source: flexpart.git/src/convmix.f90 @ c0884a8

univie
Last change on this file since c0884a8 was c0884a8, checked in by pesei <petra seibert at univie ac at>, 6 years ago

replace CTBTO code for checking type of GRIB

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