-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathdiffus.f
404 lines (401 loc) · 13.9 KB
/
diffus.f
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
subroutine diffus(igrd,ncol,nrow,nlay,nspc,nsen,deltat,dx,dy,
& idfin,vdep,rkx,rky,rkv,depth,tempk,press,mapscl,
& conc,fluxes,depfld,sens,tarray2,strz,strxy,
& ipa_xy,ipa_lay)
c
c-----CAMx v4.02 030709
c
c DIFFUS drives 3-D diffusion of concentrations. This version also
c performs diffusion of sensitivities if DDM is enabled.
c
c Copyright 1996, 1997, 1998, 1999, 2000, 2001, 2002, 2003
c ENVIRON International Corporation
c
c Modifications:
c 4/17/00 Revised diffusion equations to weight fluxes by density
c 12/07/01 added instructions for OMP
c 1/13/03 added deposited mass array
c
c Input arguments:
c igrd grid index
c ncol number of columns
c nrow number of rows
c nlay number of layers
c nspc number of species
c nsen number of species times number of DDM parameters
c or 1, whichever is larger
c deltat time step (s)
c dx cell size in x-direction (m)
c dy cell size in y-direction (m)
c idfin map of nested grids in this grid
c vdep deposition velocity (m/s)
c rkx/y horizontal diffusion coefficient (m2/s)
c rkv vertical diffusion coefficient (m2/s)
c depth layer depth (m)
c tempk temperature (K)
c press pressure (mb)
c mapscl map-scale factor at cell centroid
c conc species concentrations (umol/m3)
c sens DDM sensitivities (umol/m3/parameter unit)
c tarray2 CPU timing arguments (s)
c strz string for labeling the z diffusion process
c strxy string for labeling the x/y diffusion process
c ipa_xy 2-D gridded array to identify if cell is
c in a IPRM sub-domain
c ipa_lay 3-D gridded array to identify which IPRM sub-domain
c each layer is in
c
c Output arguments:
c conc species concentrations (umol/m3)
c sens DDM sensitivities (umol/m3/parameter unit)
c fluxes fluxes across the boundaries (umol)
c depfld 2-D array of dry deposited mass (mol/ha, g/ha)
c
c Routines Called:
c VDIFFIMP
c
c Called by:
c EMISTRNS
c
include "camx.prm"
include "filunit.com"
include "bndary.com"
include "chmstry.com"
include "App.com" !Added by Kristina 05/08/07
include "flags.com" !Added by Kristina 05/08/07
c
c=============================DDM Begin================================
c
include "tracer.com"
c
c=============================DDM End==================================
c
c
c======================== Process Analysis Begin ====================================
c
include "procan.com"
c
integer ipa_xy(ncol,nrow)
integer ipa_lay(ncol,nrow,nlay)
real fcup(MXLAYA),fcdn(MXLAYA)
logical ldoipts
c
c========================= Process Analysis End =====================================
c
dimension conc(ncol,nrow,nlay,nspc),vdep(ncol,nrow,nspc)
dimension sens(ncol,nrow,nlay,nsen)
dimension tarray2(2)
real depfld(ncol,nrow,3*nspc)
real rkx(ncol,nrow,nlay),rky(ncol,nrow,nlay),
& rkv(ncol,nrow,nlay),depth(ncol,nrow,nlay),
& tempk(ncol,nrow,nlay),press(ncol,nrow,nlay),
& mapscl(ncol,nrow),dx(nrow)
integer idfin(ncol,nrow)
dimension c1d(MXLAYA+MXLAYA*MXTRSP),d1d(MXLAYA),rk1d(MXLAYA),
& ro1d(MXLAYA)
real cnc(MXCOLA,MXROWA),sns(MXCOLA,MXROWA,MXTRSP),
& rho(MXCOLA,MXROWA)
real*8 fluxes(nspc,13),fluxbot
character*20 strz, strxy
real flxsa(8),original(MXCOL1,MXROW1,MXSOUR),Apptot(MXCOL1,MXROW1) !Added by Kristina 05/08/07
integer nx,ny,nz
c
c-----Entry point
c
c----------Added by Kristina 05/08/07------------------------------------
c
nx = ncol
ny = nrow
nz = nlay
c
c-----------End Added 05/08/07-------------------------------------------
c
c-----Vertical diffusion
c
write(*,'(a20,$)') strz
write(iout,'(a20,$)') strz
c
c$omp parallel default(shared)
c$omp& private(ispc,i1,i2,i,j,k,ro1d,d1d,c1d,rk1d,lk,
c$omp& ioff,fluxbot,isen,rokp1,ldoipts,fcup,
c$omp& fcdn,ipa_idx)
c
c$omp do schedule(dynamic)
c
do 61 ispc = 1,nspc
do 60 j = 2,nrow-1
i1 = 2
i2 = ncol - 1
if (igrd.eq.1) then
if (ibeg(j).eq.-999) goto 60
i1 = ibeg(j)
i2 = iend(j)
endif
do 50 i = i1,i2
c
c-----Skip cells occupied by child grids; load 1-D arrays
c
if (idfin(i,j).gt.igrd) goto 50
do k = 1,nlay
ro1d(k) = press(i,j,k)/tempk(i,j,k)
d1d(k) = depth(i,j,k)
c1d(k) = conc(i,j,k,ispc)
enddo
do k = 1,nlay-1
rokp1 = press(i,j,k+1)/tempk(i,j,k+1)
rk1d(k) = rkv(i,j,k)*(ro1d(k) + rokp1)/2.
enddo
rk1d(nlay) = 0.
c
c
c=============================DDM Begin================================
c
c------Load sensitivities
c
if (lddm) then
lk = nlay
do isen = 1,nddmsp
ioff = iptddm(ispc)+isen-1
do k = 1,nlay
lk = lk + 1
c1d(lk) = sens(i,j,k,ioff)
enddo
enddo
endif
c
c=============================DDM End==================================
c
c
c======================== Process Analysis Begin ====================================
c
ldoipts = .FALSE.
if( .NOT. ltrace .AND. lipr .AND. ipa_xy(i,j) .GT. 0)
& ldoipts = .TRUE.
c
c========================= Process Analysis End =====================================
c
call vdiffimp(nlay,deltat,vdep(i,j,ispc),d1d,ro1d,rk1d,c1d,
& nddmsp,fcup,fcdn,ldoipts)
c
c======================== Process Analysis Begin ====================================
c
if( ldoipts ) then
do k = 1,nlay
if( ipa_lay(i,j,k) .GT. 0) then
ipa_idx = ipa_lay(i,j,k)
if( k .GT. 1 ) then
cipr(IPR_BDIF, ipa_idx, ispc) =
& cipr(IPR_BDIF, ipa_idx, ispc) + fcdn(k)
else
cipr(IPR_DDEP, ipa_idx, ispc) =
& cipr(IPR_DDEP, ipa_idx, ispc) + fcdn(k)
endif
cipr(IPR_TDIF, ipa_idx, ispc) =
& cipr(IPR_TDIF, ipa_idx, ispc) + fcup(k)
endif
enddo
endif
c
c========================= Process Analysis End =====================================
c
do k = 1,nlay
conc(i,j,k,ispc) = c1d(k)
enddo
fluxbot = -vdep(i,j,ispc)*conc(i,j,1,ispc)
fluxes(ispc,11) = fluxes(ispc,11) + fluxbot*dx(j)*dy*deltat
do ll = 1,navspc
if (ispc .eq. lavmap(ll)) then
depfld(i,j,ll) = depfld(i,j,ll) - 1.e-2*fluxbot*deltat
goto 100
endif
enddo
100 continue
c
c=============================DDM Begin================================
c
if (lddm) then
lk = nlay
do isen = 1,nddmsp
ioff = iptddm(ispc)+isen-1
do k = 1,nlay
lk = lk + 1
sens(i,j,k,ioff) = c1d(lk)
enddo
enddo
endif
c
c=============================DDM End==================================
c
50 continue
60 continue
61 continue
c
c$omp end parallel
c$omp master
c
tcpu = dtime(tarray2)
write(*,'(a,f10.3)') ' CPU = ', tarray2(1)
write(iout,'(a,f10.3)') ' CPU = ', tarray2(1)
c
c-----Perform explicit horizontal diffusion
c
write(*,'(a20,$)') strxy
write(iout,'(a20,$)') strxy
c
c$omp end master
c
c$omp parallel default(shared)
c$omp& private(ispc,i1,i2,i,j,k,ioff,rho,cnc,sns,rhoavg,scl,
c$omp& dfxp,dfxm,dfyp,dfym,fxp,fxm,fyp,fym,isen)
c
c$omp do schedule(dynamic)
c
do 91 ispc = 1,nspc
do 90 k = 1,nlay
do 80 j = 1,nrow
do 70 i = 1,ncol
rho(i,j) = press(i,j,k)/tempk(i,j,k)
cnc(i,j) = conc(i,j,k,ispc)/rho(i,j)
70 continue
80 continue
c
c=============================DDM Begin================================
c
if (lddm) then
do isen = 1,nddmsp
ioff = iptddm(ispc)+isen-1
do j=1,nrow
do i = 1,ncol
sns(i,j,isen) = sens(i,j,k,ioff)/rho(i,j)
enddo
enddo
enddo
endif
c
c=============================DDM End==================================
c
do 85 j = 2,nrow-1
i1 = 2
i2 = ncol - 1
if (igrd.eq.1) then
if (ibeg(j).eq.-999) goto 85
i1 = ibeg(j)
i2 = iend(j)
endif
do 75 i = i1,i2
if (idfin(i,j).gt.igrd) goto 75
c
rhoavg = (rho(i,j) + rho(i+1,j))/2.
scl = (mapscl(i,j) + mapscl(i+1,j))/2.
dfxp = scl*rhoavg*rkx(i,j,k)*deltat/dx(j)/dx(j)
rhoavg = (rho(i,j) + rho(i-1,j))/2.
scl = (mapscl(i,j) + mapscl(i-1,j))/2.
dfxm = scl*rhoavg*rkx(i-1,j,k)*deltat/dx(j)/dx(j)
c
rhoavg = (rho(i,j) + rho(i,j+1))/2.
scl = (mapscl(i,j) + mapscl(i,j+1))/2.
dfyp = scl*rhoavg*rky(i,j,k)*deltat/dy/dy
rhoavg = (rho(i,j) + rho(i,j-1))/2.
scl = (mapscl(i,j) + mapscl(i,j-1))/2.
dfym = scl*rhoavg*rky(i,j-1,k)*deltat/dy/dy
c
fxp = (cnc(i+1,j) - cnc(i,j))*dfxp
fxm = (cnc(i,j) - cnc(i-1,j))*dfxm
fyp = (cnc(i,j+1) - cnc(i,j))*dfyp
fym = (cnc(i,j) - cnc(i,j-1))*dfym
conc(i,j,k,ispc) = conc(i,j,k,ispc) + mapscl(i,j)*
& ((fxp - fxm) + (fyp - fym))
c
c---------Added by Kristina 05/08/07---------------------------------------
c
c if(lApp) then
c flxsa(1) = cnc(i+1,j)*dfxp*mapscl(i,j)
c flxsa(2) = cnc(i,j)*dfxp*mapscl(i,j)
c flxsa(3) = cnc(i-1,j)*dfxm*mapscl(i,j)
c flxsa(4) = cnc(i,j)*dfxm*mapscl(i,j)
c flxsa(5) = cnc(i,j+1)*dfyp*mapscl(i,j)
c flxsa(6) = cnc(i,j)*dfyp*mapscl(i,j)
c flxsa(7) = cnc(i,j-1)*dfym*mapscl(i,j)
c flxsa(8) = cnc(i,j)*dfym*mapscl(i,j)
c
c save original values
c if (j.eq.2.and.i.eq.i1) then
c do loop1 = 1,MXCOL1
c do loop2 = 1,MXROW1
c Apptot(loop1,loop2) = 0
c do loop3 = 1,Appnum+3
c loc = loop1 + nx*(loop2-1) + nx*ny*(k-1) +
c & nx*ny*nz*(Appmap(ispc)-1) +
c & nx*ny*nz*MXTRK*(loop3-1)
c original(loop1,loop2,loop3) = Appconc(loc)
c Apptot(loop1,loop2) = Apptot(loop1,loop2)
c & + Appconc(loc)
c enddo
c enddo
c enddo
c endif
cc
c if (Appmap(ispc).ne.0)
c & call Apphdiff(flxsa,i,j,k,Appmap(ispc),original,
c & Apptot,conc(i,j,k,ispc))
c
c endif
c
c-------------End Added 05/08/07----------------------------------------------
c
c======================== Process Analysis Begin ====================================
c
if( .NOT. ltrace .AND. lipr .AND.
& ipa_lay(i,j,k) .GT. 0 ) then
ipa_idx = ipa_lay(i,j,k)
cipr(IPR_WDIF, ipa_idx, ispc) =
& cipr( IPR_WDIF, ipa_idx, ispc) -
& mapscl(i,j)*fxm
cipr(IPR_EDIF, ipa_idx, ispc) =
& cipr(IPR_EDIF, ipa_idx, ispc) +
& mapscl(i,j)*fxp
cipr(IPR_SDIF, ipa_idx, ispc) =
& cipr(IPR_SDIF, ipa_idx, ispc) -
& mapscl(i,j)*fym
cipr(IPR_NDIF, ipa_idx, ispc) =
& cipr(IPR_NDIF, ipa_idx, ispc) +
& mapscl(i,j)*fyp
endif
c
c========================= Process Analysis End =====================================
c
c
c=============================DDM Begin================================
c
if (lddm) then
do isen =1,nddmsp
fxp = (sns(i+1,j,isen) - sns(i,j,isen))*dfxp
fxm = (sns(i,j,isen) - sns(i-1,j,isen))*dfxm
fyp = (sns(i,j+1,isen) - sns(i,j,isen))*dfyp
fym = (sns(i,j,isen) - sns(i,j-1,isen))*dfym
ioff = iptddm(ispc)+isen-1
sens(i,j,k,ioff) = sens(i,j,k,ioff) + mapscl(i,j)*
& ((fxp - fxm) + (fyp - fym))
enddo
endif
c
c=============================DDM End==================================
c
75 continue
85 continue
90 continue
91 continue
c
c$omp end parallel
c$omp master
c
tcpu = dtime(tarray2)
write(*,'(a,f10.3)') ' CPU = ', tarray2(1)
write(iout,'(a,f10.3)') ' CPU = ', tarray2(1)
c
c$omp end master
c
call flush(6)
call flush(iout)
return
end