source: flexpart.git/src/init_domainfill.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: 14.3 KB
Line 
1subroutine init_domainfill
2!
3!*****************************************************************************
4!                                                                            *
5! Initializes particles equally distributed over the first release location  *
6! specified in file RELEASES. This box is assumed to be the domain for doing *
7! domain-filling trajectory calculations.                                    *
8! All particles carry the same amount of mass which alltogether comprises the*
9! mass of air within the box.                                                *
10!                                                                            *
11!     Author: A. Stohl                                                       *
12!                                                                            *
13!     15 October 2002                                                        *
14!                                                                            *
15!*****************************************************************************
16!                                                                            *
17! Variables:                                                                 *
18!                                                                            *
19! numparticlecount    consecutively counts the number of particles released  *
20! nx_we(2)       grid indices for western and eastern boundary of domain-    *
21!                filling trajectory calculations                             *
22! ny_sn(2)       grid indices for southern and northern boundary of domain-  *
23!                filling trajectory calculations                             *
24!                                                                            *
25!*****************************************************************************
26
27  use point_mod
28  use par_mod
29  use com_mod
30  use random_mod
31
32  implicit none
33
34  integer :: j,ix,jy,kz,ncolumn,numparttot
35  real :: gridarea(0:nymax-1),pp(nzmax),ylat,ylatp,ylatm,hzone
36  real :: cosfactm,cosfactp,deltacol,dz1,dz2,dz,pnew,fractus
37  real,parameter :: pih=pi/180.
38  real :: colmass(0:nxmax-1,0:nymax-1),colmasstotal,zposition
39
40  integer :: ixm,ixp,jym,jyp,indzm,indzp,in,indzh,i,jj
41  real :: pvpart,ddx,ddy,rddx,rddy,p1,p2,p3,p4,y1(2)
42
43  integer :: idummy = -11
44
45
46! Determine the release region (only full grid cells), over which particles
47! shall be initialized
48! Use 2 fields for west/east and south/north boundary
49!**************************************************************************
50
51  nx_we(1)=max(int(xpoint1(1)),0)
52  nx_we(2)=min((int(xpoint2(1))+1),nxmin1)
53  ny_sn(1)=max(int(ypoint1(1)),0)
54  ny_sn(2)=min((int(ypoint2(1))+1),nymin1)
55
56! For global simulations (both global wind data and global domain-filling),
57! set a switch, such that no boundary conditions are used
58!**************************************************************************
59  if (xglobal.and.sglobal.and.nglobal) then
60    if ((nx_we(1).eq.0).and.(nx_we(2).eq.nxmin1).and. &
61         (ny_sn(1).eq.0).and.(ny_sn(2).eq.nymin1)) then
62      gdomainfill=.true.
63    else
64      gdomainfill=.false.
65    endif
66  endif
67
68! Exit here if resuming a run from particle dump
69!***********************************************
70  if (gdomainfill.and.ipin.ne.0) return
71
72! Do not release particles twice (i.e., not at both in the leftmost and rightmost
73! grid cell) for a global domain
74!*****************************************************************************
75  if (xglobal) nx_we(2)=min(nx_we(2),nx-2)
76
77
78! Calculate area of grid cell with formula M=2*pi*R*h*dx/360,
79! see Netz, Formeln der Mathematik, 5. Auflage (1983), p.90
80!************************************************************
81
82  do jy=ny_sn(1),ny_sn(2)      ! loop about latitudes
83    ylat=ylat0+real(jy)*dy
84    ylatp=ylat+0.5*dy
85    ylatm=ylat-0.5*dy
86    if ((ylatm.lt.0).and.(ylatp.gt.0.)) then
87      hzone=1./dyconst
88    else
89      cosfactp=cos(ylatp*pih)*r_earth
90      cosfactm=cos(ylatm*pih)*r_earth
91      if (cosfactp.lt.cosfactm) then
92        hzone=sqrt(r_earth**2-cosfactp**2)- &
93             sqrt(r_earth**2-cosfactm**2)
94      else
95        hzone=sqrt(r_earth**2-cosfactm**2)- &
96             sqrt(r_earth**2-cosfactp**2)
97      endif
98    endif
99    gridarea(jy)=2.*pi*r_earth*hzone*dx/360.
100  end do
101
102! Do the same for the south pole
103
104  if (sglobal) then
105    ylat=ylat0
106    ylatp=ylat+0.5*dy
107    ylatm=ylat
108    cosfactm=0.
109    cosfactp=cos(ylatp*pih)*r_earth
110    hzone=sqrt(r_earth**2-cosfactm**2)- &
111         sqrt(r_earth**2-cosfactp**2)
112    gridarea(0)=2.*pi*r_earth*hzone*dx/360.
113  endif
114
115! Do the same for the north pole
116
117  if (nglobal) then
118    ylat=ylat0+real(nymin1)*dy
119    ylatp=ylat
120    ylatm=ylat-0.5*dy
121    cosfactp=0.
122    cosfactm=cos(ylatm*pih)*r_earth
123    hzone=sqrt(r_earth**2-cosfactp**2)- &
124         sqrt(r_earth**2-cosfactm**2)
125    gridarea(nymin1)=2.*pi*r_earth*hzone*dx/360.
126  endif
127
128
129! Calculate total mass of each grid column and of the whole atmosphere
130!*********************************************************************
131
132  colmasstotal=0.
133  do jy=ny_sn(1),ny_sn(2)          ! loop about latitudes
134    do ix=nx_we(1),nx_we(2)      ! loop about longitudes
135      pp(1)=rho(ix,jy,1,1)*r_air*tt(ix,jy,1,1)
136      pp(nz)=rho(ix,jy,nz,1)*r_air*tt(ix,jy,nz,1)
137      colmass(ix,jy)=(pp(1)-pp(nz))/ga*gridarea(jy)
138      colmasstotal=colmasstotal+colmass(ix,jy)
139    end do
140  end do
141
142  write(*,*) 'Atm. mass: ',colmasstotal
143
144
145  if (ipin.eq.0) numpart=0
146
147! Determine the particle positions
148!*********************************
149
150  numparttot=0
151  numcolumn=0
152  do jy=ny_sn(1),ny_sn(2)      ! loop about latitudes
153    ylat=ylat0+real(jy)*dy
154    do ix=nx_we(1),nx_we(2)      ! loop about longitudes
155      ncolumn=nint(0.999*real(npart(1))*colmass(ix,jy)/ &
156           colmasstotal)
157      if (ncolumn.eq.0) goto 30
158      if (ncolumn.gt.numcolumn) numcolumn=ncolumn
159
160! Calculate pressure at the altitudes of model surfaces, using the air density
161! information, which is stored as a 3-d field
162!*****************************************************************************
163
164      do kz=1,nz
165        pp(kz)=rho(ix,jy,kz,1)*r_air*tt(ix,jy,kz,1)
166      end do
167
168
169      deltacol=(pp(1)-pp(nz))/real(ncolumn)
170      pnew=pp(1)+deltacol/2.
171      jj=0
172      do j=1,ncolumn
173        jj=jj+1
174
175
176! For columns with many particles (i.e. around the equator), distribute
177! the particles equally, for columns with few particles (i.e. around the
178! poles), distribute the particles randomly
179!***********************************************************************
180
181
182        if (ncolumn.gt.20) then
183          pnew=pnew-deltacol
184        else
185          pnew=pp(1)-ran1(idummy)*(pp(1)-pp(nz))
186        endif
187
188        do kz=1,nz-1
189          if ((pp(kz).ge.pnew).and.(pp(kz+1).lt.pnew)) then
190            dz1=pp(kz)-pnew
191            dz2=pnew-pp(kz+1)
192            dz=1./(dz1+dz2)
193
194! Assign particle position
195!*************************
196! Do the following steps only if particles are not read in from previous model run
197!*****************************************************************************
198            if (ipin.eq.0) then
199              xtra1(numpart+jj)=real(ix)-0.5+ran1(idummy)
200              if (ix.eq.0) xtra1(numpart+jj)=ran1(idummy)
201              if (ix.eq.nxmin1) xtra1(numpart+jj)= &
202                   real(nxmin1)-ran1(idummy)
203              ytra1(numpart+jj)=real(jy)-0.5+ran1(idummy)
204              ztra1(numpart+jj)=(height(kz)*dz2+height(kz+1)*dz1)*dz
205              if (ztra1(numpart+jj).gt.height(nz)-0.5) &
206                   ztra1(numpart+jj)=height(nz)-0.5
207
208
209! Interpolate PV to the particle position
210!****************************************
211              ixm=int(xtra1(numpart+jj))
212              jym=int(ytra1(numpart+jj))
213              ixp=ixm+1
214              jyp=jym+1
215              ddx=xtra1(numpart+jj)-real(ixm)
216              ddy=ytra1(numpart+jj)-real(jym)
217              rddx=1.-ddx
218              rddy=1.-ddy
219              p1=rddx*rddy
220              p2=ddx*rddy
221              p3=rddx*ddy
222              p4=ddx*ddy
223              do i=2,nz
224                if (height(i).gt.ztra1(numpart+jj)) then
225                  indzm=i-1
226                  indzp=i
227                  goto 6
228                endif
229              end do
2306             continue
231              dz1=ztra1(numpart+jj)-height(indzm)
232              dz2=height(indzp)-ztra1(numpart+jj)
233              dz=1./(dz1+dz2)
234              do in=1,2
235                indzh=indzm+in-1
236                y1(in)=p1*pv(ixm,jym,indzh,1) &
237                     +p2*pv(ixp,jym,indzh,1) &
238                     +p3*pv(ixm,jyp,indzh,1) &
239                     +p4*pv(ixp,jyp,indzh,1)
240              end do
241              pvpart=(dz2*y1(1)+dz1*y1(2))*dz
242              if (ylat.lt.0.) pvpart=-1.*pvpart
243
244
245! For domain-filling option 2 (stratospheric O3), do the rest only in the stratosphere
246!*****************************************************************************
247
248              if (((ztra1(numpart+jj).gt.3000.).and. &
249                   (pvpart.gt.pvcrit)).or.(mdomainfill.eq.1)) then
250
251! Assign certain properties to the particle
252!******************************************
253                nclass(numpart+jj)=min(int(ran1(idummy)* &
254                     real(nclassunc))+1,nclassunc)
255                numparticlecount=numparticlecount+1
256                npoint(numpart+jj)=numparticlecount
257                idt(numpart+jj)=mintime
258                itra1(numpart+jj)=0
259                itramem(numpart+jj)=0
260                itrasplit(numpart+jj)=itra1(numpart+jj)+ldirect* &
261                     itsplit
262                xmass1(numpart+jj,1)=colmass(ix,jy)/real(ncolumn)
263                if (mdomainfill.eq.2) xmass1(numpart+jj,1)= &
264                     xmass1(numpart+jj,1)*pvpart*48./29.*ozonescale/10.**9
265              else
266                jj=jj-1
267              endif
268            endif
269          endif
270        end do
271      end do
272      numparttot=numparttot+ncolumn
273      if (ipin.eq.0) numpart=numpart+jj
27430    continue
275    end do
276  end do
277
278
279! Check whether numpart is really smaller than maxpart
280!*****************************************************
281
282! ESO :TODO: this warning need to be moved further up, else out-of-bounds error earlier
283  if (numpart.gt.maxpart) then
284    write(*,*) 'numpart too large: change source in init_atm_mass.f'
285    write(*,*) 'numpart: ',numpart,' maxpart: ',maxpart
286  endif
287
288
289  xmassperparticle=colmasstotal/real(numparttot)
290
291
292! Make sure that all particles are within domain
293!***********************************************
294
295  do j=1,numpart
296    if ((xtra1(j).lt.0.).or.(xtra1(j).ge.real(nxmin1)).or. &
297         (ytra1(j).lt.0.).or.(ytra1(j).ge.real(nymin1))) then
298      itra1(j)=-999999999
299    endif
300  end do
301
302
303
304
305! For boundary conditions, we need fewer particle release heights per column,
306! because otherwise it takes too long until enough mass has accumulated to
307! release a particle at the boundary (would take dx/u seconds), leading to
308! relatively large position errors of the order of one grid distance.
309! It's better to release fewer particles per column, but to do so more often.
310! Thus, use on the order of nz starting heights per column.
311! We thus repeat the above to determine fewer starting heights, that are
312! used furtheron in subroutine boundcond_domainfill.f.
313!****************************************************************************
314
315  fractus=real(numcolumn)/real(nz)
316  write(*,*) 'Total number of particles at model start: ',numpart
317  write(*,*) 'Maximum number of particles per column: ',numcolumn
318  write(*,*) 'If ',fractus,' <1, better use more particles'
319  fractus=sqrt(max(fractus,1.))/2.
320
321  do jy=ny_sn(1),ny_sn(2)      ! loop about latitudes
322    do ix=nx_we(1),nx_we(2)      ! loop about longitudes
323      ncolumn=nint(0.999/fractus*real(npart(1))*colmass(ix,jy) &
324           /colmasstotal)
325      if (ncolumn.gt.maxcolumn) stop 'maxcolumn too small'
326      if (ncolumn.eq.0) goto 80
327
328
329! Memorize how many particles per column shall be used for all boundaries
330! This is further used in subroutine boundcond_domainfill.f
331! Use 2 fields for west/east and south/north boundary
332!************************************************************************
333
334      if (ix.eq.nx_we(1)) numcolumn_we(1,jy)=ncolumn
335      if (ix.eq.nx_we(2)) numcolumn_we(2,jy)=ncolumn
336      if (jy.eq.ny_sn(1)) numcolumn_sn(1,ix)=ncolumn
337      if (jy.eq.ny_sn(2)) numcolumn_sn(2,ix)=ncolumn
338
339! Calculate pressure at the altitudes of model surfaces, using the air density
340! information, which is stored as a 3-d field
341!*****************************************************************************
342
343      do kz=1,nz
344        pp(kz)=rho(ix,jy,kz,1)*r_air*tt(ix,jy,kz,1)
345      end do
346
347! Determine the reference starting altitudes
348!*******************************************
349
350      deltacol=(pp(1)-pp(nz))/real(ncolumn)
351      pnew=pp(1)+deltacol/2.
352      do j=1,ncolumn
353        pnew=pnew-deltacol
354        do kz=1,nz-1
355          if ((pp(kz).ge.pnew).and.(pp(kz+1).lt.pnew)) then
356            dz1=pp(kz)-pnew
357            dz2=pnew-pp(kz+1)
358            dz=1./(dz1+dz2)
359            zposition=(height(kz)*dz2+height(kz+1)*dz1)*dz
360            if (zposition.gt.height(nz)-0.5) zposition=height(nz)-0.5
361
362! Memorize vertical positions where particles are introduced
363! This is further used in subroutine boundcond_domainfill.f
364!***********************************************************
365
366            if (ix.eq.nx_we(1)) zcolumn_we(1,jy,j)=zposition
367            if (ix.eq.nx_we(2)) zcolumn_we(2,jy,j)=zposition
368            if (jy.eq.ny_sn(1)) zcolumn_sn(1,ix,j)=zposition
369            if (jy.eq.ny_sn(2)) zcolumn_sn(2,ix,j)=zposition
370
371! Initialize mass that has accumulated at boundary to zero
372!*********************************************************
373
374            acc_mass_we(1,jy,j)=0.
375            acc_mass_we(2,jy,j)=0.
376            acc_mass_sn(1,jy,j)=0.
377            acc_mass_sn(2,jy,j)=0.
378          endif
379        end do
380      end do
38180    continue
382    end do
383  end do
384
385! Reduce numpart if invalid particles at end of arrays
386  do i=numpart, 1, -1
387      if (itra1(i).eq.-999999999) then
388        numpart=numpart-1
389      else
390        exit
391      end if
392    end do
393
394! If particles shall be read in to continue an existing run,
395! then the accumulated masses at the domain boundaries must be read in, too.
396! This overrides any previous calculations.
397!***************************************************************************
398
399  if ((ipin.eq.1).and.(.not.gdomainfill)) then
400    open(unitboundcond,file=path(2)(1:length(2))//'boundcond.bin', &
401         form='unformatted')
402    read(unitboundcond) numcolumn_we,numcolumn_sn, &
403         zcolumn_we,zcolumn_sn,acc_mass_we,acc_mass_sn
404    close(unitboundcond)
405  endif
406
407
408end subroutine init_domainfill
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG