source: flexpart.git/src/readpartpositions.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 4 years ago

move license from headers to a different file

  • Property mode set to 100644
File size: 4.9 KB
Line 
1subroutine readpartpositions
2
3  !*****************************************************************************
4  !                                                                            *
5  !   This routine opens the particle dump file and reads all the particle     *
6  !   positions from a previous run to initialize the current run.             *
7  !                                                                            *
8  !                                                                            *
9  !     Author: A. Stohl                                                       *
10  !                                                                            *
11  !     24 March 2000                                                          *
12  !                                                                            *
13  !*****************************************************************************
14  !                                                                            *
15  ! Variables:                                                                 *
16  !                                                                            *
17  !*****************************************************************************
18
19  use par_mod
20  use com_mod
21  use random_mod
22
23  implicit none
24
25  integer :: ibdatein,ibtimein,nspecin,itimein,numpointin,i,j,ix
26  integer :: id1,id2,it1,it2
27  real :: xlonin,ylatin,topo,hmixi,pvi,qvi,rhoi,tri,tti
28  character :: specin*7
29  real(kind=dp) :: julin,julpartin,juldate
30
31  integer :: idummy = -8
32
33  numparticlecount=0
34
35  ! Open header file of dumped particle data
36  !*****************************************
37
38  open(unitpartin,file=path(2)(1:length(2))//'header', &
39       form='unformatted',err=998)
40
41  read(unitpartin) ibdatein,ibtimein
42  read(unitpartin)
43  read(unitpartin)
44
45  read(unitpartin)
46  read(unitpartin)
47  read(unitpartin) nspecin
48  nspecin=nspecin/3
49  if ((ldirect.eq.1).and.(nspec.ne.nspecin)) goto 997
50
51  do i=1,nspecin
52    read(unitpartin)
53    read(unitpartin)
54    read(unitpartin) j,specin
55    if ((ldirect.eq.1).and.(species(i)(1:7).ne.specin)) goto 996
56  end do
57
58  read(unitpartin) numpointin
59  if (numpointin.ne.numpoint) goto 995
60999 continue
61  do i=1,numpointin
62    read(unitpartin)
63    read(unitpartin)
64    read(unitpartin)
65    read(unitpartin)
66    do j=1,nspec
67      read(unitpartin)
68      read(unitpartin)
69      read(unitpartin)
70    end do
71  end do
72  read(unitpartin)
73  read(unitpartin)
74
75  do ix=0,numxgrid-1
76    read(unitpartin)
77  end do
78
79
80  ! Open data file of dumped particle data
81  !***************************************
82
83  close(unitpartin)
84  open(unitpartin,file=path(2)(1:length(2))//'partposit_end', &
85       form='unformatted',err=998)
86 
87
88100 read(unitpartin,end=99) itimein
89  i=0
90200 i=i+1
91  read(unitpartin) npoint(i),xlonin,ylatin,ztra1(i),itramem(i), &
92       topo,pvi,qvi,rhoi,hmixi,tri,tti,(xmass1(i,j),j=1,nspec)
93 
94  if (xlonin.eq.-9999.9) goto 100
95  xtra1(i)=(xlonin-xlon0)/dx
96  ytra1(i)=(ylatin-ylat0)/dy
97  numparticlecount=max(numparticlecount,npoint(i))
98  goto 200
99
10099 numpart=i-1
101
102  close(unitpartin)
103
104  julin=juldate(ibdatein,ibtimein)+real(itimein,kind=dp)/86400._dp
105  if (abs(julin-bdate).gt.1.e-5) goto 994
106  do i=1,numpart
107    julpartin=juldate(ibdatein,ibtimein)+ &
108         real(itramem(i),kind=dp)/86400._dp
109    nclass(i)=min(int(ran1(idummy)*real(nclassunc))+1, &
110         nclassunc)
111    idt(i)=mintime
112    itra1(i)=0
113    itramem(i)=nint((julpartin-bdate)*86400.)
114    itrasplit(i)=ldirect*itsplit
115  end do
116
117  return
118
119
120994   write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### '
121  write(*,*) ' #### ENDING TIME OF PREVIOUS MODEL RUN DOES   #### '
122  write(*,*) ' #### NOT AGREE WITH STARTING TIME OF THIS RUN.#### '
123  call caldate(julin,id1,it1)
124  call caldate(bdate,id2,it2)
125  write(*,*) 'julin: ',julin,id1,it1
126  write(*,*) 'bdate: ',bdate,id2,it2
127  stop
128
129!995   write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### '
130995   write(*,*) ' #### FLEXPART MODEL WARNING IN READPARTPOSITIONS#### '
131  write(*,*) ' #### NUMBER OF RELEASE LOCATIONS DOES NOT     #### '
132  write(*,*) ' #### AGREE WITH CURRENT SETTINGS!             #### '
133!  stop
134  goto 999
135
136996   write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### '
137  write(*,*) ' #### SPECIES NAMES TO BE READ IN DO NOT       #### '
138  write(*,*) ' #### AGREE WITH CURRENT SETTINGS!             #### '
139  stop
140
141997   write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### '
142  write(*,*) ' #### THE NUMBER OF SPECIES TO BE READ IN DOES #### '
143  write(*,*) ' #### NOT AGREE WITH CURRENT SETTINGS!         #### '
144  stop
145
146998   write(*,*) ' #### FLEXPART MODEL ERROR!   THE FILE         #### '
147  write(*,*) ' #### '//path(2)(1:length(2))//'grid'//' #### '
148  write(*,*) ' #### CANNOT BE OPENED. IF A FILE WITH THIS    #### '
149  write(*,*) ' #### NAME ALREADY EXISTS, DELETE IT AND START #### '
150  write(*,*) ' #### THE PROGRAM AGAIN.                       #### '
151  stop
152
153end subroutine readpartpositions
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG