source: codes/icosagcm/devel/src/dynamics/caldyn_kernels_hevi.F90 @ 538

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

devel : macro-generated kernels for caldyn_coriolis, compute_NH_geopot, caldyn_solver, caldyn_vert_NH

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