source: trunk/LMDZ.MARS/libf/phymars/lmdz_call_atke.F90 @ 3307

Last change on this file since 3307 was 3168, checked in by llange, 2 years ago

Mars PCM
Following previous commit, add the functions for ATKE turbulence in the svn brench (yes I forgot the svn add...)
LL

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