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