[859] | 1 | MODULE compute_NH_geopot_mod |
---|
[939] | 2 | USE grid_param |
---|
[859] | 3 | IMPLICIT NONE |
---|
| 4 | PRIVATE |
---|
| 5 | |
---|
| 6 | LOGICAL, SAVE :: debug_hevi_solver = .FALSE. |
---|
| 7 | |
---|
[878] | 8 | #include "../unstructured/unstructured.h90" |
---|
[859] | 9 | |
---|
[878] | 10 | PUBLIC :: compute_NH_geopot,compute_NH_geopot_unst |
---|
| 11 | |
---|
[859] | 12 | CONTAINS |
---|
| 13 | |
---|
[878] | 14 | #ifdef BEGIN_DYSL |
---|
| 15 | |
---|
| 16 | KERNEL(compute_NH_geopot) |
---|
| 17 | tau2_g=tau*tau/g |
---|
| 18 | g2=g*g |
---|
| 19 | gm2 = 1./g2 |
---|
| 20 | vreff = Treff*cpp/preff*kappa |
---|
| 21 | gamma = 1./(1.-kappa) |
---|
| 22 | |
---|
| 23 | BARRIER |
---|
| 24 | |
---|
| 25 | ! compute Phi_star |
---|
| 26 | SEQUENCE_C1 |
---|
| 27 | BODY('1,llm+1') |
---|
| 28 | Phi_star_il(CELL) = Phi_il(CELL) + tau*g2*(W_il(CELL)/m_il(CELL)-tau) |
---|
| 29 | END_BLOCK |
---|
| 30 | END_BLOCK |
---|
| 31 | |
---|
| 32 | ! Newton-Raphson iteration : Phi_il contains current guess value |
---|
| 33 | DO iter=1,2 ! 2 iterations should be enough |
---|
| 34 | ! Compute pressure, A_ik |
---|
| 35 | SELECT CASE(caldyn_thermo) |
---|
| 36 | CASE(thermo_theta) |
---|
| 37 | {% call() compute_p_and_Aik() %} |
---|
| 38 | X_ij = (cpp/preff)*kappa*theta(CELL)*rho_ij |
---|
| 39 | p_ik(CELL) = preff*(X_ij**gamma) |
---|
| 40 | c2_mik = gamma*p_ik(CELL)/(rho_ij*m_ik(CELL)) ! c^2 = gamma*R*T = gamma*p/rho |
---|
| 41 | {% endcall %} |
---|
| 42 | CASE(thermo_entropy) |
---|
| 43 | {% call() compute_p_and_Aik() %} |
---|
| 44 | X_ij = log(vreff*rho_ij) + theta(CELL)/cpp |
---|
| 45 | p_ik(CELL) = preff*exp(X_ij*gamma) |
---|
| 46 | c2_mik = gamma*p_ik(CELL)/(rho_ij*m_ik(CELL)) ! c^2 = gamma*R*T = gamma*p/rho |
---|
| 47 | {% endcall %} |
---|
| 48 | CASE DEFAULT |
---|
| 49 | PRINT *, 'caldyn_thermo not supported by compute_NH_geopot', caldyn_thermo |
---|
| 50 | STOP |
---|
| 51 | END SELECT |
---|
| 52 | |
---|
| 53 | ! NB : A(1), A(llm), R(1), R(llm+1) = 0 => x(l)=0 at l=1,llm+1 => flat, rigid top and bottom |
---|
| 54 | ! Solve -A(l-1)x(l-1) + B(l)x(l) - A(l)x(l+1) = R(l) using Thomas algorithm |
---|
| 55 | |
---|
| 56 | SEQUENCE_C1 |
---|
| 57 | ! Compute residual R_il and B_il |
---|
| 58 | PROLOGUE(1) |
---|
| 59 | ! bottom interface l=1 |
---|
| 60 | ml_g2 = gm2*m_il(CELL) |
---|
| 61 | B_il(CELL) = A_ik(CELL) + ml_g2 + tau2_g*rho_bot |
---|
| 62 | R_il(CELL) = ml_g2*( Phi_il(CELL)-Phi_star_il(CELL)) & |
---|
| 63 | + tau2_g*( p_ik(CELL)-pbot+rho_bot*(Phi_il(CELL)-PHI_BOT(HIDX(CELL))) ) |
---|
| 64 | END_BLOCK |
---|
| 65 | BODY('2,llm') |
---|
| 66 | ! inner interfaces |
---|
| 67 | ml_g2 = gm2*m_il(CELL) |
---|
| 68 | B_il(CELL) = A_ik(CELL)+A_ik(DOWN(CELL)) + ml_g2 |
---|
| 69 | R_il(CELL) = ml_g2*( Phi_il(CELL)-Phi_star_il(CELL)) & |
---|
| 70 | + tau2_g*(p_ik(CELL)-p_ik(DOWN(CELL))) |
---|
| 71 | ! consistency check : if Wil=0 and initial state is in hydrostatic balance |
---|
| 72 | ! then Phi_star_il(CELL) = Phi_il(CELL) - tau^2*g^2 |
---|
| 73 | ! and residual = tau^2*(ml+(1/g)dl_pi)=0 |
---|
| 74 | END_BLOCK |
---|
| 75 | EPILOGUE('llm+1') |
---|
| 76 | ! top interface l=llm+1 |
---|
| 77 | ml_g2 = gm2*m_il(CELL) |
---|
| 78 | B_il(CELL) = A_ik(DOWN(CELL)) + ml_g2 |
---|
| 79 | R_il(CELL) = ml_g2*( Phi_il(CELL)-Phi_star_il(CELL)) & |
---|
| 80 | + tau2_g*( ptop-p_ik(DOWN(CELL)) ) |
---|
| 81 | END_BLOCK |
---|
| 82 | ! |
---|
| 83 | ! Forward sweep : |
---|
| 84 | ! C(0)=0, C(l) = -A(l) / (B(l)+A(l-1)C(l-1)), |
---|
| 85 | ! D(0)=0, D(l) = (R(l)+A(l-1)D(l-1)) / (B(l)+A(l-1)C(l-1)) |
---|
| 86 | PROLOGUE(1) |
---|
| 87 | X_ij = 1./B_il(CELL) |
---|
| 88 | C_ik(CELL) = -A_ik(CELL) * X_ij |
---|
| 89 | D_il(CELL) = R_il(CELL) * X_ij |
---|
| 90 | END_BLOCK |
---|
| 91 | BODY('2,llm') |
---|
| 92 | X_ij = 1./( B_il(CELL) + A_ik(DOWN(CELL))*C_ik(DOWN(CELL)) ) |
---|
| 93 | C_ik(CELL) = -A_ik(CELL) * X_ij |
---|
| 94 | D_il(CELL) = (R_il(CELL)+A_ik(DOWN(CELL))*D_il(DOWN(CELL))) * X_ij |
---|
| 95 | END_BLOCK |
---|
| 96 | EPILOGUE('llm+1') |
---|
| 97 | X_ij = 1./( B_il(CELL) + A_ik(DOWN(CELL))*C_ik(DOWN(CELL)) ) |
---|
| 98 | D_il(CELL) = (R_il(CELL)+A_ik(DOWN(CELL))*D_il(DOWN(CELL))) * X_ij |
---|
| 99 | ! Back substitution : |
---|
| 100 | ! x(i) = D(i)-C(i)x(i+1), x(llm+1)=0 |
---|
| 101 | ! + Newton-Raphson update |
---|
| 102 | ! top interface l=llm+1 |
---|
| 103 | x_il(CELL) = D_il(CELL) |
---|
| 104 | Phi_il(CELL) = Phi_il(CELL) - x_il(CELL) |
---|
| 105 | END_BLOCK |
---|
| 106 | BODY('llm,1,-1') |
---|
| 107 | ! Back substitution at lower interfaces |
---|
| 108 | x_il(CELL) = D_il(CELL) - C_ik(CELL)*x_il(UP(CELL)) |
---|
| 109 | Phi_il(CELL) = Phi_il(CELL) - x_il(CELL) |
---|
| 110 | END_BLOCK |
---|
| 111 | END_BLOCK |
---|
| 112 | |
---|
| 113 | IF(debug_hevi_solver) THEN |
---|
| 114 | PRINT *, '[hevi_solver] A,B', iter, MAXVAL(ABS(A_ik)),MAXVAL(ABS(B_il)) |
---|
| 115 | PRINT *, '[hevi_solver] C,D', iter, MAXVAL(ABS(C_ik)),MAXVAL(ABS(D_il)) |
---|
| 116 | DO l=1,llm+1 |
---|
| 117 | WRITE(*,'(A,I2.1,I3.2,E9.2)') '[hevi_solver] x_il', iter,l, MAXVAL(ABS(x_il(l,:))) |
---|
| 118 | END DO |
---|
| 119 | DO l=1,llm+1 |
---|
| 120 | WRITE(*,'(A,I2.1,I3.2,E9.2)') '[hevi_solver] R_il', iter,l, MAXVAL(ABS(R_il(l,:))) |
---|
| 121 | END DO |
---|
| 122 | END IF |
---|
| 123 | |
---|
| 124 | END DO ! Newton-Raphson |
---|
| 125 | |
---|
| 126 | BARRIER |
---|
| 127 | |
---|
| 128 | debug_hevi_solver=.FALSE. |
---|
| 129 | END_BLOCK |
---|
| 130 | |
---|
| 131 | #endif END_DYSL |
---|
| 132 | |
---|
| 133 | SUBROUTINE compute_NH_geopot_unst(tau, m_ik, m_il, theta, W_il, Phi_il) |
---|
| 134 | USE ISO_C_BINDING, only : C_DOUBLE, C_FLOAT |
---|
| 135 | USE disvert_mod, only : g,Treff,cpp,preff,kappa,caldyn_thermo,thermo_theta, & |
---|
| 136 | thermo_entropy |
---|
| 137 | USE disvert_mod, ONLY : ptop |
---|
[939] | 138 | USE data_unstructured_mod, ONLY : enter_trace, exit_trace, & |
---|
| 139 | id_NH_geopot,debug_hevi_solver_, & |
---|
| 140 | PHI_BOT,pbot,rho_bot |
---|
[878] | 141 | FIELD_MASS :: m_ik, theta ! IN*2 |
---|
| 142 | FIELD_GEOPOT :: m_il, W_il, Phi_il, Phi_star_il ! IN,INOUT*2, LOCAL |
---|
| 143 | NUM :: tau, gamma, tau2_g, tau2_g2, g2, gm2, vreff, Rd_preff |
---|
| 144 | INTEGER :: iter |
---|
| 145 | LOGICAL :: debug_hevi_solver |
---|
| 146 | DECLARE_INDICES |
---|
| 147 | NUM :: rho_ij, X_ij, Y_ij, wil, rho_c2_mik, c2_mik, ml_g2 |
---|
| 148 | #define COLUMN 0 |
---|
| 149 | #if COLUMN |
---|
| 150 | NUM1(llm) :: pk, Ak, Ck |
---|
| 151 | NUM1(llm+1):: Rl, Bl, Dl, xl |
---|
| 152 | #define p_ik(l,ij) pk(l) |
---|
| 153 | #define A_ik(l,ij) Ak(l) |
---|
| 154 | #define C_ik(l,ij) Ck(l) |
---|
| 155 | #define R_il(l,ij) Rl(l) |
---|
| 156 | #define B_il(l,ij) Bl(l) |
---|
| 157 | #define D_il(l,ij) Dl(l) |
---|
| 158 | #define x_il(l,ij) xl(l) |
---|
| 159 | #else |
---|
| 160 | FIELD_MASS :: p_ik, A_ik, C_ik |
---|
| 161 | FIELD_GEOPOT :: R_il, B_il, D_il, x_il |
---|
| 162 | #endif |
---|
| 163 | |
---|
| 164 | debug_hevi_solver=.FALSE. |
---|
| 165 | !$OMP MASTER |
---|
| 166 | debug_hevi_solver = debug_hevi_solver_ |
---|
| 167 | !$OMP END MASTER |
---|
| 168 | |
---|
| 169 | #define PHI_BOT(ij) Phi_bot |
---|
| 170 | START_TRACE(id_NH_geopot, 7,0,0) |
---|
| 171 | #include "../kernels_unst/compute_NH_geopot.k90" |
---|
| 172 | STOP_TRACE |
---|
| 173 | |
---|
| 174 | !$OMP MASTER |
---|
| 175 | debug_hevi_solver_ = debug_hevi_solver |
---|
| 176 | !$OMP END MASTER |
---|
| 177 | #undef PHI_BOT |
---|
| 178 | |
---|
| 179 | #if COLUMN |
---|
| 180 | #undef p_ik |
---|
| 181 | #undef A_ik |
---|
| 182 | #undef C_ik |
---|
| 183 | #undef R_il |
---|
| 184 | #undef B_il |
---|
| 185 | #undef D_il |
---|
| 186 | #undef x_il |
---|
| 187 | #endif |
---|
| 188 | #undef COLUMN |
---|
| 189 | END SUBROUTINE compute_NH_geopot_unst |
---|
| 190 | |
---|
[859] | 191 | SUBROUTINE compute_NH_geopot(tau, phis, m_ik, m_il, theta, W_il, Phi_il) |
---|
| 192 | USE icosa |
---|
| 193 | USE disvert_mod |
---|
| 194 | USE caldyn_vars_mod |
---|
| 195 | USE omp_para |
---|
| 196 | REAL(rstd), PARAMETER :: pbot=1e5, rho_bot=1e6 |
---|
| 197 | REAL(rstd),INTENT(IN) :: tau ! solve Phi-tau*dPhi/dt = Phi_rhs |
---|
| 198 | REAL(rstd),INTENT(IN) :: phis(iim*jjm) |
---|
| 199 | REAL(rstd),INTENT(IN) :: m_ik(iim*jjm,llm) |
---|
| 200 | REAL(rstd),INTENT(IN) :: m_il(iim*jjm,llm+1) |
---|
| 201 | REAL(rstd),INTENT(IN) :: theta(iim*jjm,llm) |
---|
| 202 | REAL(rstd),INTENT(IN) :: W_il(iim*jjm,llm+1) ! vertical momentum |
---|
| 203 | REAL(rstd),INTENT(INOUT) :: Phi_il(iim*jjm,llm+1) ! geopotential |
---|
| 204 | |
---|
| 205 | REAL(rstd) :: Phi_star_il(iim*jjm,llm+1) |
---|
| 206 | REAL(rstd) :: p_ik(iim*jjm,llm) ! pressure |
---|
| 207 | REAL(rstd) :: R_il(iim*jjm,llm+1) ! rhs of tridiag problem |
---|
| 208 | REAL(rstd) :: x_il(iim*jjm,llm+1) ! solution of tridiag problem |
---|
| 209 | REAL(rstd) :: A_ik(iim*jjm,llm) ! off-diagonal coefficients of tridiag problem |
---|
| 210 | REAL(rstd) :: B_il(iim*jjm,llm+1) ! diagonal coefficients of tridiag problem |
---|
| 211 | REAL(rstd) :: C_ik(iim*jjm,llm) ! Thomas algorithm |
---|
| 212 | REAL(rstd) :: D_il(iim*jjm,llm+1) ! Thomas algorithm |
---|
| 213 | REAL(rstd) :: gamma, rho_ij, X_ij, Y_ij |
---|
| 214 | REAL(rstd) :: wil, tau2_g, g2, gm2, ml_g2, c2_mik, vreff |
---|
| 215 | |
---|
| 216 | INTEGER :: iter, ij, l, ij_omp_begin_ext, ij_omp_end_ext |
---|
| 217 | |
---|
| 218 | CALL distrib_level(ij_begin_ext,ij_end_ext, ij_omp_begin_ext,ij_omp_end_ext) |
---|
| 219 | |
---|
| 220 | IF(dysl) THEN |
---|
| 221 | #define PHI_BOT(ij) phis(ij) |
---|
| 222 | #include "../kernels_hex/compute_NH_geopot.k90" |
---|
| 223 | #undef PHI_BOT |
---|
| 224 | ELSE |
---|
| 225 | ! FIXME : vertical OpenMP parallelism will not work |
---|
| 226 | |
---|
| 227 | tau2_g=tau*tau/g |
---|
| 228 | g2=g*g |
---|
| 229 | gm2 = g**-2 |
---|
| 230 | gamma = 1./(1.-kappa) |
---|
| 231 | |
---|
| 232 | ! compute Phi_star |
---|
| 233 | DO l=1,llm+1 |
---|
| 234 | !DIR$ SIMD |
---|
| 235 | DO ij=ij_begin_ext,ij_end_ext |
---|
| 236 | Phi_star_il(ij,l) = Phi_il(ij,l) + tau*g2*(W_il(ij,l)/m_il(ij,l)-tau) |
---|
| 237 | ENDDO |
---|
| 238 | ENDDO |
---|
| 239 | |
---|
| 240 | ! Newton-Raphson iteration : Phi_il contains current guess value |
---|
| 241 | DO iter=1,5 ! 2 iterations should be enough |
---|
| 242 | |
---|
| 243 | ! Compute pressure, A_ik |
---|
| 244 | DO l=1,llm |
---|
| 245 | !DIR$ SIMD |
---|
| 246 | DO ij=ij_begin_ext,ij_end_ext |
---|
| 247 | rho_ij = (g*m_ik(ij,l))/(Phi_il(ij,l+1)-Phi_il(ij,l)) |
---|
| 248 | X_ij = (cpp/preff)*kappa*theta(ij,l)*rho_ij |
---|
| 249 | p_ik(ij,l) = preff*(X_ij**gamma) |
---|
| 250 | c2_mik = gamma*p_ik(ij,l)/(rho_ij*m_ik(ij,l)) ! c^2 = gamma*R*T = gamma*p/rho |
---|
| 251 | A_ik(ij,l) = c2_mik*(tau/g*rho_ij)**2 |
---|
| 252 | ENDDO |
---|
| 253 | ENDDO |
---|
| 254 | |
---|
| 255 | ! Compute residual, B_il |
---|
| 256 | ! bottom interface l=1 |
---|
| 257 | !DIR$ SIMD |
---|
| 258 | DO ij=ij_begin_ext,ij_end_ext |
---|
| 259 | ml_g2 = gm2*m_il(ij,1) |
---|
| 260 | B_il(ij,1) = A_ik(ij,1) + ml_g2 + tau2_g*rho_bot |
---|
| 261 | R_il(ij,1) = ml_g2*( Phi_il(ij,1)-Phi_star_il(ij,1)) & |
---|
| 262 | + tau2_g*( p_ik(ij,1)-pbot+rho_bot*(Phi_il(ij,1)-phis(ij)) ) |
---|
| 263 | ENDDO |
---|
| 264 | ! inner interfaces |
---|
| 265 | DO l=2,llm |
---|
| 266 | !DIR$ SIMD |
---|
| 267 | DO ij=ij_begin_ext,ij_end_ext |
---|
| 268 | ml_g2 = gm2*m_il(ij,l) |
---|
| 269 | B_il(ij,l) = A_ik(ij,l)+A_ik(ij,l-1) + ml_g2 |
---|
| 270 | R_il(ij,l) = ml_g2*( Phi_il(ij,l)-Phi_star_il(ij,l)) & |
---|
| 271 | + tau2_g*(p_ik(ij,l)-p_ik(ij,l-1)) |
---|
| 272 | ! consistency check : if Wil=0 and initial state is in hydrostatic balance |
---|
| 273 | ! then Phi_star_il(ij,l) = Phi_il(ij,l) - tau^2*g^2 |
---|
| 274 | ! and residual = tau^2*(ml+(1/g)dl_pi)=0 |
---|
| 275 | ENDDO |
---|
| 276 | ENDDO |
---|
| 277 | ! top interface l=llm+1 |
---|
| 278 | !DIR$ SIMD |
---|
| 279 | DO ij=ij_begin_ext,ij_end_ext |
---|
| 280 | ml_g2 = gm2*m_il(ij,llm+1) |
---|
| 281 | B_il(ij,llm+1) = A_ik(ij,llm) + ml_g2 |
---|
| 282 | R_il(ij,llm+1) = ml_g2*( Phi_il(ij,llm+1)-Phi_star_il(ij,llm+1)) & |
---|
| 283 | + tau2_g*( ptop-p_ik(ij,llm) ) |
---|
| 284 | ENDDO |
---|
| 285 | |
---|
| 286 | ! FIXME later |
---|
| 287 | ! the lines below modify the tridiag problem |
---|
| 288 | ! for flat, rigid boundary conditions at top and bottom : |
---|
| 289 | ! zero out A(1), A(llm), R(1), R(llm+1) |
---|
| 290 | ! => x(l)=0 at l=1,llm+1 |
---|
| 291 | DO ij=ij_begin_ext,ij_end_ext |
---|
| 292 | A_ik(ij,1) = 0. |
---|
| 293 | A_ik(ij,llm) = 0. |
---|
| 294 | R_il(ij,1) = 0. |
---|
| 295 | R_il(ij,llm+1) = 0. |
---|
| 296 | ENDDO |
---|
| 297 | |
---|
| 298 | IF(debug_hevi_solver) THEN ! print Linf(residual) |
---|
| 299 | PRINT *, '[hevi_solver] R,p', iter, MAXVAL(ABS(R_il)), MAXVAL(p_ik) |
---|
| 300 | END IF |
---|
| 301 | |
---|
| 302 | ! Solve -A(l-1)x(l-1) + B(l)x(l) - A(l)x(l+1) = R(l) using Thomas algorithm |
---|
| 303 | ! Forward sweep : |
---|
| 304 | ! C(0)=0, C(l) = -A(l) / (B(l)+A(l-1)C(l-1)), |
---|
| 305 | ! D(0)=0, D(l) = (R(l)+A(l-1)D(l-1)) / (B(l)+A(l-1)C(l-1)) |
---|
| 306 | ! bottom interface l=1 |
---|
| 307 | !DIR$ SIMD |
---|
| 308 | DO ij=ij_begin_ext,ij_end_ext |
---|
| 309 | X_ij = 1./B_il(ij,1) |
---|
| 310 | C_ik(ij,1) = -A_ik(ij,1) * X_ij |
---|
| 311 | D_il(ij,1) = R_il(ij,1) * X_ij |
---|
| 312 | ENDDO |
---|
| 313 | ! inner interfaces/layers |
---|
| 314 | DO l=2,llm |
---|
| 315 | !DIR$ SIMD |
---|
| 316 | DO ij=ij_begin_ext,ij_end_ext |
---|
| 317 | X_ij = 1./(B_il(ij,l) + A_ik(ij,l-1)*C_ik(ij,l-1)) |
---|
| 318 | C_ik(ij,l) = -A_ik(ij,l) * X_ij |
---|
| 319 | D_il(ij,l) = (R_il(ij,l)+A_ik(ij,l-1)*D_il(ij,l-1)) * X_ij |
---|
| 320 | ENDDO |
---|
| 321 | ENDDO |
---|
| 322 | ! top interface l=llm+1 |
---|
| 323 | !DIR$ SIMD |
---|
| 324 | DO ij=ij_begin_ext,ij_end_ext |
---|
| 325 | X_ij = 1./(B_il(ij,llm+1) + A_ik(ij,llm)*C_ik(ij,llm)) |
---|
| 326 | D_il(ij,llm+1) = (R_il(ij,llm+1)+A_ik(ij,llm)*D_il(ij,llm)) * X_ij |
---|
| 327 | ENDDO |
---|
| 328 | |
---|
| 329 | ! Back substitution : |
---|
| 330 | ! x(i) = D(i)-C(i)x(i+1), x(N+1)=0 |
---|
| 331 | ! + Newton-Raphson update |
---|
| 332 | x_il=0. ! FIXME |
---|
| 333 | ! top interface l=llm+1 |
---|
| 334 | !DIR$ SIMD |
---|
| 335 | DO ij=ij_begin_ext,ij_end_ext |
---|
| 336 | x_il(ij,llm+1) = D_il(ij,llm+1) |
---|
| 337 | Phi_il(ij,llm+1) = Phi_il(ij,llm+1) - x_il(ij,llm+1) |
---|
| 338 | ENDDO |
---|
| 339 | ! lower interfaces |
---|
| 340 | DO l=llm,1,-1 |
---|
| 341 | !DIR$ SIMD |
---|
| 342 | DO ij=ij_begin_ext,ij_end_ext |
---|
| 343 | x_il(ij,l) = D_il(ij,l) - C_ik(ij,l)*x_il(ij,l+1) |
---|
| 344 | Phi_il(ij,l) = Phi_il(ij,l) - x_il(ij,l) |
---|
| 345 | ENDDO |
---|
| 346 | ENDDO |
---|
| 347 | |
---|
| 348 | IF(debug_hevi_solver) THEN |
---|
| 349 | PRINT *, '[hevi_solver] A,B', iter, MAXVAL(ABS(A_ik)),MAXVAL(ABS(B_il)) |
---|
| 350 | PRINT *, '[hevi_solver] C,D', iter, MAXVAL(ABS(C_ik)),MAXVAL(ABS(D_il)) |
---|
| 351 | DO l=1,llm+1 |
---|
| 352 | WRITE(*,'(A,I2.1,I3.2,E9.2)') '[hevi_solver] x', iter,l, MAXVAL(ABS(x_il(:,l))) |
---|
| 353 | END DO |
---|
| 354 | END IF |
---|
| 355 | |
---|
| 356 | END DO ! Newton-Raphson |
---|
| 357 | |
---|
| 358 | END IF ! dysl |
---|
| 359 | |
---|
| 360 | END SUBROUTINE compute_NH_geopot |
---|
| 361 | |
---|
| 362 | END MODULE compute_NH_geopot_mod |
---|