source: flexpart.git/src/conccalc_mpi.f90

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file was 92fab65, checked in by Ignacio Pisso <ip@…>, 4 years ago

add SPDX-License-Identifier to all .f90 files

  • Property mode set to 100644
File size: 18.2 KB
Line 
1! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt
2! SPDX-License-Identifier: GPL-3.0-or-later
3
4subroutine conccalc(itime,weight)
5  !                 i     i
6  !*****************************************************************************
7  !                                                                            *
8  !     Calculation of the concentrations on a regular grid using volume       *
9  !     sampling                                                               *
10  !                                                                            *
11  !     Author: A. Stohl                                                       *
12  !                                                                            *
13  !     24 May 1996                                                            *
14  !                                                                            *
15  !     April 2000: Update to calculate age spectra                            *
16  !                 Bug fix to avoid negative conc. at the domain boundaries,  *
17  !                 as suggested by Petra Seibert                              *
18  !                                                                            *
19  !     2 July 2002: re-order if-statements in order to optimize CPU time      *
20  !                                                                            *
21  !                                                                            *
22  !*****************************************************************************
23  !                                                                            *
24  ! Variables:                                                                 *
25  ! nspeciesdim     = nspec for forward runs, 1 for backward runs              *
26  !                                                                            *
27  !*****************************************************************************
28
29  use unc_mod
30  use outg_mod
31  use par_mod
32  use com_mod
33  use mpi_mod
34
35  implicit none
36
37  integer,intent(in) :: itime
38  integer :: itage,i,ix,jy,ixp,jyp,kz,ks,n,nage
39  integer :: il,ind,indz,indzp,nrelpointer
40  real,intent(in) :: weight
41  real :: rddx,rddy,p1,p2,p3,p4,dz1,dz2,dz
42  real :: hx,hy,hz,h,xd,yd,zd,xkern,r2,c(maxspec),ddx,ddy
43  real :: rhoprof(2),rhoi
44  real :: xl,yl,wx,wy,w
45  real,parameter :: factor=.596831, hxmax=6.0, hymax=4.0, hzmax=150.
46  integer :: mind2
47  ! mind2        eso: pointer to 2nd windfield in memory
48
49  ! For forward simulations, make a loop over the number of species;
50  ! for backward simulations, make an additional loop over the
51  ! releasepoints
52  !***************************************************************************
53
54  if (mp_measure_time) call mpif_mtime('conccalc',0)
55
56  mind2=memind(2)
57
58  do i=1,numpart
59    if (itra1(i).ne.itime) goto 20
60
61  ! Determine age class of the particle
62    itage=abs(itra1(i)-itramem(i))
63    do nage=1,nageclass
64      if (itage.lt.lage(nage)) goto 33
65    end do
6633   continue
67
68
69  ! For special runs, interpolate the air density to the particle position
70  !************************************************************************
71  !***********************************************************************
72  !AF IND_SOURCE switches between different units for concentrations at the source
73  !Af    NOTE that in backward simulations the release of particles takes place
74  !Af    at the receptor and the sampling at the source.
75  !Af          1="mass"
76  !Af          2="mass mixing ratio"
77  !Af IND_RECEPTOR switches between different units for concentrations at the receptor
78  !Af          1="mass"
79  !Af          2="mass mixing ratio"
80
81  !Af switches for the conccalcfile:
82  !AF IND_SAMP =  0 : xmass * 1
83  !Af IND_SAMP = -1 : xmass / rho
84
85  !Af ind_samp is defined in readcommand.f
86
87    if ( ind_samp .eq. -1 ) then
88
89      ix=int(xtra1(i))
90      jy=int(ytra1(i))
91      ixp=ix+1
92      jyp=jy+1
93      ddx=xtra1(i)-real(ix)
94      ddy=ytra1(i)-real(jy)
95      rddx=1.-ddx
96      rddy=1.-ddy
97      p1=rddx*rddy
98      p2=ddx*rddy
99      p3=rddx*ddy
100      p4=ddx*ddy
101
102! eso: Temporary fix for particle exactly at north pole
103      if (jyp >= nymax) then
104        ! write(*,*) 'WARNING: conccalc.f90 jyp >= nymax'
105        jyp=jyp-1
106      end if
107
108      do il=2,nz
109        if (height(il).gt.ztra1(i)) then
110          indz=il-1
111          indzp=il
112          goto 6
113        endif
114      end do
1156     continue
116
117      dz1=ztra1(i)-height(indz)
118      dz2=height(indzp)-ztra1(i)
119      dz=1./(dz1+dz2)
120
121  ! Take density from 2nd wind field in memory (accurate enough, no time interpolation needed)
122  !*****************************************************************************
123      do ind=indz,indzp
124        ! rhoprof(ind-indz+1)=p1*rho(ix ,jy ,ind,2) &
125        !      +p2*rho(ixp,jy ,ind,2) &
126        !      +p3*rho(ix ,jyp,ind,2) &
127        !      +p4*rho(ixp,jyp,ind,2)
128        rhoprof(ind-indz+1)=p1*rho(ix ,jy ,ind,mind2) &
129             +p2*rho(ixp,jy ,ind,mind2) &
130             +p3*rho(ix ,jyp,ind,mind2) &
131             +p4*rho(ixp,jyp,ind,mind2)
132
133      end do
134      rhoi=(dz1*rhoprof(2)+dz2*rhoprof(1))*dz
135   elseif (ind_samp.eq.0) then
136      rhoi = 1.
137   endif
138
139
140  !****************************************************************************
141  ! 1. Evaluate grid concentrations using a uniform kernel of bandwidths dx, dy
142  !****************************************************************************
143
144
145  ! For backward simulations, look from which release point the particle comes from
146  ! For domain-filling trajectory option, npoint contains a consecutive particle
147  ! number, not the release point information. Therefore, nrelpointer is set to 1
148  ! for the domain-filling option.
149  !*****************************************************************************
150
151    if ((ioutputforeachrelease.eq.0).or.(mdomainfill.eq.1)) then
152       nrelpointer=1
153    else
154       nrelpointer=npoint(i)
155    endif
156
157    do kz=1,numzgrid                ! determine height of cell
158      if (outheight(kz).gt.ztra1(i)) goto 21
159    end do
16021   continue
161    if (kz.le.numzgrid) then           ! inside output domain
162
163
164  !********************************
165  ! Do everything for mother domain
166  !********************************
167
168      xl=(xtra1(i)*dx+xoutshift)/dxout
169      yl=(ytra1(i)*dy+youtshift)/dyout
170      ix=int(xl)
171      if (xl.lt.0.) ix=ix-1
172      jy=int(yl)
173      if (yl.lt.0.) jy=jy-1
174
175  ! if (i.eq.10000) write(*,*) itime,xtra1(i),ytra1(i),ztra1(i),xl,yl
176
177
178  ! For particles aged less than 3 hours, attribute particle mass to grid cell
179  ! it resides in rather than use the kernel, in order to avoid its smoothing effect.
180  ! For older particles, use the uniform kernel.
181  ! If a particle is close to the domain boundary, do not use the kernel either.
182  !*****************************************************************************
183
184      if ((.not.lusekerneloutput).or.(itage.lt.10800).or. &
185           (xl.lt.0.5).or.(yl.lt.0.5).or. &
186           (xl.gt.real(numxgrid-1)-0.5).or. &
187           (yl.gt.real(numygrid-1)-0.5)) then             ! no kernel, direct attribution to grid cell
188        if ((ix.ge.0).and.(jy.ge.0).and.(ix.le.numxgrid-1).and. &
189             (jy.le.numygrid-1)) then
190          if (DRYBKDEP.or.WETBKDEP) then
191            do ks=1,nspec
192              gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
193                   gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
194                   xmass1(i,ks)/rhoi*weight*max(xscav_frac1(i,ks),0.0)
195            end do
196          else
197            if (lparticlecountoutput) then
198              do ks=1,nspec
199                gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
200                     gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+1
201              end do
202            else
203              do ks=1,nspec
204                gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
205                     gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
206                     xmass1(i,ks)/rhoi*weight
207              end do
208            end if
209          endif
210        endif
211
212      else                                 ! attribution via uniform kernel
213
214        ddx=xl-real(ix)                   ! distance to left cell border
215        ddy=yl-real(jy)                   ! distance to lower cell border
216        if (ddx.gt.0.5) then
217          ixp=ix+1
218          wx=1.5-ddx
219        else
220          ixp=ix-1
221          wx=0.5+ddx
222        endif
223
224        if (ddy.gt.0.5) then
225          jyp=jy+1
226          wy=1.5-ddy
227        else
228          jyp=jy-1
229          wy=0.5+ddy
230        endif
231
232
233  ! Determine mass fractions for four grid points
234  !**********************************************
235
236        if ((ix.ge.0).and.(ix.le.numxgrid-1)) then
237          if ((jy.ge.0).and.(jy.le.numygrid-1)) then
238            w=wx*wy
239            if (DRYBKDEP.or.WETBKDEP) then
240               do ks=1,nspec
241                 gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
242                   gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
243                   xmass1(i,ks)/rhoi*w*weight*max(xscav_frac1(i,ks),0.0)
244               end do
245            else
246               do ks=1,nspec
247                 gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
248                   gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
249                   xmass1(i,ks)/rhoi*weight*w
250               end do
251            endif
252          endif
253
254          if ((jyp.ge.0).and.(jyp.le.numygrid-1)) then
255            w=wx*(1.-wy)
256            if (DRYBKDEP.or.WETBKDEP) then
257              do ks=1,nspec
258                 gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
259                   gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
260                   xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
261               end do
262             else
263              do ks=1,nspec
264                 gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
265                   gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
266                   xmass1(i,ks)/rhoi*weight*w
267               end do
268             endif
269          endif
270        endif
271
272
273        if ((ixp.ge.0).and.(ixp.le.numxgrid-1)) then
274          if ((jyp.ge.0).and.(jyp.le.numygrid-1)) then
275            w=(1.-wx)*(1.-wy)
276            if (DRYBKDEP.or.WETBKDEP) then
277               do ks=1,nspec
278                 gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
279                   gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
280                   xmass1(i,ks)/rhoi*w*weight*max(xscav_frac1(i,ks),0.0)
281               end do
282            else
283               do ks=1,nspec
284                 gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
285                   gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
286                   xmass1(i,ks)/rhoi*weight*w
287               end do
288            endif
289          endif
290
291          if ((jy.ge.0).and.(jy.le.numygrid-1)) then
292            w=(1.-wx)*wy
293            if (DRYBKDEP.or.WETBKDEP) then
294               do ks=1,nspec
295                 gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
296                   gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
297                   xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
298               end do
299            else
300               do ks=1,nspec
301                 gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
302                   gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
303                   xmass1(i,ks)/rhoi*weight*w
304               end do
305            endif
306          endif
307        endif
308      endif
309
310
311
312  !************************************
313  ! Do everything for the nested domain
314  !************************************
315
316      if (nested_output.eq.1) then
317        xl=(xtra1(i)*dx+xoutshiftn)/dxoutn
318        yl=(ytra1(i)*dy+youtshiftn)/dyoutn
319        ix=int(xl)
320        if (xl.lt.0.) ix=ix-1
321        jy=int(yl)
322        if (yl.lt.0.) jy=jy-1
323
324
325  ! For particles aged less than 3 hours, attribute particle mass to grid cell
326  ! it resides in rather than use the kernel, in order to avoid its smoothing effect.
327  ! For older particles, use the uniform kernel.
328  ! If a particle is close to the domain boundary, do not use the kernel either.
329  !*****************************************************************************
330
331        if ((itage.lt.10800).or.(xl.lt.0.5).or.(yl.lt.0.5).or. &
332             (xl.gt.real(numxgridn-1)-0.5).or. &
333             (yl.gt.real(numygridn-1)-0.5).or.((.not.lusekerneloutput))) then
334! no kernel, direct attribution to grid cell
335          if ((ix.ge.0).and.(jy.ge.0).and.(ix.le.numxgridn-1).and. &
336               (jy.le.numygridn-1)) then
337            if (DRYBKDEP.or.WETBKDEP) then
338               do ks=1,nspec
339                 griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
340                   griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
341                   xmass1(i,ks)/rhoi*weight*max(xscav_frac1(i,ks),0.0)
342               end do
343            else
344              if (lparticlecountoutput) then
345                do ks=1,nspec
346                  griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
347                       griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+1
348                end do
349              else
350                do ks=1,nspec
351                  griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
352                       griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
353                       xmass1(i,ks)/rhoi*weight
354                end do
355              endif
356            endif
357          endif
358
359        else                                 ! attribution via uniform kernel
360
361          ddx=xl-real(ix)                   ! distance to left cell border
362          ddy=yl-real(jy)                   ! distance to lower cell border
363          if (ddx.gt.0.5) then
364            ixp=ix+1
365            wx=1.5-ddx
366          else
367            ixp=ix-1
368            wx=0.5+ddx
369          endif
370
371          if (ddy.gt.0.5) then
372            jyp=jy+1
373            wy=1.5-ddy
374          else
375            jyp=jy-1
376            wy=0.5+ddy
377          endif
378
379
380  ! Determine mass fractions for four grid points
381  !**********************************************
382
383          if ((ix.ge.0).and.(ix.le.numxgridn-1)) then
384            if ((jy.ge.0).and.(jy.le.numygridn-1)) then
385              w=wx*wy
386              if (DRYBKDEP.or.WETBKDEP) then
387                 do ks=1,nspec
388                   griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
389                     griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
390                     xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
391                 end do
392              else
393                do ks=1,nspec
394                   griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
395                     griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
396                     xmass1(i,ks)/rhoi*weight*w
397                 end do
398              endif
399            endif
400
401            if ((jyp.ge.0).and.(jyp.le.numygridn-1)) then
402              w=wx*(1.-wy)
403              if (DRYBKDEP.or.WETBKDEP) then
404                 do ks=1,nspec
405                   griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
406                     griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
407                     xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
408                 end do
409              else
410                 do ks=1,nspec
411                   griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
412                     griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
413                     xmass1(i,ks)/rhoi*weight*w
414                 end do
415              endif
416            endif
417          endif
418
419
420          if ((ixp.ge.0).and.(ixp.le.numxgridn-1)) then
421            if ((jyp.ge.0).and.(jyp.le.numygridn-1)) then
422              w=(1.-wx)*(1.-wy)
423              if (DRYBKDEP.or.WETBKDEP) then
424                 do ks=1,nspec
425                   griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
426                     griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
427                     xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
428                 end do
429              else
430                 do ks=1,nspec
431                   griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
432                     griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
433                     xmass1(i,ks)/rhoi*weight*w
434                 end do
435              endif
436            endif
437
438            if ((jy.ge.0).and.(jy.le.numygridn-1)) then
439              w=(1.-wx)*wy
440              if (DRYBKDEP.or.WETBKDEP) then
441                 do ks=1,nspec
442                   griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
443                     griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
444                     xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
445                 end do
446              else
447                 do ks=1,nspec
448                    griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
449                     griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
450                     xmass1(i,ks)/rhoi*weight*w
451                 end do
452              endif
453            endif
454          endif
455        endif
456
457      endif
458    endif
45920  continue
460  end do
461
462  !***********************************************************************
463  ! 2. Evaluate concentrations at receptor points, using the kernel method
464  !***********************************************************************
465
466  do n=1,numreceptor
467
468
469  ! Reset concentrations
470  !*********************
471
472    do ks=1,nspec
473      c(ks)=0.
474    end do
475
476
477  ! Estimate concentration at receptor
478  !***********************************
479
480    do i=1,numpart
481
482      if (itra1(i).ne.itime) goto 40
483      itage=abs(itra1(i)-itramem(i))
484
485      hz=min(50.+0.3*sqrt(real(itage)),hzmax)
486      zd=ztra1(i)/hz
487      if (zd.gt.1.) goto 40          ! save computing time, leave loop
488
489      hx=min((0.29+2.222e-3*sqrt(real(itage)))*dx+ &
490           real(itage)*1.2e-5,hxmax)                     ! 80 km/day
491      xd=(xtra1(i)-xreceptor(n))/hx
492      if (xd*xd.gt.1.) goto 40       ! save computing time, leave loop
493
494      hy=min((0.18+1.389e-3*sqrt(real(itage)))*dy+ &
495           real(itage)*7.5e-6,hymax)                     ! 80 km/day
496      yd=(ytra1(i)-yreceptor(n))/hy
497      if (yd*yd.gt.1.) goto 40       ! save computing time, leave loop
498      h=hx*hy*hz
499
500      r2=xd*xd+yd*yd+zd*zd
501      if (r2.lt.1.) then
502        xkern=factor*(1.-r2)
503        do ks=1,nspec
504          c(ks)=c(ks)+xmass1(i,ks)*xkern/h
505        end do
506      endif
50740    continue
508    end do
509
510    do ks=1,nspec
511      creceptor(n,ks)=creceptor(n,ks)+2.*weight*c(ks)/receptorarea(n)
512    end do
513  end do
514
515  if (mp_measure_time) call mpif_mtime('conccalc',1)
516
517end subroutine conccalc
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG