source: trunk/LMDZ.MARS/libf/phymars/calltherm_interface.F90 @ 847

Last change on this file since 847 was 709, checked in by aslmd, 13 years ago

MESOSCALE: set default r_aspect_thermals to 1 also in mesoscale. updated DEF files for reference THARSIS case (latest settings by TN).

File size: 19.7 KB
RevLine 
[161]1!
2! AC 2011-01-05
3!
[695]4      SUBROUTINE calltherm_interface (firstcall, &
[652]5     & zzlev,zzlay, &
[161]6     & ptimestep,pu,pv,pt,pq,pdu,pdv,pdt,pdq,q2, &
[185]7     & pplay,pplev,pphi,zpopsk, &
[660]8     & pdu_th,pdv_th,pdt_th,pdq_th,lmax,zmaxth,pbl_dtke, &
9     & pdhdif,hfmax,wstar,sensibFlux)
[161]10
[342]11       USE ioipsl_getincom
[161]12
13      implicit none
14#include "callkeys.h"
[185]15#include "dimensions.h"
16#include "dimphys.h"
[342]17#include "comcstfi.h"
[508]18#include "tracer.h"
[185]19
[161]20!--------------------------------------------------------
[342]21! Input Variables
[161]22!--------------------------------------------------------
23
[652]24!      REAL, INTENT(IN) :: long(ngridmx),lati(ngridmx)
[161]25      REAL, INTENT(IN) :: ptimestep
[695]26      REAL, INTENT(IN) :: pplev(ngridmx,nlayermx+1)
27      REAL, INTENT(IN) :: pplay(ngridmx,nlayermx)
[185]28      REAL, INTENT(IN) :: pphi(ngridmx,nlayermx)
29      REAL, INTENT(IN) :: pu(ngridmx,nlayermx),pv(ngridmx,nlayermx)
30      REAL, INTENT(IN) :: pt(ngridmx,nlayermx),pq(ngridmx,nlayermx,nqmx)
31      REAL, INTENT(IN) :: zzlay(ngridmx,nlayermx)
32      REAL, INTENT(IN) :: zzlev(ngridmx,nlayermx+1)
[161]33      LOGICAL, INTENT(IN) :: firstcall
[185]34      REAL, INTENT(IN) :: pdu(ngridmx,nlayermx),pdv(ngridmx,nlayermx)
[695]35      REAL, INTENT(IN) :: pdq(ngridmx,nlayermx,nqmx)
36      REAL, INTENT(IN) :: pdt(ngridmx,nlayermx)
[185]37      REAL, INTENT(IN) :: q2(ngridmx,nlayermx+1)
38      REAL, INTENT(IN) :: zpopsk(ngridmx,nlayermx)
[660]39      REAL, INTENT(IN) :: pdhdif(ngridmx,nlayermx)
40      REAL, INTENT(IN) :: sensibFlux(ngridmx)
[161]41
42!--------------------------------------------------------
[342]43! Output Variables
[161]44!--------------------------------------------------------
45
[342]46      REAL, INTENT(OUT) :: pdu_th(ngridmx,nlayermx)
47      REAL, INTENT(OUT) :: pdv_th(ngridmx,nlayermx)
48      REAL, INTENT(OUT) :: pdt_th(ngridmx,nlayermx)
49      REAL, INTENT(OUT) :: pdq_th(ngridmx,nlayermx,nqmx)
50      INTEGER, INTENT(OUT) :: lmax(ngridmx)
51      REAL, INTENT(OUT) :: zmaxth(ngridmx)
52      REAL, INTENT(OUT) :: pbl_dtke(ngridmx,nlayermx+1)
[499]53      REAL, INTENT(OUT) :: wstar(ngridmx)
[161]54
55!--------------------------------------------------------
[342]56! Thermals local variables
[161]57!--------------------------------------------------------
[342]58      REAL zu(ngridmx,nlayermx), zv(ngridmx,nlayermx)
59      REAL zt(ngridmx,nlayermx)
[185]60      REAL d_t_ajs(ngridmx,nlayermx)
61      REAL d_u_ajs(ngridmx,nlayermx), d_q_ajs(ngridmx,nlayermx,nqmx)
62      REAL d_v_ajs(ngridmx,nlayermx)
63      REAL fm_therm(ngridmx,nlayermx+1), entr_therm(ngridmx,nlayermx)
[628]64      REAL detr_therm(ngridmx,nlayermx),detrmod(ngridmx,nlayermx)
[185]65      REAL zw2(ngridmx,nlayermx+1)
[512]66      REAL fraca(ngridmx,nlayermx+1),zfraca(ngridmx,nlayermx+1)
[185]67      REAL ztla(ngridmx,nlayermx)
68      REAL q_therm(ngridmx,nlayermx), pq_therm(ngridmx,nlayermx,nqmx)
69      REAL q2_therm(ngridmx,nlayermx), dq2_therm(ngridmx,nlayermx)
[342]70      REAL lmax_real(ngridmx)
71      REAL masse(ngridmx,nlayermx)
[161]72      LOGICAL qtransport_thermals,dtke_thermals
[660]73      INTEGER l,ig,iq,ii(1),k
[628]74      CHARACTER (LEN=20) modname
[161]75
[342]76!--------------------------------------------------------
77! Local variables for sub-timestep
78!--------------------------------------------------------
[161]79
[342]80      REAL d_t_the(ngridmx,nlayermx), d_q_the(ngridmx,nlayermx,nqmx)
81      REAL d_u_the(ngridmx,nlayermx),d_v_the(ngridmx,nlayermx)
82      REAL dq2_the(ngridmx,nlayermx)
[561]83      INTEGER isplit
84      INTEGER,SAVE :: nsplit_thermals
85      REAL, SAVE :: r_aspect_thermals
[342]86      REAL fact
87      REAL zfm_therm(ngridmx,nlayermx+1),zdt
88      REAL zentr_therm(ngridmx,nlayermx),zdetr_therm(ngridmx,nlayermx)
89      REAL zheatFlux(ngridmx,nlayermx)
90      REAL zheatFlux_down(ngridmx,nlayermx)
91      REAL zbuoyancyOut(ngridmx,nlayermx)
92      REAL zbuoyancyEst(ngridmx,nlayermx)
93      REAL zzw2(ngridmx,nlayermx+1)
94      REAL zmax(ngridmx)
[628]95      INTEGER ndt,zlmax
[342]96
97!--------------------------------------------------------
98! Diagnostics
99!--------------------------------------------------------
100
[185]101      REAL heatFlux(ngridmx,nlayermx)
102      REAL heatFlux_down(ngridmx,nlayermx)
103      REAL buoyancyOut(ngridmx,nlayermx)
104      REAL buoyancyEst(ngridmx,nlayermx)
105      REAL hfmax(ngridmx),wmax(ngridmx)
[499]106      REAL pbl_teta(ngridmx),dteta(ngridmx,nlayermx)
[660]107      REAL rpdhd(ngridmx,nlayermx)
108      REAL wtdif(ngridmx,nlayermx),rho(ngridmx,nlayermx)
109      REAL wtth(ngridmx,nlayermx)
[161]110
[508]111!--------------------------------------------------------
112! Theta_m
113!--------------------------------------------------------
[342]114
[508]115      INTEGER ico2
116      SAVE ico2
[342]117
[161]118! **********************************************************************
[342]119! Initialization
[161]120! **********************************************************************
121
[621]122      lmax(:)=0
[161]123      pdu_th(:,:)=0.
124      pdv_th(:,:)=0.
125      pdt_th(:,:)=0.
126      entr_therm(:,:)=0.
127      detr_therm(:,:)=0.
128      q2_therm(:,:)=0.
129      dq2_therm(:,:)=0.
130      ztla(:,:)=0.
131      pbl_dtke(:,:)=0.
132      fm_therm(:,:)=0.
133      zw2(:,:)=0.
134      fraca(:,:)=0.
[512]135      zfraca(:,:)=0.
[161]136      if (tracer) then
137         pdq_th(:,:,:)=0.
138      end if
[342]139      d_t_ajs(:,:)=0.
140      d_u_ajs(:,:)=0.
141      d_v_ajs(:,:)=0.
142      d_q_ajs(:,:,:)=0.
143      heatFlux(:,:)=0.
144      heatFlux_down(:,:)=0.
145      buoyancyOut(:,:)=0.
146      buoyancyEst(:,:)=0.
147      zmaxth(:)=0.
148      lmax_real(:)=0.
[161]149
150
[342]151! **********************************************************************
152! Preparing inputs for the thermals
153! **********************************************************************
[161]154
[342]155       zu(:,:)=pu(:,:)+pdu(:,:)*ptimestep
156       zv(:,:)=pv(:,:)+pdv(:,:)*ptimestep
157       zt(:,:)=pt(:,:)+pdt(:,:)*ptimestep
[161]158
[342]159       pq_therm(:,:,:)=0.
160       qtransport_thermals=.true. !! default setting
161       !call getin("qtransport_thermals",qtransport_thermals)
[161]162
[342]163       if(qtransport_thermals) then
164          if(tracer) then
165                pq_therm(:,:,:)=pq(:,:,:)+pdq(:,:,:)*ptimestep
166          endif
167       endif
[161]168
[544]169       dtke_thermals=.false. !! default setting
170       call getin("dtke_thermals",dtke_thermals)
171       IF(dtke_thermals) THEN
172          DO l=1,nlayermx
173              q2_therm(:,l)=0.5*(q2(:,l)+q2(:,l+1))
174          ENDDO
175       ENDIF
[342]176
177! **********************************************************************
[508]178! Polar night mixing : theta_m
[342]179! **********************************************************************
[508]180
181      if(firstcall) then
182        ico2=0
183        if (tracer) then
184!     Prepare Special treatment if one of the tracers is CO2 gas
185           do iq=1,nqmx
186             if (noms(iq).eq."co2") then
187                ico2=iq
188             end if
189           enddo
190        endif
191      endif !of if firstcall
192
193
[342]194! **********************************************************************
[508]195! **********************************************************************
196! **********************************************************************
[342]197! CALLTHERM
198! **********************************************************************
199! **********************************************************************
200! **********************************************************************
201
[561]202!         r_aspect_thermals     ! Mainly control the shape of the temperature profile
203                                ! in the surface layer. Decreasing it goes toward
204                                ! a convective-adjustment like profile.
205!         nsplit_thermals       ! Sub-timestep for the thermals. Very dependant on the
206                                ! chosen timestep for the radiative transfer.
207                                ! It is recommended to run with 96 timestep per day and
208                                ! iradia = 1., configuration in which thermals can run
209                                ! very well with a sub-timestep of 10.
210         IF (firstcall) THEN
[709]211            r_aspect_thermals=1.  ! same value is OK for GCM and mesoscale
[342]212#ifdef MESOSCALE
[561]213            !! valid for timesteps < 200s
214            nsplit_thermals=4
[342]215#else
[592]216            IF ((ptimestep .le. 3699.*24./96.) .and. (iradia .eq. 1)) THEN
[561]217               nsplit_thermals=10
218            ELSE
219               nsplit_thermals=35
220            ENDIF
[342]221#endif
[561]222            call getin("nsplit_thermals",nsplit_thermals)
223            call getin("r_aspect_thermals",r_aspect_thermals)
224         ENDIF
[342]225
226! **********************************************************************
227! SUB-TIMESTEP LOOP
228! **********************************************************************
229
230         zdt=ptimestep/REAL(nsplit_thermals)
231
232         DO isplit=1,nsplit_thermals
233
234! Initialization of intermediary variables
235
[628]236!         zfm_therm(:,:)=0. !init is done inside
237!         zentr_therm(:,:)=0.
238!         zdetr_therm(:,:)=0.
239!         zheatFlux(:,:)=0. 
240!         zheatFlux_down(:,:)=0.
241!         zbuoyancyOut(:,:)=0.
242!         zbuoyancyEst(:,:)=0.
[342]243         zzw2(:,:)=0.
244         zmax(:)=0.
[621]245         lmax(:)=0
[628]246!         d_t_the(:,:)=0. !init is done inside
247
248!         d_u_the(:,:)=0. !transported outside
249!         d_v_the(:,:)=0.
[342]250         dq2_the(:,:)=0.
[628]251
252         if (nqmx .ne. 0 .and. ico2 .ne. 0) then
253            d_q_the(:,:,ico2)=0.
[161]254         endif
255
[342]256             CALL thermcell_main_mars(zdt  &
257     &      ,pplay,pplev,pphi,zzlev,zzlay  &
258     &      ,zu,zv,zt,pq_therm,q2_therm  &
259     &      ,d_u_the,d_v_the,d_t_the,d_q_the,dq2_the  &
260     &      ,zfm_therm,zentr_therm,zdetr_therm,lmax,zmax  &
261     &      ,r_aspect_thermals &
262     &      ,zzw2,fraca,zpopsk &
263     &      ,ztla,zheatFlux,zheatFlux_down &
264     &      ,zbuoyancyOut,zbuoyancyEst)
[161]265
[342]266      fact=1./REAL(nsplit_thermals)
[161]267
[342]268            d_t_the(:,:)=d_t_the(:,:)*ptimestep*fact
[628]269!            d_u_the(:,:)=d_u_the(:,:)*ptimestep*fact
270!            d_v_the(:,:)=d_v_the(:,:)*ptimestep*fact
271            dq2_the(:,:)=dq2_the(:,:)*fact
[508]272            if (ico2 .ne. 0) then
[624]273               d_q_the(:,:,ico2)=d_q_the(:,:,ico2)*ptimestep*fact
[508]274            endif
[161]275
[628]276            zmaxth(:)=zmaxth(:)+zmax(:)*fact
277            lmax_real(:)=lmax_real(:)+float(lmax(:))*fact
[342]278            fm_therm(:,:)=fm_therm(:,:)  &
279     &      +zfm_therm(:,:)*fact
280            entr_therm(:,:)=entr_therm(:,:)  &
281     &       +zentr_therm(:,:)*fact
282            detr_therm(:,:)=detr_therm(:,:)  &
283     &       +zdetr_therm(:,:)*fact
[512]284            zfraca(:,:)=zfraca(:,:) + fraca(:,:)*fact
[342]285
286            heatFlux(:,:)=heatFlux(:,:) &
287     &       +zheatFlux(:,:)*fact
288            heatFlux_down(:,:)=heatFlux_down(:,:) &
289     &       +zheatFlux_down(:,:)*fact
[508]290            buoyancyOut(:,:)=buoyancyOut(:,:) &
291     &       +zbuoyancyOut(:,:)*fact
292            buoyancyEst(:,:)=buoyancyEst(:,:) &
293     &       +zbuoyancyEst(:,:)*fact
[512]294 
[342]295
296            zw2(:,:)=zw2(:,:) + zzw2(:,:)*fact
297
298!  accumulation de la tendance
299
[624]300           d_t_ajs(:,:)=d_t_ajs(:,:)+d_t_the(:,:)
[628]301!           d_u_ajs(:,:)=d_u_ajs(:,:)+d_u_the(:,:)
302!           d_v_ajs(:,:)=d_v_ajs(:,:)+d_v_the(:,:)
[508]303            if (ico2 .ne. 0) then
304               d_q_ajs(:,:,ico2)=d_q_ajs(:,:,ico2)+d_q_the(:,:,ico2)
305            endif
[342]306!            dq2_therm(:,:)=dq2_therm(:,:)+dq2_the(:,:)
307!  incrementation des variables meteo
308
309            zt(:,:) = zt(:,:) + d_t_the(:,:)
[628]310!            zu(:,:) = zu(:,:) + d_u_the(:,:)
311!            zv(:,:) = zv(:,:) + d_v_the(:,:)
[508]312            if (ico2 .ne. 0) then
313             pq_therm(:,:,ico2) = &
[624]314     &          pq_therm(:,:,ico2) + d_q_the(:,:,ico2)
[508]315            endif
[342]316!            q2_therm(:,:) = q2_therm(:,:) + dq2_therm(:,:)
317
318
319         ENDDO ! isplit
320!****************************************************************
321
[621]322      lmax(:)=nint(lmax_real(:))
[628]323      zlmax=MAXVAL(lmax(:))+2
324      if (zlmax .ge. nlayermx) then
325        print*,'thermals have reached last layer of the model'
326        print*,'this is not good !'
327      endif
[621]328
[628]329
[342]330! Now that we have computed total entrainment and detrainment, we can
331! advect u, v, and q in thermals. (theta already advected). We can do
332! that separatly because u,v,and q are not used in thermcell_main for
333! any thermals-related computation : they are purely passive.
334
335! mass of cells
336      do l=1,nlayermx
337         masse(:,l)=(pplev(:,l)-pplev(:,l+1))/g
338      enddo
339
[628]340      detrmod(:,:)=0.
341      do l=1,zlmax
342         do ig=1,ngridmx
343            detrmod(ig,l)=fm_therm(ig,l)-fm_therm(ig,l+1) &
344     &      +entr_therm(ig,l)
345            if (detrmod(ig,l).lt.0.) then
346               entr_therm(ig,l)=entr_therm(ig,l)-detrmod(ig,l)
347               detrmod(ig,l)=0.
348            endif
349         enddo
350      enddo
351      ndt=10
352      call thermcell_dqup(ngridmx,nlayermx,ptimestep                &
353     &      ,fm_therm,entr_therm,detrmod,  &
354     &     masse,zu,d_u_ajs,ndt,zlmax)
[342]355
[628]356      call thermcell_dqup(ngridmx,nlayermx,ptimestep    &
357     &       ,fm_therm,entr_therm,detrmod,  &
358     &     masse,zv,d_v_ajs,ndt,zlmax)
359
[342]360      if (nqmx .ne. 0.) then
361      DO iq=1,nqmx
[508]362      if (iq .ne. ico2) then
[342]363      call thermcell_dqup(ngridmx,nlayermx,ptimestep     &
[628]364     &     ,fm_therm,entr_therm,detrmod,  &
365     &    masse,pq_therm(:,:,iq),d_q_ajs(:,:,iq),ndt,zlmax)
[508]366      endif
[342]367      ENDDO
368      endif
369
[544]370      if (dtke_thermals) then
[628]371      detrmod(:,:)=0.
[652]372      ndt=10
[628]373      do l=1,zlmax
374         do ig=1,ngridmx
375            detrmod(ig,l)=fm_therm(ig,l)-fm_therm(ig,l+1) &
376     &      +entr_therm(ig,l)
377            if (detrmod(ig,l).lt.0.) then
378               entr_therm(ig,l)=entr_therm(ig,l)-detrmod(ig,l)
379               detrmod(ig,l)=0.
380            endif
381         enddo
382      enddo
[544]383      call thermcell_dqup(ngridmx,nlayermx,ptimestep     &
[628]384     &     ,fm_therm,entr_therm,detrmod,  &
385     &    masse,q2_therm,dq2_therm,ndt,zlmax)
[544]386      endif
387
[342]388      DO ig=1,ngridmx
389         wmax(ig)=MAXVAL(zw2(ig,:))
390      ENDDO
391
392! **********************************************************************
393! **********************************************************************
394! **********************************************************************
395! CALLTHERM END
396! **********************************************************************
397! **********************************************************************
398! **********************************************************************
399
400
401! **********************************************************************
402! Preparing outputs
403! **********************************************************************
404
[628]405      do l=1,zlmax
406        pdu_th(:,l)=d_u_ajs(:,l)
407        pdv_th(:,l)=d_v_ajs(:,l)
408      enddo
[342]409
[161]410           if(qtransport_thermals) then
[342]411              if(tracer) then
[625]412               do iq=1,nqmx
413                if (iq .ne. ico2) then
[628]414                  do l=1,zlmax
415                     pdq_th(:,l,iq)=d_q_ajs(:,l,iq)
416                  enddo
[625]417                else
[628]418                  do l=1,zlmax
419                     pdq_th(:,l,iq)=d_q_ajs(:,l,iq)/ptimestep
420                  enddo
[625]421                endif
422               enddo
[342]423              endif
[161]424           endif
425
[544]426           IF(dtke_thermals) THEN
427              DO l=2,nlayermx
428                 pbl_dtke(:,l)=0.5*(dq2_therm(:,l-1)+dq2_therm(:,l))
429              ENDDO
430 
431              pbl_dtke(:,1)=0.5*dq2_therm(:,1)
432              pbl_dtke(:,nlayermx+1)=0.
433           ENDIF
[161]434
[628]435           do l=1,zlmax
436              pdt_th(:,l)=d_t_ajs(:,l)/ptimestep
437           enddo
[342]438
[499]439
[342]440! **********************************************************************
[499]441! Compute the free convection velocity scale for vdifc
442! **********************************************************************
443
444
445! Potential temperature gradient
446
447      dteta(:,nlayermx)=0.
448      DO l=1,nlayermx-1
449         DO ig=1, ngridmx
450            dteta(ig,l) = ((zt(ig,l+1)-zt(ig,l))/zpopsk(ig,l))          &
451     &              /(zzlay(ig,l+1)-zzlay(ig,l))
452         ENDDO
453      ENDDO
454
455! Computation of the pbl mixed layer temperature
456
457      DO ig=1, ngridmx
458         ii=MINLOC(abs(dteta(ig,1:lmax(ig))))
459         pbl_teta(ig) = zt(ig,ii(1))/zpopsk(ig,ii(1))
460      ENDDO
461
[660]462! we must add the heat flux from the diffusion scheme to hfmax
463
464! compute rho as it is after the diffusion
465
466      rho(:,:)=pplay(:,:)                                               &
467     & /(r*(pt(:,:)+pdhdif(:,:)*zpopsk(:,:)*ptimestep))
468
469! integrate -rho*pdhdif
470
471      rpdhd(:,:)=0.
472
473      DO ig=1,ngridmx
474       DO l=1,lmax(ig)
475        rpdhd(ig,l)=0.
476        DO k=1,l
477         rpdhd(ig,l)=rpdhd(ig,l)-rho(ig,k)*pdhdif(ig,k)*                &
478     & (zzlev(ig,k+1)-zzlev(ig,k))
479        ENDDO
480        rpdhd(ig,l)=rpdhd(ig,l)-sensibFlux(ig)/cpp
481       ENDDO
482      ENDDO
483
484! compute w'teta' from diffusion
485
486      wtdif(:,:)=rpdhd(:,:)/rho(:,:)
487
488! compute rho as it is after the thermals
489
490      rho(:,:)=pplay(:,:)/(r*(zt(:,:)))
491! integrate -rho*pdhdif
492
493      DO ig=1,ngridmx
494       DO l=1,lmax(ig)
495        rpdhd(ig,l)=0.
496        DO k=1,l
497         rpdhd(ig,l)=rpdhd(ig,l)-rho(ig,k)*(pdt_th(ig,k)/zpopsk(ig,k))* &
498     & (zzlev(ig,k+1)-zzlev(ig,k))
499        ENDDO
500        rpdhd(ig,l)=rpdhd(ig,l)+                                        &
501     &    rho(ig,1)*(heatFlux(ig,1)+heatFlux_down(ig,1))
502       ENDDO
503      ENDDO
504      rpdhd(:,nlayermx)=0.
505
506! compute w'teta' from thermals
507
508      wtth(:,:)=rpdhd(:,:)/rho(:,:)
509
510! We get the max heat flux from thermals and add the contribution from the diffusion
511
512      DO ig=1,ngridmx
513        hfmax(ig)=MAXVAL(wtth(ig,:)+wtdif(ig,:))
514      ENDDO
[499]515! We follow Spiga et. al 2010 (QJRMS)
516! ------------
517
518      DO ig=1, ngridmx
519         IF (zmax(ig) .gt. 0.) THEN
520            wstar(ig)=(g*zmaxth(ig)*hfmax(ig)/pbl_teta(ig))**(1./3.)
521         ELSE
522            wstar(ig)=0.
523         ENDIF
524      ENDDO
525
526
527
528! **********************************************************************
[342]529! Diagnostics
530! **********************************************************************
[161]531       
532        if(outptherm) then
[185]533        if (ngridmx .eq. 1) then
534        call WRITEDIAGFI(ngridmx,'entr_therm','entrainement thermique',&
[161]535     &                       'kg/m-2',1,entr_therm)
[185]536        call WRITEDIAGFI(ngridmx,'detr_therm','detrainement thermique',&
[161]537     &                       'kg/m-2',1,detr_therm)
[185]538        call WRITEDIAGFI(ngridmx,'fm_therm','flux masse thermique',&
[161]539     &                       'kg/m-2',1,fm_therm)
[185]540        call WRITEDIAGFI(ngridmx,'zw2','vitesse verticale thermique',&
[161]541     &                       'm/s',1,zw2)
[185]542        call WRITEDIAGFI(ngridmx,'heatFlux_up','heatFlux_updraft',&
[161]543     &                       'SI',1,heatFlux)
[185]544       call WRITEDIAGFI(ngridmx,'heatFlux_down','heatFlux_downdraft',&
[161]545     &                       'SI',1,heatFlux_down)
[185]546        call WRITEDIAGFI(ngridmx,'fraca','fraction coverage',&
[161]547     &                       'percent',1,fraca)
[185]548        call WRITEDIAGFI(ngridmx,'buoyancyOut','buoyancyOut',&
[161]549     &                       'm.s-2',1,buoyancyOut)
[185]550        call WRITEDIAGFI(ngridmx,'buoyancyEst','buoyancyEst',&
[161]551     &                       'm.s-2',1,buoyancyEst)
[185]552        call WRITEDIAGFI(ngridmx,'d_t_th',  &
[161]553     &         'tendance temp TH','K',1,d_t_ajs)
[619]554        call WRITEDIAGFI(ngridmx,'d_q_th',  &
555     &         'tendance traceur TH','kg/kg',1,d_q_ajs)
[185]556        call WRITEDIAGFI(ngridmx,'zmax',  &
[342]557     &         'pbl height','m',0,zmaxth)
[624]558        call WRITEDIAGFI(ngridmx,'d_u_th',  &
559     &         'tendance moment','m/s',1,pdu_th)
[660]560        call WRITEDIAGFI(ngridmx,'wtdif',  &
561     &         'heat flux from diffusion','K.m/s',1,wtdif)
562        call WRITEDIAGFI(ngridmx,'wtth',  &
563     &         'heat flux from thermals','K.m/s',1,wtth)
564        call WRITEDIAGFI(ngridmx,'wttot',  &
565     &         'heat flux PBL','K.m/s',1,wtdif(:,:)+wtth(:,:))
566
[161]567      else
568
[185]569        call WRITEDIAGFI(ngridmx,'entr_therm','entrainement thermique',&
[161]570     &                       'kg/m-2',3,entr_therm)
[185]571        call WRITEDIAGFI(ngridmx,'detr_therm','detrainement thermique',&
[161]572     &                       'kg/m-2',3,detr_therm)
[185]573        call WRITEDIAGFI(ngridmx,'fm_therm','flux masse thermique',&
[161]574     &                       'kg/m-2',3,fm_therm)
[185]575        call WRITEDIAGFI(ngridmx,'zw2','vitesse verticale thermique',&
[161]576     &                       'm/s',3,zw2)
[185]577        call WRITEDIAGFI(ngridmx,'heatFlux','heatFlux',&
[161]578     &                       'SI',3,heatFlux)
[185]579        call WRITEDIAGFI(ngridmx,'buoyancyOut','buoyancyOut',&
[161]580     &                       'SI',3,buoyancyOut)
[185]581        call WRITEDIAGFI(ngridmx,'d_t_th',  &
[161]582     &         'tendance temp TH','K',3,d_t_ajs)
583
584      endif
585      endif
586
587       END
Note: See TracBrowser for help on using the repository browser.