[358] | 1 | MODULE dynspg_flt |
---|
| 2 | !!====================================================================== |
---|
| 3 | !! *** MODULE dynspg_flt *** |
---|
| 4 | !! Ocean dynamics: surface pressure gradient trend |
---|
| 5 | !!====================================================================== |
---|
[508] | 6 | !! History 8.0 ! 98-05 (G. Roullet) free surface |
---|
| 7 | !! ! 98-10 (G. Madec, M. Imbard) release 8.2 |
---|
| 8 | !! 8.5 ! 02-08 (G. Madec) F90: Free form and module |
---|
| 9 | !! " " ! 02-11 (C. Talandier, A-M Treguier) Open boundaries |
---|
| 10 | !! 9.0 ! 04-08 (C. Talandier) New trends organization |
---|
| 11 | !! " " ! 05-11 (V. Garnier) Surface pressure gradient organization |
---|
| 12 | !! " " ! 06-07 (S. Masson) distributed restart using iom |
---|
[358] | 13 | !!---------------------------------------------------------------------- |
---|
[575] | 14 | #if defined key_dynspg_flt || defined key_esopa |
---|
[508] | 15 | !!---------------------------------------------------------------------- |
---|
[358] | 16 | !! 'key_dynspg_flt' filtered free surface |
---|
| 17 | !!---------------------------------------------------------------------- |
---|
[508] | 18 | !!---------------------------------------------------------------------- |
---|
[358] | 19 | !! dyn_spg_flt : update the momentum trend with the surface pressure |
---|
| 20 | !! gradient in the filtered free surface case with |
---|
| 21 | !! vector optimization |
---|
[508] | 22 | !! flt_rst : read/write the time-splitting restart fields in the ocean restart file |
---|
[358] | 23 | !!---------------------------------------------------------------------- |
---|
| 24 | USE oce ! ocean dynamics and tracers |
---|
| 25 | USE dom_oce ! ocean space and time domain |
---|
| 26 | USE zdf_oce ! ocean vertical physics |
---|
[719] | 27 | USE phycst ! physical constants |
---|
| 28 | USE ocesbc ! ocean surface boundary condition |
---|
| 29 | USE flxrnf ! ocean runoffs |
---|
[708] | 30 | USE sol_oce ! ocean elliptic solver |
---|
[508] | 31 | USE solver ! solver initialization |
---|
[358] | 32 | USE solpcg ! preconditionned conjugate gradient solver |
---|
| 33 | USE solsor ! Successive Over-relaxation solver |
---|
| 34 | USE solfet ! FETI solver |
---|
[719] | 35 | USE obc_oce ! Lateral open boundary condition |
---|
[358] | 36 | USE obcdyn ! ocean open boundary condition (obc_dyn routines) |
---|
| 37 | USE obcvol ! ocean open boundary condition (obc_vol routines) |
---|
| 38 | USE lib_mpp ! distributed memory computing library |
---|
| 39 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
[719] | 40 | USE cla_dynspg ! cross land advection |
---|
[358] | 41 | USE prtctl ! Print control |
---|
[413] | 42 | USE solmat ! matrix construction for elliptic solvers |
---|
[389] | 43 | USE agrif_opa_interp |
---|
[719] | 44 | USE in_out_manager ! I/O manager |
---|
[508] | 45 | USE iom |
---|
| 46 | USE restart ! only for lrst_oce |
---|
[719] | 47 | USE domvvl ! variable volume |
---|
[358] | 48 | |
---|
| 49 | IMPLICIT NONE |
---|
| 50 | PRIVATE |
---|
| 51 | |
---|
| 52 | PUBLIC dyn_spg_flt ! routine called by step.F90 |
---|
| 53 | |
---|
| 54 | !! * Substitutions |
---|
| 55 | # include "domzgr_substitute.h90" |
---|
| 56 | # include "vectopt_loop_substitute.h90" |
---|
| 57 | !!---------------------------------------------------------------------- |
---|
| 58 | !! OPA 9.0 , LOCEAN-IPSL (2005) |
---|
[719] | 59 | !! $Header: /home/opalod/NEMOCVSROOT/NEMO/OPA_SRC/DYN/dynspg_flt.F90,v 1.14 2007/06/05 10:38:27 opalod Exp $ |
---|
[508] | 60 | !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt) |
---|
[358] | 61 | !!---------------------------------------------------------------------- |
---|
| 62 | |
---|
| 63 | CONTAINS |
---|
| 64 | |
---|
| 65 | SUBROUTINE dyn_spg_flt( kt, kindic ) |
---|
| 66 | !!---------------------------------------------------------------------- |
---|
| 67 | !! *** routine dyn_spg_flt *** |
---|
| 68 | !! |
---|
| 69 | !! ** Purpose : Compute the now trend due to the surface pressure |
---|
| 70 | !! gradient in case of filtered free surface formulation and add |
---|
| 71 | !! it to the general trend of momentum equation. |
---|
| 72 | !! |
---|
| 73 | !! ** Method : Filtered free surface formulation. The surface |
---|
| 74 | !! pressure gradient is given by: |
---|
| 75 | !! spgu = 1/rau0 d/dx(ps) = 1/e1u di( sshn + rnu btda ) |
---|
| 76 | !! spgv = 1/rau0 d/dy(ps) = 1/e2v dj( sshn + rnu btda ) |
---|
| 77 | !! where sshn is the free surface elevation and btda is the after |
---|
| 78 | !! of the free surface elevation |
---|
| 79 | !! -1- compute the after sea surface elevation from the kinematic |
---|
| 80 | !! surface boundary condition: |
---|
| 81 | !! zssha = sshb + 2 rdt ( wn - emp ) |
---|
| 82 | !! Time filter applied on now sea surface elevation to avoid |
---|
| 83 | !! the divergence of two consecutive time-steps and swap of free |
---|
| 84 | !! surface arrays to start the next time step: |
---|
| 85 | !! sshb = sshn + atfp * [ sshb + zssha - 2 sshn ] |
---|
| 86 | !! sshn = zssha |
---|
| 87 | !! -2- evaluate the surface presure trend (including the addi- |
---|
| 88 | !! tional force) in three steps: |
---|
| 89 | !! a- compute the right hand side of the elliptic equation: |
---|
| 90 | !! gcb = 1/(e1t e2t) [ di(e2u spgu) + dj(e1v spgv) ] |
---|
| 91 | !! where (spgu,spgv) are given by: |
---|
| 92 | !! spgu = vertical sum[ e3u (ub+ 2 rdt ua ) ] |
---|
| 93 | !! - grav 2 rdt hu /e1u di[sshn + emp] |
---|
| 94 | !! spgv = vertical sum[ e3v (vb+ 2 rdt va) ] |
---|
| 95 | !! - grav 2 rdt hv /e2v dj[sshn + emp] |
---|
| 96 | !! and define the first guess from previous computation : |
---|
| 97 | !! zbtd = btda |
---|
| 98 | !! btda = 2 zbtd - btdb |
---|
| 99 | !! btdb = zbtd |
---|
| 100 | !! b- compute the relative accuracy to be reached by the |
---|
| 101 | !! iterative solver |
---|
| 102 | !! c- apply the solver by a call to sol... routine |
---|
| 103 | !! -3- compute and add the free surface pressure gradient inclu- |
---|
| 104 | !! ding the additional force used to stabilize the equation. |
---|
| 105 | !! |
---|
| 106 | !! ** Action : - Update (ua,va) with the surf. pressure gradient trend |
---|
| 107 | !! |
---|
[508] | 108 | !! References : Roullet and Madec 1999, JGR. |
---|
[358] | 109 | !!--------------------------------------------------------------------- |
---|
[592] | 110 | !! * Modules used |
---|
| 111 | USE oce , ONLY : zub => ta, & ! ta used as workspace |
---|
| 112 | zvb => sa ! sa " " |
---|
| 113 | |
---|
[358] | 114 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
[508] | 115 | INTEGER, INTENT( out ) :: kindic ! solver convergence flag (<0 if not converge) |
---|
| 116 | !! |
---|
| 117 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 118 | REAL(wp) :: z2dt, z2dtg, zraur, znugdt, & ! temporary scalars |
---|
[592] | 119 | & znurau, zgcb, zbtd, & ! " " |
---|
[508] | 120 | & ztdgu, ztdgv ! " " |
---|
[592] | 121 | !! Variable volume |
---|
| 122 | REAL(wp), DIMENSION(jpi,jpj) :: & |
---|
| 123 | zsshub, zsshua, zsshvb, zsshva, zssha ! 2D workspace |
---|
| 124 | REAL(wp), DIMENSION(jpi,jpj,jpk) :: & ! 3D workspace |
---|
| 125 | zfse3ub, zfse3ua, zfse3vb, zfse3va ! 3D workspace |
---|
[358] | 126 | !!---------------------------------------------------------------------- |
---|
[508] | 127 | ! |
---|
[358] | 128 | IF( kt == nit000 ) THEN |
---|
| 129 | IF(lwp) WRITE(numout,*) |
---|
| 130 | IF(lwp) WRITE(numout,*) 'dyn_spg_flt : surface pressure gradient trend' |
---|
| 131 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~~ (free surface constant volume case)' |
---|
| 132 | |
---|
| 133 | ! set to zero free surface specific arrays |
---|
| 134 | spgu(:,:) = 0.e0 ! surface pressure gradient (i-direction) |
---|
| 135 | spgv(:,:) = 0.e0 ! surface pressure gradient (j-direction) |
---|
[508] | 136 | CALL solver_init( nit000 ) ! Elliptic solver initialisation |
---|
| 137 | |
---|
| 138 | ! read filtered free surface arrays in restart file |
---|
| 139 | CALL flt_rst( nit000, 'READ' ) ! read or initialize the following fields: |
---|
| 140 | ! ! gcx, gcxb, sshb, sshn |
---|
[358] | 141 | ENDIF |
---|
| 142 | |
---|
| 143 | ! Local constant initialization |
---|
| 144 | z2dt = 2. * rdt ! time step: leap-frog |
---|
| 145 | IF( neuler == 0 .AND. kt == nit000 ) z2dt = rdt ! time step: Euler if restart from rest |
---|
[413] | 146 | IF( neuler == 0 .AND. kt == nit000+1 ) CALL sol_mat(kt) |
---|
[358] | 147 | z2dtg = grav * z2dt |
---|
| 148 | zraur = 1. / rauw |
---|
| 149 | znugdt = rnu * grav * z2dt |
---|
| 150 | znurau = znugdt * zraur |
---|
| 151 | |
---|
[592] | 152 | !! Explicit physics with thickness weighted updates |
---|
| 153 | IF( lk_vvl ) THEN ! variable volume |
---|
[358] | 154 | |
---|
[592] | 155 | DO jj = 1, jpjm1 |
---|
| 156 | DO ji = 1,jpim1 |
---|
| 157 | |
---|
| 158 | ! Sea Surface Height at u-point before |
---|
| 159 | zsshub(ji,jj) = 0.5 * umask(ji,jj,1) / ( e1u(ji,jj) * e2u(ji,jj) ) & |
---|
| 160 | & * ( e1t(ji ,jj) * e2t(ji ,jj) * sshb(ji ,jj) & |
---|
| 161 | & + e1t(ji+1,jj) * e2t(ji+1,jj) * sshb(ji+1,jj) ) |
---|
| 162 | |
---|
| 163 | ! Sea Surface Height at v-point before |
---|
| 164 | zsshvb(ji,jj) = 0.5 * vmask(ji,jj,1) / ( e1v(ji,jj) * e2v(ji,jj) ) & |
---|
| 165 | & * ( e1t(ji,jj ) * e2t(ji,jj ) * sshb(ji,jj ) & |
---|
| 166 | & + e1t(ji,jj+1) * e2t(ji,jj+1) * sshb(ji,jj+1) ) |
---|
| 167 | |
---|
| 168 | ! Sea Surface Height at u-point after |
---|
| 169 | zsshua(ji,jj) = 0.5 * umask(ji,jj,1) / ( e1u(ji,jj) * e2u(ji,jj) ) & |
---|
| 170 | & * ( e1t(ji ,jj) * e2t(ji ,jj) * ssha(ji ,jj) & |
---|
| 171 | & + e1t(ji+1,jj) * e2t(ji+1,jj) * ssha(ji+1,jj) ) |
---|
| 172 | |
---|
| 173 | ! Sea Surface Height at v-point after |
---|
| 174 | zsshva(ji,jj) = 0.5 * vmask(ji,jj,1) / ( e1v(ji,jj) * e2v(ji,jj) ) & |
---|
| 175 | & * ( e1t(ji,jj ) * e2t(ji,jj ) * ssha(ji,jj ) & |
---|
| 176 | & + e1t(ji,jj+1) * e2t(ji,jj+1) * ssha(ji,jj+1) ) |
---|
| 177 | |
---|
[358] | 178 | END DO |
---|
| 179 | END DO |
---|
| 180 | |
---|
[592] | 181 | ! Boundaries conditions |
---|
| 182 | CALL lbc_lnk( zsshub, 'U', 1. ) ; CALL lbc_lnk( zsshvb, 'V', 1. ) |
---|
| 183 | CALL lbc_lnk( zsshua, 'U', 1. ) ; CALL lbc_lnk( zsshva, 'V', 1. ) |
---|
| 184 | |
---|
| 185 | ! Scale factors at before and after time step |
---|
| 186 | ! ------------------------------------------- |
---|
[661] | 187 | CALL dom_vvl_sf( zsshub, 'U', zfse3ub ) ; CALL dom_vvl_sf( zsshua, 'U', zfse3ua ) |
---|
| 188 | CALL dom_vvl_sf( zsshvb, 'V', zfse3vb ) ; CALL dom_vvl_sf( zsshva, 'V', zfse3va ) |
---|
[592] | 189 | |
---|
[661] | 190 | |
---|
[592] | 191 | ! Thickness weighting |
---|
| 192 | ! ------------------- |
---|
[642] | 193 | DO jk = 1, jpkm1 |
---|
| 194 | DO jj = 1, jpj |
---|
| 195 | DO ji = 1, jpi |
---|
| 196 | ua(ji,jj,jk) = ua(ji,jj,jk) * fse3u(ji,jj,jk) |
---|
| 197 | va(ji,jj,jk) = va(ji,jj,jk) * fse3v(ji,jj,jk) |
---|
[592] | 198 | |
---|
[642] | 199 | zub(ji,jj,jk) = ub(ji,jj,jk) * zfse3ub(ji,jj,jk) |
---|
| 200 | zvb(ji,jj,jk) = vb(ji,jj,jk) * zfse3vb(ji,jj,jk) |
---|
| 201 | END DO |
---|
| 202 | END DO |
---|
| 203 | END DO |
---|
[592] | 204 | |
---|
| 205 | ! Evaluate the masked next velocity (effect of the additional force not included) |
---|
| 206 | DO jk = 1, jpkm1 |
---|
| 207 | DO jj = 2, jpjm1 |
---|
| 208 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 209 | ua(ji,jj,jk) = ( zub(ji,jj,jk) + z2dt * ua(ji,jj,jk) ) /zfse3ua(ji,jj,jk) * umask(ji,jj,jk) |
---|
| 210 | va(ji,jj,jk) = ( zvb(ji,jj,jk) + z2dt * va(ji,jj,jk) ) /zfse3va(ji,jj,jk) * vmask(ji,jj,jk) |
---|
| 211 | END DO |
---|
| 212 | END DO |
---|
| 213 | END DO |
---|
| 214 | |
---|
| 215 | ELSE ! fixed volume |
---|
| 216 | |
---|
| 217 | ! Surface pressure gradient (now) |
---|
[358] | 218 | DO jj = 2, jpjm1 |
---|
| 219 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[592] | 220 | spgu(ji,jj) = - grav * ( sshn(ji+1,jj) - sshn(ji,jj) ) / e1u(ji,jj) |
---|
| 221 | spgv(ji,jj) = - grav * ( sshn(ji,jj+1) - sshn(ji,jj) ) / e2v(ji,jj) |
---|
| 222 | END DO |
---|
| 223 | END DO |
---|
| 224 | |
---|
| 225 | ! Add the surface pressure trend to the general trend |
---|
| 226 | DO jk = 1, jpkm1 |
---|
| 227 | DO jj = 2, jpjm1 |
---|
| 228 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 229 | ua(ji,jj,jk) = ua(ji,jj,jk) + spgu(ji,jj) |
---|
| 230 | va(ji,jj,jk) = va(ji,jj,jk) + spgv(ji,jj) |
---|
| 231 | END DO |
---|
[358] | 232 | END DO |
---|
| 233 | END DO |
---|
[389] | 234 | |
---|
[592] | 235 | ! Evaluate the masked next velocity (effect of the additional force not included) |
---|
| 236 | DO jk = 1, jpkm1 |
---|
| 237 | DO jj = 2, jpjm1 |
---|
| 238 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 239 | ua(ji,jj,jk) = ( ub(ji,jj,jk) + z2dt * ua(ji,jj,jk) ) * umask(ji,jj,jk) |
---|
| 240 | va(ji,jj,jk) = ( vb(ji,jj,jk) + z2dt * va(ji,jj,jk) ) * vmask(ji,jj,jk) |
---|
| 241 | END DO |
---|
| 242 | END DO |
---|
| 243 | END DO |
---|
| 244 | |
---|
| 245 | ENDIF |
---|
| 246 | |
---|
[358] | 247 | #if defined key_obc |
---|
[508] | 248 | CALL obc_dyn( kt ) ! Update velocities on each open boundary with the radiation algorithm |
---|
| 249 | CALL obc_vol( kt ) ! Correction of the barotropic componant velocity to control the volume of the system |
---|
[358] | 250 | #endif |
---|
[392] | 251 | #if defined key_agrif |
---|
[508] | 252 | CALL Agrif_dyn( kt ) ! Update velocities on each coarse/fine interfaces |
---|
[389] | 253 | #endif |
---|
[358] | 254 | #if defined key_orca_r2 |
---|
| 255 | IF( n_cla == 1 ) CALL dyn_spg_cla( kt ) ! Cross Land Advection (update (ua,va)) |
---|
| 256 | #endif |
---|
| 257 | |
---|
| 258 | ! compute the next vertically averaged velocity (effect of the additional force not included) |
---|
| 259 | ! --------------------------------------------- |
---|
| 260 | DO jj = 2, jpjm1 |
---|
| 261 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 262 | spgu(ji,jj) = 0.e0 |
---|
| 263 | spgv(ji,jj) = 0.e0 |
---|
| 264 | END DO |
---|
| 265 | END DO |
---|
| 266 | |
---|
| 267 | ! vertical sum |
---|
| 268 | !CDIR NOLOOPCHG |
---|
| 269 | IF( lk_vopt_loop ) THEN ! vector opt., forced unroll |
---|
| 270 | DO jk = 1, jpkm1 |
---|
| 271 | DO ji = 1, jpij |
---|
| 272 | spgu(ji,1) = spgu(ji,1) + fse3u(ji,1,jk) * ua(ji,1,jk) |
---|
| 273 | spgv(ji,1) = spgv(ji,1) + fse3v(ji,1,jk) * va(ji,1,jk) |
---|
| 274 | END DO |
---|
| 275 | END DO |
---|
| 276 | ELSE ! No vector opt. |
---|
| 277 | DO jk = 1, jpkm1 |
---|
| 278 | DO jj = 2, jpjm1 |
---|
| 279 | DO ji = 2, jpim1 |
---|
| 280 | spgu(ji,jj) = spgu(ji,jj) + fse3u(ji,jj,jk) * ua(ji,jj,jk) |
---|
| 281 | spgv(ji,jj) = spgv(ji,jj) + fse3v(ji,jj,jk) * va(ji,jj,jk) |
---|
| 282 | END DO |
---|
| 283 | END DO |
---|
| 284 | END DO |
---|
| 285 | ENDIF |
---|
| 286 | |
---|
| 287 | ! transport: multiplied by the horizontal scale factor |
---|
| 288 | DO jj = 2, jpjm1 |
---|
| 289 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 290 | spgu(ji,jj) = spgu(ji,jj) * e2u(ji,jj) |
---|
| 291 | spgv(ji,jj) = spgv(ji,jj) * e1v(ji,jj) |
---|
| 292 | END DO |
---|
| 293 | END DO |
---|
| 294 | |
---|
| 295 | ! Boundary conditions on (spgu,spgv) |
---|
| 296 | CALL lbc_lnk( spgu, 'U', -1. ) |
---|
| 297 | CALL lbc_lnk( spgv, 'V', -1. ) |
---|
| 298 | |
---|
[592] | 299 | IF( lk_vvl ) CALL sol_mat( kt ) |
---|
| 300 | |
---|
[358] | 301 | ! Right hand side of the elliptic equation and first guess |
---|
| 302 | ! ----------------------------------------------------------- |
---|
| 303 | DO jj = 2, jpjm1 |
---|
| 304 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 305 | ! Divergence of the after vertically averaged velocity |
---|
| 306 | zgcb = spgu(ji,jj) - spgu(ji-1,jj) & |
---|
| 307 | + spgv(ji,jj) - spgv(ji,jj-1) |
---|
| 308 | gcb(ji,jj) = gcdprc(ji,jj) * zgcb |
---|
| 309 | ! First guess of the after barotropic transport divergence |
---|
| 310 | zbtd = gcx(ji,jj) |
---|
| 311 | gcx (ji,jj) = 2. * zbtd - gcxb(ji,jj) |
---|
| 312 | gcxb(ji,jj) = zbtd |
---|
| 313 | END DO |
---|
| 314 | END DO |
---|
| 315 | ! applied the lateral boundary conditions |
---|
[784] | 316 | IF( nsolv == 2 .AND. MAX( jpr2di, jpr2dj ) > 0 ) CALL lbc_lnk_e( gcb, c_solver_pt, 1. ) |
---|
[358] | 317 | |
---|
[392] | 318 | #if defined key_agrif |
---|
[413] | 319 | IF( .NOT. AGRIF_ROOT() ) THEN |
---|
[389] | 320 | ! add contribution of gradient of after barotropic transport divergence |
---|
[508] | 321 | IF( nbondi == -1 .OR. nbondi == 2 ) gcb(3 ,:) = & |
---|
| 322 | & gcb(3 ,:) - znugdt * z2dt * laplacu(2 ,:) * gcdprc(3 ,:) * hu(2 ,:) * e2u(2 ,:) |
---|
| 323 | IF( nbondi == 1 .OR. nbondi == 2 ) gcb(nlci-2,:) = & |
---|
| 324 | & gcb(nlci-2,:) + znugdt * z2dt * laplacu(nlci-2,:) * gcdprc(nlci-2,:) * hu(nlci-2,:) * e2u(nlci-2,:) |
---|
| 325 | IF( nbondj == -1 .OR. nbondj == 2 ) gcb(: ,3) = & |
---|
| 326 | & gcb(:,3 ) - znugdt * z2dt * laplacv(:,2 ) * gcdprc(:,3 ) * hv(:,2 ) * e1v(:,2 ) |
---|
| 327 | IF( nbondj == 1 .OR. nbondj == 2 ) gcb(:,nlcj-2) = & |
---|
| 328 | & gcb(:,nlcj-2) + znugdt * z2dt * laplacv(:,nlcj-2) * gcdprc(:,nlcj-2) * hv(:,nlcj-2) * e1v(:,nlcj-2) |
---|
[413] | 329 | ENDIF |
---|
[389] | 330 | #endif |
---|
| 331 | |
---|
| 332 | |
---|
[358] | 333 | ! Relative precision (computation on one processor) |
---|
| 334 | ! ------------------ |
---|
| 335 | rnorme =0. |
---|
| 336 | rnorme = SUM( gcb(1:jpi,1:jpj) * gcdmat(1:jpi,1:jpj) * gcb(1:jpi,1:jpj) * bmask(:,:) ) |
---|
| 337 | IF( lk_mpp ) CALL mpp_sum( rnorme ) ! sum over the global domain |
---|
| 338 | |
---|
| 339 | epsr = eps * eps * rnorme |
---|
| 340 | ncut = 0 |
---|
[508] | 341 | ! if rnorme is 0, the solution is 0, the solver is not called |
---|
[358] | 342 | IF( rnorme == 0.e0 ) THEN |
---|
| 343 | gcx(:,:) = 0.e0 |
---|
| 344 | res = 0.e0 |
---|
| 345 | niter = 0 |
---|
| 346 | ncut = 999 |
---|
| 347 | ENDIF |
---|
| 348 | |
---|
| 349 | ! Evaluate the next transport divergence |
---|
| 350 | ! -------------------------------------- |
---|
| 351 | ! Iterarive solver for the elliptic equation (except IF sol.=0) |
---|
| 352 | ! (output in gcx with boundary conditions applied) |
---|
| 353 | kindic = 0 |
---|
| 354 | IF( ncut == 0 ) THEN |
---|
| 355 | IF( nsolv == 1 ) THEN ! diagonal preconditioned conjuguate gradient |
---|
| 356 | CALL sol_pcg( kindic ) |
---|
| 357 | ELSEIF( nsolv == 2 ) THEN ! successive-over-relaxation |
---|
| 358 | CALL sol_sor( kindic ) |
---|
| 359 | ELSEIF( nsolv == 3 ) THEN ! FETI solver |
---|
| 360 | CALL sol_fet( kindic ) |
---|
| 361 | ELSE ! e r r o r in nsolv namelist parameter |
---|
[474] | 362 | WRITE(ctmp1,*) ' ~~~~~~~~~~~ not = ', nsolv |
---|
[784] | 363 | CALL ctl_stop( ' dyn_spg_flt : e r r o r, nsolv = 1, 2 or 3', ctmp1 ) |
---|
[358] | 364 | ENDIF |
---|
| 365 | ENDIF |
---|
| 366 | |
---|
| 367 | ! Transport divergence gradient multiplied by z2dt |
---|
| 368 | ! --------------------------------------------==== |
---|
| 369 | DO jj = 2, jpjm1 |
---|
| 370 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 371 | ! trend of Transport divergence gradient |
---|
| 372 | ztdgu = znugdt * (gcx(ji+1,jj ) - gcx(ji,jj) ) / e1u(ji,jj) |
---|
| 373 | ztdgv = znugdt * (gcx(ji ,jj+1) - gcx(ji,jj) ) / e2v(ji,jj) |
---|
| 374 | ! multiplied by z2dt |
---|
| 375 | #if defined key_obc |
---|
| 376 | ! caution : grad D = 0 along open boundaries |
---|
| 377 | spgu(ji,jj) = z2dt * ztdgu * obcumask(ji,jj) |
---|
| 378 | spgv(ji,jj) = z2dt * ztdgv * obcvmask(ji,jj) |
---|
| 379 | #else |
---|
| 380 | spgu(ji,jj) = z2dt * ztdgu |
---|
| 381 | spgv(ji,jj) = z2dt * ztdgv |
---|
| 382 | #endif |
---|
| 383 | END DO |
---|
| 384 | END DO |
---|
| 385 | |
---|
[392] | 386 | #if defined key_agrif |
---|
[413] | 387 | IF( .NOT. Agrif_Root() ) THEN |
---|
| 388 | ! caution : grad D (fine) = grad D (coarse) at coarse/fine interface |
---|
[508] | 389 | IF( nbondi == -1 .OR. nbondi == 2 ) spgu(2 ,:) = znugdt * z2dt * laplacu(2 ,:) * umask(2 ,:,1) |
---|
| 390 | IF( nbondi == 1 .OR. nbondi == 2 ) spgu(nlci-2,:) = znugdt * z2dt * laplacu(nlci-2,:) * umask(nlci-2,:,1) |
---|
| 391 | IF( nbondj == -1 .OR. nbondj == 2 ) spgv(:,2 ) = znugdt * z2dt * laplacv(:,2 ) * vmask(: ,2,1) |
---|
| 392 | IF( nbondj == 1 .OR. nbondj == 2 ) spgv(:,nlcj-2) = znugdt * z2dt * laplacv(:,nlcj-2) * vmask(:,nlcj-2,1) |
---|
[389] | 393 | ENDIF |
---|
| 394 | #endif |
---|
| 395 | ! 7. Add the trends multiplied by z2dt to the after velocity |
---|
| 396 | ! ----------------------------------------------------------- |
---|
[358] | 397 | ! ( c a u t i o n : (ua,va) here are the after velocity not the |
---|
| 398 | ! trend, the leap-frog time stepping will not |
---|
[508] | 399 | ! be done in dynnxt.F90 routine) |
---|
[358] | 400 | DO jk = 1, jpkm1 |
---|
| 401 | DO jj = 2, jpjm1 |
---|
| 402 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 403 | ua(ji,jj,jk) = (ua(ji,jj,jk) + spgu(ji,jj)) * umask(ji,jj,jk) |
---|
| 404 | va(ji,jj,jk) = (va(ji,jj,jk) + spgv(ji,jj)) * vmask(ji,jj,jk) |
---|
| 405 | END DO |
---|
| 406 | END DO |
---|
| 407 | END DO |
---|
| 408 | |
---|
| 409 | ! Sea surface elevation time stepping |
---|
| 410 | ! ----------------------------------- |
---|
[592] | 411 | IF( lk_vvl ) THEN ! after free surface elevation |
---|
| 412 | zssha(:,:) = ssha(:,:) |
---|
| 413 | ELSE |
---|
| 414 | zssha(:,:) = sshb(:,:) + z2dt * ( wn(:,:,1) - emp(:,:) * zraur ) * tmask(:,:,1) |
---|
| 415 | ENDIF |
---|
| 416 | |
---|
[358] | 417 | IF( neuler == 0 .AND. kt == nit000 ) THEN ! Euler (forward) time stepping, no time filter |
---|
[592] | 418 | ! swap of arrays |
---|
| 419 | sshb(:,:) = sshn (:,:) |
---|
| 420 | sshn(:,:) = zssha(:,:) |
---|
[358] | 421 | ELSE ! Leap-frog time stepping and time filter |
---|
[592] | 422 | ! time filter and array swap |
---|
| 423 | sshb(:,:) = atfp * ( sshb(:,:) + zssha(:,:) ) + atfp1 * sshn(:,:) |
---|
| 424 | sshn(:,:) = zssha(:,:) |
---|
[358] | 425 | ENDIF |
---|
| 426 | |
---|
[508] | 427 | ! write filtered free surface arrays in restart file |
---|
| 428 | ! -------------------------------------------------- |
---|
| 429 | IF( lrst_oce ) CALL flt_rst( kt, 'WRITE' ) |
---|
[358] | 430 | |
---|
[508] | 431 | ! print sum trends (used for debugging) |
---|
| 432 | IF(ln_ctl) CALL prt_ctl( tab2d_1=sshn, clinfo1=' spg - ssh: ', mask1=tmask ) |
---|
| 433 | ! |
---|
[358] | 434 | END SUBROUTINE dyn_spg_flt |
---|
| 435 | |
---|
[508] | 436 | |
---|
| 437 | SUBROUTINE flt_rst( kt, cdrw ) |
---|
| 438 | !!--------------------------------------------------------------------- |
---|
| 439 | !! *** ROUTINE ts_rst *** |
---|
| 440 | !! |
---|
| 441 | !! ** Purpose : Read or write filtered free surface arrays in restart file |
---|
| 442 | !!---------------------------------------------------------------------- |
---|
| 443 | INTEGER , INTENT(in) :: kt ! ocean time-step |
---|
| 444 | CHARACTER(len=*), INTENT(in) :: cdrw ! "READ"/"WRITE" flag |
---|
| 445 | !!---------------------------------------------------------------------- |
---|
| 446 | |
---|
| 447 | IF( TRIM(cdrw) == 'READ' ) THEN |
---|
[746] | 448 | IF( iom_varid( numror, 'gcx', ldstop = .FALSE. ) > 0 ) THEN |
---|
[508] | 449 | ! Caution : extra-hallow |
---|
| 450 | ! gcx and gcxb are defined as: DIMENSION(1-jpr2di:jpi+jpr2di,1-jpr2dj:jpj+jpr2dj) |
---|
[683] | 451 | CALL iom_get( numror, jpdom_autoglo, 'gcx' , gcx (1:jpi,1:jpj) ) |
---|
| 452 | CALL iom_get( numror, jpdom_autoglo, 'gcxb', gcxb(1:jpi,1:jpj) ) |
---|
| 453 | CALL iom_get( numror, jpdom_autoglo, 'sshb', sshb(:,:) ) |
---|
| 454 | CALL iom_get( numror, jpdom_autoglo, 'sshn', sshn(:,:) ) |
---|
[508] | 455 | IF( neuler == 0 ) THEN |
---|
| 456 | sshb(:,:) = sshn(:,:) |
---|
| 457 | gcxb(:,:) = gcx (:,:) |
---|
| 458 | ENDIF |
---|
| 459 | ELSE |
---|
| 460 | gcx (:,:) = 0.e0 |
---|
| 461 | gcxb(:,:) = 0.e0 |
---|
[558] | 462 | IF( nn_rstssh == 1 ) THEN |
---|
| 463 | sshb(:,:) = 0.e0 |
---|
| 464 | sshn(:,:) = 0.e0 |
---|
| 465 | ENDIF |
---|
[508] | 466 | ENDIF |
---|
| 467 | ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN |
---|
| 468 | ! Caution : extra-hallow |
---|
| 469 | ! gcx and gcxb are defined as: DIMENSION(1-jpr2di:jpi+jpr2di,1-jpr2dj:jpj+jpr2dj) |
---|
| 470 | CALL iom_rstput( kt, nitrst, numrow, 'gcx' , gcx( 1:jpi,1:jpj) ) |
---|
| 471 | CALL iom_rstput( kt, nitrst, numrow, 'gcxb', gcxb(1:jpi,1:jpj) ) |
---|
| 472 | CALL iom_rstput( kt, nitrst, numrow, 'sshb', sshb(:,:) ) |
---|
| 473 | CALL iom_rstput( kt, nitrst, numrow, 'sshn', sshn(:,:) ) |
---|
| 474 | ENDIF |
---|
| 475 | ! |
---|
| 476 | END SUBROUTINE flt_rst |
---|
| 477 | |
---|
[358] | 478 | #else |
---|
| 479 | !!---------------------------------------------------------------------- |
---|
| 480 | !! Default case : Empty module No standart free surface cst volume |
---|
| 481 | !!---------------------------------------------------------------------- |
---|
| 482 | CONTAINS |
---|
| 483 | SUBROUTINE dyn_spg_flt( kt, kindic ) ! Empty routine |
---|
| 484 | WRITE(*,*) 'dyn_spg_flt: You should not have seen this print! error?', kt, kindic |
---|
| 485 | END SUBROUTINE dyn_spg_flt |
---|
[657] | 486 | SUBROUTINE flt_rst ( kt, cdrw ) ! Empty routine |
---|
| 487 | INTEGER , INTENT(in) :: kt ! ocean time-step |
---|
| 488 | CHARACTER(len=*), INTENT(in) :: cdrw ! "READ"/"WRITE" flag |
---|
| 489 | WRITE(*,*) 'flt_rst: You should not have seen this print! error?', kt, cdrw |
---|
| 490 | END SUBROUTINE flt_rst |
---|
[358] | 491 | #endif |
---|
| 492 | |
---|
| 493 | !!====================================================================== |
---|
| 494 | END MODULE dynspg_flt |
---|