-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathboundcond_domainfill.f90
637 lines (529 loc) · 23.7 KB
/
boundcond_domainfill.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
!***********************************************************************
!* Copyright 2012,2013 *
!* Jerome Brioude, Delia Arnold, Andreas Stohl, Wayne Angevine, *
!* John Burkhart, Massimo Cassiani, Adam Dingwell, Richard C Easter, Sabine Eckhardt,*
!* Stephanie Evan, Jerome D Fast, Don Morton, Ignacio Pisso, *
!* Petra Seibert, Gerard Wotawa, Caroline Forster, Harald Sodemann, *
!* *
!* This file is part of FLEXPART WRF *
!* *
!* FLEXPART is free software: you can redistribute it and/or modify *
!* it under the terms of the GNU General Public License as published by*
!* the Free Software Foundation, either version 3 of the License, or *
!* (at your option) any later version. *
!* *
!* FLEXPART is distributed in the hope that it will be useful, *
!* but WITHOUT ANY WARRANTY; without even the implied warranty of *
!* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
!* GNU General Public License for more details. *
!* *
!* You should have received a copy of the GNU General Public License *
!* along with FLEXPART. If not, see <http://www.gnu.org/licenses/>. *
!***********************************************************************
subroutine boundcond_domainfill(itime,loutend)
! i i
!*******************************************************************************
! *
! Note: This is the FLEXPART_WRF version of subr. boundcond_domainfill. *
! The computational grid is the WRF x-y grid rather than lat-lon. *
! *
! Particles are created by this subroutine continuously throughout the *
! simulation at the boundaries of the domain-filling box. *
! All particles carry the same amount of mass which alltogether comprises the *
! mass of air within the box, which remains (more or less) constant. *
! *
! Author: A. Stohl *
! *
! 16 October 2002 *
! *
! 26 Oct 2005, R. Easter - changes to calc. of boundarea *
! associated with WRF horizontal grid. *
! Also need to get true ylat for pv calcs. *
! 11 Nov 2005, R. Easter - fixed error involving xy to latlong *
! 2012, J. Brioude: coded in fortran 90 *
!*******************************************************************************
! *
! Variables: *
! *
! nx_we(2) grid indices for western and eastern boundary of domain- *
! filling trajectory calculations *
! ny_sn(2) grid indices for southern and northern boundary of domain- *
! filling trajectory calculations *
! *
!*******************************************************************************
use point_mod
use par_mod
use com_mod
implicit none
real :: dz,dz1,dz2,ran1,dt1,dt2,dtt,xm,cosfact,accmasst
integer :: itime,in,indz,indzp,i,loutend
integer :: j,k,ix,jy,m,indzh,indexh,minpart,ipart,mmass
integer :: numactiveparticles
real :: windl(2),rhol(2),dumx,dumy,xlon,ylat
real :: windhl(2),rhohl(2)
real :: windx,rhox
real :: deltaz,boundarea,fluxofmass
integer :: ixm,ixp,jym,jyp,indzm,mm
real :: pvpart,ddx,ddy,rddx,rddy,p1,p2,p3,p4,y1(2),yh1(2)
integer :: idummy = -11
! If domain-filling is global, no boundary conditions are needed
!***************************************************************
if (gdomainfill) return
accmasst=0.
numactiveparticles=0
! Terminate trajectories that have left the domain, if domain-filling
! trajectory calculation domain is not global
!********************************************************************
do i=1,numpart
if (itra1(i).eq.itime) then
if ((ytra1(i).gt.real(ny_sn(2))).or. &
(ytra1(i).lt.real(ny_sn(1)))) itra1(i)=-999999999
if (((.not.xglobal).or.(nx_we(2).ne.(nx-2))).and. &
((xtra1(i).lt.real(nx_we(1))).or. &
(xtra1(i).gt.real(nx_we(2))))) itra1(i)=-999999999
endif
if (itra1(i).ne.-999999999) numactiveparticles= &
+ numactiveparticles+1
enddo
! Determine auxiliary variables for time interpolation
!*****************************************************
dt1=real(itime-memtime(1))
dt2=real(memtime(2)-itime)
dtt=1./(dt1+dt2)
! Initialize auxiliary variable used to search for vacant storage space
!**********************************************************************
minpart=1
!***************************************
! Western and eastern boundary condition
!***************************************
! Loop from south to north
!*************************
do jy=ny_sn(1),ny_sn(2)
! Loop over western (index 1) and eastern (index 2) boundary
!***********************************************************
do k=1,2
! for FLEXPART_WRF, x & y coords are in meters.
! In the "do 70" loop, ylat is only needed for for pv calcs,
! "if (ylat.lt.0.) pvpart=-1.*pvpart"
! Note: in the FLEXPART_ECMWF code, ylat was not defined
! in the "do 70" loop (a bug).
dumx=real(nx_we(k))
dumy=real(jy)
! Are these dumx,dumy correct ???
call xyindex_to_ll_wrf( 0, dumx, dumy, xlon, ylat )
! Loop over all release locations in a column
!********************************************
do j=1,numcolumn_we(k,jy)
! Determine, for each release location, the area of the corresponding boundary
!*****************************************************************************
if (j.eq.1) then
deltaz=(zcolumn_we(k,jy,2)+zcolumn_we(k,jy,1))/2.
else if (j.eq.numcolumn_we(k,jy)) then
! deltaz=height(nz)-(zcolumn_we(k,jy,j-1)+
! + zcolumn_we(k,jy,j))/2.
! In order to avoid taking a very high column for very many particles,
! use the deltaz from one particle below instead
deltaz=(zcolumn_we(k,jy,j)-zcolumn_we(k,jy,j-2))/2.
else
deltaz=(zcolumn_we(k,jy,j+1)-zcolumn_we(k,jy,j-1))/2.
endif
! for FLEXPART_ECMWF, dy is in degrees-lat, and 111198.5 converts
! from degrees-latitude to m
! for FLEXPART_WRF, dy is in meters
! if ((jy.eq.ny_sn(1)).or.(jy.eq.ny_sn(2))) then
! boundarea=deltaz*111198.5/2.*dy
! else
! boundarea=deltaz*111198.5*dy
! endif
if ((jy.eq.ny_sn(1)).or.(jy.eq.ny_sn(2))) then
boundarea=deltaz/2.*dy
else
boundarea=deltaz*dy
endif
! Interpolate the wind velocity and density to the release location
!******************************************************************
! Determine the model level below the release position
!*****************************************************
do i=2,nz
if (height(i).gt.zcolumn_we(k,jy,j)) then
indz=i-1
indzp=i
goto 6
endif
enddo
6 continue
! Vertical distance to the level below and above current position
!****************************************************************
dz1=zcolumn_we(k,jy,j)-height(indz)
dz2=height(indzp)-zcolumn_we(k,jy,j)
dz=1./(dz1+dz2)
! Vertical and temporal interpolation
!************************************
do m=1,2
indexh=memind(m)
do in=1,2
indzh=indz+in-1
windl(in)=uu(nx_we(k),jy,indzh,indexh)
rhol(in)=rho(nx_we(k),jy,indzh,indexh)
enddo
windhl(m)=(dz2*windl(1)+dz1*windl(2))*dz
rhohl(m)=(dz2*rhol(1)+dz1*rhol(2))*dz
enddo
windx=(windhl(1)*dt2+windhl(2)*dt1)*dtt
rhox=(rhohl(1)*dt2+rhohl(2)*dt1)*dtt
! Calculate mass flux
!********************
fluxofmass=windx*rhox*boundarea*real(lsynctime)
! If the mass flux is directed into the domain, add it to previous mass fluxes;
! if it is out of the domain, set accumulated mass flux to zero
!******************************************************************************
if (k.eq.1) then
if (fluxofmass.ge.0.) then
acc_mass_we(k,jy,j)=acc_mass_we(k,jy,j)+fluxofmass
else
acc_mass_we(k,jy,j)=0.
endif
else
if (fluxofmass.le.0.) then
acc_mass_we(k,jy,j)=acc_mass_we(k,jy,j)+abs(fluxofmass)
else
acc_mass_we(k,jy,j)=0.
endif
endif
accmasst=accmasst+acc_mass_we(k,jy,j)
! If the accumulated mass exceeds half the mass that each particle shall carry,
! one (or more) particle(s) is (are) released and the accumulated mass is
! reduced by the mass of this (these) particle(s)
!******************************************************************************
if (acc_mass_we(k,jy,j).ge.xmassperparticle/2.) then
mmass=int((acc_mass_we(k,jy,j)+xmassperparticle/2.)/ &
xmassperparticle)
acc_mass_we(k,jy,j)=acc_mass_we(k,jy,j)- &
real(mmass)*xmassperparticle
else
mmass=0
endif
do m=1,mmass
do ipart=minpart,maxpart
! If a vacant storage space is found, attribute everything to this array element
!*******************************************************************************
if (itra1(ipart).ne.itime) then
! Assign particle positions
!**************************
xtra1(ipart)=real(nx_we(k))
if (jy.eq.ny_sn(1)) then
ytra1(ipart)=real(jy)+0.5*ran1(idummy)
else if (jy.eq.ny_sn(2)) then
ytra1(ipart)=real(jy)-0.5*ran1(idummy)
else
ytra1(ipart)=real(jy)+(ran1(idummy)-.5)
endif
if (j.eq.1) then
ztra1(ipart)=zcolumn_we(k,jy,1)+(zcolumn_we(k,jy,2)- &
zcolumn_we(k,jy,1))/4.
else if (j.eq.numcolumn_we(k,jy)) then
ztra1(ipart)=(2.*zcolumn_we(k,jy,j)+ &
zcolumn_we(k,jy,j-1)+height(nz))/4.
else
ztra1(ipart)=zcolumn_we(k,jy,j-1)+ran1(idummy)* &
(zcolumn_we(k,jy,j+1)-zcolumn_we(k,jy,j-1))
endif
! Interpolate PV to the particle position
!****************************************
ixm=int(xtra1(ipart))
jym=int(ytra1(ipart))
ixp=ixm+1
jyp=jym+1
ddx=xtra1(ipart)-real(ixm)
ddy=ytra1(ipart)-real(jym)
rddx=1.-ddx
rddy=1.-ddy
p1=rddx*rddy
p2=ddx*rddy
p3=rddx*ddy
p4=ddx*ddy
do i=2,nz
if (height(i).gt.ztra1(ipart)) then
indzm=i-1
indzp=i
goto 26
endif
enddo
26 continue
dz1=ztra1(ipart)-height(indzm)
dz2=height(indzp)-ztra1(ipart)
dz=1./(dz1+dz2)
do mm=1,2
indexh=memind(mm)
do in=1,2
indzh=indzm+in-1
y1(in)=p1*pv(ixm,jym,indzh,indexh) &
+p2*pv(ixp,jym,indzh,indexh) &
+p3*pv(ixm,jyp,indzh,indexh) &
+p4*pv(ixp,jyp,indzh,indexh)
enddo
yh1(mm)=(dz2*y1(1)+dz1*y1(2))*dz
enddo
pvpart=(yh1(1)*dt2+yh1(2)*dt1)*dtt
!JB
! ylat=ylat0+ytra1(ipart)*dy
if (ylat.lt.0.) pvpart=-1.*pvpart
! For domain-filling option 2 (stratospheric O3), do the rest only in the stratosphere
!*************************************************************************************
if (((ztra1(ipart).gt.3000.).and. &
(pvpart.gt.pvcrit)).or.(mdomainfill.eq.1)) then
! if (((ztra1(ipart).lt.8000.)
! + ).or.(mdomainfill.eq.1)) then
nclass(ipart)=min(int(ran1(idummy)* &
real(nclassunc))+1,nclassunc)
numactiveparticles=numactiveparticles+1
numparticlecount=numparticlecount+1
npoint(ipart)=numparticlecount
idt(ipart)=mintime
itra1(ipart)=itime
itramem(ipart)=itra1(ipart)
itrasplit(ipart)=itra1(ipart)+ldirect*itsplit
xmass1(ipart,1)=xmassperparticle
if (mdomainfill.eq.2) xmass1(ipart,1)= &
xmass1(ipart,1)*pvpart*48./29.*ozonescale/10.**9
! + xmass1(ipart,1)*60.*48./29./10.**9
else
goto 71
endif
! Increase numpart, if necessary
!*******************************
numpart=max(numpart,ipart)
goto 73 ! Storage space has been found, stop searching
endif
enddo
if (ipart.gt.maxpart) &
stop 'boundcond_domainfill.f: too many particles required'
73 minpart=ipart+1
71 continue
enddo
enddo
enddo
enddo
!*****************************************
! Southern and northern boundary condition
!*****************************************
! Loop from west to east
!***********************
do ix=nx_we(1),nx_we(2)
! Loop over southern (index 1) and northern (index 2) boundary
!*************************************************************
do k=1,2
! for FLEXPART_WRF, x & y coords are in meters.
! ylat=ylat0+real(ny_sn(k))*dy
! cosfact=cos(ylat*pi180)
! In the "do 170" loop, ylat is only needed for for pv calcs,
! "if (ylat.lt.0.) pvpart=-1.*pvpart"
dumx=real(ix)
dumy=real(ny_sn(k))
! Are these dumx,dumy correct ???
call xyindex_to_ll_wrf( 0, dumx, dumy, xlon, ylat )
! Loop over all release locations in a column
!********************************************
do j=1,numcolumn_sn(k,ix)
! Determine, for each release location, the area of the corresponding boundary
!*****************************************************************************
if (j.eq.1) then
deltaz=(zcolumn_sn(k,ix,2)+zcolumn_sn(k,ix,1))/2.
else if (j.eq.numcolumn_sn(k,ix)) then
! deltaz=height(nz)-(zcolumn_sn(k,ix,j-1)+
! + zcolumn_sn(k,ix,j))/2.
! In order to avoid taking a very high column for very many particles,
! use the deltaz from one particle below instead
deltaz=(zcolumn_sn(k,ix,j)-zcolumn_sn(k,ix,j-2))/2.
else
deltaz=(zcolumn_sn(k,ix,j+1)-zcolumn_sn(k,ix,j-1))/2.
endif
! for FLEXPART_ECMWF, dx is in degrees-long, and 111198.5*cosfact converts
! from degrees-longitude to m
! for FLEXPART_WRF, dx is in meters
! if ((ix.eq.nx_we(1)).or.(ix.eq.nx_we(2))) then
! boundarea=deltaz*111198.5/2.*cosfact*dx
! else
! boundarea=deltaz*111198.5*cosfact*dx
! endif
if ((ix.eq.nx_we(1)).or.(ix.eq.nx_we(2))) then
boundarea=deltaz/2.*dx
else
boundarea=deltaz*dx
endif
! Interpolate the wind velocity and density to the release location
!******************************************************************
! Determine the model level below the release position
!*****************************************************
do i=2,nz
if (height(i).gt.zcolumn_sn(k,ix,j)) then
indz=i-1
indzp=i
goto 16
endif
enddo
16 continue
! Vertical distance to the level below and above current position
!****************************************************************
dz1=zcolumn_sn(k,ix,j)-height(indz)
dz2=height(indzp)-zcolumn_sn(k,ix,j)
dz=1./(dz1+dz2)
! Vertical and temporal interpolation
!************************************
do m=1,2
indexh=memind(m)
do in=1,2
indzh=indz+in-1
windl(in)=vv(ix,ny_sn(k),indzh,indexh)
rhol(in)=rho(ix,ny_sn(k),indzh,indexh)
end do
windhl(m)=(dz2*windl(1)+dz1*windl(2))*dz
rhohl(m)=(dz2*rhol(1)+dz1*rhol(2))*dz
end do
windx=(windhl(1)*dt2+windhl(2)*dt1)*dtt
rhox=(rhohl(1)*dt2+rhohl(2)*dt1)*dtt
! Calculate mass flux
!********************
fluxofmass=windx*rhox*boundarea*real(lsynctime)
! If the mass flux is directed into the domain, add it to previous mass fluxes;
! if it is out of the domain, set accumulated mass flux to zero
!******************************************************************************
if (k.eq.1) then
if (fluxofmass.ge.0.) then
acc_mass_sn(k,ix,j)=acc_mass_sn(k,ix,j)+fluxofmass
else
acc_mass_sn(k,ix,j)=0.
endif
else
if (fluxofmass.le.0.) then
acc_mass_sn(k,ix,j)=acc_mass_sn(k,ix,j)+abs(fluxofmass)
else
acc_mass_sn(k,ix,j)=0.
endif
endif
accmasst=accmasst+acc_mass_sn(k,ix,j)
! If the accumulated mass exceeds half the mass that each particle shall carry,
! one (or more) particle(s) is (are) released and the accumulated mass is
! reduced by the mass of this (these) particle(s)
!******************************************************************************
if (acc_mass_sn(k,ix,j).ge.xmassperparticle/2.) then
mmass=int((acc_mass_sn(k,ix,j)+xmassperparticle/2.)/ &
xmassperparticle)
acc_mass_sn(k,ix,j)=acc_mass_sn(k,ix,j)- &
real(mmass)*xmassperparticle
else
mmass=0
endif
do m=1,mmass
do ipart=minpart,maxpart
! If a vacant storage space is found, attribute everything to this array element
!*******************************************************************************
if (itra1(ipart).ne.itime) then
! Assign particle positions
!**************************
ytra1(ipart)=real(ny_sn(k))
if (ix.eq.nx_we(1)) then
xtra1(ipart)=real(ix)+0.5*ran1(idummy)
else if (ix.eq.nx_we(2)) then
xtra1(ipart)=real(ix)-0.5*ran1(idummy)
else
xtra1(ipart)=real(ix)+(ran1(idummy)-.5)
endif
if (j.eq.1) then
ztra1(ipart)=zcolumn_sn(k,ix,1)+(zcolumn_sn(k,ix,2)- &
zcolumn_sn(k,ix,1))/4.
else if (j.eq.numcolumn_sn(k,ix)) then
ztra1(ipart)=(2.*zcolumn_sn(k,ix,j)+ &
zcolumn_sn(k,ix,j-1)+height(nz))/4.
else
ztra1(ipart)=zcolumn_sn(k,ix,j-1)+ran1(idummy)* &
(zcolumn_sn(k,ix,j+1)-zcolumn_sn(k,ix,j-1))
endif
! Interpolate PV to the particle position
!****************************************
ixm=int(xtra1(ipart))
jym=int(ytra1(ipart))
ixp=ixm+1
jyp=jym+1
ddx=xtra1(ipart)-real(ixm)
ddy=ytra1(ipart)-real(jym)
rddx=1.-ddx
rddy=1.-ddy
p1=rddx*rddy
p2=ddx*rddy
p3=rddx*ddy
p4=ddx*ddy
do i=2,nz
if (height(i).gt.ztra1(ipart)) then
indzm=i-1
indzp=i
goto 126
endif
enddo
126 continue
dz1=ztra1(ipart)-height(indzm)
dz2=height(indzp)-ztra1(ipart)
dz=1./(dz1+dz2)
do mm=1,2
indexh=memind(mm)
do in=1,2
indzh=indzm+in-1
y1(in)=p1*pv(ixm,jym,indzh,indexh) &
+p2*pv(ixp,jym,indzh,indexh) &
+p3*pv(ixm,jyp,indzh,indexh) &
+p4*pv(ixp,jyp,indzh,indexh)
end do
yh1(mm)=(dz2*y1(1)+dz1*y1(2))*dz
end do
pvpart=(yh1(1)*dt2+yh1(2)*dt1)*dtt
if (ylat.lt.0.) pvpart=-1.*pvpart
! For domain-filling option 2 (stratospheric O3), do the rest only in the stratosphere
!*************************************************************************************
if (((ztra1(ipart).gt.3000.).and. &
(pvpart.gt.pvcrit)).or.(mdomainfill.eq.1)) then
nclass(ipart)=min(int(ran1(idummy)* &
real(nclassunc))+1,nclassunc)
numactiveparticles=numactiveparticles+1
numparticlecount=numparticlecount+1
npoint(ipart)=numparticlecount
idt(ipart)=mintime
itra1(ipart)=itime
itramem(ipart)=itra1(ipart)
itrasplit(ipart)=itra1(ipart)+ldirect*itsplit
xmass1(ipart,1)=xmassperparticle
if (mdomainfill.eq.2) xmass1(ipart,1)= &
xmass1(ipart,1)*pvpart*48./29.*ozonescale/10.**9
else
goto 171
endif
! Increase numpart, if necessary
!*******************************
numpart=max(numpart,ipart)
goto 173 ! Storage space has been found, stop searching
endif
enddo
if (ipart.gt.maxpart) &
stop 'boundcond_domainfill.f: too many particles required'
173 minpart=ipart+1
171 continue
enddo
enddo
enddo
enddo
xm=0.
do i=1,numpart
if (itra1(i).eq.itime) xm=xm+xmass1(i,1)
end do
!write(*,*) itime,numactiveparticles,numparticlecount,numpart,
! +xm,accmasst,xm+accmasst
! If particles shall be dumped, then accumulated masses at the domain boundaries
! must be dumped, too, to be used for later runs
!*****************************************************************************
if ((ipout.gt.0).and.(itime.eq.loutend)) then
open(unitboundcond,file=path(2)(1:length(2))//'boundcond.bin', &
form='unformatted')
write(unitboundcond) numcolumn_we,numcolumn_sn, &
zcolumn_we,zcolumn_sn,acc_mass_we,acc_mass_sn
close(unitboundcond)
endif
end subroutine boundcond_domainfill