source: flexpart.git/src/redist.f90 @ 3481cc1

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

move license from headers to a different file

  • Property mode set to 100644
File size: 6.7 KB
Line 
1subroutine redist (ipart,ktop,ipconv)
2
3  !**************************************************************************
4  ! Do the redistribution of particles due to convection
5  ! This subroutine is called for each particle which is assigned
6  ! a new vertical position randomly, based on the convective redistribution
7  ! matrix
8  !**************************************************************************
9
10  ! Petra Seibert, Feb 2001, Apr 2001, May 2001, Jan 2002, Nov 2002 and
11  ! Andreas Frank, Nov 2002
12
13  ! Caroline Forster:  November 2004 - February 2005
14
15  use par_mod
16  use com_mod
17  use conv_mod
18  use random_mod
19
20  implicit none
21
22  real,parameter :: const=r_air/ga
23  integer :: ipart, ktop,ipconv
24  integer :: k, kz, levnew, levold
25  real,save :: uvzlev(nuvzmax)
26  real :: wsub(nuvzmax)
27  real :: totlevmass, wsubpart
28  real :: temp_levold,temp_levold1
29  real :: sub_levold,sub_levold1
30  real :: pint, pold, rn, tv, tvold, dlevfrac
31  real :: ew,ztold,ffraction
32  real :: tv1, tv2, dlogp, dz, dz1, dz2
33  integer :: iseed = -88
34
35
36
37  ! ipart   ... number of particle to be treated
38
39  ipconv=1
40
41  ! determine height of the eta half-levels (uvzlev)
42  ! do that only once for each grid column
43  ! i.e. when ktop.eq.1
44  !**************************************************************
45
46  if (ktop .le. 1) then
47
48    tvold=tt2conv*(1.+0.378*ew(td2conv)/psconv)
49    pold=psconv
50    uvzlev(1)=0.
51
52    pint = phconv(2)
53  !  determine next virtual temperatures
54    tv1 = tconv(1)*(1.+0.608*qconv(1))
55    tv2 = tconv(2)*(1.+0.608*qconv(2))
56  !  interpolate virtual temperature to half-level
57    tv = tv1 + (tv2-tv1)*(pconv(1)-phconv(2))/(pconv(1)-pconv(2))
58    if (abs(tv-tvold).gt.0.2) then
59      uvzlev(2) = uvzlev(1) + &
60           const*log(pold/pint)* &
61           (tv-tvold)/log(tv/tvold)
62    else
63      uvzlev(2) = uvzlev(1)+ &
64           const*log(pold/pint)*tv
65    endif
66    tvold=tv
67    tv1=tv2
68    pold=pint
69
70  ! integrate profile (calculation of height agl of eta layers) as required
71    do kz = 3, nconvtop+1
72  !    note that variables defined in calcmatrix.f (pconv,tconv,qconv)
73  !    start at the first real ECMWF model level whereas kz and
74  !    thus uvzlev(kz) starts at the surface. uvzlev is defined at the
75  !    half-levels (between the tconv, qconv etc. values !)
76  !    Thus, uvzlev(kz) is the lower boundary of the tconv(kz) cell.
77      pint = phconv(kz)
78  !    determine next virtual temperatures
79      tv2 = tconv(kz)*(1.+0.608*qconv(kz))
80  !    interpolate virtual temperature to half-level
81      tv = tv1 + (tv2-tv1)*(pconv(kz-1)-phconv(kz))/ &
82           (pconv(kz-1)-pconv(kz))
83      if (abs(tv-tvold).gt.0.2) then
84        uvzlev(kz) = uvzlev(kz-1) + &
85             const*log(pold/pint)* &
86             (tv-tvold)/log(tv/tvold)
87      else
88        uvzlev(kz) = uvzlev(kz-1)+ &
89             const*log(pold/pint)*tv
90      endif
91      tvold=tv
92      tv1=tv2
93      pold=pint
94
95    end do
96
97    ktop = 2
98
99  endif ! (if ktop .le. 1) then
100
101  !  determine vertical grid position of particle in the eta system
102  !****************************************************************
103
104  ztold = ztra1(abs(ipart))
105  ! find old particle grid position
106  do kz = 2, nconvtop
107    if (uvzlev(kz) .ge. ztold ) then
108      levold = kz-1
109      goto 30
110    endif
111  end do
112
113  ! Particle is above the potentially convective domain. Skip it.
114  goto 90
115
11630   continue
117
118  ! now redistribute particles
119  !****************************
120
121  !  Choose a random number and find corresponding level of destination
122  !  Random numbers to be evenly distributed in [0,1]
123
124  rn = ran3(iseed)
125
126  ! initialize levnew
127
128  levnew = levold
129
130  ffraction = 0.
131  totlevmass=dpr(levold)/ga
132  do k = 1,nconvtop
133  ! for backward runs use the transposed matrix
134   if (ldirect.eq.1) then
135     ffraction=ffraction+fmassfrac(levold,k) &
136          /totlevmass
137   else
138     ffraction=ffraction+fmassfrac(k,levold) &
139          /totlevmass
140   endif
141   if (rn.le.ffraction) then
142     levnew=k
143  ! avoid division by zero or a too small number
144  ! if division by zero or a too small number happens the
145  ! particle is assigned to the center of the grid cell
146     if (ffraction.gt.1.e-20) then
147      if (ldirect.eq.1) then
148        dlevfrac = (ffraction-rn) / fmassfrac(levold,k) * totlevmass
149      else
150        dlevfrac = (ffraction-rn) / fmassfrac(k,levold) * totlevmass
151      endif
152     else
153       dlevfrac = 0.5
154     endif
155     goto 40
156   endif
157  end do
158
15940   continue
160
161  ! now assign new position to particle
162
163  if (levnew.le.nconvtop) then
164   if (levnew.eq.levold) then
165      ztra1(abs(ipart)) = ztold
166   else
167    dlogp = (1.-dlevfrac)* &
168         (log(phconv(levnew+1))-log(phconv(levnew)))
169    pint = log(phconv(levnew))+dlogp
170    dz1 = pint - log(phconv(levnew))
171    dz2 = log(phconv(levnew+1)) - pint
172    dz = dz1 + dz2
173    ztra1(abs(ipart)) = (uvzlev(levnew)*dz2+uvzlev(levnew+1)*dz1)/dz
174     if (ztra1(abs(ipart)).lt.0.) &
175          ztra1(abs(ipart))=-1.*ztra1(abs(ipart))
176     if (ipconv.gt.0) ipconv=-1
177   endif
178  endif
179
180  ! displace particle according to compensating subsidence
181  ! this is done to those particles, that were not redistributed
182  ! by the matrix
183  !**************************************************************
184
185  if (levnew.le.nconvtop.and.levnew.eq.levold) then
186
187  ztold = ztra1(abs(ipart))
188
189  ! determine compensating vertical velocity at the levels
190  ! above and below the particel position
191  ! increase compensating subsidence by the fraction that
192  ! is displaced by convection to this level
193
194    if (levold.gt.1) then
195     temp_levold = tconv(levold-1) + &
196          (tconv(levold)-tconv(levold-1)) &
197          *(pconv(levold-1)-phconv(levold))/ &
198          (pconv(levold-1)-pconv(levold))
199     sub_levold = sub(levold)/(1.-sub(levold)/dpr(levold)*ga)
200     wsub(levold)=-1.*sub_levold*r_air*temp_levold/(phconv(levold))
201    else
202     wsub(levold)=0.
203    endif
204
205     temp_levold1 = tconv(levold) + &
206          (tconv(levold+1)-tconv(levold)) &
207          *(pconv(levold)-phconv(levold+1))/ &
208          (pconv(levold)-pconv(levold+1))
209     sub_levold1 = sub(levold+1)/(1.-sub(levold+1)/dpr(levold+1)*ga)
210     wsub(levold+1)=-1.*sub_levold1*r_air*temp_levold1/ &
211          (phconv(levold+1))
212
213  ! interpolate wsub to the vertical particle position
214
215  dz1 = ztold - uvzlev(levold)
216  dz2 = uvzlev(levold+1) - ztold
217  dz = dz1 + dz2
218
219  wsubpart = (dz2*wsub(levold)+dz1*wsub(levold+1))/dz
220  ztra1(abs(ipart)) = ztold+wsubpart*real(lsynctime)
221  if (ztra1(abs(ipart)).lt.0.) then
222     ztra1(abs(ipart))=-1.*ztra1(abs(ipart))
223  endif
224
225  endif      !(levnew.le.nconvtop.and.levnew.eq.levold)
226
227  ! Maximum altitude .5 meter below uppermost model level
228  !*******************************************************
229
230 90   continue
231
232  if (ztra1(abs(ipart)) .gt. height(nz)-0.5) &
233       ztra1(abs(ipart)) = height(nz)-0.5
234
235end subroutine redist
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG