source: codes/icosagcm/devel/src/dynamics/compute_caldyn_solver.F90 @ 853

Last change on this file since 853 was 853, checked in by jisesh, 5 years ago

devel: separate module for compute_caldyn_solver and removed Rd duplication in it

File size: 4.6 KB
Line 
1MODULE compute_caldyn_solver_mod
2  USE grid_param, ONLY : llm
3  IMPLICIT NONE
4  PRIVATE
5
6  PUBLIC :: compute_caldyn_solver
7
8CONTAINS
9
10  SUBROUTINE compute_caldyn_solver(tau,phis, rhodz,theta,pk, geopot,W, m_il,pres, dPhi,dW,du)
11    USE icosa
12    USE caldyn_vars_mod
13    USE trace
14    USE omp_para, ONLY : ll_begin, ll_end,ll_beginp1,ll_endp1
15    USE disvert_mod, ONLY : ptop
16    USE caldyn_kernels_hevi_mod
17    REAL(rstd), PARAMETER :: pbot=1e5, rho_bot=1e6
18    REAL(rstd),INTENT(IN) :: tau ! "solve" Phi-tau*dPhi/dt = Phi_rhs
19    REAL(rstd),INTENT(IN)    :: phis(iim*jjm)
20    REAL(rstd),INTENT(IN)    :: rhodz(iim*jjm,llm)
21    REAL(rstd),INTENT(IN)    :: theta(iim*jjm,llm,nqdyn)
22    REAL(rstd),INTENT(OUT)   :: pk(iim*jjm,llm)
23    REAL(rstd),INTENT(INOUT) :: geopot(iim*jjm,llm+1)
24    REAL(rstd),INTENT(INOUT) :: W(iim*jjm,llm+1) ! OUT if tau>0
25    REAL(rstd),INTENT(OUT)   :: m_il(iim*jjm,llm+1)        ! rhodz averaged to interfaces
26    REAL(rstd),INTENT(OUT)   :: pres(iim*jjm,llm)          ! pressure
27    REAL(rstd),INTENT(OUT)   :: dW(iim*jjm,llm+1)
28    REAL(rstd),INTENT(OUT)   :: dPhi(iim*jjm,llm+1)
29    REAL(rstd),INTENT(OUT)   :: du(3*iim*jjm,llm)
30
31    REAL(rstd) :: berni(iim*jjm,llm)         ! (W/m_il)^2
32    REAL(rstd) :: berni1(iim*jjm)         ! (W/m_il)^2
33    !REAL(rstd) :: gamma, rho_ij, T_ij, X_ij, Y_ij, vreff, Rd, Cvd, Rd_preff
34    REAL(rstd) :: gamma, rho_ij, T_ij, X_ij, Y_ij, vreff, Cvd, Rd_preff
35    INTEGER    :: ij, l
36
37    CALL trace_start("compute_caldyn_solver")
38
39    !Rd=cpp*kappa
40
41    IF(dysl) THEN
42
43!$OMP BARRIER
44
45#include "../kernels_hex/caldyn_mil.k90"
46  IF(tau>0) THEN ! solve implicit problem for geopotential
47    CALL compute_NH_geopot(tau,phis, rhodz, m_il, theta, W, geopot)
48  END IF
49#define PHI_BOT(ij) phis(ij)
50#include "../kernels_hex/caldyn_solver.k90"
51#undef PHI_BOT
52!$OMP BARRIER
53
54    ELSE
55
56#define BERNI(ij) berni1(ij)
57    ! FIXME : vertical OpenMP parallelism will not work
58
59    ! average m_ik to interfaces => m_il
60    !DIR$ SIMD
61    DO ij=ij_begin_ext,ij_end_ext
62       m_il(ij,1) = .5*rhodz(ij,1)
63    ENDDO
64    DO l=2,llm
65       !DIR$ SIMD
66       DO ij=ij_begin_ext,ij_end_ext
67          m_il(ij,l) = .5*(rhodz(ij,l-1)+rhodz(ij,l))
68       ENDDO
69    ENDDO
70    !DIR$ SIMD
71    DO ij=ij_begin_ext,ij_end_ext
72       m_il(ij,llm+1) = .5*rhodz(ij,llm)
73    ENDDO
74
75    IF(tau>0) THEN ! solve implicit problem for geopotential
76       CALL compute_NH_geopot(tau, phis, rhodz, m_il, theta, W, geopot)
77    END IF
78
79    ! Compute pressure, stored temporarily in pk
80    ! kappa = R/Cp
81    ! 1-kappa = Cv/Cp
82    ! Cp/Cv = 1/(1-kappa)
83    gamma = 1./(1.-kappa)
84    DO l=1,llm
85       !DIR$ SIMD
86       DO ij=ij_begin_ext,ij_end_ext
87          rho_ij = (g*rhodz(ij,l))/(geopot(ij,l+1)-geopot(ij,l))
88          X_ij = (cpp/preff)*kappa*theta(ij,l,1)*rho_ij
89          ! kappa.theta.rho = p/exner
90          ! => X = (p/p0)/(exner/Cp)
91          !      = (p/p0)^(1-kappa)
92          pk(ij,l) = preff*(X_ij**gamma)
93       ENDDO
94    ENDDO
95
96    ! Update W, compute tendencies
97    DO l=2,llm
98       !DIR$ SIMD
99       DO ij=ij_begin_ext,ij_end_ext
100          dW(ij,l)   = (1./g)*(pk(ij,l-1)-pk(ij,l)) - m_il(ij,l)
101          W(ij,l)    = W(ij,l)+tau*dW(ij,l) ! update W
102          dPhi(ij,l) = g*g*W(ij,l)/m_il(ij,l)
103       ENDDO
104       !       PRINT *,'Max dPhi', l,ij_begin,ij_end, MAXVAL(abs(dPhi(ij_begin:ij_end,l)))
105       !       PRINT *,'Max dW', l,ij_begin,ij_end, MAXVAL(abs(dW(ij_begin:ij_end,l)))
106    ENDDO
107    ! Lower BC (FIXME : no orography yet !)
108    DO ij=ij_begin,ij_end         
109       dPhi(ij,1)=0
110       W(ij,1)=0
111       dW(ij,1)=0
112       dPhi(ij,llm+1)=0 ! rigid lid
113       W(ij,llm+1)=0
114       dW(ij,llm+1)=0
115    ENDDO
116    ! Upper BC p=ptop
117    !    DO ij=ij_omp_begin_ext,ij_omp_end_ext
118    !       dPhi(ij,llm+1) = W(ij,llm+1)/rhodz(ij,llm)
119    !       dW(ij,llm+1) = (1./g)*(pk(ij,llm)-ptop) - .5*rhodz(ij,llm)
120    !    ENDDO
121   
122    ! Compute Exner function (needed by compute_caldyn_fast) and du=-g^2.grad(w^2)
123    DO l=1,llm
124       !DIR$ SIMD
125       DO ij=ij_begin_ext,ij_end_ext
126          pk(ij,l) = cpp*((pk(ij,l)/preff)**kappa) ! other formulae possible if exponentiation is slow
127          BERNI(ij) = (-.25*g*g)*(              &
128                 (W(ij,l)/m_il(ij,l))**2       &
129               + (W(ij,l+1)/m_il(ij,l+1))**2 )
130       ENDDO
131       DO ij=ij_begin,ij_end         
132          du(ij+u_right,l) = ne_right*(BERNI(ij)-BERNI(ij+t_right))
133          du(ij+u_lup,l)   = ne_lup  *(BERNI(ij)-BERNI(ij+t_lup))
134          du(ij+u_ldown,l) = ne_ldown*(BERNI(ij)-BERNI(ij+t_ldown))
135       ENDDO
136    ENDDO
137#undef BERNI
138
139    END IF ! dysl
140
141    CALL trace_end("compute_caldyn_solver")
142   
143  END SUBROUTINE compute_caldyn_solver
144 
145END MODULE compute_caldyn_solver_mod
Note: See TracBrowser for help on using the repository browser.