[390] | 1 | ! |
---|
| 2 | Module agrif_opa_interp |
---|
[393] | 3 | #if defined key_agrif |
---|
[390] | 4 | USE par_oce |
---|
| 5 | USE oce |
---|
| 6 | USE dom_oce |
---|
| 7 | USE sol_oce |
---|
| 8 | |
---|
| 9 | CONTAINS |
---|
| 10 | SUBROUTINE Agrif_tra( kt ) |
---|
| 11 | |
---|
| 12 | Implicit none |
---|
| 13 | |
---|
| 14 | !! * Substitutions |
---|
| 15 | # include "domzgr_substitute.h90" |
---|
| 16 | # include "vectopt_loop_substitute.h90" |
---|
| 17 | ! |
---|
| 18 | INTEGER :: kt |
---|
| 19 | REAL(wp) tatemp(jpi,jpj,jpk) , satemp(jpi,jpj,jpk) |
---|
| 20 | INTEGER :: ji,jj,jk |
---|
| 21 | REAL(wp) :: rhox |
---|
| 22 | REAL(wp) :: alpha1, alpha2, alpha3, alpha4 |
---|
| 23 | REAL(wp) :: alpha5, alpha6, alpha7 |
---|
| 24 | ! |
---|
| 25 | IF (Agrif_Root()) RETURN |
---|
| 26 | |
---|
| 27 | Agrif_SpecialValue=0. |
---|
| 28 | Agrif_UseSpecialValue = .TRUE. |
---|
| 29 | tatemp = 0. |
---|
| 30 | satemp = 0. |
---|
| 31 | |
---|
| 32 | Call Agrif_Bc_variable(tatemp,tn) |
---|
| 33 | Call Agrif_Bc_variable(satemp,sn) |
---|
| 34 | Agrif_UseSpecialValue = .FALSE. |
---|
| 35 | |
---|
| 36 | rhox = Agrif_Rhox() |
---|
| 37 | |
---|
| 38 | alpha1 = (rhox-1.)/2. |
---|
| 39 | alpha2 = 1.-alpha1 |
---|
| 40 | |
---|
| 41 | alpha3 = (rhox-1)/(rhox+1) |
---|
| 42 | alpha4 = 1.-alpha3 |
---|
| 43 | |
---|
| 44 | alpha6 = 2.*(rhox-1.)/(rhox+1.) |
---|
| 45 | alpha7 = -(rhox-1)/(rhox+3) |
---|
| 46 | alpha5 = 1. - alpha6 - alpha7 |
---|
| 47 | |
---|
| 48 | ! |
---|
| 49 | If ((nbondi == 1).OR.(nbondi == 2)) THEN |
---|
| 50 | |
---|
| 51 | ta(nlci,:,:) = alpha1 * tatemp(nlci,:,:) + alpha2 * tatemp(nlci-1,:,:) |
---|
| 52 | sa(nlci,:,:) = alpha1 * satemp(nlci,:,:) + alpha2 * satemp(nlci-1,:,:) |
---|
| 53 | |
---|
| 54 | Do jk=1,jpk |
---|
| 55 | Do jj=1,jpj |
---|
| 56 | IF (umask(nlci-2,jj,jk).EQ.0.) THEN |
---|
| 57 | ta(nlci-1,jj,jk) = ta(nlci,jj,jk) * tmask(nlci-1,jj,jk) |
---|
| 58 | sa(nlci-1,jj,jk) = sa(nlci,jj,jk) * tmask(nlci-1,jj,jk) |
---|
| 59 | ELSE |
---|
| 60 | ta(nlci-1,jj,jk)=(alpha4*ta(nlci,jj,jk)+alpha3*ta(nlci-2,jj,jk))*tmask(nlci-1,jj,jk) |
---|
| 61 | sa(nlci-1,jj,jk)=(alpha4*sa(nlci,jj,jk)+alpha3*sa(nlci-2,jj,jk))*tmask(nlci-1,jj,jk) |
---|
| 62 | IF (un(nlci-2,jj,jk).GT.0.) THEN |
---|
| 63 | ta(nlci-1,jj,jk)=(alpha6*ta(nlci-2,jj,jk)+alpha5*ta(nlci,jj,jk)+alpha7*ta(nlci-3,jj,jk))*tmask(nlci-1,jj,jk) |
---|
| 64 | sa(nlci-1,jj,jk)=(alpha6*sa(nlci-2,jj,jk)+alpha5*sa(nlci,jj,jk)+alpha7*sa(nlci-3,jj,jk))*tmask(nlci-1,jj,jk) |
---|
| 65 | ENDIF |
---|
| 66 | ENDIF |
---|
| 67 | End Do |
---|
| 68 | enddo |
---|
| 69 | ENDIF |
---|
| 70 | |
---|
| 71 | If ((nbondj == 1).OR.(nbondj == 2)) THEN |
---|
| 72 | |
---|
| 73 | ta(:,nlcj,:) = alpha1 * tatemp(:,nlcj,:) + alpha2 * tatemp(:,nlcj-1,:) |
---|
| 74 | sa(:,nlcj,:) = alpha1 * satemp(:,nlcj,:) + alpha2 * satemp(:,nlcj-1,:) |
---|
| 75 | |
---|
| 76 | Do jk=1,jpk |
---|
| 77 | Do ji=1,jpi |
---|
| 78 | IF (vmask(ji,nlcj-2,jk).EQ.0.) THEN |
---|
| 79 | ta(ji,nlcj-1,jk) = ta(ji,nlcj,jk) * tmask(ji,nlcj-1,jk) |
---|
| 80 | sa(ji,nlcj-1,jk) = sa(ji,nlcj,jk) * tmask(ji,nlcj-1,jk) |
---|
| 81 | ELSE |
---|
| 82 | ta(ji,nlcj-1,jk)=(alpha4*ta(ji,nlcj,jk)+alpha3*ta(ji,nlcj-2,jk))*tmask(ji,nlcj-1,jk) |
---|
| 83 | sa(ji,nlcj-1,jk)=(alpha4*sa(ji,nlcj,jk)+alpha3*sa(ji,nlcj-2,jk))*tmask(ji,nlcj-1,jk) |
---|
| 84 | IF (vn(ji,nlcj-2,jk) .GT. 0.) THEN |
---|
| 85 | ta(ji,nlcj-1,jk)=(alpha6*ta(ji,nlcj-2,jk)+alpha5*ta(ji,nlcj,jk)+alpha7*ta(ji,nlcj-3,jk))*tmask(ji,nlcj-1,jk) |
---|
| 86 | sa(ji,nlcj-1,jk)=(alpha6*sa(ji,nlcj-2,jk)+alpha5*sa(ji,nlcj,jk)+alpha7*sa(ji,nlcj-3,jk))*tmask(ji,nlcj-1,jk) |
---|
| 87 | ENDIF |
---|
| 88 | ENDIF |
---|
| 89 | End Do |
---|
| 90 | enddo |
---|
| 91 | ENDIF |
---|
| 92 | |
---|
| 93 | IF ((nbondi == -1).OR.(nbondi == 2)) THEN |
---|
| 94 | |
---|
| 95 | ta(1,:,:) = alpha1 * tatemp(1,:,:) + alpha2 * tatemp(2,:,:) |
---|
| 96 | sa(1,:,:) = alpha1 * satemp(1,:,:) + alpha2 * satemp(2,:,:) |
---|
| 97 | |
---|
| 98 | Do jk=1,jpk |
---|
| 99 | Do jj=1,jpj |
---|
| 100 | IF (umask(2,jj,jk).EQ.0.) THEN |
---|
| 101 | ta(2,jj,jk) = ta(1,jj,jk) * tmask(2,jj,jk) |
---|
| 102 | sa(2,jj,jk) = sa(1,jj,jk) * tmask(2,jj,jk) |
---|
| 103 | ELSE |
---|
| 104 | ta(2,jj,jk)=(alpha4*ta(1,jj,jk)+alpha3*ta(3,jj,jk))*tmask(2,jj,jk) |
---|
| 105 | sa(2,jj,jk)=(alpha4*sa(1,jj,jk)+alpha3*sa(3,jj,jk))*tmask(2,jj,jk) |
---|
| 106 | IF (un(2,jj,jk).LT.0.) THEN |
---|
| 107 | ta(2,jj,jk)=(alpha6*ta(3,jj,jk)+alpha5*ta(1,jj,jk)+alpha7*ta(4,jj,jk))*tmask(2,jj,jk) |
---|
| 108 | sa(2,jj,jk)=(alpha6*sa(3,jj,jk)+alpha5*sa(1,jj,jk)+alpha7*sa(4,jj,jk))*tmask(2,jj,jk) |
---|
| 109 | ENDIF |
---|
| 110 | ENDIF |
---|
| 111 | End Do |
---|
| 112 | enddo |
---|
| 113 | ENDIF |
---|
| 114 | |
---|
| 115 | IF ((nbondj == -1).OR.(nbondj == 2)) THEN |
---|
| 116 | |
---|
| 117 | ta(:,1,:) = alpha1 * tatemp(:,1,:) + alpha2 * tatemp(:,2,:) |
---|
| 118 | sa(:,1,:) = alpha1 * satemp(:,1,:) + alpha2 * satemp(:,2,:) |
---|
| 119 | |
---|
| 120 | Do jk=1,jpk |
---|
| 121 | Do ji=1,jpi |
---|
| 122 | IF (vmask(ji,2,jk).EQ.0.) THEN |
---|
| 123 | ta(ji,2,jk)=ta(ji,1,jk) * tmask(ji,2,jk) |
---|
| 124 | sa(ji,2,jk)=sa(ji,1,jk) * tmask(ji,2,jk) |
---|
| 125 | ELSE |
---|
| 126 | ta(ji,2,jk)=(alpha4*ta(ji,1,jk)+alpha3*ta(ji,3,jk))*tmask(ji,2,jk) |
---|
| 127 | sa(ji,2,jk)=(alpha4*sa(ji,1,jk)+alpha3*sa(ji,3,jk))*tmask(ji,2,jk) |
---|
| 128 | IF (vn(ji,2,jk) .LT. 0.) THEN |
---|
| 129 | ta(ji,2,jk)=(alpha6*ta(ji,3,jk)+alpha5*ta(ji,1,jk)+alpha7*ta(ji,4,jk))*tmask(ji,2,jk) |
---|
| 130 | sa(ji,2,jk)=(alpha6*sa(ji,3,jk)+alpha5*sa(ji,1,jk)+alpha7*sa(ji,4,jk))*tmask(ji,2,jk) |
---|
| 131 | ENDIF |
---|
| 132 | ENDIF |
---|
| 133 | End Do |
---|
| 134 | enddo |
---|
| 135 | ENDIF |
---|
| 136 | |
---|
| 137 | Return |
---|
| 138 | End Subroutine Agrif_tra |
---|
| 139 | ! |
---|
| 140 | ! |
---|
| 141 | SUBROUTINE Agrif_dyn(kt) |
---|
| 142 | ! |
---|
| 143 | USE phycst |
---|
| 144 | USE sol_oce |
---|
| 145 | USE in_out_manager |
---|
| 146 | |
---|
| 147 | implicit none |
---|
| 148 | # include "domzgr_substitute.h90" |
---|
| 149 | ! |
---|
| 150 | REAL(wp) uatemp(jpi,jpj,jpk) , vatemp(jpi,jpj,jpk) |
---|
| 151 | INTEGER :: ji,jj,jk |
---|
| 152 | INTEGER kt |
---|
| 153 | REAL(wp) :: z2dt, znugdt |
---|
| 154 | REAL(wp), DIMENSION(jpi,jpj) :: uatemp2D, vatemp2D |
---|
| 155 | REAL(wp) :: timeref |
---|
| 156 | REAL(wp), DIMENSION(jpi,jpj) :: spgu1,spgv1 |
---|
| 157 | REAL(wp) :: rhox, rhoy |
---|
| 158 | |
---|
| 159 | IF (Agrif_Root()) RETURN |
---|
| 160 | |
---|
| 161 | rhox = Agrif_Rhox() |
---|
| 162 | rhoy = Agrif_Rhoy() |
---|
| 163 | |
---|
| 164 | timeref = 1. |
---|
| 165 | |
---|
| 166 | ! time step: leap-frog |
---|
| 167 | z2dt = 2. * rdt |
---|
| 168 | ! time step: Euler if restart from rest |
---|
| 169 | IF( neuler == 0 .AND. kt == nit000 ) z2dt = rdt |
---|
| 170 | ! coefficients |
---|
| 171 | znugdt = rnu * grav * z2dt |
---|
| 172 | |
---|
| 173 | Agrif_SpecialValue=0. |
---|
| 174 | Agrif_UseSpecialValue = .TRUE. |
---|
| 175 | uatemp = 0. |
---|
| 176 | vatemp = 0. |
---|
| 177 | Call Agrif_Bc_variable(uatemp,un,procname=interpu) |
---|
| 178 | Call Agrif_Bc_variable(vatemp,vn,procname=interpv) |
---|
| 179 | uatemp2d = 0. |
---|
| 180 | vatemp2d = 0. |
---|
| 181 | |
---|
| 182 | Agrif_SpecialValue=0. |
---|
| 183 | Agrif_UseSpecialValue = .TRUE. |
---|
| 184 | Call Agrif_Bc_variable(uatemp2d,e1u,calledweight=1.,procname=interpu2d) |
---|
| 185 | Call Agrif_Bc_variable(vatemp2d,e2v,calledweight=1.,procname=interpv2d) |
---|
| 186 | Agrif_UseSpecialValue = .FALSE. |
---|
| 187 | |
---|
| 188 | |
---|
| 189 | If ((nbondi == -1).OR.(nbondi == 2)) THEN |
---|
| 190 | |
---|
| 191 | DO jj=1,jpj |
---|
| 192 | laplacu(2,jj) = timeref * (uatemp2d(2,jj)/(rhoy*e2u(2,jj)))*umask(2,jj,1) |
---|
| 193 | ENDDO |
---|
| 194 | |
---|
| 195 | Do jk=1,jpkm1 |
---|
| 196 | DO jj=1,jpj |
---|
| 197 | ua(1:2,jj,jk) = (uatemp(1:2,jj,jk)/(rhoy*e2u(1:2,jj))) |
---|
[469] | 198 | #if ! defined key_zco |
---|
[390] | 199 | ua(1:2,jj,jk) = ua(1:2,jj,jk) / fse3u(1:2,jj,jk) |
---|
| 200 | #endif |
---|
| 201 | ENDDO |
---|
| 202 | ENDDO |
---|
| 203 | |
---|
| 204 | Do jk=1,jpkm1 |
---|
| 205 | DO jj=1,jpj |
---|
| 206 | ua(2,jj,jk) = (ua(2,jj,jk) - z2dt * znugdt * laplacu(2,jj))*umask(2,jj,jk) |
---|
| 207 | ENDDO |
---|
| 208 | ENDDO |
---|
| 209 | |
---|
| 210 | spgu(2,:)=0. |
---|
| 211 | |
---|
| 212 | do jk=1,jpkm1 |
---|
| 213 | do jj=1,jpj |
---|
| 214 | spgu(2,jj)=spgu(2,jj)+fse3u(2,jj,jk)*ua(2,jj,jk) |
---|
| 215 | enddo |
---|
| 216 | enddo |
---|
| 217 | |
---|
| 218 | DO jj=1,jpj |
---|
| 219 | IF (umask(2,jj,1).NE.0.) THEN |
---|
| 220 | spgu(2,jj)=spgu(2,jj)/hu(2,jj) |
---|
| 221 | ENDIF |
---|
| 222 | enddo |
---|
| 223 | |
---|
| 224 | Do jk=1,jpkm1 |
---|
| 225 | DO jj=1,jpj |
---|
| 226 | ua(2,jj,jk) = 0.25*(ua(1,jj,jk)+2.*ua(2,jj,jk)+ua(3,jj,jk)) |
---|
| 227 | ua(2,jj,jk) = ua(2,jj,jk) * umask(2,jj,jk) |
---|
| 228 | ENDDO |
---|
| 229 | ENDDO |
---|
| 230 | |
---|
| 231 | spgu1(2,:)=0. |
---|
| 232 | |
---|
| 233 | do jk=1,jpkm1 |
---|
| 234 | do jj=1,jpj |
---|
| 235 | spgu1(2,jj)=spgu1(2,jj)+fse3u(2,jj,jk)*ua(2,jj,jk) |
---|
| 236 | enddo |
---|
| 237 | enddo |
---|
| 238 | |
---|
| 239 | DO jj=1,jpj |
---|
| 240 | IF (umask(2,jj,1).NE.0.) THEN |
---|
| 241 | spgu1(2,jj)=spgu1(2,jj)/hu(2,jj) |
---|
| 242 | ENDIF |
---|
| 243 | enddo |
---|
| 244 | |
---|
| 245 | DO jk=1,jpkm1 |
---|
| 246 | DO jj=1,jpj |
---|
| 247 | ua(2,jj,jk) = (ua(2,jj,jk)+spgu(2,jj)-spgu1(2,jj))*umask(2,jj,jk) |
---|
| 248 | ENDDO |
---|
| 249 | ENDDO |
---|
| 250 | |
---|
| 251 | Do jk=1,jpkm1 |
---|
| 252 | Do jj=1,jpj |
---|
| 253 | va(2,jj,jk) = (vatemp(2,jj,jk)/(rhox*e1v(2,jj)))*vmask(2,jj,jk) |
---|
[469] | 254 | #if ! defined key_zco |
---|
[390] | 255 | va(2,jj,jk) = va(2,jj,jk) / fse3v(2,jj,jk) |
---|
| 256 | #endif |
---|
| 257 | End Do |
---|
| 258 | End Do |
---|
| 259 | |
---|
| 260 | sshn(2,:)=sshn(3,:) |
---|
| 261 | sshb(2,:)=sshb(3,:) |
---|
| 262 | |
---|
| 263 | ENDIF |
---|
| 264 | |
---|
| 265 | If ((nbondi == 1).OR.(nbondi == 2)) THEN |
---|
| 266 | |
---|
| 267 | DO jj=1,jpj |
---|
| 268 | laplacu(nlci-2,jj) = timeref * (uatemp2d(nlci-2,jj)/(rhoy*e2u(nlci-2,jj))) |
---|
| 269 | ENDDO |
---|
| 270 | |
---|
| 271 | Do jk=1,jpkm1 |
---|
| 272 | DO jj=1,jpj |
---|
| 273 | ua(nlci-2:nlci-1,jj,jk) = (uatemp(nlci-2:nlci-1,jj,jk)/(rhoy*e2u(nlci-2:nlci-1,jj))) |
---|
| 274 | |
---|
[469] | 275 | #if ! defined key_zco |
---|
[390] | 276 | ua(nlci-2:nlci-1,jj,jk) = ua(nlci-2:nlci-1,jj,jk) / fse3u(nlci-2:nlci-1,jj,jk) |
---|
| 277 | #endif |
---|
| 278 | |
---|
| 279 | ENDDO |
---|
| 280 | ENDDO |
---|
| 281 | |
---|
| 282 | Do jk=1,jpkm1 |
---|
| 283 | DO jj=1,jpj |
---|
| 284 | ua(nlci-2,jj,jk) = (ua(nlci-2,jj,jk)- z2dt * znugdt * laplacu(nlci-2,jj))*umask(nlci-2,jj,jk) |
---|
| 285 | ENDDO |
---|
| 286 | ENDDO |
---|
| 287 | |
---|
| 288 | |
---|
| 289 | spgu(nlci-2,:)=0. |
---|
| 290 | |
---|
| 291 | do jk=1,jpkm1 |
---|
| 292 | do jj=1,jpj |
---|
| 293 | spgu(nlci-2,jj)=spgu(nlci-2,jj)+fse3u(nlci-2,jj,jk)*ua(nlci-2,jj,jk) |
---|
| 294 | enddo |
---|
| 295 | enddo |
---|
| 296 | |
---|
| 297 | DO jj=1,jpj |
---|
| 298 | IF (umask(nlci-2,jj,1).NE.0.) THEN |
---|
| 299 | spgu(nlci-2,jj)=spgu(nlci-2,jj)/hu(nlci-2,jj) |
---|
| 300 | ENDIF |
---|
| 301 | enddo |
---|
| 302 | |
---|
| 303 | Do jk=1,jpkm1 |
---|
| 304 | DO jj=1,jpj |
---|
| 305 | ua(nlci-2,jj,jk) = 0.25*(ua(nlci-3,jj,jk)+2.*ua(nlci-2,jj,jk)+ua(nlci-1,jj,jk)) |
---|
| 306 | |
---|
| 307 | ua(nlci-2,jj,jk) = ua(nlci-2,jj,jk) * umask(nlci-2,jj,jk) |
---|
| 308 | |
---|
| 309 | ENDDO |
---|
| 310 | ENDDO |
---|
| 311 | |
---|
| 312 | spgu1(nlci-2,:)=0. |
---|
| 313 | |
---|
| 314 | do jk=1,jpkm1 |
---|
| 315 | do jj=1,jpj |
---|
| 316 | spgu1(nlci-2,jj)=spgu1(nlci-2,jj)+fse3u(nlci-2,jj,jk)*ua(nlci-2,jj,jk)*umask(nlci-2,jj,jk) |
---|
| 317 | enddo |
---|
| 318 | enddo |
---|
| 319 | |
---|
| 320 | DO jj=1,jpj |
---|
| 321 | IF (umask(nlci-2,jj,1).NE.0.) THEN |
---|
| 322 | spgu1(nlci-2,jj)=spgu1(nlci-2,jj)/hu(nlci-2,jj) |
---|
| 323 | ENDIF |
---|
| 324 | enddo |
---|
| 325 | |
---|
| 326 | DO jk=1,jpkm1 |
---|
| 327 | DO jj=1,jpj |
---|
| 328 | ua(nlci-2,jj,jk) = (ua(nlci-2,jj,jk)+spgu(nlci-2,jj)-spgu1(nlci-2,jj))*umask(nlci-2,jj,jk) |
---|
| 329 | ENDDO |
---|
| 330 | ENDDO |
---|
| 331 | |
---|
| 332 | Do jk=1,jpkm1 |
---|
| 333 | Do jj=1,jpj-1 |
---|
| 334 | va(nlci-1,jj,jk) = (vatemp(nlci-1,jj,jk)/(rhox*e1v(nlci-1,jj)))*vmask(nlci-1,jj,jk) |
---|
[469] | 335 | #if ! defined key_zco |
---|
[390] | 336 | va(nlci-1,jj,jk) = va(nlci-1,jj,jk) / fse3v(nlci-1,jj,jk) |
---|
| 337 | #endif |
---|
| 338 | End Do |
---|
| 339 | End Do |
---|
| 340 | |
---|
| 341 | sshn(nlci-1,:)=sshn(nlci-2,:) |
---|
| 342 | sshb(nlci-1,:)=sshb(nlci-2,:) |
---|
| 343 | ENDIF |
---|
| 344 | |
---|
| 345 | If ((nbondj == -1).OR.(nbondj == 2)) THEN |
---|
| 346 | |
---|
| 347 | DO ji=1,jpi |
---|
| 348 | laplacv(ji,2) = timeref * (vatemp2d(ji,2)/(rhox*e1v(ji,2))) |
---|
| 349 | ENDDO |
---|
| 350 | |
---|
| 351 | DO jk=1,jpkm1 |
---|
| 352 | DO ji=1,jpi |
---|
| 353 | va(ji,1:2,jk) = (vatemp(ji,1:2,jk)/(rhox*e1v(ji,1:2))) |
---|
[469] | 354 | #if ! defined key_zco |
---|
[390] | 355 | va(ji,1:2,jk) = va(ji,1:2,jk) / fse3v(ji,1:2,jk) |
---|
| 356 | #endif |
---|
| 357 | ENDDO |
---|
| 358 | ENDDO |
---|
| 359 | |
---|
| 360 | DO jk=1,jpkm1 |
---|
| 361 | DO ji=1,jpi |
---|
| 362 | va(ji,2,jk) = (va(ji,2,jk) - z2dt * znugdt * laplacv(ji,2))*vmask(ji,2,jk) |
---|
| 363 | ENDDO |
---|
| 364 | ENDDO |
---|
| 365 | |
---|
| 366 | spgv(:,2)=0. |
---|
| 367 | |
---|
| 368 | do jk=1,jpkm1 |
---|
| 369 | do ji=1,jpi |
---|
| 370 | spgv(ji,2)=spgv(ji,2)+fse3v(ji,2,jk)*va(ji,2,jk) |
---|
| 371 | enddo |
---|
| 372 | enddo |
---|
| 373 | |
---|
| 374 | DO ji=1,jpi |
---|
| 375 | IF (vmask(ji,2,1).NE.0.) THEN |
---|
| 376 | spgv(ji,2)=spgv(ji,2)/hv(ji,2) |
---|
| 377 | ENDIF |
---|
| 378 | enddo |
---|
| 379 | |
---|
| 380 | DO jk=1,jpkm1 |
---|
| 381 | DO ji=1,jpi |
---|
| 382 | va(ji,2,jk)=0.25*(va(ji,1,jk)+2.*va(ji,2,jk)+va(ji,3,jk)) |
---|
| 383 | va(ji,2,jk)=va(ji,2,jk)*vmask(ji,2,jk) |
---|
| 384 | ENDDO |
---|
| 385 | ENDDO |
---|
| 386 | |
---|
| 387 | spgv1(:,2)=0. |
---|
| 388 | |
---|
| 389 | do jk=1,jpkm1 |
---|
| 390 | do ji=1,jpi |
---|
| 391 | spgv1(ji,2)=spgv1(ji,2)+fse3v(ji,2,jk)*va(ji,2,jk)*vmask(ji,2,jk) |
---|
| 392 | enddo |
---|
| 393 | enddo |
---|
| 394 | |
---|
| 395 | DO ji=1,jpi |
---|
| 396 | IF (vmask(ji,2,1).NE.0.) THEN |
---|
| 397 | spgv1(ji,2)=spgv1(ji,2)/hv(ji,2) |
---|
| 398 | ENDIF |
---|
| 399 | enddo |
---|
| 400 | |
---|
| 401 | DO jk=1,jpkm1 |
---|
| 402 | DO ji=1,jpi |
---|
| 403 | va(ji,2,jk) = (va(ji,2,jk)+spgv(ji,2)-spgv1(ji,2))*vmask(ji,2,jk) |
---|
| 404 | ENDDO |
---|
| 405 | ENDDO |
---|
| 406 | |
---|
| 407 | DO jk=1,jpkm1 |
---|
| 408 | DO ji=1,jpi |
---|
| 409 | ua(ji,2,jk) = (uatemp(ji,2,jk)/(rhoy*e2u(ji,2)))*umask(ji,2,jk) |
---|
[469] | 410 | #if ! defined key_zco |
---|
[390] | 411 | ua(ji,2,jk) = ua(ji,2,jk) / fse3u(ji,2,jk) |
---|
| 412 | #endif |
---|
| 413 | ENDDO |
---|
| 414 | ENDDO |
---|
| 415 | |
---|
| 416 | sshn(:,2)=sshn(:,3) |
---|
| 417 | sshb(:,2)=sshb(:,3) |
---|
| 418 | ENDIF |
---|
| 419 | |
---|
| 420 | If ((nbondj == 1).OR.(nbondj == 2)) THEN |
---|
| 421 | |
---|
| 422 | DO ji=1,jpi |
---|
| 423 | laplacv(ji,nlcj-2) = timeref * (vatemp2d(ji,nlcj-2)/(rhox*e1v(ji,nlcj-2))) |
---|
| 424 | ENDDO |
---|
| 425 | |
---|
| 426 | DO jk=1,jpkm1 |
---|
| 427 | DO ji=1,jpi |
---|
| 428 | va(ji,nlcj-2:nlcj-1,jk) = (vatemp(ji,nlcj-2:nlcj-1,jk)/(rhox*e1v(ji,nlcj-2:nlcj-1))) |
---|
[469] | 429 | #if ! defined key_zco |
---|
[390] | 430 | va(ji,nlcj-2:nlcj-1,jk) = va(ji,nlcj-2:nlcj-1,jk) / fse3v(ji,nlcj-2:nlcj-1,jk) |
---|
| 431 | #endif |
---|
| 432 | ENDDO |
---|
| 433 | ENDDO |
---|
| 434 | |
---|
| 435 | DO jk=1,jpkm1 |
---|
| 436 | DO ji=1,jpi |
---|
| 437 | va(ji,nlcj-2,jk) = (va(ji,nlcj-2,jk)-z2dt * znugdt * laplacv(ji,nlcj-2))*vmask(ji,nlcj-2,jk) |
---|
| 438 | ENDDO |
---|
| 439 | ENDDO |
---|
| 440 | |
---|
| 441 | |
---|
| 442 | spgv(:,nlcj-2)=0. |
---|
| 443 | |
---|
| 444 | do jk=1,jpkm1 |
---|
| 445 | do ji=1,jpi |
---|
| 446 | spgv(ji,nlcj-2)=spgv(ji,nlcj-2)+fse3v(ji,nlcj-2,jk)*va(ji,nlcj-2,jk) |
---|
| 447 | enddo |
---|
| 448 | enddo |
---|
| 449 | |
---|
| 450 | DO ji=1,jpi |
---|
| 451 | IF (vmask(ji,nlcj-2,1).NE.0.) THEN |
---|
| 452 | spgv(ji,nlcj-2)=spgv(ji,nlcj-2)/hv(ji,nlcj-2) |
---|
| 453 | ENDIF |
---|
| 454 | enddo |
---|
| 455 | |
---|
| 456 | DO jk=1,jpkm1 |
---|
| 457 | DO ji=1,jpi |
---|
| 458 | va(ji,nlcj-2,jk)=0.25*(va(ji,nlcj-3,jk)+2.*va(ji,nlcj-2,jk)+va(ji,nlcj-1,jk)) |
---|
| 459 | va(ji,nlcj-2,jk) = va(ji,nlcj-2,jk) * vmask(ji,nlcj-2,jk) |
---|
| 460 | ENDDO |
---|
| 461 | ENDDO |
---|
| 462 | |
---|
| 463 | spgv1(:,nlcj-2)=0. |
---|
| 464 | |
---|
| 465 | do jk=1,jpkm1 |
---|
| 466 | do ji=1,jpi |
---|
| 467 | spgv1(ji,nlcj-2)=spgv1(ji,nlcj-2)+fse3v(ji,nlcj-2,jk)*va(ji,nlcj-2,jk) |
---|
| 468 | enddo |
---|
| 469 | enddo |
---|
| 470 | |
---|
| 471 | DO ji=1,jpi |
---|
| 472 | IF (vmask(ji,nlcj-2,1).NE.0.) THEN |
---|
| 473 | spgv1(ji,nlcj-2)=spgv1(ji,nlcj-2)/hv(ji,nlcj-2) |
---|
| 474 | ENDIF |
---|
| 475 | enddo |
---|
| 476 | |
---|
| 477 | DO jk=1,jpkm1 |
---|
| 478 | DO ji=1,jpi |
---|
| 479 | va(ji,nlcj-2,jk) = (va(ji,nlcj-2,jk)+spgv(ji,nlcj-2)-spgv1(ji,nlcj-2))*vmask(ji,nlcj-2,jk) |
---|
| 480 | ENDDO |
---|
| 481 | ENDDO |
---|
| 482 | |
---|
| 483 | DO jk=1,jpkm1 |
---|
| 484 | DO ji=1,jpi |
---|
| 485 | ua(ji,nlcj-1,jk) = (uatemp(ji,nlcj-1,jk)/(rhoy*e2u(ji,nlcj-1)))*umask(ji,nlcj-1,jk) |
---|
[469] | 486 | #if ! defined key_zco |
---|
[390] | 487 | ua(ji,nlcj-1,jk) = ua(ji,nlcj-1,jk) / fse3u(ji,nlcj-1,jk) |
---|
| 488 | #endif |
---|
| 489 | ENDDO |
---|
| 490 | ENDDO |
---|
| 491 | |
---|
| 492 | sshn(:,nlcj-1)=sshn(:,nlcj-2) |
---|
| 493 | sshb(:,nlcj-1)=sshb(:,nlcj-2) |
---|
| 494 | ENDIF |
---|
| 495 | |
---|
| 496 | ! |
---|
| 497 | Return |
---|
| 498 | End Subroutine Agrif_dyn |
---|
| 499 | |
---|
| 500 | |
---|
| 501 | subroutine interpu(tabres,i1,i2,j1,j2,k1,k2) |
---|
| 502 | Implicit none |
---|
| 503 | # include "domzgr_substitute.h90" |
---|
| 504 | integer i1,i2,j1,j2,k1,k2 |
---|
| 505 | integer ji,jj,jk |
---|
| 506 | real,dimension(i1:i2,j1:j2,k1:k2) :: tabres |
---|
| 507 | |
---|
| 508 | do jk=k1,k2 |
---|
| 509 | DO jj=j1,j2 |
---|
| 510 | DO ji=i1,i2 |
---|
| 511 | tabres(ji,jj,jk) = e2u(ji,jj) * un(ji,jj,jk) |
---|
[469] | 512 | #if ! defined key_zco |
---|
[390] | 513 | tabres(ji,jj,jk) = tabres(ji,jj,jk) * fse3u(ji,jj,jk) |
---|
| 514 | #endif |
---|
| 515 | ENDDO |
---|
| 516 | ENDDO |
---|
| 517 | ENDDO |
---|
| 518 | end subroutine interpu |
---|
| 519 | |
---|
| 520 | subroutine interpu2d(tabres,i1,i2,j1,j2) |
---|
| 521 | Implicit none |
---|
| 522 | integer i1,i2,j1,j2 |
---|
| 523 | integer ji,jj |
---|
| 524 | real,dimension(i1:i2,j1:j2) :: tabres |
---|
| 525 | |
---|
| 526 | DO jj=j1,j2 |
---|
| 527 | DO ji=i1,i2 |
---|
| 528 | tabres(ji,jj) = e2u(ji,jj) * ((gcx(ji+1,jj) - gcx(ji,jj))/e1u(ji,jj)) & |
---|
| 529 | *umask(ji,jj,1) |
---|
| 530 | ENDDO |
---|
| 531 | ENDDO |
---|
| 532 | end subroutine interpu2d |
---|
| 533 | |
---|
| 534 | subroutine interpv(tabres,i1,i2,j1,j2,k1,k2) |
---|
| 535 | Implicit none |
---|
| 536 | # include "domzgr_substitute.h90" |
---|
| 537 | integer i1,i2,j1,j2,k1,k2 |
---|
| 538 | integer ji,jj,jk |
---|
| 539 | real,dimension(i1:i2,j1:j2,k1:k2) :: tabres |
---|
| 540 | |
---|
| 541 | do jk=k1,k2 |
---|
| 542 | DO jj=j1,j2 |
---|
| 543 | DO ji=i1,i2 |
---|
| 544 | tabres(ji,jj,jk) = e1v(ji,jj) * vn(ji,jj,jk) |
---|
[469] | 545 | #if ! defined key_zco |
---|
[390] | 546 | tabres(ji,jj,jk) = tabres(ji,jj,jk) * fse3v(ji,jj,jk) |
---|
| 547 | #endif |
---|
| 548 | ENDDO |
---|
| 549 | ENDDO |
---|
| 550 | ENDDO |
---|
| 551 | end subroutine interpv |
---|
| 552 | |
---|
| 553 | subroutine interpv2d(tabres,i1,i2,j1,j2) |
---|
| 554 | Implicit none |
---|
| 555 | integer i1,i2,j1,j2 |
---|
| 556 | integer ji,jj |
---|
| 557 | real,dimension(i1:i2,j1:j2) :: tabres |
---|
| 558 | |
---|
| 559 | DO jj=j1,j2 |
---|
| 560 | DO ji=i1,i2 |
---|
| 561 | tabres(ji,jj) = e1v(ji,jj) * ((gcx(ji,jj+1) - gcx(ji,jj))/e2v(ji,jj)) & |
---|
| 562 | * vmask(ji,jj,1) |
---|
| 563 | ENDDO |
---|
| 564 | ENDDO |
---|
| 565 | end subroutine interpv2d |
---|
| 566 | |
---|
| 567 | #else |
---|
| 568 | CONTAINS |
---|
| 569 | subroutine Agrif_OPA_Interp_empty |
---|
| 570 | |
---|
| 571 | end subroutine Agrif_OPA_Interp_empty |
---|
| 572 | #endif |
---|
| 573 | End Module agrif_opa_interp |
---|
| 574 | |
---|