source: LMDZ6/trunk/libf/phylmd/lmdz_call_atke.F90 @ 4976

Last change on this file since 4976 was 4881, checked in by evignon, 8 months ago

extraction plus propre de la dissipation de TKE

File size: 5.6 KB
RevLine 
[4687]1module lmdz_call_atke
[4545]2
[4687]3USE lmdz_atke_exchange_coeff, ONLY :  atke_compute_km_kh
[4545]4
5implicit none
6
7
8contains
9
10subroutine call_atke(dtime,ngrid,nlay,cdrag_uv,cdrag_t,u_surf,v_surf,temp_surf, &
[4653]11                        wind_u,wind_v,temp,qvap,play,pinterf, &
[4881]12                        tke,eps,Km_out,Kh_out)
[4545]13
14
15
16
[4687]17USE lmdz_atke_turbulence_ini, ONLY : iflag_num_atke, rg, rd
[4545]18
19implicit none
20
21
22! Declarations:
23!=============
24
25REAL, INTENT(IN)    :: dtime ! timestep (s)
26INTEGER, INTENT(IN) :: ngrid ! number of horizontal index (flat grid)
27INTEGER, INTENT(IN) :: nlay ! number of vertical index 
28
29
30REAL, DIMENSION(ngrid), INTENT(IN)       :: cdrag_uv  ! drag coefficient for wind
31REAL, DIMENSION(ngrid), INTENT(IN)       :: cdrag_t   ! drag coefficient for temperature
32
33REAL, DIMENSION(ngrid), INTENT(IN)       :: u_surf    ! x wind velocity at the surface
34REAL, DIMENSION(ngrid), INTENT(IN)       :: v_surf    ! y wind velocity at the surface
35REAL, DIMENSION(ngrid), INTENT(IN)       :: temp_surf ! surface temperature
36
37REAL, DIMENSION(ngrid,nlay), INTENT(IN)       :: wind_u   ! zonal velocity (m/s)
38REAL, DIMENSION(ngrid,nlay), INTENT(IN)       :: wind_v   ! meridional velocity (m/s)
39REAL, DIMENSION(ngrid,nlay), INTENT(IN)       :: temp   ! temperature (K)
[4653]40REAL, DIMENSION(ngrid,nlay), INTENT(IN)       :: qvap   ! specific humidity (kg/kg)
[4545]41REAL, DIMENSION(ngrid,nlay), INTENT(IN)       :: play   ! pressure (Pa)
42REAL, DIMENSION(ngrid,nlay+1), INTENT(IN)     :: pinterf   ! pressure at interfaces(Pa)
43
44
[4881]45REAL, DIMENSION(ngrid,nlay+1), INTENT(INOUT)  :: tke      ! turbulent kinetic energy at interface between layers
[4545]46
[4881]47REAL, DIMENSION(ngrid,nlay+1), INTENT(OUT)    :: eps      ! output: tke dissipation rate at interface between layers
[4545]48REAL, DIMENSION(ngrid,nlay), INTENT(OUT)      :: Km_out   ! output: Exchange coefficient for momentum at interface between layers
49REAL, DIMENSION(ngrid,nlay), INTENT(OUT)      :: Kh_out   ! output: Exchange coefficient for heat flux at interface between layers
50
51
52REAL, DIMENSION(ngrid,nlay) :: wind_u_predict, wind_v_predict
53REAL, DIMENSION(ngrid) ::  wind1
54INTEGER i
55
56
[4631]57call atke_compute_km_kh(ngrid,nlay,dtime,&
[4653]58                        wind_u,wind_v,temp,qvap,play,pinterf,cdrag_uv,&
[4881]59                        tke,eps,Km_out,Kh_out)
[4545]60
[4644]61
62                               
[4545]63if (iflag_num_atke .EQ. 1) then
[4644]64   !! In this case, we make an explicit prediction of the wind shear to calculate the tke in a
65   !! forward backward way
[4545]66   !! pay attention that the treatment of the TKE
67   !! has to be adapted when solving the TKE with a prognostic equation
68
69   do i=1,ngrid
70      wind1(i)=sqrt(wind_u(i,1)**2+wind_v(i,1)**2)
71   enddo
72   call atke_explicit_prediction(ngrid,nlay,rg,rd,dtime,pinterf,play,temp,wind1,wind_u,Km_out,u_surf,cdrag_uv,wind_u_predict) 
73   call atke_explicit_prediction(ngrid,nlay,rg,rd,dtime,pinterf,play,temp,wind1,wind_v,Km_out,v_surf,cdrag_uv,wind_v_predict)
74
75
[4631]76   call atke_compute_km_kh(ngrid,nlay,dtime,&
[4653]77                        wind_u_predict,wind_v_predict,temp,qvap,play,pinterf,cdrag_uv, &
[4881]78                        tke,eps,Km_out,Kh_out)
[4545]79
80end if
81
[4644]82
83
84
85
86
[4545]87end subroutine call_atke
88
89!----------------------------------------------------------------------------------------
90
91subroutine atke_explicit_prediction(ngrid,nlay,rg,rd,dtime,pinterf,play,temp,wind1,x_in,K_in,x_surf,cdrag,x_predict)
92
93INTEGER, INTENT(IN) :: ngrid ! number of horizontal index (flat grid)
94INTEGER, INTENT(IN) :: nlay ! number of vertical index 
95REAL, INTENT(IN) :: rg,rd,dtime ! gravity, R dry air and timestep
96REAL, DIMENSION(ngrid,nlay),   INTENT(IN)  :: play      ! pressure middle of layers (Pa)
97REAL, DIMENSION(ngrid,nlay),   INTENT(IN)  :: temp      ! temperature (K)
98REAL, DIMENSION(ngrid),        INTENT(IN)  :: wind1     ! wind speed first level (m/s)
99REAL, DIMENSION(ngrid,nlay+1), INTENT(IN)  :: pinterf ! pressure at interfaces(Pa)
100REAL, DIMENSION(ngrid,nlay),   INTENT(IN)  :: x_in    ! variable at the beginning of timestep
101REAL, DIMENSION(ngrid,nlay+1), INTENT(IN)  :: K_in    ! eddy diffusivity coef at the beginning of time step
102REAL, DIMENSION(ngrid),        INTENT(IN)  :: x_surf  ! surface variable
103REAL, DIMENSION(ngrid),        INTENT(IN)  :: cdrag  ! drag coefficient
104REAL, DIMENSION(ngrid,nlay),   INTENT(OUT) :: x_predict  ! variable at the end of time step after explicit prediction
105
106
107integer i,k
108real ml,F1,rho
109real, dimension(ngrid) :: play1,temp1
110real, dimension(ngrid,nlay+1) :: K_big
111
112! computation of K_big
113
114play1(:)=play(:,1)
115temp1(:)=temp(:,1)
116
117! "big K" calculation
118do  k=2,nlay-1
119   do i=1,ngrid
120      rho=pinterf(i,k)/rd/(0.5*(temp(i,k-1)+temp(i,k)))
121      K_big(i,k)=rg*K_in(i,k)/(play(i,k)-play(i,k+1))*(rho**2)
122   enddo
123enddo
124! speficic treatment for k=nlay
125do i=1,ngrid
126   rho=pinterf(i,nlay)/rd/temp(i,nlay)
127   K_big(i,nlay)=rg*K_in(i,nlay)/(2*(play(i,nlay)-pinterf(i,nlay+1)))*(rho**2)
128enddo
129
130
131
132! x_predict calculation for 2<=k<=nlay-1
133do  k=2,nlay-1
134   do i=1,ngrid
135      ml=(pinterf(i,k)-pinterf(i,k+1))/rg
136      x_predict(i,k)=x_in(i,k)-dtime/ml*(-K_big(i,k+1)*x_in(i,k+1) &
137                  + (K_big(i,k)+K_big(i,k+1))*x_in(i,k) &
138                  - K_big(i,k)*x_in(i,k-1))
139   enddo
140enddo
141
142! Specific treatment for k=1
143do i=1,ngrid
144   ml=(pinterf(i,1)-pinterf(i,2))/rg
145   F1=-play1(i)/rd/temp1(i)*wind1(i)*cdrag(i)*(x_in(i,1)-x_surf(i)) ! attention convention sens du flux
146   x_predict(i,1)=x_in(i,1)-dtime/ml*(-K_big(i,2)*(x_in(i,2) - x_in(i,1))-F1)
147enddo
148
149! Specific treatment for k=nlay
150! flux at the top of the atmosphere=0
151do i=1,ngrid
152   ml=0.5*(pinterf(i,nlay)-pinterf(i,nlay+1))/rg
153   x_predict(i,nlay)=x_in(i,nlay)+dtime/ml*(K_big(i,nlay)*(x_in(i,nlay)-x_in(i,nlay-1)))
154enddo
155
156
157end subroutine atke_explicit_prediction
158
159
160
[4687]161end module lmdz_call_atke
Note: See TracBrowser for help on using the repository browser.