source: trunk/LMDZ.MARS/libf/phymars/calltherm_interface_mod.F90

Last change on this file was 4024, checked in by aslmd, 6 days ago

Mars GCM: a couple adding _mod to the files' name because some compiling workflows are happier with this -- progress with MESOSCALE catch-up, only a couple of efforts remaining

File size: 19.1 KB
RevLine 
[3727]1MODULE calltherm_interface_mod
2
3IMPLICIT NONE
4
5CONTAINS
[1033]6!=======================================================================
[1032]7! CALLTHERM_INTERFACE
[1033]8!=======================================================================
[1032]9! Main interface to the Martian thermal plume model
10! This interface handles sub-timesteps for this model
11! A call to this interface must be inserted in the main 'physics' routine
12!   NB: for information:
13!   In the Mars LMD-GCM, the thermal plume model is called after
14!   the vertical turbulent mixing scheme (Mellor and Yamada) 
15!   and the surface layer scheme (Richardson-based surface layer + subgrid gustiness)
16!   Other routines called before the thermals model are
17!   radiative transfer and (orographic) gravity wave drag.
18! -----------------------------------------------------------------------
[1033]19! Author : A. Colaitis 2011-01-05 (with updates 2011-2013)
20!          after C. Rio and F. Hourdin
21! Institution : Laboratoire de Meteorologie Dynamique (LMD) Paris, France
22! -----------------------------------------------------------------------
23! Corresponding author : A. Spiga aymeric.spiga_AT_upmc.fr
24! -----------------------------------------------------------------------
[1032]25! ASSOCIATED FILES
[1033]26! --> thermcell_main_mars.F90
[1032]27! --> thermcell_dqup.F90
28! --> comtherm_h.F90
29! -----------------------------------------------------------------------
30! Reference paper:
31! A. Colaïtis, A. Spiga, F. Hourdin, C. Rio, F. Forget, and E. Millour.
32! A thermal plume model for the Martian convective boundary layer.
33! Journal of Geophysical Research (Planets), 118:1468-1487, July 2013.
[1033]34! http://dx.doi.org/10.1002/jgre.20104
35! http://arxiv.org/abs/1306.6215
[1032]36! -----------------------------------------------------------------------
37! Reference paper for terrestrial plume model:
38! C. Rio and F. Hourdin.
39! A thermal plume model for the convective boundary layer : Representation of cumulus clouds.
40! Journal of the Atmospheric Sciences, 65:407-425, 2008.
41! -----------------------------------------------------------------------
42
[2823]43      SUBROUTINE calltherm_interface (ngrid,nlayer,nq, igcm_co2, &
[652]44     & zzlev,zzlay, &
[161]45     & ptimestep,pu,pv,pt,pq,pdu,pdv,pdt,pdq,q2, &
[185]46     & pplay,pplev,pphi,zpopsk, &
[660]47     & pdu_th,pdv_th,pdt_th,pdq_th,lmax,zmaxth,pbl_dtke, &
48     & pdhdif,hfmax,wstar,sensibFlux)
[161]49
[1032]50      use comtherm_h
[1036]51      use tracer_mod, only: nqmx,noms
[1226]52
53      ! SHARED VARIABLES. This needs adaptations in another climate model.
54      ! contains physical constant values such as
55      ! "g" : gravitational acceleration (m.s-2)
56      ! "r" : recuced gas constant (J.K-1.mol-1)
57      ! "cpp" : specific heat of the atmosphere (J.kg-1.K-1)
[3727]58      USE comcstfi_h, only: r, g, cpp
[1226]59
[3727]60      use thermcell_main_mars_mod, only: thermcell_main_mars
61     
62      use thermcell_dqup_mod, only: thermcell_dqup
63
[1032]64      implicit none
[161]65
66!--------------------------------------------------------
[342]67! Input Variables
[161]68!--------------------------------------------------------
69
[1032]70      INTEGER, INTENT(IN) :: ngrid ! number of horizontal grid points
71      INTEGER, INTENT(IN) :: nlayer ! number of vertical grid points
72      INTEGER, INTENT(IN) :: nq ! number of tracer species
[1028]73      REAL, INTENT(IN) :: ptimestep !timestep (s)
[1032]74      REAL, INTENT(IN) :: pplev(ngrid,nlayer+1) !intermediate pressure levels (Pa)
75      REAL, INTENT(IN) :: pplay(ngrid,nlayer) !Pressure at the middle of the layers (Pa)
76      REAL, INTENT(IN) :: pphi(ngrid,nlayer) !Geopotential at the middle of the layers (m2s-2)
77      REAL, INTENT(IN) :: pu(ngrid,nlayer),pv(ngrid,nlayer) !u,v components of the wind (ms-1)
78      REAL, INTENT(IN) :: pt(ngrid,nlayer),pq(ngrid,nlayer,nq)!temperature (K) and tracer concentration (kg/kg)
79      REAL, INTENT(IN) :: zzlay(ngrid,nlayer) ! altitude at the middle of the layers
80      REAL, INTENT(IN) :: zzlev(ngrid,nlayer+1) ! altitude at layer boundaries
81      INTEGER, INTENT(IN) :: igcm_co2 ! index of the CO2 tracer in mixing ratio array
82                                      ! --> 0 if no tracer is CO2 (or no tracer at all)
83                                      ! --> this prepares special treatment for polar night mixing
84                                      !       (see thermcell_main_mars)
85      REAL, INTENT(IN) :: pdu(ngrid,nlayer),pdv(ngrid,nlayer) ! wind velocity change from routines called
[1028]86                                                                      ! before thermals du/dt (m/s/s)
[1032]87      REAL, INTENT(IN) :: pdq(ngrid,nlayer,nq) ! tracer concentration change from routines called
[1028]88                                                     ! before thermals dq/dt (kg/kg/s)
[1032]89      REAL, INTENT(IN) :: pdt(ngrid,nlayer) ! temperature change from routines called before thermals dT/dt (K/s)
90      REAL, INTENT(IN) :: q2(ngrid,nlayer+1) ! turbulent kinetic energy
91      REAL, INTENT(IN) :: zpopsk(ngrid,nlayer) ! ratio of pressure at middle of layer to surface pressure,
[1028]92                                                   ! to the power r/cp, i.e. zpopsk=(pplay(ig,l)/pplev(ig,1))**rcp
[1032]93      REAL, INTENT(IN) :: pdhdif(ngrid,nlayer) ! potential temperature change from turbulent diffusion scheme dT/dt (K/s)
94      REAL, INTENT(IN) :: sensibFlux(ngrid) ! sensible heat flux computed from surface layer scheme
[161]95
96!--------------------------------------------------------
[342]97! Output Variables
[161]98!--------------------------------------------------------
99
[1032]100      REAL, INTENT(OUT) :: pdu_th(ngrid,nlayer) ! wind velocity change from thermals du/dt (m/s/s)
101      REAL, INTENT(OUT) :: pdv_th(ngrid,nlayer) ! wind velocity change from thermals dv/dt (m/s/s)
102      REAL, INTENT(OUT) :: pdt_th(ngrid,nlayer) ! temperature change from thermals dT/dt (K/s)
103      REAL, INTENT(OUT) :: pdq_th(ngrid,nlayer,nq) ! tracer change from thermals dq/dt (kg/kg/s)
104      INTEGER, INTENT(OUT) :: lmax(ngrid) ! layer number reacher by thermals in grid point
105      REAL, INTENT(OUT) :: zmaxth(ngrid) ! equivalent to lmax, but in (m), interpolated
106      REAL, INTENT(OUT) :: pbl_dtke(ngrid,nlayer+1) ! turbulent kinetic energy change from thermals dtke/dt
107      REAL, INTENT(OUT) :: wstar(ngrid) ! free convection velocity (m/s)
[161]108
109!--------------------------------------------------------
[342]110! Thermals local variables
[161]111!--------------------------------------------------------
[1032]112      REAL zu(ngrid,nlayer), zv(ngrid,nlayer)
113      REAL zt(ngrid,nlayer)
114      REAL d_t_ajs(ngrid,nlayer)
115      REAL d_u_ajs(ngrid,nlayer), d_q_ajs(ngrid,nlayer,nq)
116      REAL d_v_ajs(ngrid,nlayer)
117      REAL fm_therm(ngrid,nlayer+1), entr_therm(ngrid,nlayer)
118      REAL detr_therm(ngrid,nlayer),detrmod(ngrid,nlayer)
119      REAL zw2(ngrid,nlayer+1)
120      REAL fraca(ngrid,nlayer+1),zfraca(ngrid,nlayer+1)
121      REAL q_therm(ngrid,nlayer), pq_therm(ngrid,nlayer,nq)
122      REAL q2_therm(ngrid,nlayer), dq2_therm(ngrid,nlayer)
123      REAL lmax_real(ngrid)
124      REAL masse(ngrid,nlayer)
[1028]125
[660]126      INTEGER l,ig,iq,ii(1),k
[628]127      CHARACTER (LEN=20) modname
[161]128
[342]129!--------------------------------------------------------
130! Local variables for sub-timestep
131!--------------------------------------------------------
[161]132
[1032]133      REAL d_t_the(ngrid,nlayer), d_q_the(ngrid,nlayer,nq)
[561]134      INTEGER isplit
[342]135      REAL fact
[1032]136      REAL zfm_therm(ngrid,nlayer+1),zdt
137      REAL zentr_therm(ngrid,nlayer),zdetr_therm(ngrid,nlayer)
138      REAL zheatFlux(ngrid,nlayer)
139      REAL zheatFlux_down(ngrid,nlayer)
140      REAL zbuoyancyOut(ngrid,nlayer)
141      REAL zbuoyancyEst(ngrid,nlayer)
142      REAL zzw2(ngrid,nlayer+1)
143      REAL zmax(ngrid)
[1212]144      INTEGER ndt,limz
[342]145
146!--------------------------------------------------------
147! Diagnostics
148!--------------------------------------------------------
149
[1032]150      REAL heatFlux(ngrid,nlayer)
151      REAL heatFlux_down(ngrid,nlayer)
152      REAL buoyancyOut(ngrid,nlayer)
153      REAL buoyancyEst(ngrid,nlayer)
154      REAL hfmax(ngrid),wmax(ngrid)
155      REAL pbl_teta(ngrid),dteta(ngrid,nlayer)
156      REAL rpdhd(ngrid,nlayer)
157      REAL wtdif(ngrid,nlayer),rho(ngrid,nlayer)
158      REAL wtth(ngrid,nlayer)
[161]159
160! **********************************************************************
[1028]161! Initializations
162! **********************************************************************
163
[621]164      lmax(:)=0
[161]165      pdu_th(:,:)=0.
166      pdv_th(:,:)=0.
167      pdt_th(:,:)=0.
168      entr_therm(:,:)=0.
169      detr_therm(:,:)=0.
170      q2_therm(:,:)=0.
171      dq2_therm(:,:)=0.
172      pbl_dtke(:,:)=0.
173      fm_therm(:,:)=0.
174      zw2(:,:)=0.
175      fraca(:,:)=0.
[512]176      zfraca(:,:)=0.
[2823]177      pdq_th(:,:,:)=0.
[342]178      d_t_ajs(:,:)=0.
179      d_u_ajs(:,:)=0.
180      d_v_ajs(:,:)=0.
181      d_q_ajs(:,:,:)=0.
182      heatFlux(:,:)=0.
183      heatFlux_down(:,:)=0.
184      buoyancyOut(:,:)=0.
185      buoyancyEst(:,:)=0.
186      zmaxth(:)=0.
187      lmax_real(:)=0.
[161]188
189
[342]190! **********************************************************************
[1028]191! Preparing inputs for the thermals: increment tendancies
192! from other subroutines called before the thermals model
[342]193! **********************************************************************
[161]194
[1028]195       zu(:,:)=pu(:,:)+pdu(:,:)*ptimestep ! u-component of wind velocity
196       zv(:,:)=pv(:,:)+pdv(:,:)*ptimestep ! v-component of wind velocity
197       zt(:,:)=pt(:,:)+pdt(:,:)*ptimestep ! temperature
[161]198
[1028]199       pq_therm(:,:,:)=0. ! tracer concentration
[161]200
[342]201       if(qtransport_thermals) then
[2823]202         pq_therm(:,:,:)=pq(:,:,:)+pdq(:,:,:)*ptimestep ! tracer concentration
[342]203       endif
[161]204
[1028]205
[544]206       IF(dtke_thermals) THEN
[1032]207          DO l=1,nlayer
[544]208              q2_therm(:,l)=0.5*(q2(:,l)+q2(:,l+1))
209          ENDDO
210       ENDIF
[342]211
212! **********************************************************************
[1032]213! --> CALLTHERM
[342]214! SUB-TIMESTEP LOOP
215! **********************************************************************
216
[1028]217         zdt=ptimestep/REAL(nsplit_thermals) !subtimestep
[342]218
[1028]219         DO isplit=1,nsplit_thermals !temporal loop on the subtimestep
[342]220
221! Initialization of intermediary variables
222
223         zzw2(:,:)=0.
224         zmax(:)=0.
[621]225         lmax(:)=0
[628]226
[1032]227         if (nq .ne. 0 .and. igcm_co2 .ne. 0) then !initialize co2 tracer tendancy
228            d_q_the(:,:,igcm_co2)=0.
[161]229         endif
230
[1028]231! CALL to main thermal routine
[1032]232             CALL thermcell_main_mars(ngrid,nlayer,nq &
[2823]233     &      ,igcm_co2 &
[1032]234     &      ,zdt  &
[342]235     &      ,pplay,pplev,pphi,zzlev,zzlay  &
236     &      ,zu,zv,zt,pq_therm,q2_therm  &
[1028]237     &      ,d_t_the,d_q_the  &
[1212]238     &      ,zfm_therm,zentr_therm,zdetr_therm,lmax,zmax,limz  &
[342]239     &      ,zzw2,fraca,zpopsk &
[1028]240     &      ,zheatFlux,zheatFlux_down &
[342]241     &      ,zbuoyancyOut,zbuoyancyEst)
[161]242
[342]243      fact=1./REAL(nsplit_thermals)
[161]244
[1028]245! Update thermals tendancies
246            d_t_the(:,:)=d_t_the(:,:)*ptimestep*fact !temperature
[1032]247            if (igcm_co2 .ne. 0) then
248               d_q_the(:,:,igcm_co2)=d_q_the(:,:,igcm_co2)*ptimestep*fact !co2 mass mixing ratio
[508]249            endif
[1028]250            zmaxth(:)=zmaxth(:)+zmax(:)*fact !thermals height
251            lmax_real(:)=lmax_real(:)+float(lmax(:))*fact !thermals height index
252            fm_therm(:,:)=fm_therm(:,:)  & !upward mass flux
[342]253     &      +zfm_therm(:,:)*fact
[1028]254            entr_therm(:,:)=entr_therm(:,:)  & !entrainment mass flux
[342]255     &       +zentr_therm(:,:)*fact
[1028]256            detr_therm(:,:)=detr_therm(:,:)  & !detrainment mass flux
[342]257     &       +zdetr_therm(:,:)*fact
[1028]258            zfraca(:,:)=zfraca(:,:) + fraca(:,:)*fact !updraft fractional coverage
259            heatFlux(:,:)=heatFlux(:,:) & !upward heat flux
[342]260     &       +zheatFlux(:,:)*fact
[1028]261            heatFlux_down(:,:)=heatFlux_down(:,:) & !downward heat flux
[342]262     &       +zheatFlux_down(:,:)*fact
[1028]263            buoyancyOut(:,:)=buoyancyOut(:,:) & !plume final buoyancy
[508]264     &       +zbuoyancyOut(:,:)*fact
[1028]265            buoyancyEst(:,:)=buoyancyEst(:,:) & !plume estimated buoyancy used for vertical velocity computation
[508]266     &       +zbuoyancyEst(:,:)*fact
[1028]267            zw2(:,:)=zw2(:,:) + zzw2(:,:)*fact !vertical velocity
[342]268
[1028]269! Save tendancies
270           d_t_ajs(:,:)=d_t_ajs(:,:)+d_t_the(:,:) !temperature tendancy (delta T)
[1032]271            if (igcm_co2 .ne. 0) then
272               d_q_ajs(:,:,igcm_co2)=d_q_ajs(:,:,igcm_co2)+d_q_the(:,:,igcm_co2) !tracer tendancy (delta q)
[508]273            endif
[342]274
[1028]275!  Increment temperature and co2 concentration for next pass in subtimestep loop
276            zt(:,:) = zt(:,:) + d_t_the(:,:) !temperature
[1032]277            if (igcm_co2 .ne. 0) then
278             pq_therm(:,:,igcm_co2) = &
279     &          pq_therm(:,:,igcm_co2) + d_q_the(:,:,igcm_co2) !co2 tracer
[508]280            endif
[342]281
282
283         ENDDO ! isplit
284!****************************************************************
285
286! Now that we have computed total entrainment and detrainment, we can
[1028]287! advect u, v, and q in thermals. (potential temperature and co2 MMR
288! have already been advected in thermcell_main because they are coupled
289! to the determination of the thermals caracteristics). This is done
290! separatly because u,v, and q are not used in thermcell_main for
[342]291! any thermals-related computation : they are purely passive.
292
293! mass of cells
[1032]294      do l=1,nlayer
[342]295         masse(:,l)=(pplev(:,l)-pplev(:,l+1))/g
296      enddo
297
[1028]298! recompute detrainment mass flux from entrainment and upward mass flux
299! this ensure mass flux conservation
[628]300      detrmod(:,:)=0.
[1212]301      do l=1,limz
[1032]302         do ig=1,ngrid
[628]303            detrmod(ig,l)=fm_therm(ig,l)-fm_therm(ig,l+1) &
304     &      +entr_therm(ig,l)
305            if (detrmod(ig,l).lt.0.) then
306               entr_therm(ig,l)=entr_therm(ig,l)-detrmod(ig,l)
307               detrmod(ig,l)=0.
308            endif
309         enddo
310      enddo
[1028]311
312      ! u component of wind velocity advection in thermals
313      ! result is a derivative (d_u_ajs in m/s/s)
[628]314      ndt=10
[1032]315      call thermcell_dqup(ngrid,nlayer,ptimestep                &
[628]316     &      ,fm_therm,entr_therm,detrmod,  &
[1212]317     &     masse,zu,d_u_ajs,ndt,limz)
[342]318
[1028]319      ! v component of wind velocity advection in thermals
320      ! result is a derivative (d_v_ajs in m/s/s)
[1032]321      call thermcell_dqup(ngrid,nlayer,ptimestep    &
[628]322     &       ,fm_therm,entr_therm,detrmod,  &
[1212]323     &     masse,zv,d_v_ajs,ndt,limz)
[628]324
[1028]325      ! non co2 tracers advection in thermals
326      ! result is a derivative (d_q_ajs in kg/kg/s)
327
[1032]328      if (nq .ne. 0.) then
329      DO iq=1,nq
330      if (iq .ne. igcm_co2) then
331      call thermcell_dqup(ngrid,nlayer,ptimestep     &
[628]332     &     ,fm_therm,entr_therm,detrmod,  &
[1212]333     &    masse,pq_therm(:,:,iq),d_q_ajs(:,:,iq),ndt,limz)
[508]334      endif
[342]335      ENDDO
336      endif
337
[1028]338      ! tke advection in thermals
339      ! result is a tendancy (d_u_ajs in J)
[544]340      if (dtke_thermals) then
[1032]341      call thermcell_dqup(ngrid,nlayer,ptimestep     &
[628]342     &     ,fm_therm,entr_therm,detrmod,  &
[1212]343     &    masse,q2_therm,dq2_therm,ndt,limz)
[544]344      endif
345
[1028]346      ! compute wmax for diagnostics
[1032]347      DO ig=1,ngrid
[342]348         wmax(ig)=MAXVAL(zw2(ig,:))
349      ENDDO
350
351! **********************************************************************
352! **********************************************************************
353! **********************************************************************
354! CALLTHERM END
355! **********************************************************************
356! **********************************************************************
357! **********************************************************************
358
359
360! **********************************************************************
361! Preparing outputs
362! **********************************************************************
363
[1212]364      do l=1,limz
[628]365        pdu_th(:,l)=d_u_ajs(:,l)
366        pdv_th(:,l)=d_v_ajs(:,l)
367      enddo
[342]368
[1028]369           ! if tracers are transported in thermals, update output variables, else these are 0.
[161]370           if(qtransport_thermals) then
[1032]371               do iq=1,nq
372                if (iq .ne. igcm_co2) then
[1212]373                  do l=1,limz
[1028]374                     pdq_th(:,l,iq)=d_q_ajs(:,l,iq) !non-co2 tracers d_q_ajs are dq/dt (kg/kg/s)
[628]375                  enddo
[625]376                else
[1212]377                  do l=1,limz
[1028]378                     pdq_th(:,l,iq)=d_q_ajs(:,l,iq)/ptimestep !co2 tracer d_q_ajs is dq (kg/kg)
[628]379                  enddo
[625]380                endif
381               enddo
[161]382           endif
383
[1028]384           ! if tke is transported in thermals, update output variable, else this is 0.
[544]385           IF(dtke_thermals) THEN
[1032]386              DO l=2,nlayer
[544]387                 pbl_dtke(:,l)=0.5*(dq2_therm(:,l-1)+dq2_therm(:,l))
388              ENDDO
389 
390              pbl_dtke(:,1)=0.5*dq2_therm(:,1)
[1032]391              pbl_dtke(:,nlayer+1)=0.
[544]392           ENDIF
[161]393
[1028]394           ! update output variable for temperature. d_t_ajs is delta T in (K), pdt_th is dT/dt in (K/s)
[1212]395           do l=1,limz
[628]396              pdt_th(:,l)=d_t_ajs(:,l)/ptimestep
397           enddo
[342]398
[499]399
[342]400! **********************************************************************
[1028]401! SURFACE LAYER INTERFACE
402! Compute the free convection velocity w* scale for surface layer gustiness
403! speed parameterization. The computed value of w* will be used at the next
404! timestep to modify surface-atmosphere exchange fluxes, because the surface
405! layer scheme and diffusion are called BEFORE the thermals. (outside of these
406! routines)
[499]407! **********************************************************************
408
409! Potential temperature gradient
410
[1032]411      dteta(:,nlayer)=0.
412      DO l=1,nlayer-1
413         DO ig=1, ngrid
[499]414            dteta(ig,l) = ((zt(ig,l+1)-zt(ig,l))/zpopsk(ig,l))          &
415     &              /(zzlay(ig,l+1)-zzlay(ig,l))
416         ENDDO
417      ENDDO
418
[1028]419! Computation of the PBL mixed layer temperature
[499]420
[1032]421      DO ig=1, ngrid
[499]422         ii=MINLOC(abs(dteta(ig,1:lmax(ig))))
423         pbl_teta(ig) = zt(ig,ii(1))/zpopsk(ig,ii(1))
424      ENDDO
425
[1028]426! In order to have an accurate w*, we must add the heat flux from the
427! diffusion scheme to the computation of the maximum heat flux hfmax
428! Here pdhdif is the potential temperature change from the diffusion
429! scheme (Mellor and Yamada, see paper section 6, paragraph 57)
[660]430
431! compute rho as it is after the diffusion
432
433      rho(:,:)=pplay(:,:)                                               &
434     & /(r*(pt(:,:)+pdhdif(:,:)*zpopsk(:,:)*ptimestep))
435
436! integrate -rho*pdhdif
437
438      rpdhd(:,:)=0.
439
[1032]440      DO ig=1,ngrid
[660]441       DO l=1,lmax(ig)
442        rpdhd(ig,l)=0.
443        DO k=1,l
444         rpdhd(ig,l)=rpdhd(ig,l)-rho(ig,k)*pdhdif(ig,k)*                &
445     & (zzlev(ig,k+1)-zzlev(ig,k))
446        ENDDO
447        rpdhd(ig,l)=rpdhd(ig,l)-sensibFlux(ig)/cpp
448       ENDDO
449      ENDDO
450
[1028]451! compute w'theta' (vertical turbulent flux of temperature) from
452! the diffusion scheme
[660]453
454      wtdif(:,:)=rpdhd(:,:)/rho(:,:)
455
[1028]456! Now we compute the contribution of the thermals to w'theta':
[660]457! compute rho as it is after the thermals
458
459      rho(:,:)=pplay(:,:)/(r*(zt(:,:)))
[1028]460
[660]461! integrate -rho*pdhdif
462
[1032]463      DO ig=1,ngrid
[660]464       DO l=1,lmax(ig)
465        rpdhd(ig,l)=0.
466        DO k=1,l
467         rpdhd(ig,l)=rpdhd(ig,l)-rho(ig,k)*(pdt_th(ig,k)/zpopsk(ig,k))* &
468     & (zzlev(ig,k+1)-zzlev(ig,k))
469        ENDDO
470        rpdhd(ig,l)=rpdhd(ig,l)+                                        &
471     &    rho(ig,1)*(heatFlux(ig,1)+heatFlux_down(ig,1))
472       ENDDO
473      ENDDO
[1032]474      rpdhd(:,nlayer)=0.
[660]475
476! compute w'teta' from thermals
477
478      wtth(:,:)=rpdhd(:,:)/rho(:,:)
479
[1028]480! Add vertical turbulent heat fluxes from the thermals and the diffusion scheme
481! and compute the maximum
[660]482
[1032]483      DO ig=1,ngrid
[660]484        hfmax(ig)=MAXVAL(wtth(ig,:)+wtdif(ig,:))
485      ENDDO
[1028]486
487! Finally we can compute the free convection velocity scale
[499]488! We follow Spiga et. al 2010 (QJRMS)
489! ------------
490
[1032]491      DO ig=1, ngrid
[499]492         IF (zmax(ig) .gt. 0.) THEN
493            wstar(ig)=(g*zmaxth(ig)*hfmax(ig)/pbl_teta(ig))**(1./3.)
494         ELSE
495            wstar(ig)=0.
496         ENDIF
497      ENDDO
498
[3727]499       END SUBROUTINE calltherm_interface
500
501END MODULE calltherm_interface_mod
Note: See TracBrowser for help on using the repository browser.