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

Last change on this file since 4745 was 4745, checked in by evignon, 6 months ago

nettoyage et corrections dans les routines atke pour utilisation en 3D (terre + mars)

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