source: trunk/LMDZ.GENERIC/libf/phystd/convadj.F @ 1243

Last change on this file since 1243 was 787, checked in by aslmd, 12 years ago

LMDZ.GENERIC. (Sorry for long text but this is a quite major commit)

Paving the path for parallel computations. And moving towards a more flexible code.

Automatic allocation is used within all routines in phystd. No further mention to ngridmx and nqmx.

  1. ngridmx and nqmx are still used in LMDZ.GENERIC in the dyn3d part
  2. if the LMDZ4/LMDZ5 dynamical core is used, there is no more fixed dimensions ngridmx and nqmx --> a fully flexible parallel implementation is now possible (e.g. no need to recompile when changing numbers of processors)

The important stuff :

  • Compilation checked with ifort. OK with and without debug mode. No errors. Checked for: gcm, newstart, rcm1d, kcm1d
  • RUN GCM: Running an Earth test case. Comparison with previous revision --> debug mode : perfect match. bit by bit (diff command). checked with plots --> O1 mode : close match (checked with plots) --> O2 mode : sometimes up to 0.5 K departure.... BUT in this new version O2 and O1 are quite close while in previous version O1 and O2 differed by about, well, typically 0.5 K (pictures available on request)
  • RUN NEWSTART : perfect match (bit-by-bit) in either debug or normal mode.
  • RUN RCM1D : perfect match in normal mode.
  • RUN KCM1D : not tested (I don't know what is the use of kcm1d)

List of main changes :

  • Additional arguments to some subroutines (ngrid and nq)
  • F77 include strategy is obsolete and replaced by F90 module strategy In this new strategy arrays are allocatable and allocated once at first use This has to be done for all common featuring arrays defined with ngridmx or nqmx

surfdat.h >> surfdat_h.F90
tracer.h >> tracer_h.F90
comsaison.h >> comsaison_h.F90
comgeomfi.h >> comgeomfi_h.F90
comsoil.h >> comsoil_h.F90
comdiurn.h >> comdiurn_h.F90
fisice.h >> DELETED. was not used. probably a fossil.
watercap.h >> DELETED. variable put in surfdat_h.F90

  • F77 'save' strategy is obsolete and replaced by F90 'allocatable save' strategy (see previous point and e.g. new version of physiq.F90)
  • Suppressing any mention to advtrac.h which is a common in the dynamics and needs nqmx This was easily solved by adding an argument with tracer names, coming from the dynamics This is probably not a definitive solution, ... but this allows for generic physics to work easily with either LMDZ.GENERIC or LMDZ dynamical cores
  • Removing consistency tests between nq and nqmx ; and ngrid and ngridmx. No use now!
  • Adaptation of rcm1d, kcm1d, newstart given above-mentioned changes

A note on phyetat0 and soil_setting:

  • Now written so that a slice of horizontal size 'ngrid' starting at grid point 'cursor' is read in startfi.nc 'cursor' is defined in dimphys.h and initialized by inifis (or in newstart) this is useful for parallel computations. default behavior is the usual one : sequential runs, cursor is 1, size ngrid is the whole global domain

A note on an additional change :

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