1 | ! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt |
---|
2 | ! SPDX-License-Identifier: GPL-3.0-or-later |
---|
3 | |
---|
4 | subroutine boundcond_domainfill(itime,loutend) |
---|
5 | ! i i |
---|
6 | !***************************************************************************** |
---|
7 | ! * |
---|
8 | ! Particles are created by this subroutine continuously throughout the * |
---|
9 | ! simulation at the boundaries of the domain-filling box. * |
---|
10 | ! All particles carry the same amount of mass which alltogether comprises the* |
---|
11 | ! mass of air within the box, which remains (more or less) constant. * |
---|
12 | ! * |
---|
13 | ! Author: A. Stohl * |
---|
14 | ! * |
---|
15 | ! 16 October 2002 * |
---|
16 | ! * |
---|
17 | !***************************************************************************** |
---|
18 | ! * |
---|
19 | ! Variables: * |
---|
20 | ! * |
---|
21 | ! nx_we(2) grid indices for western and eastern boundary of domain- * |
---|
22 | ! filling trajectory calculations * |
---|
23 | ! ny_sn(2) grid indices for southern and northern boundary of domain- * |
---|
24 | ! filling trajectory calculations * |
---|
25 | ! * |
---|
26 | !***************************************************************************** |
---|
27 | ! CHANGES |
---|
28 | ! 08/2016 eso: MPI version: |
---|
29 | ! |
---|
30 | ! -Root process release particles and distributes to other processes. |
---|
31 | ! Temporary arrays are used, also for the non-root (receiving) processes. |
---|
32 | ! -The scheme can be improved by having all processes report numpart |
---|
33 | ! (keeping track of how many particles have left the domain), so that |
---|
34 | ! a proportional amount of new particles can be distributed (however |
---|
35 | ! we have a separate function called from timemanager that will |
---|
36 | ! redistribute particles among processes if there are imbalances) |
---|
37 | !***************************************************************************** |
---|
38 | |
---|
39 | use point_mod |
---|
40 | use par_mod |
---|
41 | use com_mod |
---|
42 | use random_mod, only: ran1 |
---|
43 | use mpi_mod |
---|
44 | |
---|
45 | implicit none |
---|
46 | |
---|
47 | real :: dz,dz1,dz2,dt1,dt2,dtt,ylat,xm,cosfact,accmasst |
---|
48 | integer :: itime,in,indz,indzp,i,loutend |
---|
49 | integer :: j,k,ix,jy,m,indzh,indexh,minpart,ipart,mmass |
---|
50 | integer :: numactiveparticles, numpart_total, rel_counter |
---|
51 | integer,allocatable,dimension(:) :: numrel_mpi !, numactiveparticles_mpi |
---|
52 | |
---|
53 | real :: windl(2),rhol(2) |
---|
54 | real :: windhl(2),rhohl(2) |
---|
55 | real :: windx,rhox |
---|
56 | real :: deltaz,boundarea,fluxofmass |
---|
57 | |
---|
58 | integer :: ixm,ixp,jym,jyp,indzm,mm |
---|
59 | real :: pvpart,ddx,ddy,rddx,rddy,p1,p2,p3,p4,y1(2),yh1(2) |
---|
60 | |
---|
61 | integer :: idummy = -11 |
---|
62 | integer :: mtag |
---|
63 | logical :: first_call=.true. |
---|
64 | ! Sizes of temporary arrays are maxpartfract*maxpart. Increase maxpartfract if |
---|
65 | ! needed. |
---|
66 | real,parameter :: maxpartfract=0.1 |
---|
67 | integer :: tmp_size = int(maxpartfract*maxpart) |
---|
68 | |
---|
69 | ! Use different seed for each process |
---|
70 | if (first_call) then |
---|
71 | idummy=idummy+mp_seed |
---|
72 | first_call=.false. |
---|
73 | end if |
---|
74 | |
---|
75 | |
---|
76 | ! If domain-filling is global, no boundary conditions are needed |
---|
77 | !*************************************************************** |
---|
78 | |
---|
79 | if (gdomainfill) return |
---|
80 | |
---|
81 | accmasst=0. |
---|
82 | numactiveparticles=0 |
---|
83 | ! Keep track of active particles on each process |
---|
84 | allocate(numrel_mpi(0:mp_partgroup_np-1)) |
---|
85 | ! numactiveparticles_mpi(0:mp_partgroup_np-1) |
---|
86 | |
---|
87 | ! New particles to be released on each process |
---|
88 | numrel_mpi(:)=0 |
---|
89 | |
---|
90 | ! Terminate trajectories that have left the domain, if domain-filling |
---|
91 | ! trajectory calculation domain is not global. Done for all processes |
---|
92 | !******************************************************************** |
---|
93 | |
---|
94 | do i=1,numpart |
---|
95 | if (itra1(i).eq.itime) then |
---|
96 | if ((ytra1(i).gt.real(ny_sn(2))).or. & |
---|
97 | (ytra1(i).lt.real(ny_sn(1)))) itra1(i)=-999999999 |
---|
98 | if (((.not.xglobal).or.(nx_we(2).ne.(nx-2))).and. & |
---|
99 | ((xtra1(i).lt.real(nx_we(1))).or. & |
---|
100 | (xtra1(i).gt.real(nx_we(2))))) itra1(i)=-999999999 |
---|
101 | endif |
---|
102 | if (itra1(i).ne.-999999999) numactiveparticles= & |
---|
103 | numactiveparticles+1 |
---|
104 | end do |
---|
105 | ! numactiveparticles_mpi(mp_partid) = numactiveparticles |
---|
106 | |
---|
107 | |
---|
108 | ! Collect number of active particles from all processes |
---|
109 | ! call MPI_Allgather(numactiveparticles, 1, MPI_INTEGER, & |
---|
110 | ! &numactiveparticles_mpi, 1, MPI_INTEGER, mp_comm_used, mp_ierr) |
---|
111 | |
---|
112 | |
---|
113 | ! Total number of new releases |
---|
114 | numpart_total = 0 |
---|
115 | |
---|
116 | |
---|
117 | ! This section only done by root process |
---|
118 | !*************************************** |
---|
119 | |
---|
120 | if (lroot) then |
---|
121 | |
---|
122 | ! Use separate arrays for newly released particles |
---|
123 | !************************************************* |
---|
124 | |
---|
125 | allocate(itra1_tmp(tmp_size),npoint_tmp(tmp_size),nclass_tmp(tmp_size),& |
---|
126 | & idt_tmp(tmp_size),itramem_tmp(tmp_size),itrasplit_tmp(tmp_size),& |
---|
127 | & xtra1_tmp(tmp_size),ytra1_tmp(tmp_size),ztra1_tmp(tmp_size),& |
---|
128 | & xmass1_tmp(tmp_size, maxspec)) |
---|
129 | |
---|
130 | ! Initialize all particles as non-existent |
---|
131 | itra1_tmp(:)=-999999999 |
---|
132 | |
---|
133 | ! Determine auxiliary variables for time interpolation |
---|
134 | !***************************************************** |
---|
135 | |
---|
136 | dt1=real(itime-memtime(1)) |
---|
137 | dt2=real(memtime(2)-itime) |
---|
138 | dtt=1./(dt1+dt2) |
---|
139 | |
---|
140 | ! Initialize auxiliary variable used to search for vacant storage space |
---|
141 | !********************************************************************** |
---|
142 | |
---|
143 | minpart=1 |
---|
144 | |
---|
145 | !*************************************** |
---|
146 | ! Western and eastern boundary condition |
---|
147 | !*************************************** |
---|
148 | |
---|
149 | ! Loop from south to north |
---|
150 | !************************* |
---|
151 | |
---|
152 | do jy=ny_sn(1),ny_sn(2) |
---|
153 | |
---|
154 | ! Loop over western (index 1) and eastern (index 2) boundary |
---|
155 | !*********************************************************** |
---|
156 | |
---|
157 | do k=1,2 |
---|
158 | |
---|
159 | ! Loop over all release locations in a column |
---|
160 | !******************************************** |
---|
161 | |
---|
162 | do j=1,numcolumn_we(k,jy) |
---|
163 | |
---|
164 | ! Determine, for each release location, the area of the corresponding boundary |
---|
165 | !***************************************************************************** |
---|
166 | |
---|
167 | if (j.eq.1) then |
---|
168 | deltaz=(zcolumn_we(k,jy,2)+zcolumn_we(k,jy,1))/2. |
---|
169 | else if (j.eq.numcolumn_we(k,jy)) then |
---|
170 | ! deltaz=height(nz)-(zcolumn_we(k,jy,j-1)+ |
---|
171 | ! + zcolumn_we(k,jy,j))/2. |
---|
172 | ! In order to avoid taking a very high column for very many particles, |
---|
173 | ! use the deltaz from one particle below instead |
---|
174 | deltaz=(zcolumn_we(k,jy,j)-zcolumn_we(k,jy,j-2))/2. |
---|
175 | else |
---|
176 | deltaz=(zcolumn_we(k,jy,j+1)-zcolumn_we(k,jy,j-1))/2. |
---|
177 | endif |
---|
178 | if ((jy.eq.ny_sn(1)).or.(jy.eq.ny_sn(2))) then |
---|
179 | boundarea=deltaz*111198.5/2.*dy |
---|
180 | else |
---|
181 | boundarea=deltaz*111198.5*dy |
---|
182 | endif |
---|
183 | |
---|
184 | |
---|
185 | ! Interpolate the wind velocity and density to the release location |
---|
186 | !****************************************************************** |
---|
187 | |
---|
188 | ! Determine the model level below the release position |
---|
189 | !***************************************************** |
---|
190 | |
---|
191 | do i=2,nz |
---|
192 | if (height(i).gt.zcolumn_we(k,jy,j)) then |
---|
193 | indz=i-1 |
---|
194 | indzp=i |
---|
195 | goto 6 |
---|
196 | endif |
---|
197 | end do |
---|
198 | 6 continue |
---|
199 | |
---|
200 | ! Vertical distance to the level below and above current position |
---|
201 | !**************************************************************** |
---|
202 | |
---|
203 | dz1=zcolumn_we(k,jy,j)-height(indz) |
---|
204 | dz2=height(indzp)-zcolumn_we(k,jy,j) |
---|
205 | dz=1./(dz1+dz2) |
---|
206 | |
---|
207 | ! Vertical and temporal interpolation |
---|
208 | !************************************ |
---|
209 | |
---|
210 | do m=1,2 |
---|
211 | indexh=memind(m) |
---|
212 | do in=1,2 |
---|
213 | indzh=indz+in-1 |
---|
214 | windl(in)=uu(nx_we(k),jy,indzh,indexh) |
---|
215 | rhol(in)=rho(nx_we(k),jy,indzh,indexh) |
---|
216 | end do |
---|
217 | |
---|
218 | windhl(m)=(dz2*windl(1)+dz1*windl(2))*dz |
---|
219 | rhohl(m)=(dz2*rhol(1)+dz1*rhol(2))*dz |
---|
220 | end do |
---|
221 | |
---|
222 | windx=(windhl(1)*dt2+windhl(2)*dt1)*dtt |
---|
223 | rhox=(rhohl(1)*dt2+rhohl(2)*dt1)*dtt |
---|
224 | |
---|
225 | ! Calculate mass flux |
---|
226 | !******************** |
---|
227 | |
---|
228 | fluxofmass=windx*rhox*boundarea*real(lsynctime) |
---|
229 | |
---|
230 | |
---|
231 | ! If the mass flux is directed into the domain, add it to previous mass fluxes; |
---|
232 | ! if it is out of the domain, set accumulated mass flux to zero |
---|
233 | !****************************************************************************** |
---|
234 | |
---|
235 | if (k.eq.1) then |
---|
236 | if (fluxofmass.ge.0.) then |
---|
237 | acc_mass_we(k,jy,j)=acc_mass_we(k,jy,j)+fluxofmass |
---|
238 | else |
---|
239 | acc_mass_we(k,jy,j)=0. |
---|
240 | endif |
---|
241 | else |
---|
242 | if (fluxofmass.le.0.) then |
---|
243 | acc_mass_we(k,jy,j)=acc_mass_we(k,jy,j)+abs(fluxofmass) |
---|
244 | else |
---|
245 | acc_mass_we(k,jy,j)=0. |
---|
246 | endif |
---|
247 | endif |
---|
248 | accmasst=accmasst+acc_mass_we(k,jy,j) |
---|
249 | |
---|
250 | ! If the accumulated mass exceeds half the mass that each particle shall carry, |
---|
251 | ! one (or more) particle(s) is (are) released and the accumulated mass is |
---|
252 | ! reduced by the mass of this (these) particle(s) |
---|
253 | !****************************************************************************** |
---|
254 | |
---|
255 | if (acc_mass_we(k,jy,j).ge.xmassperparticle/2.) then |
---|
256 | mmass=int((acc_mass_we(k,jy,j)+xmassperparticle/2.)/ & |
---|
257 | xmassperparticle) |
---|
258 | acc_mass_we(k,jy,j)=acc_mass_we(k,jy,j)- & |
---|
259 | real(mmass)*xmassperparticle |
---|
260 | else |
---|
261 | mmass=0 |
---|
262 | endif |
---|
263 | |
---|
264 | do m=1,mmass |
---|
265 | do ipart=minpart,maxpart |
---|
266 | |
---|
267 | ! If a vacant storage space is found, attribute everything to this array element |
---|
268 | ! TODO: for the MPI version this test can be removed, as all |
---|
269 | ! elements in _tmp arrays are initialized to zero |
---|
270 | !***************************************************************************** |
---|
271 | |
---|
272 | if (itra1_tmp(ipart).ne.itime) then |
---|
273 | |
---|
274 | ! Assign particle positions |
---|
275 | !************************** |
---|
276 | |
---|
277 | xtra1_tmp(ipart)=real(nx_we(k)) |
---|
278 | if (jy.eq.ny_sn(1)) then |
---|
279 | ytra1_tmp(ipart)=real(jy)+0.5*ran1(idummy) |
---|
280 | else if (jy.eq.ny_sn(2)) then |
---|
281 | ytra1_tmp(ipart)=real(jy)-0.5*ran1(idummy) |
---|
282 | else |
---|
283 | ytra1_tmp(ipart)=real(jy)+(ran1(idummy)-.5) |
---|
284 | endif |
---|
285 | if (j.eq.1) then |
---|
286 | ztra1_tmp(ipart)=zcolumn_we(k,jy,1)+(zcolumn_we(k,jy,2)- & |
---|
287 | zcolumn_we(k,jy,1))/4. |
---|
288 | else if (j.eq.numcolumn_we(k,jy)) then |
---|
289 | ztra1_tmp(ipart)=(2.*zcolumn_we(k,jy,j)+ & |
---|
290 | zcolumn_we(k,jy,j-1)+height(nz))/4. |
---|
291 | else |
---|
292 | ztra1_tmp(ipart)=zcolumn_we(k,jy,j-1)+ran1(idummy)* & |
---|
293 | (zcolumn_we(k,jy,j+1)-zcolumn_we(k,jy,j-1)) |
---|
294 | endif |
---|
295 | |
---|
296 | ! Interpolate PV to the particle position |
---|
297 | !**************************************** |
---|
298 | ixm=int(xtra1_tmp(ipart)) |
---|
299 | jym=int(ytra1_tmp(ipart)) |
---|
300 | ixp=ixm+1 |
---|
301 | jyp=jym+1 |
---|
302 | ddx=xtra1_tmp(ipart)-real(ixm) |
---|
303 | ddy=ytra1_tmp(ipart)-real(jym) |
---|
304 | rddx=1.-ddx |
---|
305 | rddy=1.-ddy |
---|
306 | p1=rddx*rddy |
---|
307 | p2=ddx*rddy |
---|
308 | p3=rddx*ddy |
---|
309 | p4=ddx*ddy |
---|
310 | do i=2,nz |
---|
311 | if (height(i).gt.ztra1_tmp(ipart)) then |
---|
312 | indzm=i-1 |
---|
313 | indzp=i |
---|
314 | goto 26 |
---|
315 | endif |
---|
316 | end do |
---|
317 | 26 continue |
---|
318 | dz1=ztra1_tmp(ipart)-height(indzm) |
---|
319 | dz2=height(indzp)-ztra1_tmp(ipart) |
---|
320 | dz=1./(dz1+dz2) |
---|
321 | do mm=1,2 |
---|
322 | indexh=memind(mm) |
---|
323 | do in=1,2 |
---|
324 | indzh=indzm+in-1 |
---|
325 | y1(in)=p1*pv(ixm,jym,indzh,indexh) & |
---|
326 | +p2*pv(ixp,jym,indzh,indexh) & |
---|
327 | +p3*pv(ixm,jyp,indzh,indexh) & |
---|
328 | +p4*pv(ixp,jyp,indzh,indexh) |
---|
329 | end do |
---|
330 | yh1(mm)=(dz2*y1(1)+dz1*y1(2))*dz |
---|
331 | end do |
---|
332 | pvpart=(yh1(1)*dt2+yh1(2)*dt1)*dtt |
---|
333 | ylat=ylat0+ytra1_tmp(ipart)*dy |
---|
334 | if (ylat.lt.0.) pvpart=-1.*pvpart |
---|
335 | |
---|
336 | |
---|
337 | ! For domain-filling option 2 (stratospheric O3), do the rest only in the stratosphere |
---|
338 | !***************************************************************************** |
---|
339 | |
---|
340 | if (((ztra1_tmp(ipart).gt.3000.).and. & |
---|
341 | (pvpart.gt.pvcrit)).or.(mdomainfill.eq.1)) then |
---|
342 | nclass_tmp(ipart)=min(int(ran1(idummy)* & |
---|
343 | real(nclassunc))+1,nclassunc) |
---|
344 | numactiveparticles=numactiveparticles+1 |
---|
345 | numparticlecount=numparticlecount+1 |
---|
346 | npoint_tmp(ipart)=numparticlecount |
---|
347 | idt_tmp(ipart)=mintime |
---|
348 | itra1_tmp(ipart)=itime |
---|
349 | itramem_tmp(ipart)=itra1_tmp(ipart) |
---|
350 | itrasplit_tmp(ipart)=itra1_tmp(ipart)+ldirect*itsplit |
---|
351 | xmass1_tmp(ipart,1)=xmassperparticle |
---|
352 | if (mdomainfill.eq.2) xmass1_tmp(ipart,1)= & |
---|
353 | xmass1_tmp(ipart,1)*pvpart*48./29.*ozonescale/10.**9 |
---|
354 | else |
---|
355 | goto 71 |
---|
356 | endif |
---|
357 | |
---|
358 | |
---|
359 | ! Increase numpart, if necessary |
---|
360 | !******************************* |
---|
361 | |
---|
362 | numpart_total=max(numpart_total,ipart) |
---|
363 | goto 73 ! Storage space has been found, stop searching |
---|
364 | endif |
---|
365 | end do |
---|
366 | if (ipart.gt.tmp_size) & |
---|
367 | stop 'boundcond_domainfill_mpi.f90: too many particles required' |
---|
368 | 73 minpart=ipart+1 |
---|
369 | 71 continue |
---|
370 | end do |
---|
371 | |
---|
372 | |
---|
373 | end do |
---|
374 | end do |
---|
375 | end do |
---|
376 | |
---|
377 | |
---|
378 | !***************************************** |
---|
379 | ! Southern and northern boundary condition |
---|
380 | !***************************************** |
---|
381 | |
---|
382 | ! Loop from west to east |
---|
383 | !*********************** |
---|
384 | |
---|
385 | do ix=nx_we(1),nx_we(2) |
---|
386 | |
---|
387 | ! Loop over southern (index 1) and northern (index 2) boundary |
---|
388 | !************************************************************* |
---|
389 | |
---|
390 | do k=1,2 |
---|
391 | ylat=ylat0+real(ny_sn(k))*dy |
---|
392 | cosfact=cos(ylat*pi180) |
---|
393 | |
---|
394 | ! Loop over all release locations in a column |
---|
395 | !******************************************** |
---|
396 | |
---|
397 | do j=1,numcolumn_sn(k,ix) |
---|
398 | |
---|
399 | ! Determine, for each release location, the area of the corresponding boundary |
---|
400 | !***************************************************************************** |
---|
401 | |
---|
402 | if (j.eq.1) then |
---|
403 | deltaz=(zcolumn_sn(k,ix,2)+zcolumn_sn(k,ix,1))/2. |
---|
404 | else if (j.eq.numcolumn_sn(k,ix)) then |
---|
405 | ! deltaz=height(nz)-(zcolumn_sn(k,ix,j-1)+ |
---|
406 | ! + zcolumn_sn(k,ix,j))/2. |
---|
407 | ! In order to avoid taking a very high column for very many particles, |
---|
408 | ! use the deltaz from one particle below instead |
---|
409 | deltaz=(zcolumn_sn(k,ix,j)-zcolumn_sn(k,ix,j-2))/2. |
---|
410 | else |
---|
411 | deltaz=(zcolumn_sn(k,ix,j+1)-zcolumn_sn(k,ix,j-1))/2. |
---|
412 | endif |
---|
413 | if ((ix.eq.nx_we(1)).or.(ix.eq.nx_we(2))) then |
---|
414 | boundarea=deltaz*111198.5/2.*cosfact*dx |
---|
415 | else |
---|
416 | boundarea=deltaz*111198.5*cosfact*dx |
---|
417 | endif |
---|
418 | |
---|
419 | |
---|
420 | ! Interpolate the wind velocity and density to the release location |
---|
421 | !****************************************************************** |
---|
422 | |
---|
423 | ! Determine the model level below the release position |
---|
424 | !***************************************************** |
---|
425 | |
---|
426 | do i=2,nz |
---|
427 | if (height(i).gt.zcolumn_sn(k,ix,j)) then |
---|
428 | indz=i-1 |
---|
429 | indzp=i |
---|
430 | goto 16 |
---|
431 | endif |
---|
432 | end do |
---|
433 | 16 continue |
---|
434 | |
---|
435 | ! Vertical distance to the level below and above current position |
---|
436 | !**************************************************************** |
---|
437 | |
---|
438 | dz1=zcolumn_sn(k,ix,j)-height(indz) |
---|
439 | dz2=height(indzp)-zcolumn_sn(k,ix,j) |
---|
440 | dz=1./(dz1+dz2) |
---|
441 | |
---|
442 | ! Vertical and temporal interpolation |
---|
443 | !************************************ |
---|
444 | |
---|
445 | do m=1,2 |
---|
446 | indexh=memind(m) |
---|
447 | do in=1,2 |
---|
448 | indzh=indz+in-1 |
---|
449 | windl(in)=vv(ix,ny_sn(k),indzh,indexh) |
---|
450 | rhol(in)=rho(ix,ny_sn(k),indzh,indexh) |
---|
451 | end do |
---|
452 | |
---|
453 | windhl(m)=(dz2*windl(1)+dz1*windl(2))*dz |
---|
454 | rhohl(m)=(dz2*rhol(1)+dz1*rhol(2))*dz |
---|
455 | end do |
---|
456 | |
---|
457 | windx=(windhl(1)*dt2+windhl(2)*dt1)*dtt |
---|
458 | rhox=(rhohl(1)*dt2+rhohl(2)*dt1)*dtt |
---|
459 | |
---|
460 | ! Calculate mass flux |
---|
461 | !******************** |
---|
462 | |
---|
463 | fluxofmass=windx*rhox*boundarea*real(lsynctime) |
---|
464 | |
---|
465 | ! If the mass flux is directed into the domain, add it to previous mass fluxes; |
---|
466 | ! if it is out of the domain, set accumulated mass flux to zero |
---|
467 | !****************************************************************************** |
---|
468 | |
---|
469 | if (k.eq.1) then |
---|
470 | if (fluxofmass.ge.0.) then |
---|
471 | acc_mass_sn(k,ix,j)=acc_mass_sn(k,ix,j)+fluxofmass |
---|
472 | else |
---|
473 | acc_mass_sn(k,ix,j)=0. |
---|
474 | endif |
---|
475 | else |
---|
476 | if (fluxofmass.le.0.) then |
---|
477 | acc_mass_sn(k,ix,j)=acc_mass_sn(k,ix,j)+abs(fluxofmass) |
---|
478 | else |
---|
479 | acc_mass_sn(k,ix,j)=0. |
---|
480 | endif |
---|
481 | endif |
---|
482 | accmasst=accmasst+acc_mass_sn(k,ix,j) |
---|
483 | |
---|
484 | ! If the accumulated mass exceeds half the mass that each particle shall carry, |
---|
485 | ! one (or more) particle(s) is (are) released and the accumulated mass is |
---|
486 | ! reduced by the mass of this (these) particle(s) |
---|
487 | !****************************************************************************** |
---|
488 | |
---|
489 | if (acc_mass_sn(k,ix,j).ge.xmassperparticle/2.) then |
---|
490 | mmass=int((acc_mass_sn(k,ix,j)+xmassperparticle/2.)/ & |
---|
491 | xmassperparticle) |
---|
492 | acc_mass_sn(k,ix,j)=acc_mass_sn(k,ix,j)- & |
---|
493 | real(mmass)*xmassperparticle |
---|
494 | else |
---|
495 | mmass=0 |
---|
496 | endif |
---|
497 | |
---|
498 | do m=1,mmass |
---|
499 | do ipart=minpart,maxpart |
---|
500 | |
---|
501 | ! If a vacant storage space is found, attribute everything to this array element |
---|
502 | !***************************************************************************** |
---|
503 | |
---|
504 | if (itra1_tmp(ipart).ne.itime) then |
---|
505 | |
---|
506 | ! Assign particle positions |
---|
507 | !************************** |
---|
508 | |
---|
509 | ytra1_tmp(ipart)=real(ny_sn(k)) |
---|
510 | if (ix.eq.nx_we(1)) then |
---|
511 | xtra1_tmp(ipart)=real(ix)+0.5*ran1(idummy) |
---|
512 | else if (ix.eq.nx_we(2)) then |
---|
513 | xtra1_tmp(ipart)=real(ix)-0.5*ran1(idummy) |
---|
514 | else |
---|
515 | xtra1_tmp(ipart)=real(ix)+(ran1(idummy)-.5) |
---|
516 | endif |
---|
517 | if (j.eq.1) then |
---|
518 | ztra1_tmp(ipart)=zcolumn_sn(k,ix,1)+(zcolumn_sn(k,ix,2)- & |
---|
519 | zcolumn_sn(k,ix,1))/4. |
---|
520 | else if (j.eq.numcolumn_sn(k,ix)) then |
---|
521 | ztra1_tmp(ipart)=(2.*zcolumn_sn(k,ix,j)+ & |
---|
522 | zcolumn_sn(k,ix,j-1)+height(nz))/4. |
---|
523 | else |
---|
524 | ztra1_tmp(ipart)=zcolumn_sn(k,ix,j-1)+ran1(idummy)* & |
---|
525 | (zcolumn_sn(k,ix,j+1)-zcolumn_sn(k,ix,j-1)) |
---|
526 | endif |
---|
527 | |
---|
528 | |
---|
529 | ! Interpolate PV to the particle position |
---|
530 | !**************************************** |
---|
531 | ixm=int(xtra1_tmp(ipart)) |
---|
532 | jym=int(ytra1_tmp(ipart)) |
---|
533 | ixp=ixm+1 |
---|
534 | jyp=jym+1 |
---|
535 | ddx=xtra1_tmp(ipart)-real(ixm) |
---|
536 | ddy=ytra1_tmp(ipart)-real(jym) |
---|
537 | rddx=1.-ddx |
---|
538 | rddy=1.-ddy |
---|
539 | p1=rddx*rddy |
---|
540 | p2=ddx*rddy |
---|
541 | p3=rddx*ddy |
---|
542 | p4=ddx*ddy |
---|
543 | do i=2,nz |
---|
544 | if (height(i).gt.ztra1_tmp(ipart)) then |
---|
545 | indzm=i-1 |
---|
546 | indzp=i |
---|
547 | goto 126 |
---|
548 | endif |
---|
549 | end do |
---|
550 | 126 continue |
---|
551 | dz1=ztra1_tmp(ipart)-height(indzm) |
---|
552 | dz2=height(indzp)-ztra1_tmp(ipart) |
---|
553 | dz=1./(dz1+dz2) |
---|
554 | do mm=1,2 |
---|
555 | indexh=memind(mm) |
---|
556 | do in=1,2 |
---|
557 | indzh=indzm+in-1 |
---|
558 | y1(in)=p1*pv(ixm,jym,indzh,indexh) & |
---|
559 | +p2*pv(ixp,jym,indzh,indexh) & |
---|
560 | +p3*pv(ixm,jyp,indzh,indexh) & |
---|
561 | +p4*pv(ixp,jyp,indzh,indexh) |
---|
562 | end do |
---|
563 | yh1(mm)=(dz2*y1(1)+dz1*y1(2))*dz |
---|
564 | end do |
---|
565 | pvpart=(yh1(1)*dt2+yh1(2)*dt1)*dtt |
---|
566 | if (ylat.lt.0.) pvpart=-1.*pvpart |
---|
567 | |
---|
568 | |
---|
569 | ! For domain-filling option 2 (stratospheric O3), do the rest only in the stratosphere |
---|
570 | !***************************************************************************** |
---|
571 | |
---|
572 | if (((ztra1_tmp(ipart).gt.3000.).and. & |
---|
573 | (pvpart.gt.pvcrit)).or.(mdomainfill.eq.1)) then |
---|
574 | nclass_tmp(ipart)=min(int(ran1(idummy)* & |
---|
575 | real(nclassunc))+1,nclassunc) |
---|
576 | numactiveparticles=numactiveparticles+1 |
---|
577 | numparticlecount=numparticlecount+1 |
---|
578 | npoint_tmp(ipart)=numparticlecount |
---|
579 | idt_tmp(ipart)=mintime |
---|
580 | itra1_tmp(ipart)=itime |
---|
581 | itramem_tmp(ipart)=itra1_tmp(ipart) |
---|
582 | itrasplit_tmp(ipart)=itra1_tmp(ipart)+ldirect*itsplit |
---|
583 | xmass1_tmp(ipart,1)=xmassperparticle |
---|
584 | if (mdomainfill.eq.2) xmass1_tmp(ipart,1)= & |
---|
585 | xmass1_tmp(ipart,1)*pvpart*48./29.*ozonescale/10.**9 |
---|
586 | else |
---|
587 | goto 171 |
---|
588 | endif |
---|
589 | |
---|
590 | |
---|
591 | ! Increase numpart, if necessary |
---|
592 | !******************************* |
---|
593 | numpart_total=max(numpart_total,ipart) |
---|
594 | goto 173 ! Storage space has been found, stop searching |
---|
595 | endif |
---|
596 | end do |
---|
597 | if (ipart.gt.tmp_size) & |
---|
598 | stop 'boundcond_domainfill.f: too many particles required' |
---|
599 | 173 minpart=ipart+1 |
---|
600 | 171 continue |
---|
601 | end do |
---|
602 | |
---|
603 | |
---|
604 | end do |
---|
605 | end do |
---|
606 | end do |
---|
607 | |
---|
608 | |
---|
609 | ! xm=0. |
---|
610 | ! do i=1,numpart_total |
---|
611 | ! if (itra1_tmp(i).eq.itime) xm=xm+xmass1(i,1) |
---|
612 | ! end do |
---|
613 | |
---|
614 | !write(*,*) itime,numactiveparticles,numparticlecount,numpart, |
---|
615 | ! +xm,accmasst,xm+accmasst |
---|
616 | |
---|
617 | end if ! if lroot |
---|
618 | |
---|
619 | ! Distribute the number of particles to be released |
---|
620 | ! ************************************************* |
---|
621 | call MPI_Bcast(numpart_total, 1, MPI_INTEGER, id_root, mp_comm_used, mp_ierr) |
---|
622 | |
---|
623 | do i=0, mp_partgroup_np-1 |
---|
624 | numrel_mpi(i) = numpart_total/mp_partgroup_np |
---|
625 | if (i.lt.mod(numpart_total,mp_partgroup_np)) numrel_mpi(i) = numrel_mpi(i) + 1 |
---|
626 | end do |
---|
627 | |
---|
628 | ! Allocate temporary arrays for receiving processes |
---|
629 | if (.not.lroot) then |
---|
630 | allocate(itra1_tmp(numrel_mpi(mp_partid)),& |
---|
631 | & npoint_tmp(numrel_mpi(mp_partid)),& |
---|
632 | & nclass_tmp(numrel_mpi(mp_partid)),& |
---|
633 | & idt_tmp(numrel_mpi(mp_partid)),& |
---|
634 | & itramem_tmp(numrel_mpi(mp_partid)),& |
---|
635 | & itrasplit_tmp(numrel_mpi(mp_partid)),& |
---|
636 | & xtra1_tmp(numrel_mpi(mp_partid)),& |
---|
637 | & ytra1_tmp(numrel_mpi(mp_partid)),& |
---|
638 | & ztra1_tmp(numrel_mpi(mp_partid)),& |
---|
639 | & xmass1_tmp(numrel_mpi(mp_partid),maxspec)) |
---|
640 | |
---|
641 | ! Initialize all particles as non-existent |
---|
642 | itra1_tmp(:)=-999999999 |
---|
643 | end if |
---|
644 | |
---|
645 | ! Distribute particles |
---|
646 | ! Keep track of released particles so far |
---|
647 | rel_counter = 0 |
---|
648 | mtag = 1000 |
---|
649 | |
---|
650 | do i=0, mp_partgroup_np-1 |
---|
651 | |
---|
652 | ! For root process, nothing to do except update release count |
---|
653 | if (i.eq.0) then |
---|
654 | rel_counter = rel_counter + numrel_mpi(i) |
---|
655 | cycle |
---|
656 | end if |
---|
657 | |
---|
658 | ! Send particles from root to non-root processes |
---|
659 | if (lroot.and.numrel_mpi(i).gt.0) then |
---|
660 | |
---|
661 | call MPI_SEND(nclass_tmp(rel_counter+1:rel_counter+numrel_mpi(i)),& |
---|
662 | &numrel_mpi(i),MPI_INTEGER,i,mtag+1*i,mp_comm_used,mp_ierr) |
---|
663 | |
---|
664 | call MPI_SEND(npoint_tmp(rel_counter+1:rel_counter+numrel_mpi(i)),& |
---|
665 | &numrel_mpi(i),MPI_INTEGER,i,mtag+2*i,mp_comm_used,mp_ierr) |
---|
666 | |
---|
667 | call MPI_SEND(itra1_tmp(rel_counter+1:rel_counter+numrel_mpi(i)),& |
---|
668 | &numrel_mpi(i),MPI_INTEGER,i,mtag+3*i,mp_comm_used,mp_ierr) |
---|
669 | |
---|
670 | call MPI_SEND(idt_tmp(rel_counter+1:rel_counter+numrel_mpi(i)),& |
---|
671 | &numrel_mpi(i),MPI_INTEGER,i,mtag+4*i,mp_comm_used,mp_ierr) |
---|
672 | |
---|
673 | call MPI_SEND(itramem_tmp(rel_counter+1:rel_counter+numrel_mpi(i)),& |
---|
674 | &numrel_mpi(i),MPI_INTEGER,i,mtag+5*i,mp_comm_used,mp_ierr) |
---|
675 | |
---|
676 | call MPI_SEND(itrasplit_tmp(rel_counter+1:rel_counter+numrel_mpi(i)),& |
---|
677 | &numrel_mpi(i),MPI_INTEGER,i,mtag+6*i,mp_comm_used,mp_ierr) |
---|
678 | |
---|
679 | call MPI_SEND(xtra1_tmp(rel_counter+1:rel_counter+numrel_mpi(i)),& |
---|
680 | &numrel_mpi(i),mp_dp,i,mtag+7*i,mp_comm_used,mp_ierr) |
---|
681 | |
---|
682 | call MPI_SEND(ytra1_tmp(rel_counter+1:rel_counter+numrel_mpi(i)),& |
---|
683 | &numrel_mpi(i),mp_dp,i,mtag+8*i,mp_comm_used,mp_ierr) |
---|
684 | |
---|
685 | call MPI_SEND(ztra1_tmp(rel_counter+1:rel_counter+numrel_mpi(i)),& |
---|
686 | &numrel_mpi(i),mp_sp,i,mtag+9*i,mp_comm_used,mp_ierr) |
---|
687 | |
---|
688 | do j=1,nspec |
---|
689 | call MPI_SEND(xmass1_tmp(rel_counter+1:rel_counter+numrel_mpi(i),j),& |
---|
690 | &numrel_mpi(i),mp_sp,i,mtag+(9+j)*i,mp_comm_used,mp_ierr) |
---|
691 | end do |
---|
692 | |
---|
693 | ! Non-root processes issue receive requests |
---|
694 | else if (i.eq.mp_partid.and.numrel_mpi(i).gt.0) then |
---|
695 | call MPI_RECV(nclass_tmp(1:numrel_mpi(i)),numrel_mpi(i),& |
---|
696 | &MPI_INTEGER,id_root,mtag+1*i,mp_comm_used,mp_status,mp_ierr) |
---|
697 | |
---|
698 | call MPI_RECV(npoint_tmp(1:numrel_mpi(i)),numrel_mpi(i),& |
---|
699 | &MPI_INTEGER,id_root,mtag+2*i,mp_comm_used,mp_status,mp_ierr) |
---|
700 | |
---|
701 | call MPI_RECV(itra1_tmp(1:numrel_mpi(i)),numrel_mpi(i),& |
---|
702 | &MPI_INTEGER,id_root,mtag+3*i,mp_comm_used,mp_status,mp_ierr) |
---|
703 | |
---|
704 | call MPI_RECV(idt_tmp(1:numrel_mpi(i)),numrel_mpi(i),& |
---|
705 | &MPI_INTEGER,id_root,mtag+4*i,mp_comm_used,mp_status,mp_ierr) |
---|
706 | |
---|
707 | call MPI_RECV(itramem_tmp(1:numrel_mpi(i)),numrel_mpi(i),& |
---|
708 | &MPI_INTEGER,id_root,mtag+5*i,mp_comm_used,mp_status,mp_ierr) |
---|
709 | |
---|
710 | call MPI_RECV(itrasplit_tmp(1:numrel_mpi(i)),numrel_mpi(i),& |
---|
711 | &MPI_INTEGER,id_root,mtag+6*i,mp_comm_used,mp_status,mp_ierr) |
---|
712 | |
---|
713 | call MPI_RECV(xtra1_tmp(1:numrel_mpi(i)),numrel_mpi(i),& |
---|
714 | &mp_dp,id_root,mtag+7*i,mp_comm_used,mp_status,mp_ierr) |
---|
715 | |
---|
716 | call MPI_RECV(ytra1_tmp(1:numrel_mpi(i)),numrel_mpi(i),& |
---|
717 | &mp_dp,id_root,mtag+8*i,mp_comm_used,mp_status,mp_ierr) |
---|
718 | |
---|
719 | call MPI_RECV(ztra1_tmp(1:numrel_mpi(i)),numrel_mpi(i),& |
---|
720 | &mp_sp,id_root,mtag+9*i,mp_comm_used,mp_status,mp_ierr) |
---|
721 | |
---|
722 | do j=1,nspec |
---|
723 | call MPI_RECV(xmass1_tmp(1:numrel_mpi(i),j),numrel_mpi(i),& |
---|
724 | &mp_sp,id_root,mtag+(9+j)*i,mp_comm_used,mp_status,mp_ierr) |
---|
725 | |
---|
726 | end do |
---|
727 | end if |
---|
728 | rel_counter = rel_counter + numrel_mpi(i) |
---|
729 | end do |
---|
730 | |
---|
731 | ! Find free storage space for the new particles. |
---|
732 | ! This section is independent of the redistribution scheme used |
---|
733 | ! ******************************************************************** |
---|
734 | |
---|
735 | ! Keep track of released particles so far |
---|
736 | minpart=1 |
---|
737 | |
---|
738 | ! The algorithm should be correct also for root process |
---|
739 | do i=1, numrel_mpi(mp_partid) |
---|
740 | do ipart=minpart, maxpart |
---|
741 | if (itra1(ipart).ne.itime) then |
---|
742 | itra1(ipart) = itra1_tmp(i) |
---|
743 | npoint(ipart) = npoint_tmp(i) |
---|
744 | nclass(ipart) = nclass_tmp(i) |
---|
745 | idt(ipart) = idt_tmp(i) |
---|
746 | itramem(ipart) = itramem_tmp(i) |
---|
747 | itrasplit(ipart) = itrasplit_tmp(i) |
---|
748 | xtra1(ipart) = xtra1_tmp(i) |
---|
749 | ytra1(ipart) = ytra1_tmp(i) |
---|
750 | ztra1(ipart) = ztra1_tmp(i) |
---|
751 | xmass1(ipart,:) = xmass1_tmp(i,:) |
---|
752 | ! Increase numpart, if necessary |
---|
753 | numpart=max(numpart,ipart) |
---|
754 | goto 200 ! Storage space has been found, stop searching |
---|
755 | end if |
---|
756 | end do |
---|
757 | 200 minpart=ipart+1 |
---|
758 | end do |
---|
759 | |
---|
760 | ! If particles shall be dumped, then accumulated masses at the domain boundaries |
---|
761 | ! must be dumped, too, to be used for later runs |
---|
762 | !***************************************************************************** |
---|
763 | |
---|
764 | if ((ipout.gt.0).and.(itime.eq.loutend)) then |
---|
765 | if (lroot) then |
---|
766 | call mpif_mtime('iotime',0) |
---|
767 | open(unitboundcond,file=path(2)(1:length(2))//'boundcond.bin', & |
---|
768 | form='unformatted') |
---|
769 | write(unitboundcond) numcolumn_we,numcolumn_sn, & |
---|
770 | zcolumn_we,zcolumn_sn,acc_mass_we,acc_mass_sn |
---|
771 | close(unitboundcond) |
---|
772 | call mpif_mtime('iotime',1) |
---|
773 | end if |
---|
774 | endif |
---|
775 | |
---|
776 | ! Deallocate temporary arrays |
---|
777 | deallocate(itra1_tmp,npoint_tmp,nclass_tmp,idt_tmp,itramem_tmp,itrasplit_tmp,& |
---|
778 | & xtra1_tmp,ytra1_tmp,ztra1_tmp,xmass1_tmp,numrel_mpi) |
---|
779 | ! numactiveparticles_mpi |
---|
780 | |
---|
781 | |
---|
782 | end subroutine boundcond_domainfill |
---|