source: branches/flexpart91_hasod/src_parallel/wetdepo.f90 @ 8

Last change on this file since 8 was 8, checked in by hasod, 11 years ago

Added parallel version of Flexpart91

File size: 13.7 KB
Line 
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
22subroutine wetdepo(itime,ltsample,loutnext)
23  !                     i      i        i
24  !*****************************************************************************
25  !                                                                            *
26  ! Calculation of wet deposition using the concept of scavenging coefficients.*
27  ! For lack of detailed information, washout and rainout are jointly treated. *
28  ! It is assumed that precipitation does not occur uniformly within the whole *
29  ! grid cell, but that only a fraction of the grid cell experiences rainfall. *
30  ! This fraction is parameterized from total cloud cover and rates of large   *
31  ! scale and convective precipitation.                                        *
32  !                                                                            *
33  !    Author: A. Stohl                                                        *
34  !                                                                            *
35  !    1 December 1996                                                         *
36  !                                                                            *
37  ! Correction by Petra Seibert, Sept 2002:                                    *
38  ! use centred precipitation data for integration                             *
39  ! Code may not be correct for decay of deposition!                           *
40  !                                                                            *
41  !*****************************************************************************
42  !                                                                            *
43  ! Variables:                                                                 *
44  ! cc [0-1]           total cloud cover                                       *
45  ! convp [mm/h]       convective precipitation rate                           *
46  ! grfraction [0-1]   fraction of grid, for which precipitation occurs        *
47  ! ix,jy              indices of output grid cell for each particle           *
48  ! itime [s]          actual simulation time [s]                              *
49  ! jpart              particle index                                          *
50  ! ldeltat [s]        interval since radioactive decay was computed           *
51  ! lfr, cfr           area fraction covered by precipitation for large scale  *
52  !                    and convective precipitation (dependent on prec. rate)  *
53  ! loutnext [s]       time for which gridded deposition is next output        *
54  ! loutstep [s]       interval at which gridded deposition is output          *
55  ! lsp [mm/h]         large scale precipitation rate                          *
56  ! ltsample [s]       interval over which mass is deposited                   *
57  ! prec [mm/h]        precipitation rate in subgrid, where precipitation occurs*
58  ! wetdeposit         mass that is wet deposited                              *
59  ! wetgrid            accumulated deposited mass on output grid               *
60  ! wetscav            scavenging coefficient                                  *
61  !                                                                            *
62  ! Constants:                                                                 *
63  !                                                                            *
64  !*****************************************************************************
65
66!  use point_mod, only:
67  use par_mod, only: nymax, nymaxn, nxmax, nxmaxn, r_air, maxspec, maxnests, &
68       nclassunc
69  use com_mod, only: nz, nyn, tccn, tcc, itra1, itramem, yresoln, xresoln,  &
70       xrn, yrn, xln, yln, xtra1, ytra1, ztra1, lage, memind, height, clouds, &
71       memtime, nx, ny, weta, cloudsnh, cloudsn, cloudsh, nxn, henry, dquer, tt, ttn, &
72       wetb, decay, npoint, xmass1, nclass, numpart, numbnests, nspec, nested_output, &
73       nageclass, lsprecn, lsprec, loutstep, ldirect, convprec, convprecn, &
74       ioutputforeachrelease
75#if (defined _OPENMP)
76  use omp_lib, only: OMP_GET_THREAD_NUM
77#endif
78       
79
80  implicit none
81
82  integer :: jpart,itime,ltsample,loutnext,ldeltat,i,j,ix,jy
83  integer :: ngrid,itage,nage,hz,il,interp_time, n, clouds_v
84  integer :: ks, kp
85  real :: S_i, act_temp, cl, cle ! in cloud scavenging
86  real :: clouds_h ! cloud height for the specific grid point
87  real :: xtn,ytn,lsp,convp,cc,grfraction,prec,wetscav
88  real :: wetdeposit(maxspec),restmass
89  integer :: griclunc
90  real :: lfr(5) = (/ 0.5,0.65,0.8,0.9,0.95/)
91  real :: cfr(5) = (/ 0.4,0.55,0.7,0.8,0.9 /)
92
93  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
94#if (defined _OPENMP)
95  integer :: thread
96#endif
97
98  ! Compute interval since radioactive decay of deposited mass was computed
99  !************************************************************************
100
101  if (itime.le.loutnext) then
102    ldeltat=itime-(loutnext-loutstep)
103  else                                  ! first half of next interval
104    ldeltat=itime-loutnext
105  endif
106
107! determine time index
108  interp_time=nint(itime-0.5*ltsample)
109  n=memind(2)
110  if (abs(memtime(1)-interp_time).lt.abs(memtime(2)-interp_time)) &
111       n=memind(1)
112
113  ! Loop over all particles
114  !************************
115
116!$OMP PARALLEL DEFAULT(none) &
117!$OMP SHARED(numpart, itra1, itime, itramem, nageclass, lage, numbnests, xtra1, xln, &
118!$OMP   xrn, yln, yrn, ytra1, xresoln, yresoln, ltsample, lsprec, convprec, tcc, &
119!$OMP   nx, ny, memind, memtime, lsprecn, convprecn, tccn, ldirect, nyn, nxn, nz, height, &
120!$OMP   ztra1, n, clouds, cloudsh, cloudsn, cloudsnh, nspec, weta, wetb, tt, ttn,   &
121!$OMP   dquer, henry, xmass1, ioutputforeachrelease, npoint, decay, ldeltat, nclass, &
122!$OMP   nested_output, interp_time, lfr, cfr) &
123!$OMP PRIVATE(jpart, itage, nage, ngrid, j, xtn, ytn, ix, jy, lsp, convp, &
124!$OMP   cc, hz, clouds_v, clouds_h, i, grfraction, prec, wetscav, wetdeposit, act_temp, &
125!$OMP   cl, s_i, cle, restmass, kp, thread, griclunc, ks, il)
126
127#if (defined _OPENMP)
128      thread = OMP_GET_THREAD_NUM()
129#endif
130
131#if (defined STATIC_SCHED)
132!$OMP DO SCHEDULE (static)
133#else
134!$OMP DO SCHEDULE (dynamic, max(1, numpart/1000))
135#endif
136
137  particleloop : do jpart=1,numpart
138
139    if (itra1(jpart).eq.-999999999) cycle
140    if(ldirect.eq.1)then
141      if (itra1(jpart).gt.itime) cycle
142    else
143      if (itra1(jpart).lt.itime) cycle
144    endif
145  ! Determine age class of the particle
146    itage=abs(itra1(jpart)-itramem(jpart))
147    do nage=1,nageclass
148      if (itage.lt.lage(nage)) exit
149    end do
150
151
152  ! Determine which nesting level to be used
153  !*****************************************
154
155    ngrid=0
156    do j=numbnests,1,-1
157      if ((xtra1(jpart).gt.xln(j)).and.(xtra1(jpart).lt.xrn(j)).and. &
158           (ytra1(jpart).gt.yln(j)).and.(ytra1(jpart).lt.yrn(j))) then
159        ngrid=j
160        exit
161      endif
162    end do
163
164
165  ! Determine nested grid coordinates
166  !**********************************
167
168    if (ngrid.gt.0) then
169      xtn=(real(xtra1(jpart))-xln(ngrid))*xresoln(ngrid)
170      ytn=(real(ytra1(jpart))-yln(ngrid))*yresoln(ngrid)
171      ix=int(xtn)
172      jy=int(ytn)
173    else
174      ix=int(xtra1(jpart))
175      jy=int(ytra1(jpart))
176    endif
177
178
179  ! Interpolate large scale precipitation, convective precipitation and
180  ! total cloud cover
181  ! Note that interpolated time refers to itime-0.5*ltsample [PS]
182  !********************************************************************
183
184    if (ngrid.eq.0) then
185      call interpol_rain(lsprec,convprec,tcc,nxmax,nymax, &
186           1,nx,ny,memind,real(xtra1(jpart)),real(ytra1(jpart)),1, &
187           memtime(1),memtime(2),interp_time,lsp,convp,cc)
188    else
189      call interpol_rain_nests(lsprecn,convprecn,tccn, &
190           nxmaxn,nymaxn,1,maxnests,ngrid,nxn,nyn,memind,xtn,ytn,1, &
191           memtime(1),memtime(2),interp_time,lsp,convp,cc)
192    endif
193
194    if ((lsp.lt.0.01).and.(convp.lt.0.01)) cycle
195
196  ! get the level were the actual particle is in
197      do il=2,nz
198        if (height(il).gt.ztra1(jpart)) then
199          hz=il-1
200          exit
201        endif
202      end do
203
204
205  ! if there is no precipitation or the particle is above the clouds no
206  ! scavenging is done
207
208  if (ngrid.eq.0) then
209     clouds_v=clouds(ix,jy,hz,n)
210     clouds_h=cloudsh(ix,jy,n)
211  else
212     clouds_v=cloudsn(ix,jy,hz,n,ngrid)
213     clouds_h=cloudsnh(ix,jy,n,ngrid)
214  endif
215  !write(*,*) 'there is
216  !    + precipitation',(clouds(ix,jy,ihz,n),ihz=1,20),lsp,convp,hz
217  if (clouds_v.le.1) cycle particleloop
218  !write (*,*) 'there is scavenging'
219
220  ! 1) Parameterization of the the area fraction of the grid cell where the
221  !    precipitation occurs: the absolute limit is the total cloud cover, but
222  !    for low precipitation rates, an even smaller fraction of the grid cell
223  !    is used. Large scale precipitation occurs over larger areas than
224  !    convective precipitation.
225  !**************************************************************************
226
227    if (lsp.gt.20.) then
228      i=5
229    else if (lsp.gt.8.) then
230      i=4
231    else if (lsp.gt.3.) then
232      i=3
233    else if (lsp.gt.1.) then
234      i=2
235    else
236      i=1
237    endif
238
239    if (convp.gt.20.) then
240      j=5
241    else if (convp.gt.8.) then
242      j=4
243    else if (convp.gt.3.) then
244      j=3
245    else if (convp.gt.1.) then
246      j=2
247    else
248      j=1
249    endif
250
251    grfraction=max(0.05,cc*(lsp*lfr(i)+convp*cfr(j))/(lsp+convp))
252
253  ! 2) Computation of precipitation rate in sub-grid cell
254  !******************************************************
255
256    prec=(lsp+convp)/grfraction
257
258  ! 3) Computation of scavenging coefficients for all species
259  !    Computation of wet deposition
260  !**********************************************************
261
262    do ks=1,nspec                                  ! loop over species
263      wetdeposit(ks)=0.
264      if (weta(ks).gt.0.) then
265        if (clouds_v.ge.4) then
266  !          BELOW CLOUD SCAVENGING
267  !          for aerosols and not highliy soluble substances weta=5E-6
268          wetscav=weta(ks)*prec**wetb(ks)                ! scavenging coeff.
269  !        write(*,*) 'bel. wetscav: ',wetscav
270        else ! below_cloud clouds_v is lt 4 and gt 1 -> in cloud scavenging
271  !        IN CLOUD SCAVENGING
272  ! BUGFIX tt for nested fields should be ttn
273  ! sec may 2008
274          if (ngrid.gt.0) then
275             act_temp=ttn(ix,jy,hz,n,ngrid)
276          else
277             act_temp=tt(ix,jy,hz,n)
278          endif
279          cl=2E-7*prec**0.36
280          if (dquer(ks).gt.0) then ! is particle
281            S_i=0.9/cl
282           else ! is gas
283            cle=(1-cl)/(henry(ks)*(r_air/3500.)*act_temp)+cl
284            S_i=1/cle
285           endif
286           wetscav=S_i*prec/3.6E6/clouds_h
287  !         write(*,*) 'in. wetscav:'
288  !    +          ,wetscav,cle,cl,act_temp,prec,clouds_h
289        endif
290
291
292  !      if (wetscav.le.0) write (*,*) 'neg, or 0 wetscav!'
293  !    +          ,wetscav,cle,cl,act_temp,prec,clouds_h,clouds_v
294        wetdeposit(ks)=xmass1(jpart,ks)* &
295             (1.-exp(-wetscav*abs(ltsample)))*grfraction  ! wet deposition
296  !      new particle mass:
297  !      if (wetdeposit(ks).gt.0) then
298  !         write(*,*) 'wetdepo: ',wetdeposit(ks),ks
299  !      endif
300        restmass = xmass1(jpart,ks)-wetdeposit(ks)
301         if (ioutputforeachrelease.eq.1) then
302           kp=npoint(jpart)
303         else
304           kp=1
305         endif
306        if (restmass .gt. smallnum) then
307          xmass1(jpart,ks)=restmass
308  !ccccccccccccccc depostatistic
309  !       wetdepo_sum(ks,kp)=wetdepo_sum(ks,kp)+wetdeposit(ks)
310  !ccccccccccccccc depostatistic
311        else
312          xmass1(jpart,ks)=0.
313        endif
314  ! Correct deposited mass to the last time step when radioactive decay of
315  ! gridded deposited mass was calculated
316       if (decay(ks).gt.0.) then
317          wetdeposit(ks)=wetdeposit(ks) &
318               *exp(abs(ldeltat)*decay(ks))
319       endif
320      else  ! weta(k)
321         wetdeposit(ks)=0.
322      endif ! weta(k)
323    end do
324
325  ! Sabine Eckhard, June 2008 create deposition runs only for forward runs
326  ! Add the wet deposition to accumulated amount on output grid and nested output grid
327  !*****************************************************************************
328
329
330    if (ldirect.eq.1) then
331
332     griclunc = mod(nclass(i)-1, nclassunc) + 1
333
334!$OMP CRITICAL
335      call wetdepokernel(griclunc,wetdeposit,real(xtra1(jpart)), &
336        real(ytra1(jpart)),nage,kp)
337      if (nested_output.eq.1) then
338         call wetdepokernel_nest(griclunc, &
339           wetdeposit,real(xtra1(jpart)),real(ytra1(jpart)), &
340           nage,kp)
341      end if
342!$OMP END CRITICAL
343    end if
344
345  end do particleloop
346
347!$OMP END DO
348!$OMP END PARALLEL
349
350end subroutine wetdepo
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG