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

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

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


Revision 1.61 - (show annotations) (download)
Mon Jun 6 14:56:13 2011 UTC (12 years, 11 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint63a, checkpoint63b, checkpoint63, checkpoint62z
Changes since 1.60: +74 -60 lines
- refine debugLevel criteria when printing messages
- print ilcall in DEBUG_STATS messages

1 C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_lsr.F,v 1.60 2011/01/13 00:13:41 heimbach Exp $
2 C $Name: $
3
4 #include "SEAICE_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: SEAICE_CHECK
8 C !INTERFACE:
9 SUBROUTINE SEAICE_LSR( myTime, myIter, myThid )
10
11 C !DESCRIPTION: \bv
12 C *==========================================================*
13 C | S/R SEAICE_LSR
14 C | o Solve ice momentum equation with an LSR dynamics solver
15 C | (see Zhang and Hibler, JGR, 102, 8691-8702, 1997
16 C | and Zhang and Rothrock, MWR, 131, 845- 861, 2003)
17 C | Written by Jinlun Zhang, PSC/UW, Feb-2001
18 C | zhang@apl.washington.edu
19 C |==========================================================*
20 C | C-grid version by Martin Losch
21 C | Since 2009/03/18: finite-Volume discretization of
22 C | stress divergence that includes all metric terms
23 C \==========================================================*
24 C \ev
25
26 C !USES:
27 IMPLICIT NONE
28
29 C === Global variables ===
30 #include "SIZE.h"
31 #include "EEPARAMS.h"
32 #include "PARAMS.h"
33 #include "DYNVARS.h"
34 #include "GRID.h"
35 #include "SEAICE.h"
36 #include "SEAICE_PARAMS.h"
37
38 #ifdef ALLOW_AUTODIFF_TAMC
39 # include "tamc.h"
40 #endif
41
42 C !INPUT/OUTPUT PARAMETERS:
43 C === Routine arguments ===
44 C myTime :: Simulation time
45 C myIter :: Simulation timestep number
46 C myThid :: my Thread Id. number
47 _RL myTime
48 INTEGER myIter
49 INTEGER myThid
50 CEOP
51
52 #ifdef SEAICE_CGRID
53 #ifdef SEAICE_ALLOW_DYNAMICS
54
55 C !LOCAL VARIABLES:
56 C === Local variables ===
57 C i,j,bi,bj :: Loop counters
58
59 INTEGER ilcall
60 INTEGER i, j, k, m, bi, bj, j1, j2, im, jm
61 INTEGER ICOUNT1, ICOUNT2
62 INTEGER kSrf
63 LOGICAL doIterate
64 CHARACTER*(MAX_LEN_MBUF) msgBuf
65
66 _RL WFAU, WFAV, WFAU1, WFAV1, WFAU2, WFAV2
67 _RL AA3, S1, S2, S1A, S2A
68 _RL hFacM, hFacP
69 _RL eplus, eminus
70
71 C coefficients of ice velocities in coefficient matrix
72 C for both U and V-equation
73 C XX: double derivative in X
74 C YY: double derivative in Y
75 C XM: metric term with derivative in X
76 C YM: metric term with derivative in Y
77 _RL UXX (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
78 _RL UYY (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
79 _RL UXM (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
80 _RL UYM (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
81 _RL VXX (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
82 _RL VYY (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
83 _RL VXM (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
84 _RL VYM (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
85 C diagonals of coefficient matrices
86 _RL AU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
87 _RL BU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
88 _RL CU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
89 _RL AV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
90 _RL BV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
91 _RL CV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
92 _RL FXY (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
93 C coefficients for lateral points (u(j+/-1),v(i+/-1))
94 _RL UVRT1(1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
95 _RL UVRT2(1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
96 C abbreviations
97 _RL etaPlusZeta (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
98 _RL zetaMinusEta(1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
99 _RL etaMeanZ (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
100 C contribution of sigma on righ hand side
101 _RL sig11(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
102 _RL sig22(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
103 _RL sig12(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
104 C auxillary fields
105 _RL URT (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
106 _RL CUU (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
107 _RL VRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
108 _RL CVV (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
109 _RL uTmp (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
110 _RL vTmp (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
111 C
112 _RL COSWAT
113 _RS SINWAT
114 _RL TEMPVAR
115 _RL UERR
116 INTEGER iMin, iMax, jMin, jMax
117 #ifdef SEAICE_VECTORIZE_LSR
118 C in this case, the copy of u(3)=u(1)/v(3)=v(1) needs to include
119 C part of the overlap, because the overlap of u/vTmp is used
120 PARAMETER ( iMin = 0, iMax = sNx+1, jMin = 0, jMax = sNy+1 )
121 #else
122 PARAMETER ( iMin = 1, iMax = sNx, jMin = 1, jMax = sNy )
123 #endif
124 #ifdef SEAICE_ALLOW_CHECK_LSR_CONVERGENCE
125 _RL resnorm, EKnorm, counter
126 #endif /* SEAICE_ALLOW_CHECK_LSR_CONVERGENCE */
127
128 #ifdef ALLOW_AUTODIFF_TAMC
129 cph break artificial dependencies
130 DO bj=myByLo(myThid),myByHi(myThid)
131 DO bi=myBxLo(myThid),myBxHi(myThid)
132 DO j=1-OLy,sNy+OLy
133 DO i=1-OLx,sNx+OLx
134 press(i,j,bi,bj)=0. _d 0
135 zeta(i,j,bi,bj)=0. _d 0
136 eta(i,j,bi,bj)=0. _d 0
137 ENDDO
138 ENDDO
139 ENDDO
140 ENDDO
141 #endif
142
143 #ifdef ALLOW_AUTODIFF_TAMC
144 DO ilcall=1,MPSEUDOTIMESTEPS
145 #else
146 DO ilcall=1,NPSEUDOTIMESTEPS
147 #endif
148 IF ( ilcall .LE. NPSEUDOTIMESTEPS ) THEN
149 #ifdef ALLOW_AUTODIFF_TAMC
150 CADJ STORE uice = comlev1_dynsol, kind=isbyte,
151 CADJ & key = ikey_dynamics + (ilcall-1)*nchklev_1
152 CADJ STORE vice = comlev1_dynsol, kind=isbyte,
153 CADJ & key = ikey_dynamics + (ilcall-1)*nchklev_1
154 CADJ STORE uicenm1 = comlev1_dynsol, kind=isbyte,
155 CADJ & key = ikey_dynamics + (ilcall-1)*nchklev_1
156 CADJ STORE vicenm1 = comlev1_dynsol, kind=isbyte,
157 CADJ & key = ikey_dynamics + (ilcall-1)*nchklev_1
158 #endif /* ALLOW_AUTODIFF_TAMC */
159
160 C surrface level
161 kSrf = 1
162 C-- introduce turning angles
163 SINWAT=SIN(SEAICE_waterTurnAngle*deg2rad)
164 COSWAT=COS(SEAICE_waterTurnAngle*deg2rad)
165
166 C SET SOME VALUES
167 WFAU1=0.95 _d 0
168 WFAV1=0.95 _d 0
169 WFAU2=ZERO
170 WFAV2=ZERO
171
172 S1A=0.80 _d 0
173 S2A=0.80 _d 0
174 WFAU=WFAU1
175 WFAV=WFAV1
176
177 ICOUNT1=SOLV_MAX_ITERS
178 ICOUNT2=SOLV_MAX_ITERS
179
180 k = 1
181
182 IF ( ilcall .EQ. 1 ) THEN
183 C NOW DO PREDICTOR TIME STEP
184 DO bj=myByLo(myThid),myByHi(myThid)
185 DO bi=myBxLo(myThid),myBxHi(myThid)
186 DO j=1-OLy,sNy+OLy
187 DO i=1-OLx,sNx+OLx
188 uIceNm1(i,j,bi,bj)=uIce(I,J,bi,bj)
189 vIceNm1(i,j,bi,bj)=vIce(I,J,bi,bj)
190 uIceC(I,J,bi,bj)=uIce(I,J,bi,bj)
191 vIceC(I,J,bi,bj)=vIce(I,J,bi,bj)
192 ENDDO
193 ENDDO
194 ENDDO
195 ENDDO
196 #ifdef ALLOW_AUTODIFF_TAMC
197 cphCADJ STORE uicec = comlev1_dynsol, kind=isbyte,
198 cphCADJ & key = ikey_dynamics + (ilcall-1)*nchklev_1
199 cphCADJ STORE vicec = comlev1_dynsol, kind=isbyte,
200 cphCADJ & key = ikey_dynamics + (ilcall-1)*nchklev_1
201 #endif
202 ELSE
203 C NOW DO MODIFIED EULER STEP
204 DO bj=myByLo(myThid),myByHi(myThid)
205 DO bi=myBxLo(myThid),myBxHi(myThid)
206 DO j=1-OLy,sNy+OLy
207 DO i=1-OLx,sNx+OLx
208 uIce(I,J,bi,bj)=HALF*(uIce(I,J,bi,bj)+uIceNm1(i,j,bi,bj))
209 vIce(I,J,bi,bj)=HALF*(vIce(I,J,bi,bj)+vIceNm1(i,j,bi,bj))
210 uIceC(I,J,bi,bj)=uIce(I,J,bi,bj)
211 vIceC(I,J,bi,bj)=vIce(I,J,bi,bj)
212 ENDDO
213 ENDDO
214 ENDDO
215 ENDDO
216 ENDIF
217 IF ( ilcall .GT. 2 ) THEN
218 C for additional (pseudo-time)steps update u/vIceNm1
219 DO bj=myByLo(myThid),myByHi(myThid)
220 DO bi=myBxLo(myThid),myBxHi(myThid)
221 DO j=1-OLy,sNy+OLy
222 DO i=1-OLx,sNx+OLx
223 uIceNm1(i,j,bi,bj)=uIce(I,J,bi,bj)
224 vIceNm1(i,j,bi,bj)=vIce(I,J,bi,bj)
225 ENDDO
226 ENDDO
227 ENDDO
228 ENDDO
229 ENDIF
230
231 #ifdef ALLOW_AUTODIFF_TAMC
232 cph That is an important one! Note, that
233 cph * lsr is called twice, thus the icall index
234 cph * this storing is still outside the iteration loop
235 CADJ STORE uice = comlev1_dynsol, kind=isbyte,
236 CADJ & key = ikey_dynamics + (ilcall-1)*nchklev_1
237 CADJ STORE vice = comlev1_dynsol, kind=isbyte,
238 CADJ & key = ikey_dynamics + (ilcall-1)*nchklev_1
239 #endif /* ALLOW_AUTODIFF_TAMC */
240
241 CALL SEAICE_CALC_STRAINRATES(
242 I uIceC, vIceC,
243 O e11, e22, e12,
244 I ilcall, myTime, myIter, myThid )
245
246 CALL SEAICE_CALC_VISCOSITIES(
247 I e11, e22, e12, zMin, zMax, hEffM, press0,
248 O eta, zeta, press,
249 I ilcall, myTime, myIter, myThid )
250
251 DO bj=myByLo(myThid),myByHi(myThid)
252 DO bi=myBxLo(myThid),myBxHi(myThid)
253 DO j=0,sNy+1
254 DO i=0,sNx+1
255 C set up non-linear water drag
256 TEMPVAR = QUART*(
257 & (uIceC(I ,J,bi,bj)-uVel(I ,J,kSrf,bi,bj)
258 & +uIceC(I+1,J,bi,bj)-uVel(I+1,J,kSrf,bi,bj))**2
259 & +(vIceC(I,J ,bi,bj)-vVel(I,J ,kSrf,bi,bj)
260 & +vIceC(I,J+1,bi,bj)-vVel(I,J+1,kSrf,bi,bj))**2)
261 IF ( YC(I,J,bi,bj) .LT. ZERO ) THEN
262 IF ( TEMPVAR .LE. (QUART/SEAICE_waterDrag_south)**2 ) THEN
263 DWATN(I,J,bi,bj)=QUART
264 ELSE
265 DWATN(I,J,bi,bj)=SEAICE_waterDrag_south*SQRT(TEMPVAR)
266 ENDIF
267 ELSE
268 IF ( TEMPVAR .LE. (QUART/SEAICE_waterDrag)**2 ) THEN
269 DWATN(I,J,bi,bj)=QUART
270 ELSE
271 DWATN(I,J,bi,bj)=SEAICE_waterDrag*SQRT(TEMPVAR)
272 ENDIF
273 ENDIF
274 DWATN(I,J,bi,bj) = DWATN(I,J,bi,bj) * HEFFM(I,J,bi,bj)
275 C set up symmettric drag
276 DRAGS(I,J,bi,bj) = DWATN(I,J,bi,bj)*COSWAT
277 ENDDO
278 ENDDO
279 C
280 DO J=1,sNy
281 DO I=1,sNx
282 C set up anti symmettric drag force and add in current force
283 C ( remember to average to correct velocity points )
284 FORCEX(I,J,bi,bj)=FORCEX0(I,J,bi,bj)+
285 & 0.5 _d 0 * ( DWATN(I,J,bi,bj)+DWATN(I-1,J,bi,bj) ) *
286 & COSWAT * uVel(I,J,kSrf,bi,bj)
287 & - SIGN(SINWAT, _fCori(I,J,bi,bj))* 0.5 _d 0 *
288 & ( DWATN(I ,J,bi,bj) * 0.5 _d 0 *
289 & (vVel(I ,J ,kSrf,bi,bj)-vIceC(I ,J ,bi,bj)
290 & +vVel(I ,J+1,kSrf,bi,bj)-vIceC(I ,J+1,bi,bj))
291 & + DWATN(I-1,J,bi,bj) * 0.5 _d 0 *
292 & (vVel(I-1,J ,kSrf,bi,bj)-vIceC(I-1,J ,bi,bj)
293 & +vVel(I-1,J+1,kSrf,bi,bj)-vIceC(I-1,J+1,bi,bj))
294 & )
295 FORCEY(I,J,bi,bj)=FORCEY0(I,J,bi,bj)+
296 & 0.5 _d 0 * ( DWATN(I,J,bi,bj)+DWATN(I,J-1,bi,bj) ) *
297 & COSWAT * vVel(I,J,kSrf,bi,bj)
298 & + SIGN(SINWAT, _fCori(I,J,bi,bj)) * 0.5 _d 0 *
299 & ( DWATN(I,J ,bi,bj) * 0.5 _d 0 *
300 & (uVel(I ,J ,kSrf,bi,bj)-uIceC(I ,J ,bi,bj)
301 & +uVel(I+1,J ,kSrf,bi,bj)-uIceC(I+1,J ,bi,bj))
302 & + DWATN(I,J-1,bi,bj) * 0.5 _d 0 *
303 & (uVel(I ,J-1,kSrf,bi,bj)-uIceC(I ,J-1,bi,bj)
304 & +uVel(I+1,J-1,kSrf,bi,bj)-uIceC(I+1,J-1,bi,bj))
305 & )
306 ENDDO
307 ENDDO
308 C this is the rhs contribution of the time derivative
309 DO j=1,sNy
310 DO i=1,sNx
311 FORCEX(I,J,bi,bj)=FORCEX(I,J,bi,bj)
312 & +seaiceMassU(I,J,bi,bj)/SEAICE_deltaTdyn
313 & *uIceNm1(i,j,bi,bj)
314 FORCEY(I,J,bi,bj)=FORCEY(I,J,bi,bj)
315 & +seaiceMassV(I,J,bi,bj)/SEAICE_deltaTdyn
316 & *vIceNm1(i,j,bi,bj)
317 FORCEX(I,J,bi,bj)=FORCEX(I,J,bi,bj)* seaiceMaskU(I,J,bi,bj)
318 FORCEY(I,J,bi,bj)=FORCEY(I,J,bi,bj)* seaiceMaskV(I,J,bi,bj)
319 ENDDO
320 ENDDO
321 ENDDO
322 ENDDO
323 C
324 C some abbreviations
325 C
326 DO bj=myByLo(myThid),myByHi(myThid)
327 DO bi=myBxLo(myThid),myBxHi(myThid)
328 DO J=0,sNy
329 DO I=0,sNx
330 etaPlusZeta (I,J,bi,bj) = ETA (I,J,bi,bj)+ZETA(I,J,bi,bj)
331 zetaMinusEta(I,J,bi,bj) = ZETA(I,J,bi,bj)-ETA (I,J,bi,bj)
332 ENDDO
333 ENDDO
334 DO J=1,sNy+1
335 DO I=1,sNx+1
336 etaMeanZ (I,J,bi,bj) =
337 & ( ETA (I,J ,bi,bj) + ETA (I-1,J ,bi,bj)
338 & + ETA (I,J-1,bi,bj) + ETA (I-1,J-1,bi,bj) )
339 & / MAX(1.D0,maskC(I,J, k,bi,bj)+maskC(I-1,J, k,bi,bj)
340 & + maskC(I,J-1,k,bi,bj)+maskC(I-1,J-1,k,bi,bj) )
341 ENDDO
342 ENDDO
343 C free-slip means no lateral stress, which is best achieved masking
344 C eta on vorticity(=Z)-points; from now on we only need to worry
345 C about the no-slip boundary conditions
346 IF (.NOT.SEAICE_no_slip) THEN
347 DO J=1,sNy+1
348 DO I=1,sNx+1
349 etaMeanZ (I,J,bi,bj) = etaMeanZ(I,J,bi,bj)
350 & *maskC(I,J, k,bi,bj)*maskC(I-1,J, k,bi,bj)
351 & *maskC(I,J-1,k,bi,bj)*maskC(I-1,J-1,k,bi,bj)
352 ENDDO
353 ENDDO
354 ENDIF
355 C coefficients of uIce(I,J) and vIce(I,J) belonging to ...
356 DO J=1,sNy
357 DO I=0,sNx
358 C ... d/dx (eta+zeta) d/dx u
359 UXX(I,J,bi,bj) = _dyF(I,J,bi,bj) * etaPlusZeta(I,J,bi,bj)
360 & * _recip_dxF(I,J,bi,bj)
361 C ... d/dx (zeta-eta) k1 u
362 UXM(I,J,bi,bj) = _dyF(I,J,bi,bj) * zetaMinusEta(I,J,bi,bj)
363 & * k1AtC(I,J,bi,bj) * 0.5 _d 0
364 ENDDO
365 ENDDO
366 DO J=1,sNy+1
367 DO I=1,sNx
368 C ... d/dy eta d/dy u
369 UYY(I,J,bi,bj) = _dxV(I,J,bi,bj) * etaMeanZ(I,J,bi,bj)
370 & * _recip_dyU(I,J,bi,bj)
371 C ... d/dy eta k2 u
372 UYM(I,J,bi,bj) = _dxV(I,J,bi,bj) * etaMeanZ(I,J,bi,bj)
373 & * k2AtZ(I,J,bi,bj) * 0.5 _d 0
374 ENDDO
375 ENDDO
376 DO J=1,sNy
377 DO I=1,sNx+1
378 C ... d/dx eta dv/dx
379 VXX(I,J,bi,bj) = _dyU(I,J,bi,bj) * etaMeanZ(I,J,bi,bj)
380 & * _recip_dxV(I,J,bi,bj)
381 C ... d/dx eta k1 v
382 VXM(I,J,bi,bj) = _dyU(I,J,bi,bj) * etaMeanZ(I,J,bi,bj)
383 & * k1AtZ(I,J,bi,bj) * 0.5 _d 0
384 ENDDO
385 ENDDO
386 DO J=0,sNy
387 DO I=1,sNx
388 C ... d/dy eta+zeta dv/dy
389 VYY(I,J,bi,bj) = _dxF(I,J,bi,bj) * etaPlusZeta(I,J,bi,bj)
390 & * _recip_dyF(I,J,bi,bj)
391 C ... d/dy (zeta-eta) k2 v
392 VYM(I,J,bi,bj) = _dxF(I,J,bi,bj) * zetaMinusEta(I,J,bi,bj)
393 & * k2AtC(I,J,bi,bj) * 0.5 _d 0
394 ENDDO
395 ENDDO
396 ENDDO
397 ENDDO
398
399 C SOLVE FOR uIce
400 DO bj=myByLo(myThid),myByHi(myThid)
401 DO bi=myBxLo(myThid),myBxHi(myThid)
402 C assemble coefficient matrix, beware of sign convention: because this
403 C is the left hand side we calculate -grad(sigma), but the coefficients
404 C of U(I,J+/-1) are counted on the right hand side
405 DO J=1,sNy
406 DO I=1,sNx
407 C coefficients for UICE(I-1,J)
408 AU(I,J,bi,bj)= ( - UXX(I-1,J,bi,bj) + UXM(I-1,J,bi,bj) )
409 & * seaiceMaskU(I,J,bi,bj)
410 C coefficients for UICE(I+1,J)
411 CU(I,J,bi,bj)= ( - UXX(I ,J,bi,bj) - UXM(I ,J,bi,bj) )
412 & * seaiceMaskU(I,J,bi,bj)
413 C coefficients for UICE(I,J)
414 BU(I,J,bi,bj)=(ONE - seaiceMaskU(I,J,bi,bj)) +
415 & ( UXX(I-1,J ,bi,bj) + UXX(I,J,bi,bj)
416 & + UYY(I ,J+1,bi,bj) + UYY(I,J,bi,bj)
417 & + UXM(I-1,J ,bi,bj) - UXM(I,J,bi,bj)
418 & + UYM(I ,J+1,bi,bj) - UYM(I,J,bi,bj)
419 & ) * seaiceMaskU(I,J,bi,bj)
420 C coefficients of uIce(I,J-1)
421 UVRT1(I,J,bi,bj)= UYY(I,J ,bi,bj) + UYM(I,J ,bi,bj)
422 C coefficients of uIce(I,J+1)
423 UVRT2(I,J,bi,bj)= UYY(I,J+1,bi,bj) - UYM(I,J+1,bi,bj)
424 ENDDO
425 ENDDO
426
427 C apply boundary conditions according to slip factor
428 C for no slip, set u on boundary to zero: u(j+/-1)=-u(j)
429 C for the free slip case sigma_12 = 0
430 DO J=1,sNy
431 DO I=1,sNx
432 hFacM = seaiceMaskU(I,J-1,bi,bj)
433 hFacP = seaiceMaskU(I,J+1,bi,bj)
434 C copy contributions to coefficient of U(I,J)
435 C beware of sign convection: UVRT1/2 have the opposite sign convention
436 C than BU, hence the minus sign
437 BU(I,J,bi,bj)=BU(I,J,bi,bj) + seaiceMaskU(I,J,bi,bj) *
438 & ( ( 1. _d 0 - hFacM )
439 & * ( UYY(I ,J ,bi,bj) + UYM(I ,J ,bi,bj) )
440 & + ( 1. _d 0 - hFacP )
441 & * ( UYY(I ,J+1,bi,bj) - UYM(I ,J+1,bi,bj) ) )
442 C reset coefficients of U(I,J-1) and U(I,J+1)
443 UVRT1(I,J,bi,bj)=UVRT1(I,J,bi,bj) * hFacM
444 UVRT2(I,J,bi,bj)=UVRT2(I,J,bi,bj) * hFacP
445 ENDDO
446 ENDDO
447
448 C now we need to normalize everything by the grid cell area
449 DO J=1,sNy
450 DO I=1,sNx
451 AU(I,J,bi,bj) = AU(I,J,bi,bj) * recip_rAw(I,J,bi,bj)
452 CU(I,J,bi,bj) = CU(I,J,bi,bj) * recip_rAw(I,J,bi,bj)
453 C here we need ad in the contribution from the time derivative
454 C and the symmetric drag term; must be done after normalizing
455 BU(I,J,bi,bj) = BU(I,J,bi,bj) * recip_rAw(I,J,bi,bj)
456 & + seaiceMaskU(I,J,bi,bj) *
457 & ( seaiceMassU(I,J,bi,bj)/SEAICE_deltaTdyn
458 & + 0.5 _d 0*( DRAGS(I,J,bi,bj) + DRAGS(I-1,J,bi,bj) ) )
459 UVRT1(I,J,bi,bj) = UVRT1(I,J,bi,bj) * recip_rAw(I,J,bi,bj)
460 UVRT2(I,J,bi,bj) = UVRT2(I,J,bi,bj) * recip_rAw(I,J,bi,bj)
461 ENDDO
462 ENDDO
463
464 DO J=1,sNy
465 AU(1,J,bi,bj)=ZERO
466 CU(sNx,J,bi,bj)=ZERO
467 CU(1,J,bi,bj)=CU(1,J,bi,bj)/BU(1,J,bi,bj)
468 ENDDO
469
470 C now set up right-hand side
471 C contribution of sigma11 to rhs
472 DO J=1,sNy
473 DO I=0,sNx
474 sig11(I,J) = zetaMinusEta(I,J,bi,bj)
475 & * ( vIceC(I,J+1,bi,bj) - vIceC(I,J,bi,bj) )
476 & * _recip_dyF(I,J,bi,bj)
477 & + etaPlusZeta(I,J,bi,bj) * k2AtC(I,J,bi,bj)
478 & * 0.5 _d 0 * ( vIceC(I,J+1,bi,bj) + vIceC(I,J,bi,bj) )
479 & - 0.5 _d 0 * PRESS(I,J,bi,bj)
480 ENDDO
481 ENDDO
482 C contribution of sigma12 to rhs of u-equation
483 DO J=1,sNy+1
484 DO I=1,sNx
485 hFacM = seaiceMaskV(I,J,bi,bj) - seaiceMaskV(I-1,J,bi,bj)
486 sig12(I,J) = etaMeanZ(I,J,bi,bj) * (
487 & ( vIceC(I,J,bi,bj) - vIceC(I-1,J,bi,bj) )
488 & * _recip_dxV(I,J,bi,bj)
489 & - k1AtZ(I,J,bi,bj)
490 & * 0.5 _d 0 * ( vIceC(I,J,bi,bj) + vIceC(I-1,J,bi,bj) )
491 & )
492 C free slip conditions (sig12=0) are taken care of by masking sig12
493 & *maskC(I ,J ,k,bi,bj)*maskC(I-1,J ,k,bi,bj)
494 & *maskC(I ,J-1,k,bi,bj)*maskC(I-1,J-1,k,bi,bj)
495 C no slip boundary conditions (v(i-1)=-v(i))
496 C v(i)+v(i-1) = 0 is also taken care of by masking sig12, so that we
497 C only need to deal with v(i)-v(i-1)
498 & + etaMeanZ(I,J,bi,bj) * _recip_dxV(I,J,bi,bj)
499 & * ( vIceC(I,J,bi,bj) + vIceC(I-1,J,bi,bj) )
500 & * hFacM * 2. _d 0
501 ENDDO
502 ENDDO
503
504 DO J=1,sNy
505 DO I=1,sNx
506 C coriolis and other forcing
507 FXY(I,J,bi,bj)=0.5 _d 0 *
508 & ( seaiceMassC(I ,J,bi,bj) * _fCori(I ,J,bi,bj)
509 & *0.5 _d 0*( vIceC( i ,j,bi,bj)+vIceC( i ,j+1,bi,bj) )
510 & + seaiceMassC(I-1,J,bi,bj) * _fCori(I-1,J,bi,bj)
511 & *0.5 _d 0*( vIceC(i-1,j,bi,bj)+vIceC(i-1,j+1,bi,bj) ) )
512 & +FORCEX(I,J,bi,bj)
513 C contribution to the rhs part of grad(sigma)_x
514 & + recip_rAw(I,J,bi,bj) * seaiceMaskU(I,J,bi,bj) *
515 & ( _dyF(I ,J ,bi,bj)*sig11(I ,J )
516 & - _dyF(I-1,J ,bi,bj)*sig11(I-1,J )
517 & + _dxV(I ,J+1,bi,bj)*sig12(I ,J+1)
518 & - _dxV(I ,J ,bi,bj)*sig12(I ,J ) )
519 ENDDO
520 ENDDO
521
522 ENDDO
523 ENDDO
524
525 #ifdef ALLOW_DEBUG
526 IF ( debugLevel .GE. debLevD ) THEN
527 WRITE(msgBuf,'(A,I3,A)')
528 & 'Uice pre iter (SEAICE_LSR', MOD(ilcall,1000), ')'
529 CALL DEBUG_STATS_RL( 1, UICE, msgBuf, myThid )
530 ENDIF
531 #endif /* ALLOW_DEBUG */
532 C NOW DO ITERATION
533
534 cph--- iteration starts here
535 cph--- need to kick out goto
536 doIterate = .TRUE.
537
538 C ITERATION START -----------------------------------------------------
539 #ifdef ALLOW_AUTODIFF_TAMC
540 CADJ LOOP = iteration uice
541 #endif /* ALLOW_AUTODIFF_TAMC */
542 DO m = 1,solv_max_iters
543 IF ( doIterate ) THEN
544
545 DO bj=myByLo(myThid),myByHi(myThid)
546 DO bi=myBxLo(myThid),myBxHi(myThid)
547 C NOW SET U(3)=U(1)
548 C save uIce prior to iteration
549 DO J=jMin,jMax
550 DO I=1,sNx
551 uTmp(I,J,bi,bj)=uIce(I,J,bi,bj)
552 ENDDO
553 ENDDO
554
555 DO J=1,sNy
556 DO I=1,sNx
557 IF(I.EQ.1) THEN
558 AA3=( UXX(I-1,J,bi,bj) - UXM(I-1,J,bi,bj) )
559 & * uIce(I-1,J,bi,bj) * seaiceMaskU(I,J,bi,bj)
560 ELSEIF(I.EQ.sNx) THEN
561 AA3=( UXX(I ,J,bi,bj) + UXM(I ,J,bi,bj) )
562 & * uIce(I+1,J,bi,bj) * seaiceMaskU(I,J,bi,bj)
563 ELSE
564 AA3=ZERO
565 ENDIF
566 URT(I,J)=FXY(I,J,bi,bj)
567 & + AA3 * recip_rAw(I,J,bi,bj)
568 #ifdef SEAICE_VECTORIZE_LSR
569 & + UVRT1(I,J,bi,bj)*uTmp(I,J-1,bi,bj)
570 & + UVRT2(I,J,bi,bj)*uTmp(I,J+1,bi,bj)
571 #else
572 & + UVRT1(I,J,bi,bj)*uIce(I,J-1,bi,bj)
573 & + UVRT2(I,J,bi,bj)*uIce(I,J+1,bi,bj)
574 #endif /* SEAICE_VECTORIZE_LSR */
575 URT(I,J)=URT(I,J)* seaiceMaskU(I,J,bi,bj)
576 ENDDO
577
578 DO I=1,sNx
579 CUU(I,J)=CU(I,J,bi,bj)
580 ENDDO
581 URT(1,J)=URT(1,J)/BU(1,J,bi,bj)
582 #ifdef SEAICE_VECTORIZE_LSR
583 ENDDO
584 C start a new loop with reversed order to support automatic vectorization
585 DO I=2,sNx
586 IM=I-1
587 DO J=1,sNy
588 #else /* do not SEAICE_VECTORIZE_LSR */
589 DO I=2,sNx
590 IM=I-1
591 #endif /* SEAICE_VECTORIZE_LSR */
592 CUU(I,J)=CUU(I,J)/(BU(I,J,bi,bj)-AU(I,J,bi,bj)*CUU(IM,J))
593 URT(I,J)=(URT(I,J)-AU(I,J,bi,bj)*URT(IM,J))
594 & /(BU(I,J,bi,bj)-AU(I,J,bi,bj)*CUU(IM,J))
595 ENDDO
596 #ifdef SEAICE_VECTORIZE_LSR
597 ENDDO
598 C go back to original order
599 DO J=1,sNy
600 #endif /* SEAICE_VECTORIZE_LSR */
601 DO I=1,sNx-1
602 J1=sNx-I
603 J2=J1+1
604 URT(J1,J)=URT(J1,J)-CUU(J1,J)*URT(J2,J)
605 ENDDO
606 DO I=1,sNx
607 uIce(I,J,bi,bj)=uTmp(I,J,bi,bj)
608 & +WFAU*(URT(I,J)-uTmp(I,J,bi,bj))
609 ENDDO
610 ENDDO
611 C end bi,bj-loops
612 ENDDO
613 ENDDO
614
615 IF(MOD(m,SOLV_NCHECK).EQ.0) THEN
616 S1=ZERO
617 DO bj=myByLo(myThid),myByHi(myThid)
618 DO bi=myBxLo(myThid),myBxHi(myThid)
619 DO J=1,sNy
620 DO I=1,sNx
621 UERR=(uIce(I,J,bi,bj)-uTmp(I,J,bi,bj))
622 & * seaiceMaskU(I,J,bi,bj)
623 S1=MAX(ABS(UERR),S1)
624 ENDDO
625 ENDDO
626 ENDDO
627 ENDDO
628 _GLOBAL_MAX_RL( S1, myThid )
629 c WRITE(standardMessageUnit,'(A,2I6,1P4E16.9)')
630 c & ' U iters,error,WF = ',ilcall,M,S1,S1A,WFAU
631 C SAFEGUARD AGAINST BAD FORCING ETC
632 IF(m.GT.1.AND.S1.GT.S1A) WFAU=WFAU2
633 S1A=S1
634 IF(S1.LT.LSR_ERROR) THEN
635 ICOUNT1=m
636 doIterate = .FALSE.
637 ENDIF
638 ENDIF
639 CALL EXCH_UV_XY_RL( uIce, vIce,.TRUE.,myThid)
640
641 ENDIF
642 ENDDO
643 C ITERATION END -----------------------------------------------------
644
645 #ifdef ALLOW_DEBUG
646 IF ( debugLevel .GE. debLevC ) THEN
647 _BEGIN_MASTER( myThid )
648 WRITE(standardMessageUnit,'(A,I6,1P2E22.14)')
649 & ' U lsr iters, error = ',ICOUNT1,S1
650 _END_MASTER( myThid )
651 IF ( debugLevel .GE. debLevD ) THEN
652 WRITE(msgBuf,'(A,I3,A)')
653 & 'Uice post iter (SEAICE_LSR', MOD(ilcall,1000), ')'
654 CALL DEBUG_STATS_RL( 1, UICE, msgBuf, myThid )
655 ENDIF
656 IF (S1.EQ.0.D0.AND.ICOUNT1.GT.SOLV_NCHECK)
657 & STOP 'ABNORMAL END: S/R SEAICE_LSR did not converge (uIce)'
658 ENDIF
659 #endif /* ALLOW_DEBUG */
660
661 C NOW FOR vIce
662 DO bj=myByLo(myThid),myByHi(myThid)
663 DO bi=myBxLo(myThid),myBxHi(myThid)
664 C assemble coefficient matrix, beware of sign convention: because this
665 C is the left hand side we calculate -grad(sigma), but the coefficients
666 C of U(I,J+/-1) are counted on the right hand side
667 DO J=1,sNy
668 DO I=1,sNx
669 C coefficients for VICE(I,J-1)
670 AV(I,J,bi,bj)=( - VYY(I,J-1,bi,bj) + VYM(I,J-1,bi,bj)
671 & ) * seaiceMaskV(I,J,bi,bj)
672 C coefficients for VICE(I,J+1)
673 CV(I,J,bi,bj)=( - VYY(I,J ,bi,bj) - VYM(I,J ,bi,bj)
674 & ) * seaiceMaskV(I,J,bi,bj)
675 C coefficients for VICE(I,J)
676 BV(I,J,bi,bj)= (ONE - seaiceMaskV(I,J,bi,bj)) +
677 & ( VXX(I,J,bi,bj) + VXX(I+1,J ,bi,bj)
678 & + VYY(I,J,bi,bj) + VYY(I ,J-1,bi,bj)
679 & - VXM(I,J,bi,bj) + VXM(I+1,J ,bi,bj)
680 & - VYM(I,J,bi,bj) + VYM(I ,J-1,bi,bj)
681 & ) * seaiceMaskV(I,J,bi,bj)
682 C coefficients for V(I-1,J)
683 UVRT1(I,J,bi,bj) = VXX(I ,J,bi,bj) + VXM(I ,J,bi,bj)
684 C coefficients for V(I+1,J)
685 UVRT2(I,J,bi,bj) = VXX(I+1,J,bi,bj) - VXM(I+1,J,bi,bj)
686 ENDDO
687 ENDDO
688
689 C apply boundary conditions according to slip factor
690 C for no slip, set u on boundary to zero: v(i+/-1)=-v(i)
691 C for the free slip case sigma_12 = 0
692 DO J=1,sNy
693 DO I=1,sNx
694 hFacM = seaiceMaskV(i-1,j,bi,bj)
695 hFacP = seaiceMaskV(i+1,j,bi,bj)
696 C copy contributions to coefficient of V(I,J)
697 C beware of sign convection: UVRT1/2 have the opposite sign convention
698 C than BV, hence the minus sign
699 BV(I,J,bi,bj)=BV(I,J,bi,bj) + seaiceMaskV(I,J,bi,bj) *
700 & ( ( 1. _d 0 - hFacM )
701 & * ( VXX(I ,J,bi,bj) + VXM(I ,J,bi,bj) )
702 & + ( 1. _d 0 - hFacP )
703 & * ( VXX(I+1,J,bi,bj) - VXM(I+1,J,bi,bj) ) )
704 C reset coefficients of V(I-1,J) and V(I+1,J)
705 UVRT1(I,J,bi,bj)=UVRT1(I,J,bi,bj) * hFacM
706 UVRT2(I,J,bi,bj)=UVRT2(I,J,bi,bj) * hFacP
707 ENDDO
708 ENDDO
709
710 C now we need to normalize everything by the grid cell area
711 DO J=1,sNy
712 DO I=1,sNx
713 AV(I,J,bi,bj) = AV(I,J,bi,bj) * recip_rAs(I,J,bi,bj)
714 CV(I,J,bi,bj) = CV(I,J,bi,bj) * recip_rAs(I,J,bi,bj)
715 C here we need ad in the contribution from the time derivative
716 C and the symmetric drag term; must be done after normalizing
717 BV(I,J,bi,bj) = BV(I,J,bi,bj) * recip_rAs(I,J,bi,bj)
718 & + seaiceMaskV(I,J,bi,bj) *
719 & ( seaiceMassV(I,J,bi,bj)/SEAICE_deltaTdyn
720 & + 0.5 _d 0 * ( DRAGS(I,J,bi,bj) + DRAGS(I,J-1,bi,bj) ) )
721 UVRT1(I,J,bi,bj) = UVRT1(I,J,bi,bj) * recip_rAs(I,J,bi,bj)
722 UVRT2(I,J,bi,bj) = UVRT2(I,J,bi,bj) * recip_rAs(I,J,bi,bj)
723 ENDDO
724 ENDDO
725
726 DO I=1,sNx
727 AV(I,1,bi,bj)=ZERO
728 CV(I,sNy,bi,bj)=ZERO
729 CV(I,1,bi,bj)=CV(I,1,bi,bj)/BV(I,1,bi,bj)
730 ENDDO
731
732 C now set up right-hand-side
733 C contribution of sigma22 to rhs
734 DO J=0,sNy
735 DO I=1,sNx
736 sig22(I,J) = zetaMinusEta(I,J,bi,bj)
737 & * ( uIceC(I+1,J,bi,bj) - uIceC(I,J,bi,bj) )
738 & * _recip_dxF(I,J,bi,bj)
739 & + etaPlusZeta(I,J,bi,bj) * k1AtC(I,J,bi,bj)
740 & * 0.5 _d 0 * ( uIceC(I+1,J,bi,bj) + uIceC(I,J,bi,bj) )
741 & - 0.5 _d 0 * PRESS(I,J,bi,bj)
742 ENDDO
743 ENDDO
744 C contribution of sigma12 to rhs of v-equation
745 DO J=1,sNy
746 DO I=1,sNx+1
747 hFacM = seaiceMaskU(i,j,bi,bj) - seaiceMaskU(i,j-1,bi,bj)
748 sig12(I,J) = etaMeanZ(I,J,bi,bj) * (
749 & ( uIceC(I,J,bi,bj) - uIceC(I,J-1,bi,bj) )
750 & * _recip_dyU(I,J,bi,bj)
751 & - k2AtZ(I,J,bi,bj)
752 & * 0.5 _d 0 * ( uIceC(I,J,bi,bj) + uIceC(I,J-1,bi,bj) )
753 & )
754 C free slip conditions (sig12=0) are taken care of by masking sig12,
755 & *maskC(I ,J ,k,bi,bj)*maskC(I-1,J ,k,bi,bj)
756 & *maskC(I ,J-1,k,bi,bj)*maskC(I-1,J-1,k,bi,bj)
757 C no slip boundary conditions (u(j-1)=-u(j))
758 C u(j)+u(j-1) = 0 is also taken care of by masking sig12, so that we
759 C only need to deal with u(j)-u(j-1)
760 & + etaMeanZ(I,J,bi,bj) * _recip_dyU(I,J,bi,bj)
761 & * ( uIceC(I,J,bi,bj) + uIceC(I,J-1,bi,bj) )
762 & * hFacM * 2. _d 0
763 ENDDO
764 ENDDO
765
766 DO J=1,sNy
767 DO I=1,sNx
768 C coriols and other foring
769 FXY(I,J,bi,bj)= - 0.5 _d 0 *
770 & ( seaiceMassC(I,J ,bi,bj) * _fCori(I,J ,bi,bj)
771 & *0.5 _d 0*( uIceC(i ,j ,bi,bj)+uIceC(i+1, j,bi,bj) )
772 & + seaiceMassC(I,J-1,bi,bj) * _fCori(I,J-1,bi,bj)
773 & *0.5 _d 0*( uIceC(i ,j-1,bi,bj)+uIceC(i+1,j-1,bi,bj) ) )
774 & + FORCEY(I,J,bi,bj)
775 C contribution to the rhs part of grad(sigma)_y
776 & + recip_rAs(I,J,bi,bj) * seaiceMaskV(I,J,bi,bj) *
777 & ( _dyU(I+1,J ,bi,bj) * sig12(I+1,J )
778 & - _dyU(I ,J ,bi,bj) * sig12(I ,J )
779 & + _dxF(I ,J ,bi,bj) * sig22(I ,J )
780 & - _dxF(I ,J-1,bi,bj) * sig22(I ,J-1) )
781 ENDDO
782 ENDDO
783
784 ENDDO
785 ENDDO
786
787 #ifdef ALLOW_DEBUG
788 IF ( debugLevel .GE. debLevD ) THEN
789 WRITE(msgBuf,'(A,I3,A)')
790 & 'Vice pre iter (SEAICE_LSR', MOD(ilcall,1000), ')'
791 CALL DEBUG_STATS_RL( 1, VICE, msgBuf, myThid )
792 ENDIF
793 #endif /* ALLOW_DEBUG */
794
795 C NOW DO ITERATION
796
797 cph--- iteration starts here
798 cph--- need to kick out goto
799 doIterate = .TRUE.
800
801 C ITERATION START -----------------------------------------------------
802 #ifdef ALLOW_AUTODIFF_TAMC
803 CADJ LOOP = iteration vice
804 #endif /* ALLOW_AUTODIFF_TAMC */
805 DO m = 1,solv_max_iters
806 IF ( doIterate ) THEN
807
808 DO bj=myByLo(myThid),myByHi(myThid)
809 DO bi=myBxLo(myThid),myBxHi(myThid)
810 C NOW SET V(3)=V(1)
811 C save vIce prior to iteration
812 DO J=1,sNy
813 DO I=iMin,iMax
814 vTmp(I,J,bi,bj)=vIce(I,J,bi,bj)
815 ENDDO
816 ENDDO
817
818 DO I=1,sNx
819 DO J=1,sNy
820 IF(J.EQ.1) THEN
821 AA3=( VYY(I,J-1,bi,bj) - VYM(I,J-1,bi,bj)
822 & ) * vIce(I,J-1,bi,bj) * seaiceMaskV(I,J,bi,bj)
823 ELSEIF(J.EQ.sNy) THEN
824 AA3=( VYY(I,J ,bi,bj) + VYM(I,J ,bi,bj)
825 & ) * vIce(I,J+1,bi,bj) * seaiceMaskV(I,J,bi,bj)
826 ELSE
827 AA3=ZERO
828 ENDIF
829
830 VRT(I,J)=FXY(I,J,bi,bj)
831 & + AA3 * recip_rAs(I,J,bi,bj)
832 #ifdef SEAICE_VECTORIZE_LSR
833 & + UVRT1(I,J,bi,bj)*vTmp(I-1,J,bi,bj)
834 & + UVRT2(I,J,bi,bj)*vTmp(I+1,J,bi,bj)
835 #else
836 & + UVRT1(I,J,bi,bj)*vIce(I-1,J,bi,bj)
837 & + UVRT2(I,J,bi,bj)*vIce(I+1,J,bi,bj)
838 #endif /* SEAICE_VECTORIZE_LSR */
839 VRT(I,J)=VRT(I,J)* seaiceMaskV(I,J,bi,bj)
840 ENDDO
841
842 DO J=1,sNy
843 CVV(I,J)=CV(I,J,bi,bj)
844 ENDDO
845 VRT(I,1)=VRT(I,1)/BV(I,1,bi,bj)
846 DO J=2,sNy
847 JM=J-1
848 CVV(I,J)=CVV(I,J)/(BV(I,J,bi,bj)-AV(I,J,bi,bj)*CVV(I,JM))
849 VRT(I,J)=(VRT(I,J)-AV(I,J,bi,bj)*VRT(I,JM))
850 & /(BV(I,J,bi,bj)-AV(I,J,bi,bj)*CVV(I,JM))
851 ENDDO
852 DO J=1,sNy-1
853 J1=sNy-J
854 J2=J1+1
855 VRT(I,J1)=VRT(I,J1)-CVV(I,J1)*VRT(I,J2)
856 ENDDO
857 DO J=1,sNy
858 vIce(I,J,bi,bj)=vTmp(I,J,bi,bj)
859 & +WFAV*(VRT(I,J)-vTmp(I,J,bi,bj))
860 ENDDO
861 ENDDO
862 C end bi,bj-loops
863 ENDDO
864 ENDDO
865
866 IF(MOD(m,SOLV_NCHECK).EQ.0) THEN
867 S2=ZERO
868 DO bj=myByLo(myThid),myByHi(myThid)
869 DO bi=myBxLo(myThid),myBxHi(myThid)
870 DO J=1,sNy
871 DO I=1,sNx
872 UERR=(vIce(I,J,bi,bj)-vTmp(I,J,bi,bj))
873 & * seaiceMaskV(I,J,bi,bj)
874 S2=MAX(ABS(UERR),S2)
875 ENDDO
876 ENDDO
877 ENDDO
878 ENDDO
879 _GLOBAL_MAX_RL( S2, myThid )
880 C SAFEGUARD AGAINST BAD FORCING ETC
881 IF(m.GT.1.AND.S2.GT.S2A) WFAV=WFAV2
882 S2A=S2
883 IF(S2.LT.LSR_ERROR) THEN
884 ICOUNT2=m
885 doIterate = .FALSE.
886 ENDIF
887 ENDIF
888 CALL EXCH_UV_XY_RL( uIce, vIce,.TRUE.,myThid)
889
890 ENDIF
891 ENDDO
892 C ITERATION END -----------------------------------------------------
893
894 #ifdef ALLOW_DEBUG
895 IF ( debugLevel .GE. debLevC ) THEN
896 _BEGIN_MASTER( myThid )
897 WRITE(standardMessageUnit,'(A,I6,1P2E22.14)')
898 & ' V lsr iters, error = ',ICOUNT2,S2
899 _END_MASTER( myThid )
900 IF ( debugLevel .GE. debLevD ) THEN
901 WRITE(msgBuf,'(A,I3,A)')
902 & 'Vice post iter (SEAICE_LSR', MOD(ilcall,1000), ')'
903 CALL DEBUG_STATS_RL( 1, VICE, msgBuf, myThid )
904 ENDIF
905 IF (S2.EQ.0.D0.AND.ICOUNT2.GT.SOLV_NCHECK)
906 & STOP 'ABNORMAL END: S/R SEAICE_LSR did not converge (vIce)'
907 ENDIF
908 #endif /* ALLOW_DEBUG */
909
910 C APPLY MASKS
911 DO bj=myByLo(myThid),myByHi(myThid)
912 DO bi=myBxLo(myThid),myBxHi(myThid)
913 DO J=1-Oly,sNy+Oly
914 DO I=1-Olx,sNx+Olx
915 uIce(I,J,bi,bj)=uIce(I,J,bi,bj)* seaiceMaskU(I,J,bi,bj)
916 vIce(I,J,bi,bj)=vIce(I,J,bi,bj)* seaiceMaskV(I,J,bi,bj)
917 ENDDO
918 ENDDO
919 ENDDO
920 ENDDO
921 CML CALL EXCH_UV_XY_RL( uIce, vIce,.TRUE.,myThid)
922
923 #ifdef SEAICE_ALLOW_CHECK_LSR_CONVERGENCE
924 IF ( debugLevel .GE. debLevC ) THEN
925 resnorm = 0. _d 0
926 EKnorm = 0. _d 0
927 counter = 0. _d 0
928 DO bj=myByLo(myThid),myByHi(myThid)
929 DO bi=myBxLo(myThid),myBxHi(myThid)
930 C-- Compute residual norm and print
931 DO j=1,sNy
932 DO i=1,sNx
933 resnorm = resnorm + 0.5 _d 0 *
934 & ( ( (uIceNm1(i,j,bi,bj)+uIceNm1(i+1,j,bi,bj))
935 & - (uIce(i,j,bi,bj)+uIce(i+1,j,bi,bj)) )**2
936 & + ( (vIceNm1(i,j,bi,bj)+vIceNm1(i,j+1,bi,bj))
937 & - (vIce(i,j,bi,bj)+vIce(i,j+1,bi,bj)) )**2 )
938 IF ( area(i,j,bi,bj) .gt. 0.5 _d 0 ) THEN
939 EKnorm = EKnorm + 0.5 _d 0 * heff(i,j,bi,bj) *
940 & ( ( (uIce(i,j,bi,bj)+uIce(i+1,j,bi,bj)) )**2
941 & + ( (vIce(i,j,bi,bj)+vIce(i,j+1,bi,bj)) )**2 )
942 counter = counter + 1. _d 0
943 ENDIF
944 ENDDO
945 ENDDO
946 ENDDO
947 ENDDO
948 _GLOBAL_SUM_RL( resnorm, myThid )
949 _GLOBAL_SUM_RL( EKnorm, myThid )
950 _GLOBAL_SUM_RL( counter, myThid )
951 IF ( counter .gt. 0. _d 0 ) EKnorm = EKnorm/counter
952 _BEGIN_MASTER( myThid )
953 WRITE(standardMessageUnit,'(A,I7,1X,2E22.14)')
954 & 'S/R SEAICE_LSR: IPSEUDO, RESNORM, EKNORM = ',
955 & ilcall, sqrt(resnorm), EKnorm
956 _END_MASTER( myThid )
957 ENDIF
958 #endif /* SEAICE_ALLOW_CHECK_LSR_CONVERGENCE */
959
960 ENDIF
961 C end outer pseudo-time stepping loop
962 ENDDO
963
964 IF ( useHB87StressCoupling ) THEN
965 C compute the divergence of stress here to be used later
966 C
967 C compute strain rate from latest velocities
968 CALL SEAICE_CALC_STRAINRATES(
969 I uIce, vIce,
970 O e11, e22, e12,
971 I 3, myTime, myIter, myThid )
972 C compute internal stresses with updated ice velocities
973 DO bj=myByLo(myThid),myByHi(myThid)
974 DO bi=myBxLo(myThid),myBxHi(myThid)
975 DO j=1-Oly,sNy+Oly
976 DO i=1-Olx,sNx+Olx
977 sig11(I,J) = 0. _d 0
978 sig22(I,J) = 0. _d 0
979 sig12(I,J) = 0. _d 0
980 ENDDO
981 ENDDO
982
983 DO j=0,sNy
984 DO i=0,sNx
985 eplus = e11(I,J,bi,bj) + e22(I,J,bi,bj)
986 eminus= e11(I,J,bi,bj) - e22(I,J,bi,bj)
987 sig11(I,J) = zeta(I,J,bi,bj)*eplus + eta(I,J,bi,bj)*eminus
988 & - 0.5 _d 0 * PRESS(I,J,bi,bj)
989 sig22(I,J) = zeta(I,J,bi,bj)*eplus - eta(I,J,bi,bj)*eminus
990 & - 0.5 _d 0 * PRESS(I,J,bi,bj)
991 ENDDO
992 ENDDO
993
994 DO j=1,sNy+1
995 DO i=1,sNx+1
996 sig12(I,J) = 2. _d 0 * e12(I,J,bi,bj) * etaMeanZ(I,J,bi,bj)
997 ENDDO
998 ENDDO
999 C evaluate divergence of stress and apply to forcing
1000 DO J=1,sNy
1001 DO I=1,sNx
1002 stressDivergenceX(I,J,bi,bj) =
1003 & ( sig11(I ,J ) * _dyF(I ,J ,bi,bj)
1004 & - sig11(I-1,J ) * _dyF(I-1,J ,bi,bj)
1005 & + sig12(I ,J+1) * _dxV(I ,J+1,bi,bj)
1006 & - sig12(I ,J ) * _dxV(I ,J ,bi,bj)
1007 & ) * recip_rAw(I,J,bi,bj)
1008 stressDivergenceY(I,J,bi,bj) =
1009 & ( sig22(I ,J ) * _dxF(I ,J ,bi,bj)
1010 & - sig22(I ,J-1) * _dxF(I ,J-1,bi,bj)
1011 & + sig12(I+1,J ) * _dyU(I+1,J ,bi,bj)
1012 & - sig12(I ,J ) * _dyU(I ,J ,bi,bj)
1013 & ) * recip_rAs(I,J,bi,bj)
1014 ENDDO
1015 ENDDO
1016 ENDDO
1017 ENDDO
1018 C endif useHB87StressCoupling
1019 ENDIF
1020
1021 #endif /* SEAICE_ALLOW_DYNAMICS */
1022 #endif /* SEAICE_CGRID */
1023
1024 RETURN
1025 END

  ViewVC Help
Powered by ViewVC 1.1.22