1 | module lmdz_call_atke |
---|
2 | |
---|
3 | USE lmdz_atke_exchange_coeff, ONLY : atke_compute_km_kh |
---|
4 | |
---|
5 | implicit none |
---|
6 | |
---|
7 | |
---|
8 | contains |
---|
9 | |
---|
10 | subroutine 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 | |
---|
17 | USE lmdz_atke_turbulence_ini, ONLY : iflag_num_atke, rg, rd |
---|
18 | |
---|
19 | implicit none |
---|
20 | |
---|
21 | |
---|
22 | ! Declarations: |
---|
23 | !============= |
---|
24 | |
---|
25 | REAL, INTENT(IN) :: dtime ! timestep (s) |
---|
26 | INTEGER, INTENT(IN) :: ngrid ! number of horizontal index (flat grid) |
---|
27 | INTEGER, INTENT(IN) :: nlay ! number of vertical index |
---|
28 | |
---|
29 | |
---|
30 | REAL, DIMENSION(ngrid), INTENT(IN) :: cdrag_uv ! drag coefficient for wind |
---|
31 | REAL, DIMENSION(ngrid), INTENT(IN) :: cdrag_t ! drag coefficient for temperature |
---|
32 | |
---|
33 | REAL, DIMENSION(ngrid), INTENT(IN) :: u_surf ! x wind velocity at the surface |
---|
34 | REAL, DIMENSION(ngrid), INTENT(IN) :: v_surf ! y wind velocity at the surface |
---|
35 | REAL, DIMENSION(ngrid), INTENT(IN) :: temp_surf ! surface temperature |
---|
36 | |
---|
37 | REAL, DIMENSION(ngrid,nlay), INTENT(IN) :: wind_u ! zonal velocity (m/s) |
---|
38 | REAL, DIMENSION(ngrid,nlay), INTENT(IN) :: wind_v ! meridional velocity (m/s) |
---|
39 | REAL, DIMENSION(ngrid,nlay), INTENT(IN) :: temp ! temperature (K) |
---|
40 | REAL, DIMENSION(ngrid,nlay), INTENT(IN) :: qvap ! specific humidity (kg/kg) |
---|
41 | REAL, DIMENSION(ngrid,nlay), INTENT(IN) :: play ! pressure (Pa) |
---|
42 | REAL, DIMENSION(ngrid,nlay+1), INTENT(IN) :: pinterf ! pressure at interfaces(Pa) |
---|
43 | |
---|
44 | |
---|
45 | REAL, DIMENSION(ngrid,nlay+1), INTENT(INOUT) :: tke ! turbulent kinetic energy at interface between layers |
---|
46 | |
---|
47 | REAL, DIMENSION(ngrid,nlay), INTENT(OUT) :: Km_out ! output: Exchange coefficient for momentum at interface between layers |
---|
48 | REAL, DIMENSION(ngrid,nlay), INTENT(OUT) :: Kh_out ! output: Exchange coefficient for heat flux at interface between layers |
---|
49 | |
---|
50 | |
---|
51 | REAL, DIMENSION(ngrid,nlay) :: wind_u_predict, wind_v_predict |
---|
52 | REAL, DIMENSION(ngrid) :: wind1 |
---|
53 | INTEGER i |
---|
54 | |
---|
55 | |
---|
56 | call 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 | |
---|
62 | if (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 | |
---|
79 | end if |
---|
80 | |
---|
81 | |
---|
82 | |
---|
83 | |
---|
84 | |
---|
85 | |
---|
86 | end subroutine call_atke |
---|
87 | |
---|
88 | !---------------------------------------------------------------------------------------- |
---|
89 | |
---|
90 | subroutine atke_explicit_prediction(ngrid,nlay,rg,rd,dtime,pinterf,play,temp,wind1,x_in,K_in,x_surf,cdrag,x_predict) |
---|
91 | |
---|
92 | INTEGER, INTENT(IN) :: ngrid ! number of horizontal index (flat grid) |
---|
93 | INTEGER, INTENT(IN) :: nlay ! number of vertical index |
---|
94 | REAL, INTENT(IN) :: rg,rd,dtime ! gravity, R dry air and timestep |
---|
95 | REAL, DIMENSION(ngrid,nlay), INTENT(IN) :: play ! pressure middle of layers (Pa) |
---|
96 | REAL, DIMENSION(ngrid,nlay), INTENT(IN) :: temp ! temperature (K) |
---|
97 | REAL, DIMENSION(ngrid), INTENT(IN) :: wind1 ! wind speed first level (m/s) |
---|
98 | REAL, DIMENSION(ngrid,nlay+1), INTENT(IN) :: pinterf ! pressure at interfaces(Pa) |
---|
99 | REAL, DIMENSION(ngrid,nlay), INTENT(IN) :: x_in ! variable at the beginning of timestep |
---|
100 | REAL, DIMENSION(ngrid,nlay+1), INTENT(IN) :: K_in ! eddy diffusivity coef at the beginning of time step |
---|
101 | REAL, DIMENSION(ngrid), INTENT(IN) :: x_surf ! surface variable |
---|
102 | REAL, DIMENSION(ngrid), INTENT(IN) :: cdrag ! drag coefficient |
---|
103 | REAL, DIMENSION(ngrid,nlay), INTENT(OUT) :: x_predict ! variable at the end of time step after explicit prediction |
---|
104 | |
---|
105 | |
---|
106 | integer i,k |
---|
107 | real ml,F1,rho |
---|
108 | real, dimension(ngrid) :: play1,temp1 |
---|
109 | real, dimension(ngrid,nlay+1) :: K_big |
---|
110 | |
---|
111 | ! computation of K_big |
---|
112 | |
---|
113 | play1(:)=play(:,1) |
---|
114 | temp1(:)=temp(:,1) |
---|
115 | |
---|
116 | ! "big K" calculation |
---|
117 | do 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 |
---|
122 | enddo |
---|
123 | ! speficic treatment for k=nlay |
---|
124 | do 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) |
---|
127 | enddo |
---|
128 | |
---|
129 | |
---|
130 | |
---|
131 | ! x_predict calculation for 2<=k<=nlay-1 |
---|
132 | do 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 |
---|
139 | enddo |
---|
140 | |
---|
141 | ! Specific treatment for k=1 |
---|
142 | do 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) |
---|
146 | enddo |
---|
147 | |
---|
148 | ! Specific treatment for k=nlay |
---|
149 | ! flux at the top of the atmosphere=0 |
---|
150 | do 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))) |
---|
153 | enddo |
---|
154 | |
---|
155 | |
---|
156 | end subroutine atke_explicit_prediction |
---|
157 | |
---|
158 | |
---|
159 | |
---|
160 | end module lmdz_call_atke |
---|