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