source: flexpart.git/src/readpartpositions_mpi.f90 @ 766d157

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 766d157 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: 7.0 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  !   CHANGES                                                                  *
14  !     12/2014 eso: MPI version                                               *
15  !                  Root process reads positions and distributes the data     *
16  !                  :TODO: The above..                                        *
17  !*****************************************************************************
18  !                                                                            *
19  ! Variables:                                                                 *
20  !                                                                            *
21  !*****************************************************************************
22
23  use par_mod
24  use com_mod
25  use random_mod, only: ran1
26  use mpi_mod !, only: mp_seed
27
28  implicit none
29
30  integer :: ibdatein,ibtimein,nspecin,itimein,numpointin,i,j,ix,ip
31  integer :: id1,id2,it1,it2
32  integer :: addone,numparticlecount_all,numpart_all,lbnd,ubnd
33  real :: xlonin,ylatin,topo,hmixi,pvi,qvi,rhoi,tri,tti
34  character :: specin*7
35  real(kind=dp) :: julin,julpartin,juldate
36
37  integer :: idummy = -8
38
39  ! These variables are allocated at the root process for all particles in file.
40  real,dimension(maxpart) :: xtra1_all,ytra1_all,ztra1_all
41  real,dimension(maxpart,maxspec) :: xmass1_all
42  integer,dimension(maxpart) :: npoint_all,itramem_all
43
44  ! Different seed for each process
45  idummy=idummy+mp_seed
46
47  numparticlecount=0
48  numparticlecount_all=0
49  numpart_all=0
50
51  ! Open header file of dumped particle data
52  ! Each MPI process sequentially access file (just in case)
53  !*********************************************************
54
55  do ip=0, mp_partgroup_np-1
56    call mpif_mpi_barrier
57    open(unitpartin,file=path(2)(1:length(2))//'header', &
58         form='unformatted',err=998)
59
60    read(unitpartin) ibdatein,ibtimein
61    read(unitpartin)
62    read(unitpartin)
63
64    read(unitpartin)
65    read(unitpartin)
66    read(unitpartin) nspecin
67    nspecin=nspecin/3
68    if ((ldirect.eq.1).and.(nspec.ne.nspecin)) goto 997
69
70    do i=1,nspecin
71      read(unitpartin)
72      read(unitpartin)
73      read(unitpartin) j,specin
74      if ((ldirect.eq.1).and.(species(i)(1:7).ne.specin)) goto 996
75    end do
76   
77    read(unitpartin) numpointin
78    if (numpointin.ne.numpoint) then ! goto 995
79
80! eso 2016: moved this warning here to avoid out-of-block goto
81!995   write(*,*) ' #### FLEXPART MODEL WARNING IN READPARTPOSITIONS#### '
82      write(*,*) ' #### FLEXPART MODEL WARNING IN READPARTPOSITIONS#### '
83      write(*,*) ' #### NUMBER OF RELEASE LOCATIONS DOES NOT     #### '
84      write(*,*) ' #### AGREE WITH CURRENT SETTINGS!             #### '
85!  stop
86      goto 999
87    end if
88
89999 continue
90    do i=1,numpointin
91      read(unitpartin)
92      read(unitpartin)
93      read(unitpartin)
94      read(unitpartin)
95      do j=1,nspec
96        read(unitpartin)
97        read(unitpartin)
98        read(unitpartin)
99      end do
100    end do
101    read(unitpartin)
102    read(unitpartin)
103
104    do ix=0,numxgrid-1
105      read(unitpartin)
106    end do
107    close(unitpartin)
108
109
110  ! Open data file of dumped particle data
111  ! All processes read the whole file
112  !***************************************
113
114    open(unitpartin,file=path(2)(1:length(2))//'partposit_end', &
115         form='unformatted',err=998)
116   
117100 read(unitpartin,end=99) itimein
118    i=0
119200 i=i+1
120    read(unitpartin) npoint_all(i),xlonin,ylatin,ztra1_all(i),itramem_all(i), &
121         topo,pvi,qvi,rhoi,hmixi,tri,tti,(xmass1_all(i,j),j=1,nspec)
122   
123    if (xlonin.eq.-9999.9) goto 100
124    xtra1_all(i)=(xlonin-xlon0)/dx
125    ytra1_all(i)=(ylatin-ylat0)/dy
126    numparticlecount_all=max(numparticlecount_all,npoint(i))
127    goto 200
128
12999  numpart_all=i-1
130   
131    close(unitpartin)
132
133  end do
134
135
136  ! Distribute particles among processes
137  !**************************************
138  lbnd=1
139  ubnd=0
140
141  do ip=0, mp_partgroup_np-1
142   
143! Extra particle distributed in case remainder in the division
144    if (ip.lt.mod(numpart_all,mp_partgroup_np)) then
145      addone=1
146    else
147      addone=0
148    end if
149
150    if (ip==0) then
151      ubnd=ubnd + numpart_all/mp_partgroup_np + addone
152    else
153      ubnd=lbnd + numpart_all/mp_partgroup_np + addone - 1
154    end if
155   
156    if (ip==mp_pid) then
157     
158      numparticlecount=numparticlecount_all/mp_partgroup_np+addone
159      numpart=numpart_all/mp_partgroup_np+addone
160
161      xtra1(1:numpart) = xtra1_all(lbnd:ubnd)
162      ytra1(1:numpart) = ytra1_all(lbnd:ubnd)
163      ztra1(1:numpart) = ztra1_all(lbnd:ubnd)
164      xmass1(1:numpart,:) = xmass1_all(lbnd:ubnd,:)
165      npoint(1:numpart) = npoint_all(lbnd:ubnd)
166      itramem(1:numpart) = itramem_all(lbnd:ubnd)
167
168    end if
169
170    lbnd=ubnd+1
171
172  end do
173
174  julin=juldate(ibdatein,ibtimein)+real(itimein,kind=dp)/86400._dp
175  if (abs(julin-bdate).gt.1.e-5) goto 994
176  do i=1,numpart
177    julpartin=juldate(ibdatein,ibtimein)+ &
178         real(itramem(i),kind=dp)/86400._dp
179    nclass(i)=min(int(ran1(idummy)*real(nclassunc))+1, &
180         nclassunc)
181    idt(i)=mintime
182    itra1(i)=0
183    itramem(i)=nint((julpartin-bdate)*86400.)
184    itrasplit(i)=ldirect*itsplit
185  end do
186
187  return
188
189
190994   write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### '
191  write(*,*) ' #### ENDING TIME OF PREVIOUS MODEL RUN DOES   #### '
192  write(*,*) ' #### NOT AGREE WITH STARTING TIME OF THIS RUN.#### '
193  call caldate(julin,id1,it1)
194  call caldate(bdate,id2,it2)
195  write(*,*) 'julin: ',julin,id1,it1
196  write(*,*) 'bdate: ',bdate,id2,it2
197  stop
198
199996   write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### '
200  write(*,*) ' #### SPECIES NAMES TO BE READ IN DO NOT       #### '
201  write(*,*) ' #### AGREE WITH CURRENT SETTINGS!             #### '
202  stop
203
204997   write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### '
205  write(*,*) ' #### THE NUMBER OF SPECIES TO BE READ IN DOES #### '
206  write(*,*) ' #### NOT AGREE WITH CURRENT SETTINGS!         #### '
207  stop
208
209998   write(*,*) ' #### FLEXPART MODEL ERROR!   THE FILE         #### '
210  write(*,*) ' #### '//path(2)(1:length(2))//'grid'//' #### '
211  write(*,*) ' #### CANNOT BE OPENED. IF A FILE WITH THIS    #### '
212  write(*,*) ' #### NAME ALREADY EXISTS, DELETE IT AND START #### '
213  write(*,*) ' #### THE PROGRAM AGAIN.                       #### '
214  stop
215
216end subroutine readpartpositions
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG