source: trunk/LMDZ.MARS/libf/phymars/convadj.F @ 3026

Last change on this file since 3026 was 2823, checked in by emillour, 2 years ago

Mars GCM:
Remove the "tracer" (logical) flag as we now always run with at least
one tracer.
EM

File size: 11.2 KB
RevLine 
[38]1      subroutine convadj(ngrid,nlay,nq,ptimestep,
[161]2     &                   pplay,pplev,ppopsk,lmax_th,
[38]3     &                   pu,pv,ph,pq,
4     &                   pdufi,pdvfi,pdhfi,pdqfi,
5     &                   pduadj,pdvadj,pdhadj,
6     &                   pdqadj)
7
[1036]8      use tracer_mod, only: noms, ! tracer names
9     &                      igcm_h2o_vap ! index of water vapor tracer
[2311]10      use comcstfi_h, only: g
[38]11      implicit none
12
13!==================================================================
14!     
15!     Purpose
16!     -------
17!     Calculates dry convective adjustment. If one tracer is CO2,
18!     we take into account the molecular mass variation (e.g. when
19!     CO2 condenses) to trigger convection (F. Forget 01/2005)
20!
21!     Authors
22!     -------
23!     Original author unknown.
24!     Modif. 2005 by F. Forget.
25!     Modif. 2010 by R. Wordsworth
[161]26!     Modif. 2011 by A. Colaitis
27!
[38]28!==================================================================
29
30!     ------------
31!     Declarations
32!     ------------
33
[2311]34      include "callkeys.h"
[38]35
36
37!     Arguments
38!     ---------
39
[161]40      INTEGER,intent(in) :: ngrid,nlay,lmax_th(ngrid)
[38]41      REAL,intent(in) :: ptimestep
42      REAL,intent(in) :: ph(ngrid,nlay)
43      REAL,intent(in) :: pdhfi(ngrid,nlay)
44      REAL,intent(out) :: pdhadj(ngrid,nlay)
45      REAL,intent(in) :: pplay(ngrid,nlay),pplev(ngrid,nlay+1)
46      REAL,intent(in) :: ppopsk(ngrid,nlay)
47      REAL,intent(in) :: pu(ngrid,nlay),pdufi(ngrid,nlay)
48      REAL,intent(out) :: pduadj(ngrid,nlay)
49      REAL,intent(in) :: pv(ngrid,nlay),pdvfi(ngrid,nlay)
50      REAL,intent(out) :: pdvadj(ngrid,nlay)
51
52!     Tracers
53      integer,intent(in) :: nq
54      real,intent(in) :: pq(ngrid,nlay,nq), pdqfi(ngrid,nlay,nq)
55      real,intent(out) :: pdqadj(ngrid,nlay,nq)
56
57
58!     Local
59!     -----
60
61      INTEGER ig,i,l,l1,l2,jj
[1047]62      INTEGER jcnt, jadrs(ngrid)
[38]63
[1047]64      REAL sig(nlay+1),sdsig(nlay),dsig(nlay)
65      REAL zu(ngrid,nlay),zv(ngrid,nlay)
66      REAL zh(ngrid,nlay)
67      REAL zu2(ngrid,nlay),zv2(ngrid,nlay)
68      REAL zh2(ngrid,nlay), zhc(ngrid,nlay)
[38]69      REAL zhm,zsm,zdsm,zum,zvm,zalpha,zhmc
70
71!     Tracers
[2311]72      INTEGER iq
73      integer,save :: ico2
[1047]74      REAL zq(ngrid,nlay,nq), zq2(ngrid,nlay,nq)
[1036]75      REAL zqm(nq),zqco2m
[2311]76      real m_co2, m_noco2
77      real,save :: A, B
[38]78
79      real mtot1, mtot2 , mm1, mm2
[2311]80      integer l1ref, l2ref
81      LOGICAL vtest(ngrid),down
82      logical,save :: firstcall=.true.
[38]83
[2586]84
85!$OMP THREADPRIVATE(ico2,A,B,firstcall)
86
[38]87!     for conservation test
88      real masse,cadjncons
89
90!     --------------
91!     Initialisation
92!     --------------
93
[1779]94      ! AS: OK firstcall absolute
[38]95      IF (firstcall) THEN
[1047]96
[38]97        ico2=0
98!     Prepare Special treatment if one of the tracers is CO2 gas
[1036]99           do iq=1,nq
[38]100             if (noms(iq).eq."co2") then
101                ico2=iq
102                m_co2 = 44.01E-3  ! CO2 molecular mass (kg/mol)   
103                m_noco2 = 33.37E-3  ! Non condensible mol mass (kg/mol)   
104!               Compute A and B coefficient use to compute
105!               mean molecular mass Mair defined by
106!               1/Mair = q(ico2)/m_co2 + (1-q(ico2))/m_noco2
107!               1/Mair = A*q(ico2) + B
108                A =(1/m_co2 - 1/m_noco2)
109                B=1/m_noco2
110             end if
111           enddo
112        firstcall=.false.
113      ENDIF ! of IF (firstcall)
114
115      DO l=1,nlay
116         DO ig=1,ngrid
117            zh(ig,l)=ph(ig,l)+pdhfi(ig,l)*ptimestep
118            zu(ig,l)=pu(ig,l)+pdufi(ig,l)*ptimestep
119            zv(ig,l)=pv(ig,l)+pdvfi(ig,l)*ptimestep
120         ENDDO
121      ENDDO
122
[2823]123      zq(1:ngrid,1:nlay,1:nq)=pq(1:ngrid,1:nlay,1:nq)+
124     &                        pdqfi(1:ngrid,1:nlay,1:nq)*ptimestep
[38]125
126      zh2(:,:)=zh(:,:)
127      zu2(:,:)=zu(:,:)
128      zv2(:,:)=zv(:,:)
129      zq2(:,:,:)=zq(:,:,:)
130
131!     -----------------------------
132!     Detection of unstable columns
133!     -----------------------------
134!     If ph(above) < ph(below) we set vtest=.true.
135
136      DO ig=1,ngrid
137        vtest(ig)=.false.
138      ENDDO
139
140      if (ico2.ne.0) then
141!     Special case if one of the tracers is CO2 gas
142         DO l=1,nlay
143           DO ig=1,ngrid
144             zhc(ig,l) = zh2(ig,l)*(A*zq2(ig,l,ico2)+B)
145           ENDDO
146         ENDDO
147       else
148          zhc(:,:)=zh2(:,:)
149       end if
150
151!     Find out which grid points are convectively unstable
152      DO l=2,nlay
153        DO ig=1,ngrid
[161]154          IF((zhc(ig,l).LT.zhc(ig,l-1))                                 $
155     $  .AND. (l .GT. lmax_th(ig))) vtest(ig)=.true.
156       ENDDO
[38]157      ENDDO
158
159!     Make a list of them
160      jcnt=0
161      DO ig=1,ngrid
162         IF(vtest(ig)) THEN
163            jcnt=jcnt+1
164            jadrs(jcnt)=ig
165         ENDIF
166      ENDDO
167
168
169!     ---------------------------------------------------------------
170!     Adjustment of the "jcnt" unstable profiles indicated by "jadrs"
171!     ---------------------------------------------------------------
172
173      DO jj = 1, jcnt   ! loop on every convective grid point
174
175          i = jadrs(jj)
176 
177!     Calculate sigma in this column
178          DO l=1,nlay+1
179            sig(l)=pplev(i,l)/pplev(i,1)
180       
181          ENDDO
182         DO l=1,nlay
183            dsig(l)=sig(l)-sig(l+1)
184            sdsig(l)=ppopsk(i,l)*dsig(l)
185         ENDDO
186          l2 = 1
187
188!     Test loop upwards on l2
189
[161]190      DO
191       l2 = l2 + 1
192       IF (l2 .GT. nlay) EXIT
193       IF ((zhc(i, l2) .LT. zhc(i, l2-1)).AND.(l2 .GT. lmax_th(i))) THEN
194
[38]195!     l2 is the highest level of the unstable column
196 
197              l1 = l2 - 1
198              l  = l1
199              zsm = sdsig(l2)
200              zdsm = dsig(l2)
201              zhm = zh2(i, l2)
202              if(ico2.ne.0) zqco2m = zq2(i,l2,ico2)
203
204!     Test loop downwards
205
206              DO
207                zsm = zsm + sdsig(l)
208                zdsm = zdsm + dsig(l)
209                zhm = zhm + sdsig(l) * (zh2(i, l) - zhm) / zsm
210                if(ico2.ne.0) then
211                  zqco2m =
212     &            zqco2m + dsig(l) * (zq2(i,l,ico2) - zqco2m) / zdsm
213                  zhmc = zhm*(A*zqco2m+B)
214                else
215                  zhmc = zhm
216                end if
217 
218!     do we have to extend the column downwards?
219 
[161]220            down = .false.
221            IF (l1 .ne. 1) then    !--  and then
222              IF ((zhmc .lt. zhc(i, l1-1)).and.(l1.gt.lmax_th(i))) then
223                down = .true.
224              END IF
225            END IF
[38]226 
227                ! this could be a problem...
228
229                if (down) then
230 
231                  l1 = l1 - 1
232                  l  = l1
233 
234                else
235 
236!     can we extend the column upwards?
237 
238                  if (l2 .eq. nlay) exit
239 
240                  if (zhc(i, l2+1) .ge. zhmc) exit
241
242                  l2 = l2 + 1
243                  l  = l2
244
245                end if ! of if (down)
246
247              enddo ! of do ! Test loop forward
248
249!     New constant profile (average value)
250
251
252              zalpha=0.
253              zum=0.
254              zvm=0.
255              do iq=1,nq
256                zqm(iq) = 0.
257              end do
258              DO l = l1, l2
259                if(ico2.ne.0) then
260                  zalpha=zalpha+
261     &            ABS(zhc(i,l)/(A+B*zqco2m) -zhm)*dsig(l)
262                else
263                  zalpha=zalpha+ABS(zh2(i,l)-zhm)*dsig(l)
264                endif
265                zh2(i, l) = zhm
266!     modifs by RDW !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
267                zum=zum+dsig(l)*zu2(i,l)
268                zvm=zvm+dsig(l)*zv2(i,l)
269!                zum=zum+dsig(l)*zu(i,l)
270!                zvm=zvm+dsig(l)*zv(i,l)
271                do iq=1,nq
272                   zqm(iq) = zqm(iq)+dsig(l)*zq2(i,l,iq)
273!                   zqm(iq) = zqm(iq)+dsig(l)*zq(i,l,iq)
274!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
275
276!     to conserve tracers/ KE, we must calculate zum, zvm and zqm using
277!     the up-to-date column values. If we do not do this, there are cases
278!     where convection stops at one level and starts at the next where we
279!     can break conservation of stuff (particularly tracers) significantly.
280
281                end do
282              ENDDO ! of DO l = l1, l2
283              zalpha=zalpha/(zhm*(sig(l1)-sig(l2+1)))
284              zum=zum/(sig(l1)-sig(l2+1))
285              zvm=zvm/(sig(l1)-sig(l2+1))
286              do iq=1,nq
287                 zqm(iq) = zqm(iq)/(sig(l1)-sig(l2+1))
288              end do
289
290              IF(zalpha.GT.1.) THEN
291                 zalpha=1.
292              ELSE
293!                IF(zalpha.LT.0.) STOP
294                 IF(zalpha.LT.1.e-4) zalpha=1.e-4
295              ENDIF
296
297              DO l=l1,l2
298                 zu2(i,l)=zu2(i,l)+zalpha*(zum-zu2(i,l))
299                 zv2(i,l)=zv2(i,l)+zalpha*(zvm-zv2(i,l))
300                 do iq=1,nq
301!                  zq2(i,l,iq)=zq2(i,l,iq)+zalpha*(zqm(iq)-zq2(i,l,iq))
302                   zq2(i,l,iq)=zqm(iq)
303                 end do
304              ENDDO
305              if (ico2.ne.0) then
306                DO l=l1, l2
307                  zhc(i,l) = zh2(i,l)*(A*zq2(i,l,ico2)+B)
308                ENDDO
309              end if
310
311
312              l2 = l2 + 1
313
314            END IF   ! End of l1 to l2 instability treatment
315                     ! We now continue to test from l2 upwards
316
317          ENDDO   ! End of upwards loop on l2
318
319
320!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
321!     check conservation
322         cadjncons=0.0
323         if(water)then
324          iq    = igcm_h2o_vap
325          do l = 1, nlay
326            masse = (pplev(i,l) - pplev(i,l+1))/g
327            cadjncons = cadjncons +
328     &           masse*(zq2(i,l,iq)-zq(i,l,iq))/ptimestep
329          end do
330         endif
331
332         if(cadjncons.lt.-1.e-6)then
333            print*,'convadj has just crashed...'
334            print*,'i  = ',i
335            print*,'l1 = ',l1
336            print*,'l2 = ',l2
337            print*,'cadjncons        = ',cadjncons
338         do l = 1, nlay
339            print*,'dsig         = ',dsig(l)
340         end do         
341         do l = 1, nlay
342            print*,'dsig         = ',dsig(l)
343         end do
344         do l = 1, nlay
345            print*,'sig         = ',sig(l)
346         end do
347         do l = 1, nlay
348            print*,'pplay(ig,:)         = ',pplay(i,l)
349         end do
350         do l = 1, nlay+1
351            print*,'pplev(ig,:)         = ',pplev(i,l)
352         end do
353         do l = 1, nlay
354            print*,'ph(ig,:)         = ',ph(i,l)
355         end do
356         do l = 1, nlay
357            print*,'zh(ig,:)         = ',zh(i,l)
358         end do
359         do l = 1, nlay
360            print*,'zh2(ig,:)        = ',zh2(i,l)
361         end do
362         do l = 1, nlay
363            print*,'zq(ig,:,vap)     = ',zq(i,l,igcm_h2o_vap)
364         end do
365         do l = 1, nlay
366            print*,'zq2(ig,:,vap)    = ',zq2(i,l,igcm_h2o_vap)
367         end do
368            print*,'zqm(vap)         = ',zqm(igcm_h2o_vap)
369            print*,'jadrs=',jadrs
370
[2311]371            call abort_physic("convadj","crashed",1)
[38]372         endif
373!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
374
375      ENDDO ! of DO jj = 1, jcnt   ! loop on every convective grid point
376
377
378      DO l=1,nlay
379        DO ig=1,ngrid
380          pdhadj(ig,l)=(zh2(ig,l)-zh(ig,l))/ptimestep
381          pduadj(ig,l)=(zu2(ig,l)-zu(ig,l))/ptimestep
382          pdvadj(ig,l)=(zv2(ig,l)-zv(ig,l))/ptimestep
383        ENDDO
384      ENDDO
385
[2823]386      pdqadj(1:ngrid,1:nlay,1:nq)=(zq2(1:ngrid,1:nlay,1:nq)
387     &                             -zq(1:ngrid,1:nlay,1:nq))/ptimestep
[38]388
389      end
390     
Note: See TracBrowser for help on using the repository browser.