source: flexpart.git/src/wetdepo.f90 @ 02095e3

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 02095e3 was db375cc, checked in by Espen Sollum ATMOS <eso@…>, 7 years ago

Removed some unused variables

  • Property mode set to 100644
File size: 7.6 KB
RevLine 
[e200b7a]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)
[db712a8]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! ix,jy              indices of output grid cell for each particle           *
45! itime [s]          actual simulation time [s]                              *
46! jpart              particle index                                          *
47! ldeltat [s]        interval since radioactive decay was computed           *
48! loutnext [s]       time for which gridded deposition is next output        *
49! loutstep [s]       interval at which gridded deposition is output          *
50! ltsample [s]       interval over which mass is deposited                   *
51! wetdeposit         mass that is wet deposited                              *
52! wetgrid            accumulated deposited mass on output grid               *
53! wetscav            scavenging coefficient                                  *
54!                                                                            *
55! Constants:                                                                 *
56!                                                                            *
57!*****************************************************************************
[e200b7a]58
59  use point_mod
60  use par_mod
61  use com_mod
62
63  implicit none
64
[db375cc]65  integer :: jpart,itime,ltsample,loutnext,ldeltat
66  integer :: itage,nage
[4fbe7a5]67  integer :: ks, kp
[6985a98]68  integer(selected_int_kind(16)), dimension(nspec) :: blc_count, inc_count
[c9cf570]69  real :: grfraction(3),wetscav
[e200b7a]70  real :: wetdeposit(maxspec),restmass
71  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
[8a65cb0]72
[db712a8]73! Compute interval since radioactive decay of deposited mass was computed
74!************************************************************************
[e200b7a]75
76  if (itime.le.loutnext) then
77    ldeltat=itime-(loutnext-loutstep)
78  else                                  ! first half of next interval
79    ldeltat=itime-loutnext
80  endif
81
[db712a8]82! Loop over all particles
83!************************
[e200b7a]84
[c8fc724]85  blc_count(:)=0
86  inc_count(:)=0
[8a65cb0]87
[e200b7a]88  do jpart=1,numpart
[4fbe7a5]89
[e200b7a]90    if (itra1(jpart).eq.-999999999) goto 20
91    if(ldirect.eq.1)then
92      if (itra1(jpart).gt.itime) goto 20
93    else
94      if (itra1(jpart).lt.itime) goto 20
95    endif
[8a65cb0]96
[c9cf570]97! Determine age class of the particle - nage is used for the kernel
[341f4b7]98!******************************************************************
[c9cf570]99     itage=abs(itra1(jpart)-itramem(jpart))
100     do nage=1,nageclass
101       if (itage.lt.lage(nage)) goto 33
102     end do
103 33  continue
[d6a0977]104
[92a74b2]105    do ks=1,nspec      ! loop over species
[8a65cb0]106
[6985a98]107      if (WETDEPSPEC(ks).eqv..false.) cycle
[db712a8]108
109!**************************************************
110! CALCULATE DEPOSITION
111!**************************************************
[0539b8f]112!       wetscav=0.
[0e085a8]113       
[7919911]114!        write(*,*) ks,dquer(ks), crain_aero(ks),csnow_aero(ks)
[0e085a8]115!       if (((dquer(ks).le.0.).and.(weta_gas(ks).gt.0..or.wetb_gas(ks).gt.0.)) &
116!          .or. &
117!          ((dquer(ks).gt.0.).and.(crain_aero(ks).gt.0..or.csnow_aero(ks).gt.0.).or. &
118!            (ccn_aero(ks).gt0) .or. (in_aero(ks).gt.0) .or. (henry(ks).gt.0)))  then
[585a533]119
[92a74b2]120      call get_wetscav(itime,ltsample,loutnext,jpart,ks,grfraction,inc_count,blc_count,wetscav)
[0e085a8]121     
[e200b7a]122
[4fbe7a5]123      if (wetscav.gt.0.) then
[e200b7a]124        wetdeposit(ks)=xmass1(jpart,ks)* &
[5f9d14a]125             (1.-exp(-wetscav*abs(ltsample)))*grfraction(1)  ! wet deposition
[4fbe7a5]126      else ! if no scavenging
127        wetdeposit(ks)=0.
128      endif
[c9cf570]129 
[4fbe7a5]130      restmass = xmass1(jpart,ks)-wetdeposit(ks)
131      if (ioutputforeachrelease.eq.1) then
132        kp=npoint(jpart)
133      else
134        kp=1
135      endif
136      if (restmass .gt. smallnum) then
137        xmass1(jpart,ks)=restmass
[db712a8]138!   depostatistic
139!   wetdepo_sum(ks,kp)=wetdepo_sum(ks,kp)+wetdeposit(ks)
140!   depostatistic
[4fbe7a5]141      else
142        xmass1(jpart,ks)=0.
143      endif
[db712a8]144!   Correct deposited mass to the last time step when radioactive decay of
145!   gridded deposited mass was calculated
[4fbe7a5]146      if (decay(ks).gt.0.) then
147        wetdeposit(ks)=wetdeposit(ks)*exp(abs(ldeltat)*decay(ks))
148      endif
[f13406c]149
[7919911]150!    endif ! no deposition
[92a74b2]151    end do ! loop over species
[e200b7a]152
[db712a8]153! Sabine Eckhardt, June 2008 create deposition runs only for forward runs
154! Add the wet deposition to accumulated amount on output grid and nested output grid
155!*****************************************************************************
[e200b7a]156
[4fbe7a5]157    if (ldirect.eq.1) then
158      call wetdepokernel(nclass(jpart),wetdeposit,real(xtra1(jpart)), &
[5f9d14a]159           real(ytra1(jpart)),nage,kp)
[4fbe7a5]160      if (nested_output.eq.1) call wetdepokernel_nest(nclass(jpart), &
[5f9d14a]161           wetdeposit,real(xtra1(jpart)),real(ytra1(jpart)),nage,kp)
[4fbe7a5]162    endif
[e200b7a]163
16420  continue
[4fbe7a5]165  end do ! all particles
[e200b7a]166
[db712a8]167! count the total number of below-cloud and in-cloud occurences:
[c8fc724]168  tot_blc_count(1:nspec)=tot_blc_count(1:nspec)+blc_count(1:nspec)
169  tot_inc_count(1:nspec)=tot_inc_count(1:nspec)+inc_count(1:nspec)
[8a65cb0]170
[e200b7a]171end subroutine wetdepo
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG