source: flexpart.git/src/convmix.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: 9.8 KB
Line 
1subroutine convmix(itime,metdata_format)
2  !                     i
3  !**************************************************************
4  !handles all the calculations related to convective mixing
5  !Petra Seibert, Bernd C. Krueger, Feb 2001
6  !nested grids included, Bernd C. Krueger, May 2001
7  !
8  !Changes by Caroline Forster, April 2004 - February 2005:
9  !  convmix called every lsynctime seconds
10  !CHANGES by A. Stohl:
11  !  various run-time optimizations - February 2005
12  ! CHANGES by C. Forster, November 2005, NCEP GFS version
13  !      in the ECMWF version convection is calculated on the
14  !      original eta-levels
15  !      in the GFS version convection is calculated on the
16  !      FLEXPART levels
17  !
18  !   Unified ECMWF and GFS builds                                             
19  !   Marian Harustak, 12.5.2017                                             
20  !     - Merged convmix and convmix_gfs into one routine using if-then           
21  !       for meteo-type dependent code                                       
22  !**************************************************************
23
24  use flux_mod
25  use par_mod
26  use com_mod
27  use conv_mod
28  use class_gribfile
29
30  implicit none
31
32  integer :: igr,igrold, ipart, itime, ix, j, inest
33  integer :: ipconv
34  integer :: jy, kpart, ktop, ngrid,kz
35  integer :: igrid(maxpart), ipoint(maxpart), igridn(maxpart,maxnests)
36  integer :: metdata_format
37
38  ! itime [s]                 current time
39  ! igrid(maxpart)            horizontal grid position of each particle
40  ! igridn(maxpart,maxnests)  dto. for nested grids
41  ! ipoint(maxpart)           pointer to access particles according to grid position
42
43  logical :: lconv
44  real :: x, y, xtn,ytn, ztold, delt
45  real :: dt1,dt2,dtt
46  integer :: mind1,mind2
47  ! dt1,dt2,dtt,mind1,mind2       variables used for time interpolation
48  integer :: itage,nage
49  real,parameter :: eps=nxmax/3.e5
50
51
52  !monitoring variables
53  !real sumconv,sumall
54
55
56  ! Calculate auxiliary variables for time interpolation
57  !*****************************************************
58
59  dt1=real(itime-memtime(1))
60  dt2=real(memtime(2)-itime)
61  dtt=1./(dt1+dt2)
62  mind1=memind(1)
63  mind2=memind(2)
64  delt=real(abs(lsynctime))
65
66
67  lconv = .false.
68
69  ! if no particles are present return after initialization
70  !********************************************************
71
72  if (numpart.le.0) return
73
74  ! Assign igrid and igridn, which are pseudo grid numbers indicating particles
75  ! that are outside the part of the grid under consideration
76  ! (e.g. particles near the poles or particles in other nests).
77  ! Do this for all nests but use only the innermost nest; for all others
78  ! igrid shall be -1
79  ! Also, initialize index vector ipoint
80  !************************************************************************
81
82  do ipart=1,numpart
83    igrid(ipart)=-1
84    do j=numbnests,1,-1
85      igridn(ipart,j)=-1
86    end do
87    ipoint(ipart)=ipart
88  ! do not consider particles that are (yet) not part of simulation
89    if (itra1(ipart).ne.itime) goto 20
90    x = xtra1(ipart)
91    y = ytra1(ipart)
92
93  ! Determine which nesting level to be used
94  !**********************************************************
95
96    ngrid=0
97    if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then
98    do j=numbnests,1,-1
99      if ( x.gt.xln(j)+eps .and. x.lt.xrn(j)-eps .and. &
100           y.gt.yln(j)+eps .and. y.lt.yrn(j)-eps ) then
101        ngrid=j
102        goto 23
103      endif
104    end do
105    else
106      do j=numbnests,1,-1
107        if ( x.gt.xln(j) .and. x.lt.xrn(j) .and. &
108             y.gt.yln(j) .and. y.lt.yrn(j) ) then
109          ngrid=j
110          goto 23
111        endif
112      end do
113    endif
114 23   continue
115
116  ! Determine nested grid coordinates
117  !**********************************
118
119    if (ngrid.gt.0) then
120  ! nested grids
121      xtn=(x-xln(ngrid))*xresoln(ngrid)
122      ytn=(y-yln(ngrid))*yresoln(ngrid)
123      ix=nint(xtn)
124      jy=nint(ytn)
125      igridn(ipart,ngrid) = 1 + jy*nxn(ngrid) + ix
126    else if(ngrid.eq.0) then
127  ! mother grid
128      ix=nint(x)
129      jy=nint(y)
130      igrid(ipart) = 1 + jy*nx + ix
131    endif
132
133 20 continue
134  end do
135
136  ! sumall = 0.
137  ! sumconv = 0.
138
139  !*****************************************************************************
140  ! 1. Now, do everything for the mother domain and, later, for all of the nested domains
141  ! While all particles have to be considered for redistribution, the Emanuel convection
142  ! scheme only needs to be called once for every grid column where particles are present.
143  ! Therefore, particles are sorted according to their grid position. Whenever a new grid
144  ! cell is encountered by looping through the sorted particles, the convection scheme is called.
145  !*****************************************************************************
146
147  ! sort particles according to horizontal position and calculate index vector IPOINT
148
149  call sort2(numpart,igrid,ipoint)
150
151  ! Now visit all grid columns where particles are present
152  ! by going through the sorted particles
153
154  igrold = -1
155  do kpart=1,numpart
156    igr = igrid(kpart)
157    if (igr .eq. -1) goto 50
158    ipart = ipoint(kpart)
159  ! sumall = sumall + 1
160    if (igr .ne. igrold) then
161  ! we are in a new grid column
162      jy = (igr-1)/nx
163      ix = igr - jy*nx - 1
164
165  ! Interpolate all meteorological data needed for the convection scheme
166      psconv=(ps(ix,jy,1,mind1)*dt2+ps(ix,jy,1,mind2)*dt1)*dtt
167      tt2conv=(tt2(ix,jy,1,mind1)*dt2+tt2(ix,jy,1,mind2)*dt1)*dtt
168      td2conv=(td2(ix,jy,1,mind1)*dt2+td2(ix,jy,1,mind2)*dt1)*dtt
169!!$      do kz=1,nconvlev+1      !old
170      if (metdata_format.eq.GRIBFILE_CENTRE_ECMWF) then
171        do kz=1,nuvz-1           !bugfix
172        tconv(kz)=(tth(ix,jy,kz+1,mind1)*dt2+ &
173             tth(ix,jy,kz+1,mind2)*dt1)*dtt
174        qconv(kz)=(qvh(ix,jy,kz+1,mind1)*dt2+ &
175             qvh(ix,jy,kz+1,mind2)*dt1)*dtt
176      end do
177      else
178        do kz=1,nuvz-1           !bugfix
179          pconv(kz)=(pplev(ix,jy,kz,mind1)*dt2+ &
180              pplev(ix,jy,kz,mind2)*dt1)*dtt
181          tconv(kz)=(tt(ix,jy,kz,mind1)*dt2+ &
182              tt(ix,jy,kz,mind2)*dt1)*dtt
183          qconv(kz)=(qv(ix,jy,kz,mind1)*dt2+ &
184              qv(ix,jy,kz,mind2)*dt1)*dtt
185        end do
186      end if
187
188  ! Calculate translocation matrix
189      call calcmatrix(lconv,delt,cbaseflux(ix,jy),metdata_format)
190      igrold = igr
191      ktop = 0
192    endif
193
194  ! treat particle only if column has convection
195    if (lconv .eqv. .true.) then
196  ! assign new vertical position to particle
197
198      ztold=ztra1(ipart)
199      call redist(ipart,ktop,ipconv)
200  !    if (ipconv.le.0) sumconv = sumconv+1
201
202  ! Calculate the gross fluxes across layer interfaces
203  !***************************************************
204
205      if (iflux.eq.1) then
206        itage=abs(itra1(ipart)-itramem(ipart))
207        do nage=1,nageclass
208          if (itage.lt.lage(nage)) goto 37
209        end do
210 37     continue
211
212        if (nage.le.nageclass) &
213             call calcfluxes(nage,ipart,real(xtra1(ipart)), &
214             real(ytra1(ipart)),ztold)
215      endif
216
217    endif   !(lconv .eqv. .true)
218 50 continue
219  end do
220
221
222  !*****************************************************************************
223  ! 2. Nested domains
224  !*****************************************************************************
225
226  ! sort particles according to horizontal position and calculate index vector IPOINT
227
228  do inest=1,numbnests
229    do ipart=1,numpart
230      ipoint(ipart)=ipart
231      igrid(ipart) = igridn(ipart,inest)
232    enddo
233    call sort2(numpart,igrid,ipoint)
234
235  ! Now visit all grid columns where particles are present
236  ! by going through the sorted particles
237
238    igrold = -1
239    do kpart=1,numpart
240      igr = igrid(kpart)
241      if (igr .eq. -1) goto 60
242      ipart = ipoint(kpart)
243      ! sumall = sumall + 1
244      if (igr .ne. igrold) then
245  ! we are in a new grid column
246        jy = (igr-1)/nxn(inest)
247        ix = igr - jy*nxn(inest) - 1
248
249  ! Interpolate all meteorological data needed for the convection scheme
250        psconv=(psn(ix,jy,1,mind1,inest)*dt2+ &
251             psn(ix,jy,1,mind2,inest)*dt1)*dtt
252        tt2conv=(tt2n(ix,jy,1,mind1,inest)*dt2+ &
253             tt2n(ix,jy,1,mind2,inest)*dt1)*dtt
254        td2conv=(td2n(ix,jy,1,mind1,inest)*dt2+ &
255             td2n(ix,jy,1,mind2,inest)*dt1)*dtt
256!!$        do kz=1,nconvlev+1    !old
257        do kz=1,nuvz-1           !bugfix
258          tconv(kz)=(tthn(ix,jy,kz+1,mind1,inest)*dt2+ &
259               tthn(ix,jy,kz+1,mind2,inest)*dt1)*dtt
260          qconv(kz)=(qvhn(ix,jy,kz+1,mind1,inest)*dt2+ &
261               qvhn(ix,jy,kz+1,mind2,inest)*dt1)*dtt
262        end do
263
264  ! calculate translocation matrix
265  !*******************************
266        call calcmatrix(lconv,delt,cbasefluxn(ix,jy,inest),metdata_format)
267        igrold = igr
268        ktop = 0
269      endif
270
271  ! treat particle only if column has convection
272      if (lconv .eqv. .true.) then
273  ! assign new vertical position to particle
274        ztold=ztra1(ipart)
275        call redist(ipart,ktop,ipconv)
276  !      if (ipconv.le.0) sumconv = sumconv+1
277
278  ! Calculate the gross fluxes across layer interfaces
279  !***************************************************
280
281        if (iflux.eq.1) then
282          itage=abs(itra1(ipart)-itramem(ipart))
283          do nage=1,nageclass
284            if (itage.lt.lage(nage)) goto 47
285          end do
286 47       continue
287
288          if (nage.le.nageclass) &
289               call calcfluxes(nage,ipart,real(xtra1(ipart)), &
290               real(ytra1(ipart)),ztold)
291        endif
292
293      endif !(lconv .eqv. .true.)
294
295
29660    continue
297    end do
298  end do
299  !--------------------------------------------------------------------------
300  ! write(*,*)'############################################'
301  ! write(*,*)'TIME=',&
302  !    &  itime
303  ! write(*,*)'fraction of particles under convection',&
304  !    &  sumconv/(sumall+0.001)
305  ! write(*,*)'total number of particles',&
306  !    &  sumall
307  ! write(*,*)'number of particles under convection',&
308  !    &  sumconv
309  ! write(*,*)'############################################'
310
311  return
312end subroutine convmix
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG