source: codes/icosagcm/trunk/src/caldyn_kernels.f90 @ 408

Last change on this file since 408 was 362, checked in by dubos, 9 years ago

Introduced DEC convention for velocity - HEVI only - cleanup to follow

File size: 16.0 KB
Line 
1MODULE caldyn_kernels_mod
2  USE icosa
3  USE transfert_mod
4  USE caldyn_kernels_base_mod
5  IMPLICIT NONE
6  PRIVATE
7
8  PUBLIC :: compute_planetvel, compute_pvort, compute_geopot, &
9       compute_caldyn_horiz, compute_caldyn_vert
10CONTAINS
11
12  SUBROUTINE compute_planetvel(planetvel)
13    USE wind_mod
14    REAL(rstd),INTENT(OUT)  :: planetvel(iim*3*jjm)
15    REAL(rstd) :: ulon(iim*3*jjm)
16    REAL(rstd) :: ulat(iim*3*jjm)
17    REAL(rstd) :: lon,lat
18    INTEGER :: ij
19    DO ij=ij_begin_ext,ij_end_ext
20       ulon(ij+u_right)=radius*omega*cos(lat_e(ij+u_right))
21       ulat(ij+u_right)=0
22
23       ulon(ij+u_lup)=radius*omega*cos(lat_e(ij+u_lup))
24       ulat(ij+u_lup)=0
25
26       ulon(ij+u_ldown)=radius*omega*cos(lat_e(ij+u_ldown))
27       ulat(ij+u_ldown)=0
28    END DO
29    CALL compute_wind2D_perp_from_lonlat_compound(ulon, ulat, planetvel)
30  END SUBROUTINE compute_planetvel
31
32  SUBROUTINE compute_pvort(ps,u,theta_rhodz, rhodz,theta,qu,qv)
33    USE icosa
34    USE disvert_mod, ONLY : mass_dak, mass_dbk, caldyn_eta, eta_mass
35    USE exner_mod
36    USE trace
37    USE omp_para
38    IMPLICIT NONE
39    REAL(rstd),INTENT(IN)  :: u(iim*3*jjm,llm)
40    REAL(rstd),INTENT(IN)  :: ps(iim*jjm)
41    REAL(rstd),INTENT(IN)  :: theta_rhodz(iim*jjm,llm)
42    REAL(rstd),INTENT(INOUT) :: rhodz(iim*jjm,llm)
43    REAL(rstd),INTENT(OUT) :: theta(iim*jjm,llm)
44    REAL(rstd),INTENT(OUT) :: qu(iim*3*jjm,llm)
45    REAL(rstd),INTENT(OUT) :: qv(iim*2*jjm,llm)
46
47    INTEGER :: i,j,ij,l
48    REAL(rstd) :: etav,hv, m
49    CALL trace_start("compute_pvort") 
50
51    IF(caldyn_eta==eta_mass) THEN
52       CALL wait_message(req_ps)  ! COM00
53    ELSE
54       CALL wait_message(req_mass) ! COM00
55    END IF
56    CALL wait_message(req_theta_rhodz) ! COM01
57
58    IF(caldyn_eta==eta_mass) THEN ! Compute mass & theta
59       DO l = ll_begin,ll_end
60          CALL test_message(req_u) 
61          !DIR$ SIMD
62          DO ij=ij_begin_ext,ij_end_ext
63             m = ( mass_dak(l)+ps(ij)*mass_dbk(l) )/g
64             rhodz(ij,l) = m
65             theta(ij,l) = theta_rhodz(ij,l)/rhodz(ij,l)
66          ENDDO
67       ENDDO
68    ELSE ! Compute only theta
69       DO l = ll_begin,ll_end
70          CALL test_message(req_u) 
71          !DIR$ SIMD
72          DO ij=ij_begin_ext,ij_end_ext
73             theta(ij,l) = theta_rhodz(ij,l)/rhodz(ij,l)
74          ENDDO
75       ENDDO
76    END IF
77
78    CALL wait_message(req_u) ! COM02
79
80!!! Compute shallow-water potential vorticity
81    DO l = ll_begin,ll_end
82       !DIR$ SIMD
83       DO ij=ij_begin_ext,ij_end_ext
84          etav= 1./Av(ij+z_up)*(  ne_rup        * u(ij+u_rup,l)        * de(ij+u_rup)         &
85               + ne_left * u(ij+t_rup+u_left,l) * de(ij+t_rup+u_left)  &
86               - ne_lup        * u(ij+u_lup,l)        * de(ij+u_lup) )                               
87
88          hv =  Riv2(ij,vup)          * rhodz(ij,l)            &
89               + Riv2(ij+t_rup,vldown) * rhodz(ij+t_rup,l)     &
90               + Riv2(ij+t_lup,vrdown) * rhodz(ij+t_lup,l)
91
92          qv(ij+z_up,l) = ( etav+fv(ij+z_up) )/hv
93
94          etav = 1./Av(ij+z_down)*(  ne_ldown         * u(ij+u_ldown,l)          * de(ij+u_ldown)          &
95               + ne_right * u(ij+t_ldown+u_right,l)  * de(ij+t_ldown+u_right)  &
96               - ne_rdown         * u(ij+u_rdown,l)          * de(ij+u_rdown) )
97
98          hv =  Riv2(ij,vdown)        * rhodz(ij,l)          &
99               + Riv2(ij+t_ldown,vrup) * rhodz(ij+t_ldown,l)  &
100               + Riv2(ij+t_rdown,vlup) * rhodz(ij+t_rdown,l)
101
102          qv(ij+z_down,l) =( etav+fv(ij+z_down) )/hv
103
104       ENDDO
105
106       !DIR$ SIMD
107       DO ij=ij_begin,ij_end
108          qu(ij+u_right,l) = 0.5*(qv(ij+z_rdown,l)+qv(ij+z_rup,l)) 
109          qu(ij+u_lup,l) = 0.5*(qv(ij+z_up,l)+qv(ij+z_lup,l)) 
110          qu(ij+u_ldown,l) = 0.5*(qv(ij+z_ldown,l)+qv(ij+z_down,l)) 
111       END DO
112
113    ENDDO
114
115    CALL trace_end("compute_pvort")
116  END SUBROUTINE compute_pvort
117
118  !************************* caldyn_horiz = caldyn_fast + caldyn_slow **********************
119
120  SUBROUTINE compute_caldyn_horiz(u,rhodz,qu,theta,pk,geopot, hflux,convm, dtheta_rhodz, du)
121    USE icosa
122    USE disvert_mod
123    USE exner_mod
124    USE trace
125    USE omp_para
126    IMPLICIT NONE
127    REAL(rstd),INTENT(IN)  :: u(iim*3*jjm,llm)    ! prognostic "velocity"
128    REAL(rstd),INTENT(IN)  :: rhodz(iim*jjm,llm)
129    REAL(rstd),INTENT(IN)  :: qu(iim*3*jjm,llm)
130    REAL(rstd),INTENT(IN)  :: theta(iim*jjm,llm)  ! potential temperature
131    REAL(rstd),INTENT(INOUT) :: pk(iim*jjm,llm) ! Exner function
132    REAL(rstd),INTENT(IN)  :: geopot(iim*jjm,llm+1)    ! geopotential
133
134    REAL(rstd),INTENT(OUT) :: hflux(iim*3*jjm,llm) ! hflux in kg/s
135    REAL(rstd),INTENT(OUT) :: convm(iim*jjm,llm)  ! mass flux convergence
136    REAL(rstd),INTENT(OUT) :: dtheta_rhodz(iim*jjm,llm)
137    REAL(rstd),INTENT(OUT) :: du(iim*3*jjm,llm)
138
139    REAL(rstd) :: cor_NT(iim*jjm,llm)  ! NT coriolis force u.(du/dPhi)
140    REAL(rstd) :: urel(3*iim*jjm,llm)  ! relative velocity
141    REAL(rstd) :: Ftheta(3*iim*jjm,llm) ! theta flux
142    REAL(rstd) :: berni(iim*jjm,llm)  ! Bernoulli function
143
144    INTEGER :: i,j,ij,l
145    REAL(rstd) :: ww,uu
146
147    CALL trace_start("compute_caldyn_horiz")
148
149    !    CALL wait_message(req_theta_rhodz)
150
151    DO l = ll_begin, ll_end
152!!!  Compute mass and theta fluxes
153       IF (caldyn_conserv==energy) CALL test_message(req_qu) 
154       !DIR$ SIMD
155       DO ij=ij_begin_ext,ij_end_ext         
156          hflux(ij+u_right,l)=0.5*(rhodz(ij,l)+rhodz(ij+t_right,l))*u(ij+u_right,l)*le(ij+u_right)
157          hflux(ij+u_lup,l)=0.5*(rhodz(ij,l)+rhodz(ij+t_lup,l))*u(ij+u_lup,l)*le(ij+u_lup)
158          hflux(ij+u_ldown,l)=0.5*(rhodz(ij,l)+rhodz(ij+t_ldown,l))*u(ij+u_ldown,l)*le(ij+u_ldown)
159
160          Ftheta(ij+u_right,l)=0.5*(theta(ij,l)+theta(ij+t_right,l))*hflux(ij+u_right,l)
161          Ftheta(ij+u_lup,l)=0.5*(theta(ij,l)+theta(ij+t_lup,l))*hflux(ij+u_lup,l)
162          Ftheta(ij+u_ldown,l)=0.5*(theta(ij,l)+theta(ij+t_ldown,l))*hflux(ij+u_ldown,l)
163       ENDDO
164
165!!! compute horizontal divergence of fluxes
166       !DIR$ SIMD
167       DO ij=ij_begin,ij_end         
168          ! convm = -div(mass flux), sign convention as in Ringler et al. 2012, eq. 21
169          convm(ij,l)= -1./Ai(ij)*(ne_right*hflux(ij+u_right,l)   +  &
170               ne_rup*hflux(ij+u_rup,l)       +  & 
171               ne_lup*hflux(ij+u_lup,l)       +  & 
172               ne_left*hflux(ij+u_left,l)     +  & 
173               ne_ldown*hflux(ij+u_ldown,l)   +  & 
174               ne_rdown*hflux(ij+u_rdown,l))   
175
176          ! signe ? attention d (rho theta dz)
177          ! dtheta_rhodz =  -div(flux.theta)
178          dtheta_rhodz(ij,l)=-1./Ai(ij)*(ne_right*Ftheta(ij+u_right,l)    +  &
179               ne_rup*Ftheta(ij+u_rup,l)        +  & 
180               ne_lup*Ftheta(ij+u_lup,l)        +  & 
181               ne_left*Ftheta(ij+u_left,l)      +  & 
182               ne_ldown*Ftheta(ij+u_ldown,l)    +  & 
183               ne_rdown*Ftheta(ij+u_rdown,l))
184       ENDDO
185
186    END DO
187
188!!! Compute potential vorticity (Coriolis) contribution to du
189
190    SELECT CASE(caldyn_conserv)
191    CASE(energy) ! energy-conserving TRiSK
192
193       CALL wait_message(req_qu) ! COM03
194
195       DO l=ll_begin,ll_end
196          !DIR$ SIMD
197          DO ij=ij_begin,ij_end         
198
199             uu = wee(ij+u_right,1,1)*hflux(ij+u_rup,l)*(qu(ij+u_right,l)+qu(ij+u_rup,l))+                             &
200                  wee(ij+u_right,2,1)*hflux(ij+u_lup,l)*(qu(ij+u_right,l)+qu(ij+u_lup,l))+                             &
201                  wee(ij+u_right,3,1)*hflux(ij+u_left,l)*(qu(ij+u_right,l)+qu(ij+u_left,l))+                           &
202                  wee(ij+u_right,4,1)*hflux(ij+u_ldown,l)*(qu(ij+u_right,l)+qu(ij+u_ldown,l))+                         &
203                  wee(ij+u_right,5,1)*hflux(ij+u_rdown,l)*(qu(ij+u_right,l)+qu(ij+u_rdown,l))+                         & 
204                  wee(ij+u_right,1,2)*hflux(ij+t_right+u_ldown,l)*(qu(ij+u_right,l)+qu(ij+t_right+u_ldown,l))+         &
205                  wee(ij+u_right,2,2)*hflux(ij+t_right+u_rdown,l)*(qu(ij+u_right,l)+qu(ij+t_right+u_rdown,l))+         &
206                  wee(ij+u_right,3,2)*hflux(ij+t_right+u_right,l)*(qu(ij+u_right,l)+qu(ij+t_right+u_right,l))+         &
207                  wee(ij+u_right,4,2)*hflux(ij+t_right+u_rup,l)*(qu(ij+u_right,l)+qu(ij+t_right+u_rup,l))+             &
208                  wee(ij+u_right,5,2)*hflux(ij+t_right+u_lup,l)*(qu(ij+u_right,l)+qu(ij+t_right+u_lup,l))
209             du(ij+u_right,l) = .5*uu/de(ij+u_right)
210
211             uu = wee(ij+u_lup,1,1)*hflux(ij+u_left,l)*(qu(ij+u_lup,l)+qu(ij+u_left,l)) +                        &
212                  wee(ij+u_lup,2,1)*hflux(ij+u_ldown,l)*(qu(ij+u_lup,l)+qu(ij+u_ldown,l)) +                      &
213                  wee(ij+u_lup,3,1)*hflux(ij+u_rdown,l)*(qu(ij+u_lup,l)+qu(ij+u_rdown,l)) +                      &
214                  wee(ij+u_lup,4,1)*hflux(ij+u_right,l)*(qu(ij+u_lup,l)+qu(ij+u_right,l)) +                      &
215                  wee(ij+u_lup,5,1)*hflux(ij+u_rup,l)*(qu(ij+u_lup,l)+qu(ij+u_rup,l)) +                          & 
216                  wee(ij+u_lup,1,2)*hflux(ij+t_lup+u_right,l)*(qu(ij+u_lup,l)+qu(ij+t_lup+u_right,l)) +          &
217                  wee(ij+u_lup,2,2)*hflux(ij+t_lup+u_rup,l)*(qu(ij+u_lup,l)+qu(ij+t_lup+u_rup,l)) +              &
218                  wee(ij+u_lup,3,2)*hflux(ij+t_lup+u_lup,l)*(qu(ij+u_lup,l)+qu(ij+t_lup+u_lup,l)) +              &
219                  wee(ij+u_lup,4,2)*hflux(ij+t_lup+u_left,l)*(qu(ij+u_lup,l)+qu(ij+t_lup+u_left,l)) +            &
220                  wee(ij+u_lup,5,2)*hflux(ij+t_lup+u_ldown,l)*(qu(ij+u_lup,l)+qu(ij+t_lup+u_ldown,l))
221             du(ij+u_lup,l) = .5*uu/de(ij+u_lup)
222
223
224             uu = wee(ij+u_ldown,1,1)*hflux(ij+u_rdown,l)*(qu(ij+u_ldown,l)+qu(ij+u_rdown,l)) +                      &
225                  wee(ij+u_ldown,2,1)*hflux(ij+u_right,l)*(qu(ij+u_ldown,l)+qu(ij+u_right,l)) +                      &
226                  wee(ij+u_ldown,3,1)*hflux(ij+u_rup,l)*(qu(ij+u_ldown,l)+qu(ij+u_rup,l)) +                          &
227                  wee(ij+u_ldown,4,1)*hflux(ij+u_lup,l)*(qu(ij+u_ldown,l)+qu(ij+u_lup,l)) +                          &
228                  wee(ij+u_ldown,5,1)*hflux(ij+u_left,l)*(qu(ij+u_ldown,l)+qu(ij+u_left,l)) +                        & 
229                  wee(ij+u_ldown,1,2)*hflux(ij+t_ldown+u_lup,l)*(qu(ij+u_ldown,l)+qu(ij+t_ldown+u_lup,l)) +          &
230                  wee(ij+u_ldown,2,2)*hflux(ij+t_ldown+u_left,l)*(qu(ij+u_ldown,l)+qu(ij+t_ldown+u_left,l)) +        &
231                  wee(ij+u_ldown,3,2)*hflux(ij+t_ldown+u_ldown,l)*(qu(ij+u_ldown,l)+qu(ij+t_ldown+u_ldown,l)) +      &
232                  wee(ij+u_ldown,4,2)*hflux(ij+t_ldown+u_rdown,l)*(qu(ij+u_ldown,l)+qu(ij+t_ldown+u_rdown,l)) +      &
233                  wee(ij+u_ldown,5,2)*hflux(ij+t_ldown+u_right,l)*(qu(ij+u_ldown,l)+qu(ij+t_ldown+u_right,l))
234             du(ij+u_ldown,l) = .5*uu/de(ij+u_ldown)
235
236          ENDDO
237       ENDDO
238
239    CASE(enstrophy) ! enstrophy-conserving TRiSK
240
241       DO l=ll_begin,ll_end
242          !DIR$ SIMD
243          DO ij=ij_begin,ij_end         
244
245             uu = wee(ij+u_right,1,1)*hflux(ij+u_rup,l)+                           &
246                  wee(ij+u_right,2,1)*hflux(ij+u_lup,l)+                           &
247                  wee(ij+u_right,3,1)*hflux(ij+u_left,l)+                          &
248                  wee(ij+u_right,4,1)*hflux(ij+u_ldown,l)+                         &
249                  wee(ij+u_right,5,1)*hflux(ij+u_rdown,l)+                         & 
250                  wee(ij+u_right,1,2)*hflux(ij+t_right+u_ldown,l)+                 &
251                  wee(ij+u_right,2,2)*hflux(ij+t_right+u_rdown,l)+                 &
252                  wee(ij+u_right,3,2)*hflux(ij+t_right+u_right,l)+                 &
253                  wee(ij+u_right,4,2)*hflux(ij+t_right+u_rup,l)+                   &
254                  wee(ij+u_right,5,2)*hflux(ij+t_right+u_lup,l)
255             du(ij+u_right,l) = qu(ij+u_right,l)*uu/de(ij+u_right)
256
257
258             uu = wee(ij+u_lup,1,1)*hflux(ij+u_left,l)+                        &
259                  wee(ij+u_lup,2,1)*hflux(ij+u_ldown,l)+                       &
260                  wee(ij+u_lup,3,1)*hflux(ij+u_rdown,l)+                       &
261                  wee(ij+u_lup,4,1)*hflux(ij+u_right,l)+                       &
262                  wee(ij+u_lup,5,1)*hflux(ij+u_rup,l)+                         & 
263                  wee(ij+u_lup,1,2)*hflux(ij+t_lup+u_right,l)+                 &
264                  wee(ij+u_lup,2,2)*hflux(ij+t_lup+u_rup,l)+                   &
265                  wee(ij+u_lup,3,2)*hflux(ij+t_lup+u_lup,l)+                   &
266                  wee(ij+u_lup,4,2)*hflux(ij+t_lup+u_left,l)+                  &
267                  wee(ij+u_lup,5,2)*hflux(ij+t_lup+u_ldown,l)
268             du(ij+u_lup,l) = qu(ij+u_lup,l)*uu/de(ij+u_lup)
269
270             uu = wee(ij+u_ldown,1,1)*hflux(ij+u_rdown,l)+                      &
271                  wee(ij+u_ldown,2,1)*hflux(ij+u_right,l)+                      &
272                  wee(ij+u_ldown,3,1)*hflux(ij+u_rup,l)+                        &
273                  wee(ij+u_ldown,4,1)*hflux(ij+u_lup,l)+                        &
274                  wee(ij+u_ldown,5,1)*hflux(ij+u_left,l)+                       & 
275                  wee(ij+u_ldown,1,2)*hflux(ij+t_ldown+u_lup,l)+                &
276                  wee(ij+u_ldown,2,2)*hflux(ij+t_ldown+u_left,l)+               &
277                  wee(ij+u_ldown,3,2)*hflux(ij+t_ldown+u_ldown,l)+              &
278                  wee(ij+u_ldown,4,2)*hflux(ij+t_ldown+u_rdown,l)+              &
279                  wee(ij+u_ldown,5,2)*hflux(ij+t_ldown+u_right,l)
280             du(ij+u_ldown,l) = qu(ij+u_ldown,l)*uu/de(ij+u_ldown)
281
282          ENDDO
283       ENDDO
284
285    CASE DEFAULT
286       STOP
287    END SELECT
288
289!!! Compute bernouilli term = Kinetic Energy + geopotential
290    IF(boussinesq) THEN
291       DO l=ll_begin,ll_end
292          !DIR$ SIMD
293          DO ij=ij_begin,ij_end         
294
295             berni(ij,l) = pk(ij,l) + &
296                  1/(4*Ai(ij))*(le(ij+u_right)*de(ij+u_right)*u(ij+u_right,l)**2 +    &
297                  le(ij+u_rup)*de(ij+u_rup)*u(ij+u_rup,l)**2 +          &
298                  le(ij+u_lup)*de(ij+u_lup)*u(ij+u_lup,l)**2 +          &
299                  le(ij+u_left)*de(ij+u_left)*u(ij+u_left,l)**2 +       &
300                  le(ij+u_ldown)*de(ij+u_ldown)*u(ij+u_ldown,l)**2 +    &
301                  le(ij+u_rdown)*de(ij+u_rdown)*u(ij+u_rdown,l)**2 ) 
302             ! from now on pk contains the vertically-averaged geopotential
303             pk(ij,l) = .5*(geopot(ij,l)+geopot(ij,l+1))
304          ENDDO
305       ENDDO
306
307    ELSE ! compressible
308
309       DO l=ll_begin,ll_end
310          !DIR$ SIMD
311          DO ij=ij_begin,ij_end         
312
313             berni(ij,l) = .5*(geopot(ij,l)+geopot(ij,l+1)) &
314                  + 1/(4*Ai(ij))*(le(ij+u_right)*de(ij+u_right)*u(ij+u_right,l)**2 +    &
315                  le(ij+u_rup)*de(ij+u_rup)*u(ij+u_rup,l)**2 +          &
316                  le(ij+u_lup)*de(ij+u_lup)*u(ij+u_lup,l)**2 +          &
317                  le(ij+u_left)*de(ij+u_left)*u(ij+u_left,l)**2 +       &
318                  le(ij+u_ldown)*de(ij+u_ldown)*u(ij+u_ldown,l)**2 +    &
319                  le(ij+u_rdown)*de(ij+u_rdown)*u(ij+u_rdown,l)**2 ) 
320          ENDDO
321       ENDDO
322
323    END IF ! Boussinesq/compressible
324
325!!! Add gradients of Bernoulli and Exner functions to du
326    DO l=ll_begin,ll_end
327       !DIR$ SIMD
328       DO ij=ij_begin,ij_end         
329
330          du(ij+u_right,l) = du(ij+u_right,l) + 1/de(ij+u_right) * (             &
331               0.5*(theta(ij,l)+theta(ij+t_right,l))                &
332               *( ne_right*pk(ij,l)+ne_left*pk(ij+t_right,l))    &
333               + ne_right*berni(ij,l)+ne_left*berni(ij+t_right,l) )
334
335
336          du(ij+u_lup,l) = du(ij+u_lup,l) +  1/de(ij+u_lup) * (                  &
337               0.5*(theta(ij,l)+theta(ij+t_lup,l))                  &
338               *( ne_lup*pk(ij,l)+ne_rdown*pk(ij+t_lup,l))       &
339               + ne_lup*berni(ij,l)+ne_rdown*berni(ij+t_lup,l) )
340
341          du(ij+u_ldown,l) = du(ij+u_ldown,l) + 1/de(ij+u_ldown) * (             &
342               0.5*(theta(ij,l)+theta(ij+t_ldown,l))                &
343               *( ne_ldown*pk(ij,l)+ne_rup*pk(ij+t_ldown,l))     &
344               + ne_ldown*berni(ij,l)+ne_rup*berni(ij+t_ldown,l) )
345
346       ENDDO
347    ENDDO
348
349    CALL trace_end("compute_caldyn_horiz")
350
351  END SUBROUTINE compute_caldyn_horiz
352
353END MODULE caldyn_kernels_mod
Note: See TracBrowser for help on using the repository browser.