source: LMDZ5/branches/LF-private/libf/phylmd/ener_conserv.F90 @ 5080

Last change on this file since 5080 was 1761, checked in by Laurent Fairhead, 11 years ago

New version of Mellor et Yamada pronostic TKE

... based on energy transfer from the mean state.

The new version is yamada_c.
It must be called after vertical diffusion rather than just before
since the source terms u_z w'u' + v_z w'theta' and g/theta w'theta'
are diagnosed from the vertical diffusion (as energy loss from the mean
state) rather than computed as K (u_z2+v_z2) or g/\theta K theta_z
(where _z means vertical derivative).
The call to this version is controled by iflag_pbl.

iflag_pbl = 20 : with TKE computed at interfaces between layers
iflag_pbl = 25 : with TKE computed within the layer
In both cases, the dissipation of turbulence is translated into heat, and
passed to the physics as dtdiss (temperature tendency due to dissipation

of TKE).

The diffusion coefficient being computed after dissipation, it must be
kept for diffusion at the next time step, and thus be stored in the
restartphy file.
This coefficient must be computed and stored for each sub-surface.

A new way of managing subsurface variables is introduced.
For arrays of the form X(:,nbsrf) are extented to X(:,nbsrf+1), where
is_ave=nbsrf+1, is an additional sub-surface which contains the averaged values.

coef_diff_turb_mod.F90 : change of flags.
ener_conserv.F90 : energy conservation must not be applied in those

cases

indicesol.h : definition of is_ave
pbl_surface_mod.F90 : call to yamada_c and changes in the management of

coefh/coefm

physiq.F : Change in the initialisation of pmflxr/s and

modified calls to pbl_surface_mod (introduction
of dtdiss, initialisation
of pbl_tke and coefh in 1D).

phys_local_var_mod.F90 : declaration of d_t_diss
phys_output_mod.F90 : new outputs (bils_tke and bils_diss) and

coefh->coefh(:,:,is_ave)

phys_state_var_mod.F90 : modified declaration for coefh and coefm

(nbsrf -> nbsrf+1)

FH

File size: 5.3 KB
Line 
1subroutine ener_conserv(klon,klev,pdtphys, &
2 &    puo,pvo,pto,pqo,pun,pvn,ptn,pqn,dtke,masse,exner,d_t_ec)
3
4!=============================================================
5! Energy conservation
6! Based on the TKE equation
7! The M2 and N2 terms at the origin of TKE production are
8! concerted into heating in the d_t_ec term
9! Option 1 is the standard
10!        101 is for M2 term only
11!        101 for N2 term only
12!         -1 is a previours treatment for kinetic energy only
13!  FH (hourdin@lmd.jussieu.fr), 2013/04/25
14!=============================================================
15
16!=============================================================
17! Declarations
18!=============================================================
19
20! From module
21USE phys_local_var_mod, ONLY : d_u_vdf,d_v_vdf,d_t_vdf,d_u_ajs,d_v_ajs,d_t_ajs,d_u_con,d_v_con,d_t_con,d_t_diss
22USE phys_output_var_mod, ONLY : bils_ec,bils_tke,bils_kinetic,bils_enthalp,bils_latent,bils_diss
23
24IMPLICIT none
25#include "YOMCST.h"
26#include "YOETHF.h"
27#include "clesphys.h"
28#include "compbl.h"
29
30! Arguments
31INTEGER, INTENT(IN) :: klon,klev
32REAL, INTENT(IN) :: pdtphys
33REAL, DIMENSION(klon,klev),INTENT(IN) :: puo,pvo,pto,pqo
34REAL, DIMENSION(klon,klev),INTENT(IN) :: pun,pvn,ptn,pqn
35REAL, DIMENSION(klon,klev),INTENT(IN) :: masse,exner
36REAL, DIMENSION(klon,klev+1),INTENT(IN) :: dtke
37REAL, DIMENSION(klon,klev),INTENT(OUT) :: d_t_ec
38      integer k,i
39
40! Local
41REAL, DIMENSION(klon,klev+1) :: fluxu,fluxv,fluxt
42REAL, DIMENSION(klon,klev+1) :: dddu,dddv,dddt
43REAL, DIMENSION(klon,klev) :: d_u,d_v,d_t,zv,zu
44REAL ZRCPD
45
46character*80 abort_message
47character*20 :: modname
48
49
50modname='ener_conser'
51d_t_ec(:,:)=0.
52
53IF (iflag_ener_conserv==-1) THEN
54!+jld ec_conser
55   DO k = 1, klev
56   DO i = 1, klon
57      ZRCPD = RCPD*(1.0+RVTMP2*pqn(i,k))
58      d_t_ec(i,k)=0.5/ZRCPD &
59 &      *(puo(i,k)**2+pvo(i,k)**2-pun(i,k)**2-pvn(i,k)**2)
60      ENDDO
61      ENDDO
62!-jld ec_conser
63
64
65
66ELSEIF (iflag_ener_conserv>=1) THEN
67
68   IF (iflag_ener_conserv<=2) THEN
69!     print*,'ener_conserv pbl=',iflag_pbl
70      IF (iflag_pbl>=20 .AND. iflag_pbl<=27) THEN !d_t_diss accounts for conserv
71         d_t(:,:)=d_t_ajs(:,:)   ! d_t_ajs = adjust + thermals
72         d_u(:,:)=d_u_ajs(:,:)+d_u_con(:,:)
73         d_v(:,:)=d_v_ajs(:,:)+d_v_con(:,:)
74      ELSE
75         d_t(:,:)=d_t_vdf(:,:)+d_t_ajs(:,:)   ! d_t_ajs = adjust + thermals
76         d_u(:,:)=d_u_vdf(:,:)+d_u_ajs(:,:)+d_u_con(:,:)
77         d_v(:,:)=d_v_vdf(:,:)+d_v_ajs(:,:)+d_v_con(:,:)
78      ENDIF
79   ELSEIF (iflag_ener_conserv==101) THEN
80      d_t(:,:)=0.
81      d_u(:,:)=d_u_vdf(:,:)+d_u_ajs(:,:)+d_u_con(:,:)
82      d_v(:,:)=d_v_vdf(:,:)+d_v_ajs(:,:)+d_v_con(:,:)
83   ELSEIF (iflag_ener_conserv==110) THEN
84      d_t(:,:)=d_t_vdf(:,:)+d_t_ajs(:,:)
85      d_u(:,:)=0.
86      d_v(:,:)=0.
87   ELSE
88      abort_message = 'iflag_ener_conserv non prevu'
89      CALL abort_gcm (modname,abort_message,1)
90   ENDIF
91
92!----------------------------------------------------------------------------
93! Two options wether we consider time integration in the energy conservation
94!----------------------------------------------------------------------------
95
96   if (iflag_ener_conserv==2) then
97      zu(:,:)=puo(:,:)
98      zv(:,:)=pvo(:,:)
99   else
100      IF (iflag_pbl>=20 .AND. iflag_pbl<=27) THEN
101         zu(:,:)=puo(:,:)+d_u_vdf(:,:)+0.5*d_u(:,:)
102         zv(:,:)=pvo(:,:)+d_v_vdf(:,:)+0.5*d_v(:,:)
103      ELSE
104         zu(:,:)=puo(:,:)+0.5*d_u(:,:)
105         zv(:,:)=pvo(:,:)+0.5*d_v(:,:)
106      ENDIF
107   endif
108
109   fluxu(:,klev+1)=0.
110   fluxv(:,klev+1)=0.
111   fluxt(:,klev+1)=0.
112
113   do k=klev,1,-1
114      fluxu(:,k)=fluxu(:,k+1)+masse(:,k)*d_u(:,k)
115      fluxv(:,k)=fluxv(:,k+1)+masse(:,k)*d_v(:,k)
116      fluxt(:,k)=fluxt(:,k+1)+masse(:,k)*d_t(:,k)/exner(:,k)
117   enddo
118
119   dddu(:,1)=2*zu(:,1)*fluxu(:,1)
120   dddv(:,1)=2*zv(:,1)*fluxv(:,1)
121   dddt(:,1)=(exner(:,1)-1.)*fluxt(:,1)
122
123   do k=2,klev
124      dddu(:,k)=(zu(:,k)-zu(:,k-1))*fluxu(:,k)
125      dddv(:,k)=(zv(:,k)-zv(:,k-1))*fluxv(:,k)
126      dddt(:,k)=(exner(:,k)-exner(:,k-1))*fluxt(:,k)
127   enddo
128   dddu(:,klev+1)=0.
129   dddv(:,klev+1)=0.
130   dddt(:,klev+1)=0.
131
132   do k=1,klev
133      d_t_ec(:,k)=-(dddu(:,k)+dddu(:,k+1)+dddv(:,k)+dddv(:,k+1) &
134   &  +rcpd*(dddt(:,k)+dddt(:,k+1)))/(2.*rcpd*masse(:,k))
135   enddo
136! d_t_ec=0.
137
138ENDIF
139
140!================================================================
141!  Computation of integrated enthalpie and kinetic energy variation
142!  FH (hourdin@lmd.jussieu.fr), 2013/04/25
143!================================================================
144
145      bils_ec(:)=0.
146      bils_tke(:)=0.
147      bils_diss(:)=0.
148      bils_kinetic(:)=0.
149      bils_enthalp(:)=0.
150      bils_latent(:)=0.
151      DO k=1,klev
152        bils_ec(:)=bils_ec(:)-d_t_ec(:,k)*masse(:,k)
153        bils_tke(:)=bils_tke(:)+0.5*(dtke(:,k)+dtke(:,k+1))*masse(:,k)
154        bils_diss(:)=bils_diss(:)-d_t_diss(:,k)*masse(:,k)
155        bils_kinetic(:)=bils_kinetic(:)+masse(:,k)* &
156     &           (pun(:,k)*pun(:,k)+pvn(:,k)*pvn(:,k) &
157     &            -puo(:,k)*puo(:,k)-pvo(:,k)*pvo(:,k))
158        bils_enthalp(:)= &
159     &  bils_enthalp(:)+masse(:,k)*(ptn(:,k)-pto(:,k)+d_t_ec(:,k))
160        bils_latent(:)=bils_latent(:)+masse(:,k)* &
161     &             (pqn(:,k)-pqo(:,k))
162      ENDDO
163      bils_ec(:)=rcpd*bils_ec(:)/pdtphys
164      bils_tke(:)=bils_tke(:)/pdtphys
165      bils_diss(:)=rcpd*bils_diss(:)/pdtphys
166      bils_kinetic(:)= 0.5*bils_kinetic(:)/pdtphys
167      bils_enthalp(:)=rcpd*bils_enthalp(:)/pdtphys
168      bils_latent(:)=rlvtt*bils_latent(:)/pdtphys
169RETURN
170
171END
Note: See TracBrowser for help on using the repository browser.