source: flexpart.git/src/convmix.f90 @ 54cbd6c

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 54cbd6c was 8a65cb0, checked in by Espen Sollum ATMOS <espen@…>, 9 years ago

Added code, makefile for dev branch

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