source: codes/icosagcm/trunk/src/caldyn_kernels_hevi.f90 @ 521

Last change on this file since 521 was 521, checked in by dubos, 7 years ago

Cleanup after fixing RK2.5

File size: 35.5 KB
Line 
1MODULE caldyn_kernels_hevi_mod
2  USE icosa
3  USE trace
4  USE omp_para
5  USE disvert_mod
6  USE transfert_mod
7  USE caldyn_kernels_base_mod
8  IMPLICIT NONE
9  PRIVATE
10
11  REAL(rstd), PARAMETER :: pbot=1e5, Phi_bot=0., rho_bot=1e6 ! FIXME
12
13  LOGICAL, PARAMETER :: debug_hevi_solver = .FALSE.
14
15  PUBLIC :: compute_theta, compute_pvort_only, compute_caldyn_Coriolis, &
16       compute_caldyn_slow_hydro, compute_caldyn_slow_NH, &
17       compute_caldyn_solver, compute_caldyn_fast
18
19CONTAINS
20
21  SUBROUTINE compute_theta(ps,theta_rhodz, rhodz,theta)
22    REAL(rstd),INTENT(IN)    :: ps(iim*jjm)
23    REAL(rstd),INTENT(IN)    :: theta_rhodz(iim*jjm,llm,nqdyn)
24    REAL(rstd),INTENT(INOUT) :: rhodz(iim*jjm,llm)
25    REAL(rstd),INTENT(OUT)   :: theta(iim*jjm,llm,nqdyn)
26    INTEGER :: ij,l,iq
27    REAL(rstd) :: m
28    CALL trace_start("compute_theta")
29
30    IF(caldyn_eta==eta_mass) THEN ! Compute mass
31       DO l = ll_begin,ll_end
32          !DIR$ SIMD
33          DO ij=ij_begin_ext,ij_end_ext
34             m = mass_dak(l)+(ps(ij)*g+ptop)*mass_dbk(l) ! ps is actually Ms
35             rhodz(ij,l) = m/g
36          END DO
37       END DO
38    END IF
39
40    DO l = ll_begin,ll_end
41       DO iq=1,nqdyn
42          !DIR$ SIMD
43          DO ij=ij_begin_ext,ij_end_ext
44             theta(ij,l,iq) = theta_rhodz(ij,l,iq)/rhodz(ij,l)
45          END DO
46       END DO
47    END DO
48
49    CALL trace_end("compute_theta")
50  END SUBROUTINE compute_theta
51
52  SUBROUTINE compute_pvort_only(u,rhodz,qu,qv)
53    REAL(rstd),INTENT(IN)  :: u(iim*3*jjm,llm)
54    REAL(rstd),INTENT(INOUT) :: rhodz(iim*jjm,llm)
55    REAL(rstd),INTENT(OUT) :: qu(iim*3*jjm,llm)
56    REAL(rstd),INTENT(OUT) :: qv(iim*2*jjm,llm)
57
58    INTEGER :: ij,l
59    REAL(rstd) :: etav,hv,radius_m2
60
61    CALL trace_start("compute_pvort_only") 
62!!! Compute shallow-water potential vorticity
63    radius_m2=radius**(-2)
64    DO l = ll_begin,ll_end
65       !DIR$ SIMD
66       DO ij=ij_begin_ext,ij_end_ext
67          etav= 1./Av(ij+z_up)*(  ne_rup  * u(ij+u_rup,l)   &
68                  + ne_left * u(ij+t_rup+u_left,l)          &
69                  - ne_lup  * u(ij+u_lup,l) )                               
70          hv =   Riv2(ij,vup)          * rhodz(ij,l)        &
71               + Riv2(ij+t_rup,vldown) * rhodz(ij+t_rup,l)  &
72               + Riv2(ij+t_lup,vrdown) * rhodz(ij+t_lup,l)
73          qv(ij+z_up,l) = ( etav+fv(ij+z_up) )/hv
74         
75          etav = 1./Av(ij+z_down)*(  ne_ldown * u(ij+u_ldown,l)   &
76               + ne_right * u(ij+t_ldown+u_right,l)               &
77               - ne_rdown * u(ij+u_rdown,l) )
78          hv =   Riv2(ij,vdown)        * rhodz(ij,l)              &
79               + Riv2(ij+t_ldown,vrup) * rhodz(ij+t_ldown,l)      &
80               + Riv2(ij+t_rdown,vlup) * rhodz(ij+t_rdown,l)
81          qv(ij+z_down,l) =( etav+fv(ij+z_down) )/hv
82       ENDDO
83
84       !DIR$ SIMD
85       DO ij=ij_begin,ij_end
86          qu(ij+u_right,l) = 0.5*(qv(ij+z_rdown,l)+qv(ij+z_rup,l)) 
87          qu(ij+u_lup,l) = 0.5*(qv(ij+z_up,l)+qv(ij+z_lup,l)) 
88          qu(ij+u_ldown,l) = 0.5*(qv(ij+z_ldown,l)+qv(ij+z_down,l)) 
89       END DO
90
91    ENDDO
92
93    CALL trace_end("compute_pvort_only")
94
95  END SUBROUTINE compute_pvort_only
96
97  SUBROUTINE compute_NH_geopot(tau, m_ik, m_il, theta, W_il, Phi_il)
98    REAL(rstd),INTENT(IN)    :: tau ! solve Phi-tau*dPhi/dt = Phi_rhs
99    REAL(rstd),INTENT(IN)    :: m_ik(iim*jjm,llm)
100    REAL(rstd),INTENT(IN)    :: m_il(iim*jjm,llm+1)
101    REAL(rstd),INTENT(IN)    :: theta(iim*jjm,llm)
102    REAL(rstd),INTENT(IN)    :: W_il(iim*jjm,llm+1)   ! vertical momentum
103    REAL(rstd),INTENT(INOUT) :: Phi_il(iim*jjm,llm+1) ! geopotential
104
105    REAL(rstd) :: Phi_star_il(iim*jjm,llm+1) 
106    REAL(rstd) :: p_ik(iim*jjm,llm)      ! pressure
107    REAL(rstd) :: R_il(iim*jjm,llm+1)    ! rhs of tridiag problem
108    REAL(rstd) :: x_il(iim*jjm,llm+1)    ! solution of tridiag problem
109    REAL(rstd) :: A_ik(iim*jjm,llm)      ! off-diagonal coefficients of tridiag problem
110    REAL(rstd) :: B_il(iim*jjm,llm+1)    ! diagonal coefficients of tridiag problem
111    REAL(rstd) :: C_ik(iim*jjm,llm)      ! Thomas algorithm
112    REAL(rstd) :: D_il(iim*jjm,llm+1)    ! Thomas algorithm
113    REAL(rstd) :: gamma, rho_ij, X_ij, Y_ij
114    REAL(rstd) :: wil, tau2_g, g2, gm2, ml_g2, c2_mik
115
116    INTEGER    :: iter, ij, l
117
118!    FIXME : vertical OpenMP parallelism will not work
119   
120    tau2_g=tau*tau/g
121    g2=g*g
122    gm2 = g**-2
123    gamma = 1./(1.-kappa)
124   
125    ! compute Phi_star
126    DO l=1,llm+1
127       !DIR$ SIMD
128       DO ij=ij_begin_ext,ij_end_ext
129          Phi_star_il(ij,l) = Phi_il(ij,l) + tau*g2*(W_il(ij,l)/m_il(ij,l)-tau)
130       ENDDO
131    ENDDO
132
133    ! Newton-Raphson iteration : Phi_il contains current guess value
134    DO iter=1,5 ! 2 iterations should be enough
135
136       ! Compute pressure, A_ik
137       DO l=1,llm
138          !DIR$ SIMD
139          DO ij=ij_begin_ext,ij_end_ext
140             rho_ij = (g*m_ik(ij,l))/(Phi_il(ij,l+1)-Phi_il(ij,l))
141             X_ij = (cpp/preff)*kappa*theta(ij,l)*rho_ij
142             p_ik(ij,l) = preff*(X_ij**gamma)
143             c2_mik = gamma*p_ik(ij,l)/(rho_ij*m_ik(ij,l)) ! c^2 = gamma*R*T = gamma*p/rho
144             A_ik(ij,l) = c2_mik*(tau/g*rho_ij)**2
145          ENDDO
146       ENDDO
147
148       ! Compute residual, B_il
149       ! bottom interface l=1
150       !DIR$ SIMD
151       DO ij=ij_begin_ext,ij_end_ext
152          ml_g2 = gm2*m_il(ij,1) 
153          B_il(ij,1) = A_ik(ij,1) + ml_g2 + tau2_g*rho_bot
154          R_il(ij,1) = ml_g2*( Phi_il(ij,1)-Phi_star_il(ij,1))  &
155               + tau2_g*( p_ik(ij,1)-pbot+rho_bot*(Phi_il(ij,1)-Phi_bot) )
156       ENDDO
157       ! inner interfaces
158       DO l=2,llm
159          !DIR$ SIMD
160          DO ij=ij_begin_ext,ij_end_ext
161             ml_g2 = gm2*m_il(ij,l) 
162             B_il(ij,l) = A_ik(ij,l)+A_ik(ij,l-1) + ml_g2
163             R_il(ij,l) = ml_g2*( Phi_il(ij,l)-Phi_star_il(ij,l))  &
164                     + tau2_g*(p_ik(ij,l)-p_ik(ij,l-1))
165             ! consistency check : if Wil=0 and initial state is in hydrostatic balance
166             ! then Phi_star_il(ij,l) = Phi_il(ij,l) - tau^2*g^2
167             ! and residual = tau^2*(ml+(1/g)dl_pi)=0
168          ENDDO
169       ENDDO
170       ! top interface l=llm+1
171       !DIR$ SIMD
172       DO ij=ij_begin_ext,ij_end_ext
173          ml_g2 = gm2*m_il(ij,llm+1) 
174          B_il(ij,llm+1) = A_ik(ij,llm) + ml_g2
175          R_il(ij,llm+1) = ml_g2*( Phi_il(ij,llm+1)-Phi_star_il(ij,llm+1))  &
176               + tau2_g*( ptop-p_ik(ij,llm) )
177       ENDDO
178
179       ! FIXME later
180       ! the lines below modify the tridiag problem
181       ! for flat, rigid boundary conditions at top and bottom :
182       ! zero out A(1), A(llm), R(1), R(llm+1)
183       ! => x(l)=0 at l=1,llm+1
184       DO ij=ij_begin_ext,ij_end_ext
185          A_ik(ij,1) = 0.
186          A_ik(ij,llm) = 0.
187          R_il(ij,1) = 0.
188          R_il(ij,llm+1) = 0.
189       ENDDO
190
191       IF(debug_hevi_solver) THEN ! print Linf(residual)
192          PRINT *, '[hevi_solver] R,p', iter, MAXVAL(ABS(R_il)), MAXVAL(p_ik)
193       END IF
194
195       ! Solve -A(l-1)x(l-1) + B(l)x(l) - A(l)x(l+1) = R(l) using Thomas algorithm
196       ! Forward sweep :
197       ! C(0)=0, C(l) = -A(l) / (B(l)+A(l-1)C(l-1)),
198       ! D(0)=0, D(l) = (R(l)+A(l-1)D(l-1)) / (B(l)+A(l-1)C(l-1))
199       ! bottom interface l=1
200       !DIR$ SIMD
201       DO ij=ij_begin_ext,ij_end_ext
202          X_ij = 1./B_il(ij,1)
203          C_ik(ij,1) = -A_ik(ij,1) * X_ij 
204          D_il(ij,1) =  R_il(ij,1) * X_ij
205       ENDDO
206       ! inner interfaces/layers
207       DO l=2,llm
208          !DIR$ SIMD
209          DO ij=ij_begin_ext,ij_end_ext
210             X_ij = 1./(B_il(ij,l) + A_ik(ij,l-1)*C_ik(ij,l-1))
211             C_ik(ij,l) = -A_ik(ij,l) * X_ij 
212             D_il(ij,l) =  (R_il(ij,l)+A_ik(ij,l-1)*D_il(ij,l-1)) * X_ij
213          ENDDO
214       ENDDO
215       ! top interface l=llm+1
216       !DIR$ SIMD
217       DO ij=ij_begin_ext,ij_end_ext
218          X_ij = 1./(B_il(ij,llm+1) + A_ik(ij,llm)*C_ik(ij,llm))
219          D_il(ij,llm+1) =  (R_il(ij,llm+1)+A_ik(ij,llm)*D_il(ij,llm)) * X_ij
220       ENDDO
221       
222       ! Back substitution :
223       ! x(i) = D(i)-C(i)x(i+1), x(N+1)=0
224       ! + Newton-Raphson update
225       x_il=0. ! FIXME
226       ! top interface l=llm+1
227       !DIR$ SIMD
228       DO ij=ij_begin_ext,ij_end_ext
229          x_il(ij,llm+1) = D_il(ij,llm+1)
230          Phi_il(ij,llm+1) = Phi_il(ij,llm+1) - x_il(ij,llm+1)
231       ENDDO
232       ! lower interfaces
233       DO l=llm,1,-1
234          !DIR$ SIMD
235          DO ij=ij_begin_ext,ij_end_ext
236             x_il(ij,l) = D_il(ij,l) - C_ik(ij,l)*x_il(ij,l+1)
237             Phi_il(ij,l) = Phi_il(ij,l) - x_il(ij,l)
238          ENDDO
239       ENDDO
240
241       IF(debug_hevi_solver) THEN
242          PRINT *, '[hevi_solver] A,B', iter, MAXVAL(ABS(A_ik)),MAXVAL(ABS(B_il))
243          PRINT *, '[hevi_solver] C,D', iter, MAXVAL(ABS(C_ik)),MAXVAL(ABS(D_il))
244          DO l=1,llm+1
245             WRITE(*,'(A,I2.1,I3.2,E9.2)'), '[hevi_solver] x', iter,l, MAXVAL(ABS(x_il(:,l)))
246          END DO
247       END IF
248
249    END DO ! Newton-Raphson
250   
251  END SUBROUTINE compute_NH_geopot
252
253  SUBROUTINE compute_caldyn_solver(tau,rhodz,theta,pk, geopot,W, dPhi,dW,du)
254    REAL(rstd),INTENT(IN) :: tau ! "solve" Phi-tau*dPhi/dt = Phi_rhs
255    REAL(rstd),INTENT(IN)    :: rhodz(iim*jjm,llm)
256    REAL(rstd),INTENT(IN)    :: theta(iim*jjm,llm)
257    REAL(rstd),INTENT(OUT)   :: pk(iim*jjm,llm)
258    REAL(rstd),INTENT(INOUT) :: geopot(iim*jjm,llm+1)
259    REAL(rstd),INTENT(INOUT) :: W(iim*jjm,llm+1) ! OUT if tau>0
260    REAL(rstd),INTENT(OUT)   :: dW(iim*jjm,llm+1)
261    REAL(rstd),INTENT(OUT)   :: dPhi(iim*jjm,llm+1)
262    REAL(rstd),INTENT(OUT)   :: du(3*iim*jjm,llm)
263
264    REAL(rstd) :: m_il(iim*jjm,llm+1)        ! rhodz averaged to interfaces
265    REAL(rstd) :: berni(iim*jjm)             ! (W/m_il)^2
266    REAL(rstd) :: gamma, rho_ij, X_ij, Y_ij
267    INTEGER    :: ij, l
268
269    CALL trace_start("compute_caldyn_solver")
270
271    ! FIXME : vertical OpenMP parallelism will not work
272
273    ! average m_ik to interfaces => m_il
274    !DIR$ SIMD
275    DO ij=ij_begin_ext,ij_end_ext
276       m_il(ij,1) = .5*rhodz(ij,1)
277    ENDDO
278    DO l=2,llm
279       !DIR$ SIMD
280       DO ij=ij_begin_ext,ij_end_ext
281          m_il(ij,l) = .5*(rhodz(ij,l-1)+rhodz(ij,l))
282       ENDDO
283    ENDDO
284    !DIR$ SIMD
285    DO ij=ij_begin_ext,ij_end_ext
286       m_il(ij,llm+1) = .5*rhodz(ij,llm)
287    ENDDO
288
289    IF(tau>0) THEN ! solve implicit problem for geopotential
290       CALL compute_NH_geopot(tau, rhodz, m_il, theta, W, geopot)
291    END IF
292
293    ! Compute pressure, stored temporarily in pk
294    ! kappa = R/Cp
295    ! 1-kappa = Cv/Cp
296    ! Cp/Cv = 1/(1-kappa)
297    gamma = 1./(1.-kappa)
298    DO l=1,llm
299       !DIR$ SIMD
300       DO ij=ij_begin_ext,ij_end_ext
301          rho_ij = (g*rhodz(ij,l))/(geopot(ij,l+1)-geopot(ij,l))
302          X_ij = (cpp/preff)*kappa*theta(ij,l)*rho_ij
303          ! kappa.theta.rho = p/exner
304          ! => X = (p/p0)/(exner/Cp)
305          !      = (p/p0)^(1-kappa)
306          pk(ij,l) = preff*(X_ij**gamma)
307       ENDDO
308    ENDDO
309
310    ! Update W, compute tendencies
311    DO l=2,llm
312       !DIR$ SIMD
313       DO ij=ij_begin_ext,ij_end_ext
314          dW(ij,l)   = (1./g)*(pk(ij,l-1)-pk(ij,l)) - m_il(ij,l)
315          W(ij,l)    = W(ij,l)+tau*dW(ij,l) ! update W
316          dPhi(ij,l) = g*g*W(ij,l)/m_il(ij,l)
317       ENDDO
318       !       PRINT *,'Max dPhi', l,ij_begin,ij_end, MAXVAL(abs(dPhi(ij_begin:ij_end,l)))
319       !       PRINT *,'Max dW', l,ij_begin,ij_end, MAXVAL(abs(dW(ij_begin:ij_end,l)))
320    ENDDO
321    ! Lower BC (FIXME : no orography yet !)
322    DO ij=ij_begin,ij_end         
323       dPhi(ij,1)=0
324       W(ij,1)=0
325       dW(ij,1)=0
326       dPhi(ij,llm+1)=0 ! rigid lid
327       W(ij,llm+1)=0
328       dW(ij,llm+1)=0
329    ENDDO
330    ! Upper BC p=ptop
331    !    DO ij=ij_omp_begin_ext,ij_omp_end_ext
332    !       dPhi(ij,llm+1) = W(ij,llm+1)/rhodz(ij,llm)
333    !       dW(ij,llm+1) = (1./g)*(pk(ij,llm)-ptop) - .5*rhodz(ij,llm)
334    !    ENDDO
335   
336    ! Compute Exner function (needed by compute_caldyn_fast) and du=-g^2.grad(w^2)
337    DO l=1,llm
338       !DIR$ SIMD
339       DO ij=ij_begin_ext,ij_end_ext
340          pk(ij,l) = cpp*((pk(ij,l)/preff)**kappa) ! other formulae possible if exponentiation is slow
341          berni(ij) = (-.25*g*g)*(              &
342                 (W(ij,l)/m_il(ij,l))**2       &
343               + (W(ij,l+1)/m_il(ij,l+1))**2 )
344       ENDDO
345       DO ij=ij_begin,ij_end         
346          du(ij+u_right,l) = ne_right*(berni(ij)-berni(ij+t_right))
347          du(ij+u_lup,l)   = ne_lup  *(berni(ij)-berni(ij+t_lup))
348          du(ij+u_ldown,l) = ne_ldown*(berni(ij)-berni(ij+t_ldown))
349       ENDDO
350    ENDDO
351   
352    CALL trace_end("compute_caldyn_solver")
353   
354  END SUBROUTINE compute_caldyn_solver
355 
356  SUBROUTINE compute_caldyn_fast(tau,u,rhodz,theta,pk,geopot,du)
357    REAL(rstd),INTENT(IN)    :: tau                ! "solve" u-tau*du/dt = rhs
358    REAL(rstd),INTENT(INOUT) :: u(iim*3*jjm,llm)   ! OUT if tau>0
359    REAL(rstd),INTENT(IN)    :: rhodz(iim*jjm,llm)
360    REAL(rstd),INTENT(IN)    :: theta(iim*jjm,llm,nqdyn)
361    REAL(rstd),INTENT(INOUT) :: pk(iim*jjm,llm)
362    REAL(rstd),INTENT(INOUT) :: geopot(iim*jjm,llm+1)
363    REAL(rstd),INTENT(INOUT)   :: du(iim*3*jjm,llm)
364    REAL(rstd) :: berni(iim*jjm,llm)  ! Bernoulli function
365    REAL(rstd) :: berniv(iim*jjm,llm)  ! moist Bernoulli function
366
367    INTEGER :: i,j,ij,l
368    REAL(rstd) :: Rd, qv, temp, chi, nu, due_right, due_lup, due_ldown
369
370    CALL trace_start("compute_caldyn_fast")
371
372    Rd=cpp*kappa
373
374    ! Compute Bernoulli term
375    IF(boussinesq) THEN
376       DO l=ll_begin,ll_end
377          !DIR$ SIMD
378          DO ij=ij_begin,ij_end         
379             berni(ij,l) = pk(ij,l)
380             ! from now on pk contains the vertically-averaged geopotential
381             pk(ij,l) = .5*(geopot(ij,l)+geopot(ij,l+1))
382          END DO
383       END DO
384    ELSE ! compressible
385
386       DO l=ll_begin,ll_end
387          SELECT CASE(caldyn_thermo)
388          CASE(thermo_theta) ! vdp = theta.dpi => B = Phi
389             !DIR$ SIMD
390             DO ij=ij_begin,ij_end         
391                berni(ij,l) = .5*(geopot(ij,l)+geopot(ij,l+1))
392             END DO
393          CASE(thermo_entropy) ! vdp = dG + sdT => B = Phi + G, G=h-Ts=T*(cpp-s)
394             !DIR$ SIMD
395             DO ij=ij_begin,ij_end         
396                berni(ij,l) = .5*(geopot(ij,l)+geopot(ij,l+1)) &
397                     + pk(ij,l)*(cpp-theta(ij,l,1)) ! pk=temperature, theta=entropy
398             END DO
399          CASE(thermo_moist) 
400             !DIR$ SIMD
401             DO ij=ij_begin,ij_end         
402                ! du/dt = grad(Bd)+rv.grad(Bv)+s.grad(T)
403                ! Bd = Phi + gibbs_d
404                ! Bv = Phi + gibbs_v
405                ! pk=temperature, theta=entropy
406                qv = theta(ij,l,2)
407                temp = pk(ij,l)
408                chi = log(temp/Treff)
409                nu = (chi*(cpp+qv*cppv)-theta(ij,l,1))/(Rd+qv*Rv) ! log(p/preff)
410                berni(ij,l) = .5*(geopot(ij,l)+geopot(ij,l+1)) &
411                     + temp*(cpp*(1.-chi)+Rd*nu)
412                berniv(ij,l) = .5*(geopot(ij,l)+geopot(ij,l+1)) &
413                     + temp*(cppv*(1.-chi)+Rv*nu)
414            END DO
415          END SELECT
416       END DO
417
418    END IF ! Boussinesq/compressible
419
420!!! u:=u+tau*du, du = -grad(B)-theta.grad(pi)
421    DO l=ll_begin,ll_end
422       IF(caldyn_thermo == thermo_moist) THEN
423          !DIR$ SIMD
424          DO ij=ij_begin,ij_end         
425             due_right =  berni(ij+t_right,l)-berni(ij,l)       &
426                  + 0.5*(theta(ij,l,1)+theta(ij+t_right,l,1))   &
427                       *(pk(ij+t_right,l)-pk(ij,l))             &
428                  + 0.5*(theta(ij,l,2)+theta(ij+t_right,l,2))   &
429                       *(berniv(ij+t_right,l)-berniv(ij,l))
430             
431             due_lup = berni(ij+t_lup,l)-berni(ij,l)            &
432                  + 0.5*(theta(ij,l,1)+theta(ij+t_lup,l,1))     &
433                       *(pk(ij+t_lup,l)-pk(ij,l))               &
434                  + 0.5*(theta(ij,l,2)+theta(ij+t_lup,l,2))     &
435                       *(berniv(ij+t_lup,l)-berniv(ij,l))
436             
437             due_ldown = berni(ij+t_ldown,l)-berni(ij,l)        &
438                  + 0.5*(theta(ij,l,1)+theta(ij+t_ldown,l,1))   &
439                       *(pk(ij+t_ldown,l)-pk(ij,l))             &
440                  + 0.5*(theta(ij,l,2)+theta(ij+t_ldown,l,2))   &
441                       *(berniv(ij+t_ldown,l)-berniv(ij,l))
442             
443             du(ij+u_right,l) = du(ij+u_right,l) - ne_right*due_right
444             du(ij+u_lup,l)   = du(ij+u_lup,l)   - ne_lup*due_lup
445             du(ij+u_ldown,l) = du(ij+u_ldown,l) - ne_ldown*due_ldown
446             u(ij+u_right,l) = u(ij+u_right,l) + tau*du(ij+u_right,l)
447             u(ij+u_lup,l)   = u(ij+u_lup,l)   + tau*du(ij+u_lup,l)
448             u(ij+u_ldown,l) = u(ij+u_ldown,l) + tau*du(ij+u_ldown,l)
449          END DO
450       ELSE
451          !DIR$ SIMD
452          DO ij=ij_begin,ij_end         
453             due_right =  0.5*(theta(ij,l,1)+theta(ij+t_right,l,1))  &
454                  *(pk(ij+t_right,l)-pk(ij,l))        &
455                  +  berni(ij+t_right,l)-berni(ij,l)
456             due_lup = 0.5*(theta(ij,l,1)+theta(ij+t_lup,l,1))    &
457                  *(pk(ij+t_lup,l)-pk(ij,l))          &
458                  +  berni(ij+t_lup,l)-berni(ij,l)
459             due_ldown = 0.5*(theta(ij,l,1)+theta(ij+t_ldown,l,1)) &
460                  *(pk(ij+t_ldown,l)-pk(ij,l))      &
461                  +   berni(ij+t_ldown,l)-berni(ij,l)
462             du(ij+u_right,l) = du(ij+u_right,l) - ne_right*due_right
463             du(ij+u_lup,l)   = du(ij+u_lup,l)   - ne_lup*due_lup
464             du(ij+u_ldown,l) = du(ij+u_ldown,l) - ne_ldown*due_ldown
465             u(ij+u_right,l) = u(ij+u_right,l) + tau*du(ij+u_right,l)
466             u(ij+u_lup,l)   = u(ij+u_lup,l)   + tau*du(ij+u_lup,l)
467             u(ij+u_ldown,l) = u(ij+u_ldown,l) + tau*du(ij+u_ldown,l)
468          END DO
469       END IF
470    END DO
471       
472    CALL trace_end("compute_caldyn_fast")
473
474  END SUBROUTINE compute_caldyn_fast
475
476  SUBROUTINE compute_caldyn_Coriolis(hflux,theta,qu, convm,dtheta_rhodz,du)
477    REAL(rstd),INTENT(IN)  :: hflux(3*iim*jjm,llm)  ! hflux in kg/s
478    REAL(rstd),INTENT(IN)  :: theta(iim*jjm,llm,nqdyn) ! active scalars
479    REAL(rstd),INTENT(IN)  :: qu(3*iim*jjm,llm)
480    REAL(rstd),INTENT(OUT) :: convm(iim*jjm,llm)  ! mass flux convergence
481    REAL(rstd),INTENT(OUT) :: dtheta_rhodz(iim*jjm,llm,nqdyn)
482    REAL(rstd),INTENT(INOUT) :: du(3*iim*jjm,llm)
483   
484    REAL(rstd) :: Ftheta(3*iim*jjm)  ! potential temperature flux
485    REAL(rstd) :: uu_right, uu_lup, uu_ldown
486    INTEGER :: ij,iq,l,kdown
487
488    CALL trace_start("compute_caldyn_Coriolis")
489
490    DO l=ll_begin, ll_end
491       ! compute theta flux
492       DO iq=1,nqdyn
493       !DIR$ SIMD
494          DO ij=ij_begin_ext,ij_end_ext     
495             Ftheta(ij+u_right) = 0.5*(theta(ij,l,iq)+theta(ij+t_right,l,iq)) &
496                                  * hflux(ij+u_right,l)
497             Ftheta(ij+u_lup)   = 0.5*(theta(ij,l,iq)+theta(ij+t_lup,l,iq)) &
498                  * hflux(ij+u_lup,l)
499             Ftheta(ij+u_ldown) = 0.5*(theta(ij,l,iq)+theta(ij+t_ldown,l,iq)) &
500                  * hflux(ij+u_ldown,l)
501          END DO
502          ! horizontal divergence of fluxes
503       !DIR$ SIMD
504          DO ij=ij_begin,ij_end         
505             ! dtheta_rhodz =  -div(flux.theta)
506             dtheta_rhodz(ij,l,iq)= &
507                  -1./Ai(ij)*(ne_right*Ftheta(ij+u_right)    +  &
508                  ne_rup*Ftheta(ij+u_rup)        +  & 
509                  ne_lup*Ftheta(ij+u_lup)        +  & 
510                  ne_left*Ftheta(ij+u_left)      +  & 
511                  ne_ldown*Ftheta(ij+u_ldown)    +  &
512                  ne_rdown*Ftheta(ij+u_rdown) )
513          END DO
514       END DO
515
516       !DIR$ SIMD
517       DO ij=ij_begin,ij_end         
518          ! convm = -div(mass flux), sign convention as in Ringler et al. 2012, eq. 21
519          convm(ij,l)= -1./Ai(ij)*(ne_right*hflux(ij+u_right,l)   +  &
520               ne_rup*hflux(ij+u_rup,l)       +  & 
521               ne_lup*hflux(ij+u_lup,l)       +  & 
522               ne_left*hflux(ij+u_left,l)     +  & 
523               ne_ldown*hflux(ij+u_ldown,l)   +  & 
524               ne_rdown*hflux(ij+u_rdown,l))
525       END DO ! ij
526    END DO ! llm
527
528!!! Compute potential vorticity (Coriolis) contribution to du
529    SELECT CASE(caldyn_conserv)
530
531    CASE(energy) ! energy-conserving TRiSK
532
533       DO l=ll_begin,ll_end
534          !DIR$ SIMD
535          DO ij=ij_begin,ij_end         
536             uu_right = &
537                  wee(ij+u_right,1,1)*hflux(ij+u_rup,l)*(qu(ij+u_right,l)+qu(ij+u_rup,l))+                             &
538                  wee(ij+u_right,2,1)*hflux(ij+u_lup,l)*(qu(ij+u_right,l)+qu(ij+u_lup,l))+                             &
539                  wee(ij+u_right,3,1)*hflux(ij+u_left,l)*(qu(ij+u_right,l)+qu(ij+u_left,l))+                           &
540                  wee(ij+u_right,4,1)*hflux(ij+u_ldown,l)*(qu(ij+u_right,l)+qu(ij+u_ldown,l))+                         &
541                  wee(ij+u_right,5,1)*hflux(ij+u_rdown,l)*(qu(ij+u_right,l)+qu(ij+u_rdown,l))+                         & 
542                  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))+         &
543                  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))+         &
544                  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))+         &
545                  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))+             &
546                  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))
547             uu_lup = &
548                  wee(ij+u_lup,1,1)*hflux(ij+u_left,l)*(qu(ij+u_lup,l)+qu(ij+u_left,l)) +                        &
549                  wee(ij+u_lup,2,1)*hflux(ij+u_ldown,l)*(qu(ij+u_lup,l)+qu(ij+u_ldown,l)) +                      &
550                  wee(ij+u_lup,3,1)*hflux(ij+u_rdown,l)*(qu(ij+u_lup,l)+qu(ij+u_rdown,l)) +                      &
551                  wee(ij+u_lup,4,1)*hflux(ij+u_right,l)*(qu(ij+u_lup,l)+qu(ij+u_right,l)) +                      &
552                  wee(ij+u_lup,5,1)*hflux(ij+u_rup,l)*(qu(ij+u_lup,l)+qu(ij+u_rup,l)) +                          & 
553                  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)) +          &
554                  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)) +              &
555                  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)) +              &
556                  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)) +            &
557                  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))
558             uu_ldown = &
559                  wee(ij+u_ldown,1,1)*hflux(ij+u_rdown,l)*(qu(ij+u_ldown,l)+qu(ij+u_rdown,l)) +                      &
560                  wee(ij+u_ldown,2,1)*hflux(ij+u_right,l)*(qu(ij+u_ldown,l)+qu(ij+u_right,l)) +                      &
561                  wee(ij+u_ldown,3,1)*hflux(ij+u_rup,l)*(qu(ij+u_ldown,l)+qu(ij+u_rup,l)) +                          &
562                  wee(ij+u_ldown,4,1)*hflux(ij+u_lup,l)*(qu(ij+u_ldown,l)+qu(ij+u_lup,l)) +                          &
563                  wee(ij+u_ldown,5,1)*hflux(ij+u_left,l)*(qu(ij+u_ldown,l)+qu(ij+u_left,l)) +                        & 
564                  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)) +          &
565                  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)) +        &
566                  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)) +      &
567                  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)) +      &
568                  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))
569             du(ij+u_right,l) = du(ij+u_right,l) + .5*uu_right
570             du(ij+u_lup,l)   = du(ij+u_lup,l)   + .5*uu_lup
571             du(ij+u_ldown,l) = du(ij+u_ldown,l)   + .5*uu_ldown
572          ENDDO
573       ENDDO
574
575    CASE(enstrophy) ! enstrophy-conserving TRiSK
576
577       DO l=ll_begin,ll_end
578          !DIR$ SIMD
579          DO ij=ij_begin,ij_end         
580             uu_right = &
581                  wee(ij+u_right,1,1)*hflux(ij+u_rup,l)+                           &
582                  wee(ij+u_right,2,1)*hflux(ij+u_lup,l)+                           &
583                  wee(ij+u_right,3,1)*hflux(ij+u_left,l)+                          &
584                  wee(ij+u_right,4,1)*hflux(ij+u_ldown,l)+                         &
585                  wee(ij+u_right,5,1)*hflux(ij+u_rdown,l)+                         & 
586                  wee(ij+u_right,1,2)*hflux(ij+t_right+u_ldown,l)+                 &
587                  wee(ij+u_right,2,2)*hflux(ij+t_right+u_rdown,l)+                 &
588                  wee(ij+u_right,3,2)*hflux(ij+t_right+u_right,l)+                 &
589                  wee(ij+u_right,4,2)*hflux(ij+t_right+u_rup,l)+                   &
590                  wee(ij+u_right,5,2)*hflux(ij+t_right+u_lup,l)
591             uu_lup = &
592                  wee(ij+u_lup,1,1)*hflux(ij+u_left,l)+                        &
593                  wee(ij+u_lup,2,1)*hflux(ij+u_ldown,l)+                       &
594                  wee(ij+u_lup,3,1)*hflux(ij+u_rdown,l)+                       &
595                  wee(ij+u_lup,4,1)*hflux(ij+u_right,l)+                       &
596                  wee(ij+u_lup,5,1)*hflux(ij+u_rup,l)+                         & 
597                  wee(ij+u_lup,1,2)*hflux(ij+t_lup+u_right,l)+                 &
598                  wee(ij+u_lup,2,2)*hflux(ij+t_lup+u_rup,l)+                   &
599                  wee(ij+u_lup,3,2)*hflux(ij+t_lup+u_lup,l)+                   &
600                  wee(ij+u_lup,4,2)*hflux(ij+t_lup+u_left,l)+                  &
601                  wee(ij+u_lup,5,2)*hflux(ij+t_lup+u_ldown,l)
602             uu_ldown = &
603                  wee(ij+u_ldown,1,1)*hflux(ij+u_rdown,l)+                      &
604                  wee(ij+u_ldown,2,1)*hflux(ij+u_right,l)+                      &
605                  wee(ij+u_ldown,3,1)*hflux(ij+u_rup,l)+                        &
606                  wee(ij+u_ldown,4,1)*hflux(ij+u_lup,l)+                        &
607                  wee(ij+u_ldown,5,1)*hflux(ij+u_left,l)+                       & 
608                  wee(ij+u_ldown,1,2)*hflux(ij+t_ldown+u_lup,l)+                &
609                  wee(ij+u_ldown,2,2)*hflux(ij+t_ldown+u_left,l)+               &
610                  wee(ij+u_ldown,3,2)*hflux(ij+t_ldown+u_ldown,l)+              &
611                  wee(ij+u_ldown,4,2)*hflux(ij+t_ldown+u_rdown,l)+              &
612                  wee(ij+u_ldown,5,2)*hflux(ij+t_ldown+u_right,l)
613
614             du(ij+u_right,l) = du(ij+u_right,l) + .5*uu_right
615             du(ij+u_lup,l)   = du(ij+u_lup,l)   + .5*uu_lup
616             du(ij+u_ldown,l) = du(ij+u_ldown,l)   + .5*uu_ldown
617          END DO
618       END DO
619    CASE DEFAULT
620       STOP
621    END SELECT
622
623    CALL trace_end("compute_caldyn_Coriolis")
624
625  END SUBROUTINE compute_caldyn_Coriolis
626
627  SUBROUTINE compute_caldyn_slow_hydro(u,rhodz,hflux,du, zero)
628    LOGICAL, INTENT(IN) :: zero
629    REAL(rstd),INTENT(IN)  :: u(3*iim*jjm,llm)    ! prognostic "velocity"
630    REAL(rstd),INTENT(IN)  :: rhodz(iim*jjm,llm)
631    REAL(rstd),INTENT(OUT) :: hflux(3*iim*jjm,llm) ! hflux in kg/s
632    REAL(rstd),INTENT(INOUT) :: du(3*iim*jjm,llm)
633   
634    REAL(rstd) :: berni(iim*jjm)  ! Bernoulli function
635    REAL(rstd) :: uu_right, uu_lup, uu_ldown
636    INTEGER :: ij,l
637
638    CALL trace_start("compute_caldyn_slow_hydro")
639
640    le_de(:) = le(:)/de(:) ! FIXME - make sure le_de is what we expect
641
642    DO l = ll_begin, ll_end
643       !  Compute mass fluxes
644       IF (caldyn_conserv==energy) CALL test_message(req_qu) 
645       !DIR$ SIMD
646       DO ij=ij_begin_ext,ij_end_ext
647          uu_right=0.5*(rhodz(ij,l)+rhodz(ij+t_right,l))*u(ij+u_right,l)
648          uu_lup=0.5*(rhodz(ij,l)+rhodz(ij+t_lup,l))*u(ij+u_lup,l)
649          uu_ldown=0.5*(rhodz(ij,l)+rhodz(ij+t_ldown,l))*u(ij+u_ldown,l)
650          uu_right= uu_right*le_de(ij+u_right)
651          uu_lup  = uu_lup  *le_de(ij+u_lup)
652          uu_ldown= uu_ldown*le_de(ij+u_ldown)
653          hflux(ij+u_right,l)=uu_right
654          hflux(ij+u_lup,l)  =uu_lup
655          hflux(ij+u_ldown,l)=uu_ldown
656       ENDDO
657       ! Compute Bernoulli=kinetic energy
658       !DIR$ SIMD
659       DO ij=ij_begin,ij_end         
660          berni(ij) = &
661               1/(4*Ai(ij))*(le_de(ij+u_right)*u(ij+u_right,l)**2 +    &
662               le_de(ij+u_rup)*u(ij+u_rup,l)**2 +          &
663               le_de(ij+u_lup)*u(ij+u_lup,l)**2 +          &
664               le_de(ij+u_left)*u(ij+u_left,l)**2 +       &
665               le_de(ij+u_ldown)*u(ij+u_ldown,l)**2 +    &
666               le_de(ij+u_rdown)*u(ij+u_rdown,l)**2 ) 
667       ENDDO
668       ! Compute du=-grad(Bernoulli)
669       IF(zero) THEN
670          !DIR$ SIMD
671          DO ij=ij_begin,ij_end
672             du(ij+u_right,l) = ne_right*(berni(ij)-berni(ij+t_right))
673             du(ij+u_lup,l)   = ne_lup*(berni(ij)-berni(ij+t_lup))
674             du(ij+u_ldown,l) = ne_ldown*(berni(ij)-berni(ij+t_ldown))
675          END DO
676       ELSE
677          !DIR$ SIMD
678          DO ij=ij_begin,ij_end
679             du(ij+u_right,l) = du(ij+u_right,l) + &
680                  ne_right*(berni(ij)-berni(ij+t_right))
681             du(ij+u_lup,l)   = du(ij+u_lup,l) + &
682                  ne_lup*(berni(ij)-berni(ij+t_lup))
683             du(ij+u_ldown,l) = du(ij+u_ldown,l) + &
684                  ne_ldown*(berni(ij)-berni(ij+t_ldown))
685          END DO
686       END IF
687    END DO
688
689    CALL trace_end("compute_caldyn_slow_hydro")   
690  END SUBROUTINE compute_caldyn_slow_hydro
691
692  SUBROUTINE compute_caldyn_slow_NH(u,rhodz,Phi,W, hflux,du,dPhi,dW)
693    REAL(rstd),INTENT(IN)  :: u(3*iim*jjm,llm)    ! prognostic "velocity"
694    REAL(rstd),INTENT(IN)  :: rhodz(iim*jjm,llm)  ! rho*dz
695    REAL(rstd),INTENT(IN)  :: Phi(iim*jjm,llm+1)  ! prognostic geopotential
696    REAL(rstd),INTENT(IN)  :: W(iim*jjm,llm+1)    ! prognostic vertical momentum
697
698    REAL(rstd),INTENT(OUT) :: hflux(3*iim*jjm,llm) ! hflux in kg/s
699    REAL(rstd),INTENT(OUT) :: du(3*iim*jjm,llm)
700    REAL(rstd),INTENT(OUT) :: dW(iim*jjm,llm+1)
701    REAL(rstd),INTENT(OUT) :: dPhi(iim*jjm,llm+1)
702
703    REAL(rstd) :: w_il(3*iim*jjm,llm+1) ! Wil/mil
704    REAL(rstd) :: F_el(3*iim*jjm,llm+1) ! NH mass flux
705    REAL(rstd) :: GradPhi2(3*iim*jjm,llm+1) ! grad_Phi**2
706    REAL(rstd) :: DePhil(3*iim*jjm,llm+1) ! grad(Phi)
707    REAL(rstd) :: berni(iim*jjm)  ! Bernoulli function
708    REAL(rstd) :: G_el(3*iim*jjm) ! horizontal flux of W
709    REAL(rstd) :: v_el(3*iim*jjm)
710   
711    INTEGER :: ij,l,kdown,kup
712    REAL(rstd) :: uu_right, uu_lup, uu_ldown, W_el, W2_el
713
714    CALL trace_start("compute_caldyn_slow_NH")
715
716    le_de(:) = le(:)/de(:) ! FIXME - make sure le_de is what we expect
717
718    DO l=ll_begin, ll_endp1 ! compute on l levels (interfaces)
719       IF(l==1) THEN
720          kdown=1
721       ELSE
722          kdown=l-1
723       END IF
724       IF(l==llm+1) THEN
725          kup=llm
726       ELSE
727          kup=l
728       END IF
729       ! below : "checked" means "formula also valid when kup=kdown (top/bottom)"
730       ! compute mil, wil=Wil/mil
731       DO ij=ij_begin_ext, ij_end_ext
732          w_il(ij,l) = 2.*W(ij,l)/(rhodz(ij,kdown)+rhodz(ij,kup)) ! checked
733       END DO
734       ! compute DePhi, v_el, G_el, F_el
735       ! v_el, W2_el and therefore G_el incorporate metric factor le_de
736       ! while DePhil, W_el and F_el don't
737       DO ij=ij_begin_ext, ij_end_ext
738          ! Compute on edge 'right'
739          W_el  = .5*( W(ij,l)+W(ij+t_right,l) )
740          DePhil(ij+u_right,l) = ne_right*(Phi(ij+t_right,l)-Phi(ij,l))
741          F_el(ij+u_right,l)   = DePhil(ij+u_right,l)*W_el
742          W2_el = .5*le_de(ij+u_right) * &
743               ( W(ij,l)*w_il(ij,l) + W(ij+t_right,l)*w_il(ij+t_right,l) )
744          v_el(ij+u_right) = .5*le_de(ij+u_right)*(u(ij+u_right,kup)+u(ij+u_right,kdown)) ! checked
745          G_el(ij+u_right) = v_el(ij+u_right)*W_el - DePhil(ij+u_right,l)*W2_el
746          ! Compute on edge 'lup'
747          W_el  = .5*( W(ij,l)+W(ij+t_lup,l) )
748          DePhil(ij+u_lup,l) = ne_lup*(Phi(ij+t_lup,l)-Phi(ij,l))
749          F_el(ij+u_lup,l)   = DePhil(ij+u_lup,l)*W_el
750          W2_el = .5*le_de(ij+u_lup) * &
751               ( W(ij,l)*w_il(ij,l) + W(ij+t_lup,l)*w_il(ij+t_lup,l) )
752          v_el(ij+u_lup) = .5*le_de(ij+u_lup)*( u(ij+u_lup,kup) + u(ij+u_lup,kdown)) ! checked
753          G_el(ij+u_lup) = v_el(ij+u_lup)*W_el - DePhil(ij+u_lup,l)*W2_el
754          ! Compute on edge 'ldown'
755          W_el  = .5*( W(ij,l)+W(ij+t_ldown,l) )
756          DePhil(ij+u_ldown,l) = ne_ldown*(Phi(ij+t_ldown,l)-Phi(ij,l))
757          F_el(ij+u_ldown,l) = DePhil(ij+u_ldown,l)*W_el
758          W2_el = .5*le_de(ij+u_ldown) * &
759               ( W(ij,l)*w_il(ij,l) + W(ij+t_ldown,l)*w_il(ij+t_ldown,l) )
760          v_el(ij+u_ldown) = .5*le_de(ij+u_ldown)*( u(ij+u_ldown,kup) + u(ij+u_ldown,kdown)) ! checked
761          G_el(ij+u_ldown) = v_el(ij+u_ldown)*W_el - DePhil(ij+u_ldown,l)*W2_el
762       END DO
763       ! compute GradPhi2, dPhi, dW
764       DO ij=ij_begin_ext, ij_end_ext
765          gradPhi2(ij,l) = &
766               1/(2*Ai(ij))*(le_de(ij+u_right)*DePhil(ij+u_right,l)**2 + &
767               le_de(ij+u_rup)*DePhil(ij+u_rup,l)**2 +      &
768               le_de(ij+u_lup)*DePhil(ij+u_lup,l)**2 +      &
769               le_de(ij+u_left)*DePhil(ij+u_left,l)**2 +    &
770               le_de(ij+u_ldown)*DePhil(ij+u_ldown,l)**2 +  &
771               le_de(ij+u_rdown)*DePhil(ij+u_rdown,l)**2 )
772!          gradPhi2(ij,l) = 0. ! FIXME !!
773
774          dPhi(ij,l) = gradPhi2(ij,l)*w_il(ij,l) -1/(2*Ai(ij))* & 
775               ( DePhil(ij+u_right,l)*v_el(ij+u_right) + & ! -v.gradPhi,
776                 DePhil(ij+u_rup,l)*v_el(ij+u_rup) +     & ! v_el already has le_de
777                 DePhil(ij+u_lup,l)*v_el(ij+u_lup) +     &
778                 DePhil(ij+u_left,l)*v_el(ij+u_left) +   &
779                 DePhil(ij+u_ldown,l)*v_el(ij+u_ldown) + &
780                 DePhil(ij+u_rdown,l)*v_el(ij+u_rdown) )
781
782          dW(ij,l) = -1./Ai(ij)*(           & ! -div(G_el),
783               ne_right*G_el(ij+u_right) +  & ! G_el already has le_de
784               ne_rup*G_el(ij+u_rup) +      &
785               ne_lup*G_el(ij+u_lup) +      & 
786               ne_left*G_el(ij+u_left) +    &
787               ne_ldown*G_el(ij+u_ldown) +  &
788               ne_rdown*G_el(ij+u_rdown))
789       END DO
790    END DO
791    ! FIXME !!
792 !   F_el(:,:)=0.
793 !   dPhi(:,:)=0.
794 !   dW(:,:)=0.
795
796    DO l=ll_begin, ll_end ! compute on k levels (layers)
797       ! Compute berni at scalar points
798       DO ij=ij_begin_ext, ij_end_ext
799          berni(ij) = &
800               1/(4*Ai(ij))*( &
801                    le_de(ij+u_right)*u(ij+u_right,l)**2 +    &
802                    le_de(ij+u_rup)*u(ij+u_rup,l)**2 +          &
803                    le_de(ij+u_lup)*u(ij+u_lup,l)**2 +          &
804                    le_de(ij+u_left)*u(ij+u_left,l)**2 +       &
805                    le_de(ij+u_ldown)*u(ij+u_ldown,l)**2 +    &
806                    le_de(ij+u_rdown)*u(ij+u_rdown,l)**2 ) &
807               - .25*( gradPhi2(ij,l)  *w_il(ij,l)**2 +   &
808                      gradPhi2(ij,l+1)*w_il(ij,l+1)**2 )
809       END DO
810       ! Compute mass flux and grad(berni) at edges
811       DO ij=ij_begin_ext, ij_end_ext
812          ! Compute on edge 'right'
813          uu_right = 0.5*(rhodz(ij,l)+rhodz(ij+t_right,l))*u(ij+u_right,l) &
814                    -0.5*(F_el(ij+u_right,l)+F_el(ij+u_right,l+1))
815          hflux(ij+u_right,l) = uu_right*le_de(ij+u_right)
816          du(ij+u_right,l) = ne_right*(berni(ij)-berni(ij+t_right))
817          ! Compute on edge 'lup'
818          uu_lup = 0.5*(rhodz(ij,l)+rhodz(ij+t_lup,l))*u(ij+u_lup,l) &
819                  -0.5*(F_el(ij+u_lup,l)+F_el(ij+u_lup,l+1))
820          hflux(ij+u_lup,l) = uu_lup*le_de(ij+u_lup)
821          du(ij+u_lup,l) = ne_lup*(berni(ij)-berni(ij+t_lup))
822          ! Compute on edge 'ldown'
823          uu_ldown = 0.5*(rhodz(ij,l)+rhodz(ij+t_ldown,l))*u(ij+u_ldown,l) &
824                    -0.5*(F_el(ij+u_ldown,l)+F_el(ij+u_ldown,l+1))
825          hflux(ij+u_ldown,l) = uu_ldown*le_de(ij+u_ldown)
826          du(ij+u_ldown,l) = ne_ldown*(berni(ij)-berni(ij+t_ldown))
827       END DO
828    END DO
829    ! FIXME !!
830    ! du(:,:)=0.
831    ! hflux(:,:)=0.
832
833    CALL trace_end("compute_caldyn_slow_NH")
834
835  END SUBROUTINE compute_caldyn_slow_NH
836
837END MODULE caldyn_kernels_hevi_mod
Note: See TracBrowser for help on using the repository browser.