/[MITgcm]/MITgcm/pkg/seaice/seaice_evp.F
ViewVC logotype

Contents of /MITgcm/pkg/seaice/seaice_evp.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.48 - (show annotations) (download)
Wed Mar 13 09:57:20 2013 UTC (11 years, 2 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64g, checkpoint64f
Changes since 1.47: +20 -7 lines
mask the asymmetric contribution of the ice-ocean stress (when the
turningAngle is non-zero) over open water to avoid stripes in the
EVP-solutions after many EVP-substeps. This does not change the
verification experiment (because turningAngle=0 in that experiment),
but will affect other EVP solutions with turningAngle.ne.0

1 C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_evp.F,v 1.47 2013/01/04 14:33:56 mlosch Exp $
2 C $Name: $
3
4 #include "SEAICE_OPTIONS.h"
5 #ifdef ALLOW_OBCS
6 # include "OBCS_OPTIONS.h"
7 #else
8 # define OBCS_UVICE_OLD
9 #endif
10
11 CBOP
12 C !ROUTINE: SEAICE_EVP
13 C !INTERFACE:
14 SUBROUTINE SEAICE_EVP( myTime, myIter, myThid )
15
16 C !DESCRIPTION: \bv
17 C *==========================================================*
18 C | SUBROUTINE SEAICE_EVP
19 C | o Ice dynamics using an EVP solver following
20 C | E. C. Hunke and J. K. Dukowicz. An
21 C | Elastic-Viscous-Plastic Model for Sea Ice Dynamics,
22 C | J. Phys. Oceanogr., 27, 1849-1867 (1997).
23 C *==========================================================*
24 C | written by Martin Losch, March 2006
25 C *==========================================================*
26 C \ev
27
28 C !USES:
29 IMPLICIT NONE
30
31 C === Global variables ===
32 #include "SIZE.h"
33 #include "EEPARAMS.h"
34 #include "PARAMS.h"
35 #include "DYNVARS.h"
36 #include "GRID.h"
37 #include "SEAICE_SIZE.h"
38 #include "SEAICE_PARAMS.h"
39 #include "SEAICE.h"
40
41 #ifdef ALLOW_AUTODIFF_TAMC
42 # include "tamc.h"
43 #endif
44
45 C !INPUT/OUTPUT PARAMETERS:
46 C === Routine arguments ===
47 C myTime :: Simulation time
48 C myIter :: Simulation timestep number
49 C myThid :: my Thread Id. number
50 _RL myTime
51 INTEGER myIter
52 INTEGER myThid
53 CEOP
54
55 #if ( defined (SEAICE_CGRID) && \
56 defined (SEAICE_ALLOW_EVP) && \
57 defined (SEAICE_ALLOW_DYNAMICS) )
58
59 C === Local variables ===
60 C i,j,bi,bj :: Loop counters
61 C kSrf :: vertical index of surface layer
62 C nEVPstep :: number of timesteps within the EVP solver
63 C iEVPstep :: Loop counter
64 C SIN/COSWAT :: sine/cosine of turning angle
65 C (recip_)ecc2 :: (one over) eccentricity squared
66 C recip_evpAlpha :: 1/SEAICE_evpAlpha
67 C recip_deltaT :: 1/SEAICE_deltaTdyn
68 C evpStarFac :: 1 if SEAICEuseEVPstar = .true., 0 otherwise
69 C betaFac :: SEAICE_evpBeta/SEAICE_deltaTdyn=1/SEAICE_deltaTevp
70 C betaFacP1 :: betaFac + evpStarFac/SEAICE_deltaTdyn
71 C e11,e12,e22 :: components of strain rate tensor
72 C seaice_div :: divergence strain rates at C-points times P
73 C /divided by Delta minus 1
74 C seaice_tension :: tension strain rates at C-points times P
75 C /divided by Delta
76 C seaice_shear :: shear strain rates, defined at Z-points times P
77 C /divided by Delta
78 C sig11, sig22 :: sum and difference of diagonal terms of stress tensor
79
80 INTEGER i, j, bi, bj
81 INTEGER kSrf
82 INTEGER nEVPstep, iEVPstep
83 #ifdef ALLOW_AUTODIFF_TAMC
84 INTEGER ikeyloc
85 #else
86 INTEGER nEVPstepMax
87 #endif
88
89 _RL COSWAT
90 _RS SINWAT
91 _RL ecc2, recip_ecc2, recip_evpAlpha, recip_deltaT
92 _RL betaFacP1, betaFac, evpStarFac
93
94 _RL seaice_div (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95 _RL seaice_tension(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96 _RL seaice_shear (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97 _RL sig11 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98 _RL sig22 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99 C auxilliary variables
100 _RL ep (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
101 _RL em (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
102 _RL pressC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
103 _RL zetaC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104 _RL deltaZ (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105 _RL zetaZ (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106 _RL deltaC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
107 _RL deltaCreg, deltaZreg, pressZ
108 #ifdef SEAICE_ALLOW_TEM
109 _RL etaDenC, zetaMaxC, etaDenZ, zetaMaxZ
110 #endif /* SEAICE_ALLOW_TEM */
111 #ifdef SEAICE_ALLOW_CLIPZETA
112 _RL zMaxZ, zMinZ, fac
113 #endif /* SEAICE_ALLOW_CLIPZETA */
114 _RL denom1, denom2, facZ
115 _RL locMaskU, locMaskV
116
117 C surface level
118 kSrf = 1
119 C-- introduce turning angles
120 SINWAT=SIN(SEAICE_waterTurnAngle*deg2rad)
121 COSWAT=COS(SEAICE_waterTurnAngle*deg2rad)
122
123 C abbreviation eccentricity squared
124 ecc2=SEAICE_eccen**2
125 recip_ecc2 = 0. _d 0
126 IF ( ecc2 .NE. 0. _d 0 ) recip_ecc2=ONE/ecc2
127 C determine number of internal time steps
128 nEVPstep = INT(SEAICE_deltaTdyn/SEAICE_deltaTevp)
129 IF ( SEAICEnEVPstarSteps.NE.UNSET_I ) nEVPstep=SEAICEnEVPstarSteps
130 C SEAICE_evpAlpha = 2. * SEAICE_evpTauRelax/SEAICE_deltaTevp
131 denom1 = SEAICE_evpAlpha / ( SEAICE_evpAlpha + 1. _d 0 )
132 denom2 = SEAICE_evpAlpha / ( SEAICE_evpAlpha + ecc2 )
133
134 recip_deltaT = 1. _d 0 / SEAICE_deltaTdyn
135 recip_evpAlpha = 0. _d 0
136 IF ( SEAICE_evpAlpha .GT. 0. _d 0 )
137 & recip_evpAlpha = 1. _d 0 / SEAICE_evpAlpha
138 evpStarFac = 0. _d 0
139 IF ( SEAICEuseEVPstar ) evpStarFac = 1. _d 0
140
141 betaFac = SEAICE_evpBeta*recip_deltaT
142 betaFacP1 = betaFac + evpStarFac*recip_deltaT
143 #ifndef ALLOW_AUTODIFF_TAMC
144 nEVPstepMax = nEVPstep
145 #endif
146
147 DO bj=myByLo(myThid),myByHi(myThid)
148 DO bi=myBxLo(myThid),myBxHi(myThid)
149 DO j=1-OLy,sNy+OLy
150 DO i=1-OLx,sNx+OLx
151 C use u/vIce as work arrays: copy previous time step to u/vIceNm1
152 uIceNm1(I,J,bi,bj) = uIce(I,J,bi,bj)
153 vIceNm1(I,J,bi,bj) = vIce(I,J,bi,bj)
154 C initialise strain rates
155 e11 (I,J,bi,bj) = 0. _d 0
156 e22 (I,J,bi,bj) = 0. _d 0
157 e12 (I,J,bi,bj) = 0. _d 0
158 ENDDO
159 ENDDO
160 ENDDO
161 ENDDO
162 #ifdef SEAICE_ALLOW_CLIPZETA
163 C damping constraint (Hunke, J.Comp.Phys.,2001)
164 IF ( SEAICE_evpDampC .GT. 0. _d 0 ) THEN
165 CML fac = HALF * SEAICE_evpDampC * SEAICE_evpTauRelax
166 CML & /SEAICE_deltaTevp**2
167 fac = 0.25 _d 0 * SEAICE_evpDampC * SEAICE_evpAlpha
168 & /SEAICE_deltaTevp
169 DO bj=myByLo(myThid),myByHi(myThid)
170 DO bi=myBxLo(myThid),myBxHi(myThid)
171 DO j=1-OLy,sNy+OLy
172 DO i=1-OLx,sNx+OLx
173 zMax (I,J,bi,bj) = _rA(I,J,bi,bj) * fac
174 ENDDO
175 ENDDO
176 ENDDO
177 ENDDO
178 ENDIF
179 #endif /* SEAICE_ALLOW_CLIPZETA */
180 C
181 C start of the main time loop
182 DO iEVPstep = 1, nEVPstepMax
183 IF (iEVPstep.LE.nEVPstep) THEN
184 C
185 #ifdef ALLOW_AUTODIFF_TAMC
186 ikeyloc = iEVPstep +
187 & (ikey_dynamics-1)*nEVPstepMax
188 CADJ STORE uice = comlev1_evp,
189 CADJ & key = ikeyloc, byte = isbyte
190 CADJ STORE vice = comlev1_evp,
191 CADJ & key = ikeyloc, byte = isbyte
192 CADJ STORE seaice_sigma1 = comlev1_evp,
193 CADJ & key = ikeyloc, byte = isbyte
194 CADJ STORE seaice_sigma2 = comlev1_evp,
195 CADJ & key = ikeyloc, byte = isbyte
196 CADJ STORE seaice_sigma12 = comlev1_evp,
197 CADJ & key = ikeyloc, byte = isbyte
198 #endif /* ALLOW_AUTODIFF_TAMC */
199 C
200 C first calculate strain rates and bulk moduli/viscosities
201 C
202 CALL SEAICE_CALC_STRAINRATES(
203 I uIce, vIce,
204 O e11, e22, e12,
205 I iEVPstep, myTime, myIter, myThid )
206
207 #ifdef ALLOW_AUTODIFF_TAMC
208 CADJ STORE e11 = comlev1_evp,key = ikeyloc, byte = isbyte
209 CADJ STORE e12 = comlev1_evp,key = ikeyloc, byte = isbyte
210 CADJ STORE e22 = comlev1_evp,key = ikeyloc, byte = isbyte
211 #endif /* ALLOW_AUTODIFF_TAMC */
212
213 DO bj=myByLo(myThid),myByHi(myThid)
214 DO bi=myBxLo(myThid),myBxHi(myThid)
215 DO j=1-OLy,sNy+OLy
216 DO i=1-OLx,sNx+OLx
217 seaice_div (I,J) = 0. _d 0
218 seaice_tension(I,J) = 0. _d 0
219 seaice_shear (I,J) = 0. _d 0
220 pressC (I,J) = 0. _d 0
221 zetaC (I,J) = 0. _d 0
222 deltaZ (I,J) = 0. _d 0
223 zetaZ (I,J) = 0. _d 0
224 deltaC (I,J) = 0. _d 0
225 ENDDO
226 ENDDO
227 DO j=1-OLy,sNy+OLy
228 DO i=1-OLx,sNx+OLx
229 ep(i,j) = e11(i,j,bi,bj) + e22(i,j,bi,bj)
230 em(i,j) = e11(i,j,bi,bj) - e22(i,j,bi,bj)
231 ENDDO
232 ENDDO
233 DO j=0,sNy+1
234 DO i=0,sNx+1
235 C average strain rates to Z and C points
236 facZ = 0.25 _d 0
237 CML facZ = 1.0 _d 0
238 CML & / MAX(1.D0,maskC(I, J, kSrf,bi,bj)
239 CML & + maskC(I-1,J, kSrf,bi,bj)
240 CML & + maskC(I, J-1,kSrf,bi,bj)
241 CML & + maskC(I-1,J-1,kSrf,bi,bj) )
242 C Averaging the squares gives more accurate viscous-plastic behavior
243 C than squaring the averages
244 deltaC(I,J) =
245 & ep(I,J)**2 + recip_ecc2*em(I,J)**2
246 & + recip_ecc2*
247 & ( e12(I, J,bi,bj)**2 + e12(I+1,J+1,bi,bj)**2
248 & + e12(I+1,J,bi,bj)**2 + e12(I, J+1,bi,bj)**2 )
249 deltaZ(I,J) = facZ *
250 & ( ep(I,J )**2 + ep(I-1,J )**2
251 & + ep(I,J-1)**2 + ep(I-1,J-1)**2
252 & +(em(I,J )**2 + em(I-1,J )**2
253 & +em(I,J-1)**2 + em(I-1,J-1)**2)*recip_ecc2
254 & )
255 & + 4. _d 0*recip_ecc2*e12(I,J,bi,bj)**2
256 #ifdef ALLOW_AUTODIFF_TAMC
257 C avoid sqrt of 0
258 IF ( deltaC(I,J) .GT. 0. _d 0 )
259 & deltaC(I,J) = SQRT(deltaC(I,J))
260 IF ( deltaZ(I,J) .GT. 0. _d 0 )
261 & deltaZ(I,J) = SQRT(deltaZ(I,J))
262 #else
263 deltaC(I,J) = SQRT(deltaC(I,J))
264 deltaZ(I,J) = SQRT(deltaZ(I,J))
265 #endif /* ALLOW_AUTODIFF_TAMC */
266 deltaCreg = MAX(deltaC(I,J),SEAICE_EPS)
267 deltaZreg = MAX(deltaZ(I,J),SEAICE_EPS)
268 C modify pressure (copied from seaice_calc_viscosities)
269 zetaC(I,J) = HALF*press0(I,J,bi,bj)/deltaCreg
270 pressZ = (deltaZ(I,J)/deltaZreg) * facZ
271 & * ( PRESS0(I,J ,bi,bj) + PRESS0(I-1,J ,bi,bj)
272 & + PRESS0(I,J-1,bi,bj) + PRESS0(I-1,J-1,bi,bj) )
273 zetaZ(I,J) = HALF/deltaZreg * pressZ
274 ENDDO
275 ENDDO
276 #ifdef SEAICE_ALLOW_CLIPZETA
277 #ifdef ALLOW_AUTODIFF_TAMC
278 CADJ STORE zetac = comlev1_evp,key = ikeyloc, byte = isbyte
279 CADJ STORE zetaz = comlev1_evp,key = ikeyloc, byte = isbyte
280 #endif /* ALLOW_AUTODIFF_TAMC */
281 C regularize zeta if necessary
282 DO j=0,sNy+1
283 DO i=0,sNx+1
284 zetaC(I,J) = MAX(zMin(I,J,bi,bj),MIN(zMax(I,J,bi,bj)
285 & ,zetaC(I,J)))
286 CML zetaC(I,J) = zetaC(I,J)*hEffM(I,J,bi,bj)
287 C
288 C zMin, zMax are defined at C-points, make sure that they are not
289 C masked by boundaries/land points
290 zMaxZ = MAX(
291 & MAX(zMax(I, J,bi,bj),zMax(I, J-1,bi,bj)),
292 & MAX(zMax(I-1,J,bi,bj),zMax(I-1,J-1,bi,bj)) )
293 zMinZ = MAX(
294 & MAX(zMin(I, J,bi,bj),zMin(I, J-1,bi,bj)),
295 & MAX(zMin(I-1,J,bi,bj),zMin(I-1,J-1,bi,bj)) )
296 zetaZ(I,J) = MAX(zMinZ,MIN(zMaxZ,zetaZ(I,J)))
297 ENDDO
298 ENDDO
299 #endif /* SEAICE_ALLOW_CLIPZETA */
300 C recompute pressure
301 DO j=0,sNy+1
302 DO i=0,sNx+1
303 pressC(I,J) = TWO*zetaC(I,J)*deltaC(I,J)
304 ENDDO
305 ENDDO
306 #ifdef ALLOW_DIAGNOSTICS
307 IF ( useDiagnostics ) THEN
308 C save eta, zeta, and pressure for diagnostics
309 DO j=1,sNy
310 DO i=1,sNx
311 press(I,J,bi,bj) = pressC(I,J)
312 zeta (I,J,bi,bj) = zetaC(I,J)
313 eta (I,J,bi,bj) = zetaC(I,J)*recip_ecc2
314 ENDDO
315 ENDDO
316 ENDIF
317 #endif /* ALLOW_DIAGNOSTICS */
318 C Calculate the RHS of the stress equations. Do this now in order to
319 C avoid multiple divisions by Delta
320 C P * ( D_d/Delta - 1 ) = 2*zeta*D_d - P = 2*zeta*D_d - 2*zeta*Delta
321 C P * ( D_t/Delta ) = 2*zeta*D_t
322 C P * ( D_s/Delta ) = 2*zeta*D_s
323 #ifdef SEAICE_ALLOW_TEM
324 IF ( SEAICEuseTEM ) THEN
325 DO j=0,sNy+1
326 DO i=0,sNx+1
327 facZ = 0.25 _d 0
328 CML facZ = 1.0 _d 0
329 CML & / MAX(1.D0,maskC(I, J, kSrf,bi,bj)
330 CML & + maskC(I-1,J, kSrf,bi,bj)
331 CML & + maskC(I, J-1,kSrf,bi,bj)
332 CML & + maskC(I-1,J-1,kSrf,bi,bj) )
333 C Averaging the squares gives more accurate viscous-plastic behavior
334 C than squaring the averages
335 etaDenC = em(I,J)**2 +
336 & ( e12(I, J,bi,bj)**2 + e12(I+1,J+1,bi,bj)**2
337 & + e12(I+1,J,bi,bj)**2 + e12(I, J+1,bi,bj)**2 )
338 etaDenC = SQRT(MAX(SEAICE_EPS_SQ,etaDenC))
339 zetaMaxC = ecc2*zetaC(I,J)*(deltaC(I,J)-ep(I,J))/etaDenC
340 etaDenZ =
341 & facZ * ( em(I, J )**2 + em(I-1,J-1)**2
342 & + em(I-1,J )**2 + em(I, J-1)**2 )
343 & + 4. _d 0*e12(I,J,bi,bj)**2
344 etaDenZ = SQRT(MAX(SEAICE_EPS_SQ,etaDenZ))
345 zetaMaxZ = ecc2*zetaZ(I,J) * ( deltaZ(I,J)
346 & - facZ * ( ep(I,J ) + ep(I-1,J )
347 & + ep(I,J-1) + ep(I-1,J-1) )
348 & )/etaDenZ
349 #ifdef ALLOW_DIAGNOSTICS
350 C save new eta for diagnostics
351 eta(I,J,bi,bj) = MIN(zetaC(I,J),zetaMaxC)*recip_ecc2
352 #endif /* ALLOW_DIAGNOSTICS */
353 seaice_div (I,J) =
354 & ( 2. _d 0 *zetaC(I,J)*ep(I,J) - pressC(I,J)
355 & ) * hEffM(I,J,bi,bj)
356 seaice_tension(I,J) = 2. _d 0*MIN(zetaC(I,J),zetaMaxC)
357 & * em(I,J) * hEffM(I,J,bi,bj)
358 seaice_shear (I,J) = 2. _d 0*MIN(zetaZ(I,J),zetaMaxZ)
359 & * 2. _d 0*e12(I,J,bi,bj)
360 ENDDO
361 ENDDO
362 ELSE
363 #else
364 IF (.TRUE. ) THEN
365 #endif /* SEAICE_ALLOW_TEM */
366 DO j=0,sNy+1
367 DO i=0,sNx+1
368 seaice_div (I,J) =
369 & ( 2. _d 0 *zetaC(I,J)*ep(I,J) - pressC(I,J)
370 & ) * hEffM(I,J,bi,bj)
371 seaice_tension(I,J) = 2. _d 0*zetaC(I,J)
372 & * em(I,J) * hEffM(I,J,bi,bj)
373 seaice_shear (I,J) =
374 & 2. _d 0*zetaZ(I,J)*2. _d 0*e12(I,J,bi,bj)
375 ENDDO
376 ENDDO
377 ENDIF
378 C
379 C first step stress equations
380 C
381 DO j=0,sNy
382 DO i=0,sNx
383 C sigma1 and sigma2 are computed on C points
384 seaice_sigma1 (I,J,bi,bj) = ( seaice_sigma1 (I,J,bi,bj)
385 & + recip_evpAlpha * seaice_div(I,J)
386 & ) * denom1
387 & *hEffM(I,J,bi,bj)
388 seaice_sigma2 (I,J,bi,bj) = ( seaice_sigma2 (I,J,bi,bj)
389 & + recip_evpAlpha * seaice_tension(I,J)
390 & ) * denom2
391 & *hEffM(I,J,bi,bj)
392 #ifdef SEAICE_EVP_ELIMINATE_UNDERFLOWS
393 C Code to avoid very small numbers that can degrade performance.
394 C Many compilers can handle this more efficiently with the help of
395 C a flag (copied from CICE after correspondence with Elizabeth Hunke)
396 seaice_sigma1(I,J,bi,bj) = SIGN(MAX(
397 & ABS( seaice_sigma1(I,J,bi,bj) ), SEAICE_EPS ),
398 & seaice_sigma1(I,J,bi,bj) )
399 seaice_sigma2(I,J,bi,bj) = SIGN(MAX(
400 & ABS( seaice_sigma2(I,J,bi,bj) ), SEAICE_EPS ),
401 & seaice_sigma2(I,J,bi,bj) )
402 #endif /* SEAICE_EVP_ELIMINATE_UNDERFLOWS */
403 C recover sigma11 and sigma22
404 sig11(I,J) = 0.5 _d 0 *
405 & ( seaice_sigma1(I,J,bi,bj)+seaice_sigma2(I,J,bi,bj) )
406 sig22(I,J) = 0.5 _d 0 *
407 & ( seaice_sigma1(I,J,bi,bj)-seaice_sigma2(I,J,bi,bj) )
408 ENDDO
409 ENDDO
410
411 C sigma12 is computed on Z points
412 DO j=1,sNy+1
413 DO i=1,sNx+1
414 seaice_sigma12(I,J,bi,bj) = ( seaice_sigma12(I,J,bi,bj)
415 & + 0.5 _d 0 * recip_evpAlpha * seaice_shear(I,J)
416 & ) * denom2
417 #ifdef SEAICE_EVP_ELIMINATE_UNDERFLOWS
418 seaice_sigma12(I,J,bi,bj) = SIGN(MAX(
419 & ABS( seaice_sigma12(I,J,bi,bj) ), SEAICE_EPS ),
420 & seaice_sigma12(I,J,bi,bj) )
421 #endif /* SEAICE_EVP_ELIMINATE_UNDERFLOWS */
422 ENDDO
423 ENDDO
424 C
425 C compute divergence of stress tensor
426 C
427 DO J=1,sNy
428 DO I=1,sNx
429 stressDivergenceX(I,J,bi,bj) =
430 & ( sig11(I ,J ) * _dyF(I ,J,bi,bj)
431 & - sig11(I-1,J ) * _dyF(I-1,J,bi,bj)
432 & + seaice_sigma12(I,J+1,bi,bj) * _dxV(I,J+1,bi,bj)
433 & - seaice_sigma12(I,J ,bi,bj) * _dxV(I,J ,bi,bj)
434 & ) * recip_rAw(I,J,bi,bj)
435 stressDivergenceY(I,J,bi,bj) =
436 & ( sig22(I,J ) * _dxF(I,J ,bi,bj)
437 & - sig22(I,J-1) * _dxF(I,J-1,bi,bj)
438 & + seaice_sigma12(I+1,J,bi,bj) * _dyU(I+1,J,bi,bj)
439 & - seaice_sigma12(I ,J,bi,bj) * _dyU(I ,J,bi,bj)
440 & ) * recip_rAs(I,J,bi,bj)
441 ENDDO
442 ENDDO
443 ENDDO
444 ENDDO
445
446 #ifdef ALLOW_AUTODIFF_TAMC
447 CADJ STORE stressDivergenceX = comlev1_evp,
448 CADJ & key = ikeyloc, byte = isbyte
449 CADJ STORE stressDivergenceY = comlev1_evp,
450 CADJ & key = ikeyloc, byte = isbyte
451 #endif /* ALLOW_AUTODIFF_TAMC */
452
453
454 #ifdef ALLOW_AUTODIFF_TAMC
455 #ifdef SEAICE_DYN_STABLE_ADJOINT
456 Cgf zero out adjoint fields to stabilize pkg/seaice dyna. adjoint
457 CALL ZERO_ADJ( 1, stressDivergenceX, myThid)
458 CALL ZERO_ADJ( 1, stressDivergenceY, myThid)
459 #endif
460 #endif /* ALLOW_AUTODIFF_TAMC */
461
462 C
463 C set up rhs for stepping the velocity field
464 C
465 CALL SEAICE_OCEANDRAG_COEFFS(
466 I uIce, vIce,
467 O DWATN,
468 I iEVPstep, myTime, myIter, myThid )
469
470 DO bj=myByLo(myThid),myByHi(myThid)
471 DO bi=myBxLo(myThid),myBxHi(myThid)
472 DO J=1,sNy
473 DO I=1,sNx
474 C over open water, all terms that contain sea ice mass drop out and
475 C the balance is determined by the atmosphere-ice and ice-ocean stress;
476 C the staggering of uIce and vIce can cause stripes in the open ocean
477 C solution when the turning angle is non-zero (SINWAT.NE.0);
478 C we mask this term here in order to avoid the stripes and because
479 C over open ocean, u/vIce do not advect anything, so that the associated
480 C error is small and most likely only confinded to the ice edge but may
481 C propagate into the ice covered regions.
482 locMaskU = seaiceMassU(I,J,bi,bj)
483 locMaskV = seaiceMassV(I,J,bi,bj)
484 IF ( locMaskU .NE. 0. _d 0 ) locMaskU = 1. _d 0
485 IF ( locMaskV .NE. 0. _d 0 ) locMaskV = 1. _d 0
486 C to recover old results replace the above lines with the line below
487 C locMaskU = 1. _d 0
488 C locMaskV = 1. _d 0
489 C set up anti symmetric drag force and add in ice ocean stress
490 C ( remember to average to correct velocity points )
491 FORCEX(I,J,bi,bj)=FORCEX0(I,J,bi,bj)+
492 & 0.5 _d 0 * ( DWATN(I,J,bi,bj)+DWATN(I-1,J,bi,bj) ) *
493 & COSWAT * uVel(I,J,kSrf,bi,bj)
494 & - SIGN(SINWAT, _fCori(I,J,bi,bj))* 0.5 _d 0 *
495 & ( DWATN(I ,J,bi,bj) * 0.5 _d 0 *
496 & (vVel(I ,J ,kSrf,bi,bj)-vIce(I ,J ,bi,bj)
497 & +vVel(I ,J+1,kSrf,bi,bj)-vIce(I ,J+1,bi,bj))
498 & + DWATN(I-1,J,bi,bj) * 0.5 _d 0 *
499 & (vVel(I-1,J ,kSrf,bi,bj)-vIce(I-1,J ,bi,bj)
500 & +vVel(I-1,J+1,kSrf,bi,bj)-vIce(I-1,J+1,bi,bj))
501 & )*locMaskU
502 FORCEY(I,J,bi,bj)=FORCEY0(I,J,bi,bj)+
503 & 0.5 _d 0 * ( DWATN(I,J,bi,bj)+DWATN(I,J-1,bi,bj) ) *
504 & COSWAT * vVel(I,J,kSrf,bi,bj)
505 & + SIGN(SINWAT, _fCori(I,J,bi,bj)) * 0.5 _d 0 *
506 & ( DWATN(I,J ,bi,bj) * 0.5 _d 0 *
507 & (uVel(I ,J ,kSrf,bi,bj)-uIce(I ,J ,bi,bj)
508 & +uVel(I+1,J ,kSrf,bi,bj)-uIce(I+1,J ,bi,bj))
509 & + DWATN(I,J-1,bi,bj) * 0.5 _d 0 *
510 & (uVel(I ,J-1,kSrf,bi,bj)-uIce(I ,J-1,bi,bj)
511 & +uVel(I+1,J-1,kSrf,bi,bj)-uIce(I+1,J-1,bi,bj))
512 & )*locMaskV
513 C coriols terms
514 FORCEX(I,J,bi,bj)=FORCEX(I,J,bi,bj) + 0.5 _d 0*(
515 & seaiceMassC(I ,J,bi,bj) * _fCori(I ,J,bi,bj)
516 & * 0.5 _d 0*( vIce(I ,J,bi,bj)+vIce(I ,J+1,bi,bj) )
517 & + seaiceMassC(I-1,J,bi,bj) * _fCori(I-1,J,bi,bj)
518 & * 0.5 _d 0*( vIce(I-1,J,bi,bj)+vIce(I-1,J+1,bi,bj) )
519 & )
520 FORCEY(I,J,bi,bj)=FORCEY(I,J,bi,bj) - 0.5 _d 0*(
521 & seaiceMassC(I,J ,bi,bj) * _fCori(I,J ,bi,bj)
522 & * 0.5 _d 0*( uIce(I,J ,bi,bj)+uIce(I+1, J,bi,bj) )
523 & + seaiceMassC(I,J-1,bi,bj) * _fCori(I,J-1,bi,bj)
524 & * 0.5 _d 0*( uIce(I,J-1,bi,bj)+uIce(I+1,J-1,bi,bj) )
525 & )
526 ENDDO
527 ENDDO
528 C
529 C step momentum equations with ice-ocean stress term treated implicitly
530 C
531 DO J=1,sNy
532 DO I=1,sNx
533 uIce(I,J,bi,bj) = seaiceMaskU(I,J,bi,bj) *
534 & ( seaiceMassU(I,J,bi,bj)*betaFac
535 & * uIce(I,J,bi,bj)
536 & + seaiceMassU(I,J,bi,bj)*recip_deltaT*evpStarFac
537 & * uIceNm1(I,J,bi,bj)
538 & + FORCEX(I,J,bi,bj)
539 & + stressDivergenceX(I,J,bi,bj) )
540 & /( 1. _d 0 - seaiceMaskU(I,J,bi,bj)
541 & + seaiceMassU(I,J,bi,bj)*betaFacP1
542 & + 0.5 _d 0*( DWATN(I,J,bi,bj) + DWATN(I-1,J,bi,bj) )
543 & * COSWAT )
544 vIce(I,J,bi,bj) = seaiceMaskV(I,J,bi,bj) *
545 & ( seaiceMassV(I,J,bi,bj)*betaFac
546 & * vIce(I,J,bi,bj)
547 & + seaiceMassV(I,J,bi,bj)*recip_deltaT*evpStarFac
548 & * vIceNm1(I,J,bi,bj)
549 & + FORCEY(I,J,bi,bj)
550 & + stressDivergenceY(I,J,bi,bj) )
551 & /( 1. _d 0 - seaiceMaskV(I,J,bi,bj)
552 & + seaiceMassV(I,J,bi,bj)*betaFacP1
553 & + 0.5 _d 0*( DWATN(I,J,bi,bj) + DWATN(I,J-1,bi,bj) )
554 & * COSWAT )
555 C-- to change results of lab_sea.hb87 test exp. (only preserve 2 digits for cg2d)
556 c uIce(i,j,bi,bj) = uIceNm1(i,j,bi,bj)
557 c & +( uIce(i,j,bi,bj) - uIceNm1(i,j,bi,bj) )
558 c vIce(i,j,bi,bj) = vIceNm1(i,j,bi,bj)
559 c & +( vIce(i,j,bi,bj) - vIceNm1(i,j,bi,bj) )
560 ENDDO
561 ENDDO
562 #ifndef OBCS_UVICE_OLD
563 DO j=1,sNy
564 DO i=1,sNx
565 locMaskU = maskInC(i,j,bi,bj)*maskInC(i-1,j,bi,bj)
566 locMaskV = maskInC(i,j,bi,bj)*maskInC(i,j-1,bi,bj)
567 uIce(i,j,bi,bj) = uIce(i,j,bi,bj)*locMaskU
568 & + uIceNm1(i,j,bi,bj)*(ONE-locMaskU)
569 vIce(i,j,bi,bj) = vIce(i,j,bi,bj)*locMaskV
570 & + vIceNm1(i,j,bi,bj)*(ONE-locMaskV)
571 ENDDO
572 ENDDO
573 #endif /* OBCS_UVICE_OLD */
574 ENDDO
575 ENDDO
576
577 #ifdef ALLOW_AUTODIFF_TAMC
578 CADJ STORE uIce = comlev1_evp, key = ikeyloc, byte = isbyte
579 CADJ STORE vIce = comlev1_evp, key = ikeyloc, byte = isbyte
580 #endif /* ALLOW_AUTODIFF_TAMC */
581
582 CALL EXCH_UV_XY_RL(uIce,vIce,.TRUE.,myThid)
583
584 ENDIF
585 C end of the main time loop
586 ENDDO
587
588 #endif /* SEAICE_ALLOW_DYNAMICS and SEAICE_CGRID and SEAICE_ALLOW_EVP */
589
590 RETURN
591 END

  ViewVC Help
Powered by ViewVC 1.1.22