MODULE compute_caldyn_solver_mod USE grid_param USE compute_NH_geopot_mod, ONLY : compute_NH_geopot, pbot, rho_bot IMPLICIT NONE PRIVATE #include "../unstructured/unstructured.h90" PUBLIC :: compute_caldyn_solver CONTAINS #ifdef BEGIN_DYSL KERNEL(caldyn_solver) ! ! Compute pressure (pres) and Exner function (pk) ! kappa = R/Cp ! 1-kappa = Cv/Cp ! Cp/Cv = 1/(1-kappa) gamma = 1./(1.-kappa) vreff = Rd*Treff/preff ! reference specific volume Cvd = 1./(cpp-Rd) Rd_preff = kappa*cpp/preff FORALL_CELLS_EXT() SELECT CASE(caldyn_thermo) CASE(thermo_theta) ON_PRIMAL rho_ij = 1./(geopot(UP(CELL))-geopot(CELL)) rho_ij = rho_ij*g*rhodz(CELL) X_ij = Rd_preff*theta(CELL,1)*rho_ij ! kappa.theta.rho = p/exner ! => X = (p/p0)/(exner/Cp) ! = (p/p0)^(1-kappa) pres(CELL) = preff*(X_ij**gamma) ! pressure ! Compute Exner function (needed by compute_caldyn_fast) ! other formulae possible if exponentiation is slow pk(CELL) = cpp*((pres(CELL)/preff)**kappa) ! Exner END_BLOCK CASE(thermo_entropy) ON_PRIMAL rho_ij = 1./(geopot(UP(CELL))-geopot(CELL)) rho_ij = rho_ij*g*rhodz(CELL) T_ij = Treff*exp( (theta(CELL,1)+Rd*log(vreff*rho_ij))*Cvd ) pres(CELL) = rho_ij*Rd*T_ij pk(CELL) = T_ij END_BLOCK CASE DEFAULT STOP END SELECT END_BLOCK ! We need a barrier here because we compute pres above and do a vertical difference below BARRIER FORALL_CELLS_EXT('1', 'llm+1') ON_PRIMAL CST_IFTHEN(IS_BOTTOM_LEVEL) ! Lower BC dW(CELL) = (1./g)*(pbot-rho_bot*(geopot(CELL)-PHI_BOT(HIDX(CELL)))-pres(CELL)) - m_il(CELL) CST_ELSEIF(IS_TOP_INTERFACE) ! Top BC dW(CELL) = (1./g)*(pres(DOWN(CELL))-ptop) - m_il(CELL) CST_ELSE dW(CELL) = (1./g)*(pres(DOWN(CELL))-pres(CELL)) - m_il(CELL) CST_ENDIF W(CELL) = W(CELL)+tau*dW(CELL) ! update W dPhi(CELL) = g*g*W(CELL)/m_il(CELL) END_BLOCK END_BLOCK ! We need a barrier here because we update W above and do a vertical average below BARRIER FORALL_CELLS_EXT() ON_PRIMAL ! compute du = -0.5*g^2.grad(w^2) berni(CELL) = (-.25*g*g)*((W(CELL)/m_il(CELL))**2 + (W(UP(CELL))/m_il(UP(CELL)))**2 ) END_BLOCK END_BLOCK FORALL_CELLS_EXT() ON_EDGES du(EDGE) = SIGN*(berni(CELL1)-berni(CELL2)) END_BLOCK END_BLOCK END_BLOCK #endif END_DYSL SUBROUTINE compute_caldyn_solver_unst(tau,rhodz,theta, berni,pres,m_il, pk,geopot,W,dPhi,dW,du) USE ISO_C_BINDING, only : C_DOUBLE, C_FLOAT USE earth_const USE trace USE grid_param, ONLY : nqdyn USE disvert_mod, ONLY : ptop USE data_unstructured_mod, ONLY : enter_trace, exit_trace, & id_solver, left, right, PHI_BOT USE compute_NH_geopot_mod, ONLY : compute_NH_geopot_unst NUM, INTENT(IN) :: tau FIELD_MASS :: rhodz,pk,berni,pres ! IN, OUT, BUF*2 FIELD_THETA :: theta ! IN FIELD_GEOPOT :: geopot,W,dPhi,dW, m_il ! INOUT,INOUT, OUT,OUT, BUF FIELD_U :: du ! OUT DECLARE_INDICES NUM :: X_ij, rho_ij, T_ij, gamma, Cvd, vreff, Rd_preff #define PHI_BOT(ij) Phi_bot #include "../kernels_unst/caldyn_mil.k90" IF(tau>0) THEN ! solve implicit problem for geopotential CALL compute_NH_geopot_unst(tau, rhodz, m_il, theta, W, geopot) END IF START_TRACE(id_solver, 7,0,1) #include "../kernels_unst/caldyn_solver.k90" STOP_TRACE #undef PHI_BOT END SUBROUTINE compute_caldyn_solver_unst SUBROUTINE compute_caldyn_solver(tau,phis, rhodz,theta,pk, geopot,W, m_il,pres, dPhi,dW,du) USE icosa USE caldyn_vars_mod USE trace USE omp_para, ONLY : ll_begin, ll_end,ll_beginp1,ll_endp1 USE disvert_mod, ONLY : ptop REAL(rstd),INTENT(IN) :: tau ! "solve" Phi-tau*dPhi/dt = Phi_rhs REAL(rstd),INTENT(IN) :: phis(iim*jjm) REAL(rstd),INTENT(IN) :: rhodz(iim*jjm,llm) REAL(rstd),INTENT(IN) :: theta(iim*jjm,llm,nqdyn) REAL(rstd),INTENT(OUT) :: pk(iim*jjm,llm) REAL(rstd),INTENT(INOUT) :: geopot(iim*jjm,llm+1) REAL(rstd),INTENT(INOUT) :: W(iim*jjm,llm+1) ! OUT if tau>0 REAL(rstd),INTENT(OUT) :: m_il(iim*jjm,llm+1) ! rhodz averaged to interfaces REAL(rstd),INTENT(OUT) :: pres(iim*jjm,llm) ! pressure REAL(rstd),INTENT(OUT) :: dW(iim*jjm,llm+1) REAL(rstd),INTENT(OUT) :: dPhi(iim*jjm,llm+1) REAL(rstd),INTENT(OUT) :: du(3*iim*jjm,llm) REAL(rstd) :: berni(iim*jjm,llm) ! (W/m_il)^2 REAL(rstd) :: berni1(iim*jjm) ! (W/m_il)^2 !REAL(rstd) :: gamma, rho_ij, T_ij, X_ij, Y_ij, vreff, Rd, Cvd, Rd_preff REAL(rstd) :: gamma, rho_ij, T_ij, X_ij, Y_ij, vreff, Cvd, Rd_preff INTEGER :: ij, l CALL trace_start("compute_caldyn_solver") !Rd=cpp*kappa IF(dysl) THEN !$OMP BARRIER #include "../kernels_hex/caldyn_mil.k90" IF(tau>0) THEN ! solve implicit problem for geopotential CALL compute_NH_geopot(tau,phis, rhodz, m_il, theta, W, geopot) END IF #define PHI_BOT(ij) phis(ij) #include "../kernels_hex/caldyn_solver.k90" #undef PHI_BOT !$OMP BARRIER ELSE #define BERNI(ij) berni1(ij) ! FIXME : vertical OpenMP parallelism will not work ! average m_ik to interfaces => m_il !DIR$ SIMD DO ij=ij_begin_ext,ij_end_ext m_il(ij,1) = .5*rhodz(ij,1) ENDDO DO l=2,llm !DIR$ SIMD DO ij=ij_begin_ext,ij_end_ext m_il(ij,l) = .5*(rhodz(ij,l-1)+rhodz(ij,l)) ENDDO ENDDO !DIR$ SIMD DO ij=ij_begin_ext,ij_end_ext m_il(ij,llm+1) = .5*rhodz(ij,llm) ENDDO IF(tau>0) THEN ! solve implicit problem for geopotential CALL compute_NH_geopot(tau, phis, rhodz, m_il, theta, W, geopot) END IF ! Compute pressure, stored temporarily in pk ! kappa = R/Cp ! 1-kappa = Cv/Cp ! Cp/Cv = 1/(1-kappa) gamma = 1./(1.-kappa) DO l=1,llm !DIR$ SIMD DO ij=ij_begin_ext,ij_end_ext rho_ij = (g*rhodz(ij,l))/(geopot(ij,l+1)-geopot(ij,l)) X_ij = (cpp/preff)*kappa*theta(ij,l,1)*rho_ij ! kappa.theta.rho = p/exner ! => X = (p/p0)/(exner/Cp) ! = (p/p0)^(1-kappa) pk(ij,l) = preff*(X_ij**gamma) ENDDO ENDDO ! Update W, compute tendencies DO l=2,llm !DIR$ SIMD DO ij=ij_begin_ext,ij_end_ext dW(ij,l) = (1./g)*(pk(ij,l-1)-pk(ij,l)) - m_il(ij,l) W(ij,l) = W(ij,l)+tau*dW(ij,l) ! update W dPhi(ij,l) = g*g*W(ij,l)/m_il(ij,l) ENDDO ! PRINT *,'Max dPhi', l,ij_begin,ij_end, MAXVAL(abs(dPhi(ij_begin:ij_end,l))) ! PRINT *,'Max dW', l,ij_begin,ij_end, MAXVAL(abs(dW(ij_begin:ij_end,l))) ENDDO ! Lower BC (FIXME : no orography yet !) DO ij=ij_begin,ij_end dPhi(ij,1)=0 W(ij,1)=0 dW(ij,1)=0 dPhi(ij,llm+1)=0 ! rigid lid W(ij,llm+1)=0 dW(ij,llm+1)=0 ENDDO ! Upper BC p=ptop ! DO ij=ij_omp_begin_ext,ij_omp_end_ext ! dPhi(ij,llm+1) = W(ij,llm+1)/rhodz(ij,llm) ! dW(ij,llm+1) = (1./g)*(pk(ij,llm)-ptop) - .5*rhodz(ij,llm) ! ENDDO ! Compute Exner function (needed by compute_caldyn_fast) and du=-g^2.grad(w^2) DO l=1,llm !DIR$ SIMD DO ij=ij_begin_ext,ij_end_ext pk(ij,l) = cpp*((pk(ij,l)/preff)**kappa) ! other formulae possible if exponentiation is slow BERNI(ij) = (-.25*g*g)*( & (W(ij,l)/m_il(ij,l))**2 & + (W(ij,l+1)/m_il(ij,l+1))**2 ) ENDDO DO ij=ij_begin,ij_end du(ij+u_right,l) = ne_right*(BERNI(ij)-BERNI(ij+t_right)) du(ij+u_lup,l) = ne_lup *(BERNI(ij)-BERNI(ij+t_lup)) du(ij+u_ldown,l) = ne_ldown*(BERNI(ij)-BERNI(ij+t_ldown)) ENDDO ENDDO #undef BERNI END IF ! dysl CALL trace_end("compute_caldyn_solver") END SUBROUTINE compute_caldyn_solver END MODULE compute_caldyn_solver_mod