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

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

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


Revision 1.21 - (show annotations) (download)
Thu Apr 4 07:02:51 2013 UTC (12 years, 3 months ago) by mlosch
Branch: MAIN
Changes since 1.20: +3 -5 lines
simplify the use of CPP flags (compile when SEAICE_ALLOW_JFNK is defined)

1 C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_jfnk.F,v 1.20 2013/03/02 04:35:05 jmc Exp $
2 C $Name: $
3
4 #include "SEAICE_OPTIONS.h"
5
6 C-- File seaice_jfnk.F: seaice jfnk dynamical solver S/R:
7 C-- Contents
8 C-- o SEAICE_JFNK
9 C-- o SEAICE_JFNK_UPDATE
10
11 CBOP
12 C !ROUTINE: SEAICE_JFNK
13 C !INTERFACE:
14 SUBROUTINE SEAICE_JFNK( myTime, myIter, myThid )
15
16 C !DESCRIPTION: \bv
17 C *==========================================================*
18 C | SUBROUTINE SEAICE_JFNK
19 C | o Ice dynamics using a Jacobian-free Newton-Krylov solver
20 C | following J.-F. Lemieux et al. Improving the numerical
21 C | convergence of viscous-plastic sea ice models with the
22 C | Jacobian-free Newton-Krylov method. J. Comp. Phys. 229,
23 C | 2840-2852 (2010).
24 C | o The logic follows JFs code.
25 C *==========================================================*
26 C | written by Martin Losch, Oct 2012
27 C *==========================================================*
28 C \ev
29
30 C !USES:
31 IMPLICIT NONE
32
33 C === Global variables ===
34 #include "SIZE.h"
35 #include "EEPARAMS.h"
36 #include "PARAMS.h"
37 #include "DYNVARS.h"
38 #include "GRID.h"
39 #include "SEAICE_SIZE.h"
40 #include "SEAICE_PARAMS.h"
41 #include "SEAICE.h"
42
43 #ifdef ALLOW_AUTODIFF_TAMC
44 # include "tamc.h"
45 #endif
46
47 C !INPUT/OUTPUT PARAMETERS:
48 C === Routine arguments ===
49 C myTime :: Simulation time
50 C myIter :: Simulation timestep number
51 C myThid :: my Thread Id. number
52 _RL myTime
53 INTEGER myIter
54 INTEGER myThid
55
56 #ifdef SEAICE_ALLOW_JFNK
57 C !FUNCTIONS:
58 LOGICAL DIFFERENT_MULTIPLE
59 EXTERNAL DIFFERENT_MULTIPLE
60
61 C !LOCAL VARIABLES:
62 C === Local variables ===
63 C i,j,bi,bj :: loop indices
64 INTEGER i,j,bi,bj
65 C loop indices
66 INTEGER newtonIter
67 INTEGER krylovIter, krylovFails
68 INTEGER totalKrylovItersLoc, totalNewtonItersLoc
69 C FGMRES flag that determines amount of output messages of fgmres
70 INTEGER iOutFGMRES
71 C FGMRES flag that indicates what fgmres wants us to do next
72 INTEGER iCode
73 _RL JFNKresidual
74 _RL JFNKresidualKm1
75 C parameters to compute convergence criterion
76 _RL phi_e, alp_e, JFNKgamma_lin
77 _RL FGMRESeps
78 _RL JFNKtol
79
80 _RL recip_deltaT
81 LOGICAL JFNKconverged, krylovConverged
82 LOGICAL writeNow
83 CHARACTER*(MAX_LEN_MBUF) msgBuf
84
85 C u/vIceRes :: residual of sea-ice momentum equations
86 _RL uIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
87 _RL vIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
88 C du/vIce :: ice velocity increment to be added to u/vIce
89 _RL duIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
90 _RL dvIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
91 C precomputed (= constant per Newton iteration) versions of
92 C zeta, eta, and DWATN, press
93 _RL zetaPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
94 _RL etaPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
95 _RL etaZPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
96 _RL dwatPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
97 CEOP
98
99 C Initialise
100 newtonIter = 0
101 krylovFails = 0
102 totalKrylovItersLoc = 0
103 JFNKconverged = .FALSE.
104 JFNKtol = 0. _d 0
105 JFNKresidual = 0. _d 0
106 JFNKresidualKm1 = 0. _d 0
107 FGMRESeps = 0. _d 0
108 recip_deltaT = 1. _d 0 / SEAICE_deltaTdyn
109
110 iOutFGMRES=0
111 C with iOutFgmres=1, seaice_fgmres prints the residual at each iteration
112 IF ( debugLevel.GE.debLevC .AND.
113 & DIFFERENT_MULTIPLE( SEAICE_monFreq, myTime, deltaTClock ) )
114 & iOutFGMRES=1
115
116 DO bj=myByLo(myThid),myByHi(myThid)
117 DO bi=myBxLo(myThid),myBxHi(myThid)
118 DO J=1-OLy,sNy+OLy
119 DO I=1-OLx,sNx+OLx
120 uIceRes(I,J,bi,bj) = 0. _d 0
121 vIceRes(I,J,bi,bj) = 0. _d 0
122 duIce (I,J,bi,bj) = 0. _d 0
123 dvIce (I,J,bi,bj) = 0. _d 0
124 uIceNm1(I,J,bi,bj) = uIce(I,J,bi,bj)
125 vIceNm1(I,J,bi,bj) = vIce(I,J,bi,bj)
126 ENDDO
127 ENDDO
128 C Compute things that do no change during the Newton iteration:
129 C sea-surface tilt and wind stress:
130 C FORCEX/Y0 - mass*(u/vIceNm1)/deltaT
131 DO J=1-OLy,sNy+OLy
132 DO I=1-OLx,sNx+OLx
133 FORCEX(I,J,bi,bj) = FORCEX0(I,J,bi,bj)
134 & + seaiceMassU(I,J,bi,bj)*uIceNm1(I,J,bi,bj)*recip_deltaT
135 FORCEY(I,J,bi,bj) = FORCEY0(I,J,bi,bj)
136 & + seaiceMassV(I,J,bi,bj)*vIceNm1(I,J,bi,bj)*recip_deltaT
137 ENDDO
138 ENDDO
139 ENDDO
140 ENDDO
141 C Start nonlinear Newton iteration: outer loop iteration
142 DO WHILE ( newtonIter.LT.SEAICEnewtonIterMax .AND.
143 & .NOT.JFNKconverged )
144 newtonIter = newtonIter + 1
145 C Compute initial residual F(u), (includes computation of global
146 C variables DWATN, zeta, and eta)
147 IF ( newtonIter .EQ. 1 ) CALL SEAICE_JFNK_UPDATE(
148 I duIce, dvIce,
149 U uIce, vIce, JFNKresidual,
150 O uIceRes, vIceRes,
151 I newtonIter, myTime, myIter, myThid )
152 C local copies of precomputed coefficients that are to stay
153 C constant for the preconditioner
154 DO bj=myByLo(myThid),myByHi(myThid)
155 DO bi=myBxLo(myThid),myBxHi(myThid)
156 DO j=1-OLy,sNy+OLy
157 DO i=1-OLx,sNx+OLx
158 zetaPre(I,J,bi,bj) = zeta(I,J,bi,bj)
159 etaPre(I,J,bi,bj) = eta(I,J,bi,bj)
160 etaZPre(I,J,bi,bj) = etaZ(I,J,bi,bj)
161 dwatPre(I,J,bi,bj) = DWATN(I,J,bi,bj)
162 ENDDO
163 ENDDO
164 ENDDO
165 ENDDO
166 C compute convergence criterion for linear preconditioned FGMRES
167 JFNKgamma_lin = JFNKgamma_lin_max
168 IF ( newtonIter.GT.1.AND.newtonIter.LE.SEAICE_JFNK_tolIter
169 & .AND.JFNKresidual.LT.JFNKres_t ) THEN
170 C Eisenstat, 1996, equ.(2.6)
171 phi_e = 1. _d 0
172 alp_e = 1. _d 0
173 JFNKgamma_lin = phi_e*( JFNKresidual/JFNKresidualKm1 )**alp_e
174 JFNKgamma_lin = min(JFNKgamma_lin_max, JFNKgamma_lin)
175 JFNKgamma_lin = max(JFNKgamma_lin_min, JFNKgamma_lin)
176 ENDIF
177 C save the residual for the next iteration
178 JFNKresidualKm1 = JFNKresidual
179
180 C The Krylov iteration using FGMRES, the preconditioner is LSOR
181 C for now. The code is adapted from SEAICE_LSR, but heavily stripped
182 C down.
183 C krylovIter is mapped into "its" in seaice_fgmres and is incremented
184 C in that routine
185 krylovIter = 0
186 iCode = 0
187
188 JFNKconverged = JFNKresidual.LT.JFNKtol
189
190 C do Krylov loop only if convergence is not reached
191
192 IF ( .NOT.JFNKconverged ) THEN
193
194 C start Krylov iteration (FGMRES)
195
196 krylovConverged = .FALSE.
197 FGMRESeps = JFNKgamma_lin * JFNKresidual
198 DO WHILE ( .NOT.krylovConverged )
199 C solution vector sol = du/vIce
200 C residual vector (rhs) Fu = u/vIceRes
201 C output work vectors wk1, -> input work vector wk2
202
203 CALL SEAICE_FGMRES_DRIVER(
204 I uIceRes, vIceRes,
205 U duIce, dvIce, iCode,
206 I FGMRESeps, iOutFGMRES,
207 I newtonIter, krylovIter, myTime, myIter, myThid )
208 C FGMRES returns iCode either asking for an new preconditioned vector
209 C or product of matrix (Jacobian) times vector. For iCode = 0, terminate
210 C iteration
211 IF (iCode.EQ.1) THEN
212 C Call preconditioner
213 IF ( SOLV_MAX_ITERS .GT. 0 )
214 & CALL SEAICE_PRECONDITIONER(
215 U duIce, dvIce,
216 I zetaPre, etaPre, etaZpre, dwatPre,
217 I newtonIter, krylovIter, myTime, myIter, myThid )
218 ELSEIF (iCode.GE.2) THEN
219 C Compute Jacobian times vector
220 CALL SEAICE_JACVEC(
221 I uIce, vIce, uIceRes, vIceRes,
222 U duIce, dvIce,
223 I newtonIter, krylovIter, myTime, myIter, myThid )
224 ENDIF
225 krylovConverged = iCode.EQ.0
226 C End of Krylov iterate
227 ENDDO
228 totalKrylovItersLoc = totalKrylovItersLoc + krylovIter
229 C some output diagnostics
230 IF ( debugLevel.GE.debLevA ) THEN
231 _BEGIN_MASTER( myThid )
232 totalNewtonItersLoc =
233 & SEAICEnewtonIterMax*(myIter-nIter0)+newtonIter
234 WRITE(msgBuf,'(2A,2(1XI6),2E12.5)')
235 & ' S/R SEAICE_JFNK: Newton iterate / total, ',
236 & 'JFNKgamma_lin, initial norm = ',
237 & newtonIter, totalNewtonItersLoc,
238 & JFNKgamma_lin,JFNKresidual
239 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
240 & SQUEEZE_RIGHT, myThid )
241 WRITE(msgBuf,'(3(A,I6))')
242 & ' S/R SEAICE_JFNK: Newton iterate / total = ',newtonIter,
243 & ' / ', totalNewtonItersLoc,
244 & ', Nb. of FGMRES iterations = ', krylovIter
245 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
246 & SQUEEZE_RIGHT, myThid )
247 _END_MASTER( myThid )
248 ENDIF
249 IF ( krylovIter.EQ.SEAICEkrylovIterMax ) THEN
250 krylovFails = krylovFails + 1
251 ENDIF
252 C Set the stopping criterion for the Newton iteration and the
253 C criterion for the transition from accurate to approximate FGMRES
254 IF ( newtonIter .EQ. 1 ) THEN
255 JFNKtol=JFNKgamma_nonlin*JFNKresidual
256 IF ( JFNKres_tFac .NE. UNSET_RL )
257 & JFNKres_t = JFNKresidual * JFNKres_tFac
258 ENDIF
259 C Update linear solution vector and return to Newton iteration
260 C Do a linesearch if necessary, and compute a new residual.
261 C Note that it should be possible to do the following operations
262 C at the beginning of the Newton iteration, thereby saving us from
263 C the extra call of seaice_jfnk_update, but unfortunately that
264 C changes the results, so we leave the stuff here for now.
265 CALL SEAICE_JFNK_UPDATE(
266 I duIce, dvIce,
267 U uIce, vIce, JFNKresidual,
268 O uIceRes, vIceRes,
269 I newtonIter, myTime, myIter, myThid )
270 C reset du/vIce here instead of setting sol = 0 in seaice_fgmres_driver
271 DO bj=myByLo(myThid),myByHi(myThid)
272 DO bi=myBxLo(myThid),myBxHi(myThid)
273 DO J=1-OLy,sNy+OLy
274 DO I=1-OLx,sNx+OLx
275 duIce(I,J,bi,bj)= 0. _d 0
276 dvIce(I,J,bi,bj)= 0. _d 0
277 ENDDO
278 ENDDO
279 ENDDO
280 ENDDO
281 ENDIF
282 C end of Newton iterate
283 ENDDO
284
285 C-- Output diagnostics
286
287 IF ( SEAICE_monFreq .GT. 0. _d 0 ) THEN
288 C Count iterations
289 totalJFNKtimeSteps = totalJFNKtimeSteps + 1
290 totalNewtonIters = totalNewtonIters + newtonIter
291 totalKrylovIters = totalKrylovIters + totalKrylovItersLoc
292 C Record failure
293 totalKrylovFails = totalKrylovFails + krylovFails
294 IF ( newtonIter .EQ. SEAICEnewtonIterMax ) THEN
295 totalNewtonFails = totalNewtonFails + 1
296 ENDIF
297 ENDIF
298 C Decide whether it is time to dump and reset the counter
299 writeNow = DIFFERENT_MULTIPLE(SEAICE_monFreq,
300 & myTime+deltaTClock, deltaTClock)
301 #ifdef ALLOW_CAL
302 IF ( useCAL ) THEN
303 CALL CAL_TIME2DUMP(
304 I zeroRL, SEAICE_monFreq, deltaTClock,
305 U writeNow,
306 I myTime+deltaTclock, myIter+1, myThid )
307 ENDIF
308 #endif
309 IF ( writeNow ) THEN
310 _BEGIN_MASTER( myThid )
311 WRITE(msgBuf,'(A)')
312 &' // ======================================================='
313 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
314 & SQUEEZE_RIGHT, myThid )
315 WRITE(msgBuf,'(A)') ' // Begin JFNK statistics'
316 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
317 & SQUEEZE_RIGHT, myThid )
318 WRITE(msgBuf,'(A)')
319 &' // ======================================================='
320 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
321 & SQUEEZE_RIGHT, myThid )
322 WRITE(msgBuf,'(A,I10)')
323 & ' %JFNK_MON: time step = ', myIter+1
324 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
325 & SQUEEZE_RIGHT, myThid )
326 WRITE(msgBuf,'(A,I10)')
327 & ' %JFNK_MON: Nb. of time steps = ', totalJFNKtimeSteps
328 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
329 & SQUEEZE_RIGHT, myThid )
330 WRITE(msgBuf,'(A,I10)')
331 & ' %JFNK_MON: Nb. of Newton steps = ', totalNewtonIters
332 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
333 & SQUEEZE_RIGHT, myThid )
334 WRITE(msgBuf,'(A,I10)')
335 & ' %JFNK_MON: Nb. of Krylov steps = ', totalKrylovIters
336 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
337 & SQUEEZE_RIGHT, myThid )
338 WRITE(msgBuf,'(A,I10)')
339 & ' %JFNK_MON: Nb. of Newton failures = ', totalNewtonFails
340 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
341 & SQUEEZE_RIGHT, myThid )
342 WRITE(msgBuf,'(A,I10)')
343 & ' %JFNK_MON: Nb. of Krylov failures = ', totalKrylovFails
344 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
345 & SQUEEZE_RIGHT, myThid )
346 WRITE(msgBuf,'(A)')
347 &' // ======================================================='
348 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
349 & SQUEEZE_RIGHT, myThid )
350 WRITE(msgBuf,'(A)') ' // End JFNK statistics'
351 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
352 & SQUEEZE_RIGHT, myThid )
353 WRITE(msgBuf,'(A)')
354 &' // ======================================================='
355 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
356 & SQUEEZE_RIGHT, myThid )
357 _END_MASTER( myThid )
358 C reset and start again
359 totalJFNKtimeSteps = 0
360 totalNewtonIters = 0
361 totalKrylovIters = 0
362 totalKrylovFails = 0
363 totalNewtonFails = 0
364 ENDIF
365
366 C Print more debugging information
367 IF ( debugLevel.GE.debLevA ) THEN
368 IF ( newtonIter .EQ. SEAICEnewtonIterMax ) THEN
369 _BEGIN_MASTER( myThid )
370 WRITE(msgBuf,'(A,I10)')
371 & ' S/R SEAICE_JFNK: JFNK did not converge in timestep ',
372 & myIter+1
373 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
374 & SQUEEZE_RIGHT, myThid )
375 _END_MASTER( myThid )
376 ENDIF
377 IF ( krylovFails .GT. 0 ) THEN
378 _BEGIN_MASTER( myThid )
379 WRITE(msgBuf,'(A,I4,A,I10)')
380 & ' S/R SEAICE_JFNK: FGMRES did not converge ',
381 & krylovFails, ' times in timestep ', myIter+1
382 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
383 & SQUEEZE_RIGHT, myThid )
384 _END_MASTER( myThid )
385 ENDIF
386 _BEGIN_MASTER( myThid )
387 WRITE(msgBuf,'(A,I6,A,I10)')
388 & ' S/R SEAICE_JFNK: Total number FGMRES iterations = ',
389 & totalKrylovItersLoc, ' in timestep ', myIter+1
390 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
391 & SQUEEZE_RIGHT, myThid )
392 _END_MASTER( myThid )
393 ENDIF
394
395 RETURN
396 END
397
398 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
399 CBOP
400 C !ROUTINE: SEAICE_JFNK_UPDATE
401 C !INTERFACE:
402
403 SUBROUTINE SEAICE_JFNK_UPDATE(
404 I duIce, dvIce,
405 U uIce, vIce, JFNKresidual,
406 O uIceRes, vIceRes,
407 I newtonIter, myTime, myIter, myThid )
408
409 C !DESCRIPTION: \bv
410 C *==========================================================*
411 C | SUBROUTINE SEAICE_JFNK_UPDATE
412 C | o Update velocities with incremental solutions of FGMRES
413 C | o compute residual of updated solutions and do
414 C | o linesearch:
415 C | reduce update until residual is smaller than previous
416 C | one (input)
417 C *==========================================================*
418 C | written by Martin Losch, Jan 2013
419 C *==========================================================*
420 C \ev
421
422 C !USES:
423 IMPLICIT NONE
424
425 C === Global variables ===
426 #include "SIZE.h"
427 #include "EEPARAMS.h"
428 #include "PARAMS.h"
429 #include "SEAICE_SIZE.h"
430 #include "SEAICE_PARAMS.h"
431
432 C !INPUT/OUTPUT PARAMETERS:
433 C === Routine arguments ===
434 C myTime :: Simulation time
435 C myIter :: Simulation timestep number
436 C myThid :: my Thread Id. number
437 C newtonIter :: current iterate of Newton iteration
438 _RL myTime
439 INTEGER myIter
440 INTEGER myThid
441 INTEGER newtonIter
442 C JFNKresidual :: Residual at the beginning of the FGMRES iteration,
443 C changes with newtonIter (updated)
444 _RL JFNKresidual
445 C du/vIce :: ice velocity increment to be added to u/vIce (input)
446 _RL duIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
447 _RL dvIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
448 C u/vIce :: ice velocity increment to be added to u/vIce (updated)
449 _RL uIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
450 _RL vIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
451 C u/vIceRes :: residual of sea-ice momentum equations (output)
452 _RL uIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
453 _RL vIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
454
455 C !LOCAL VARIABLES:
456 C === Local variables ===
457 C i,j,bi,bj :: loop indices
458 INTEGER i,j,bi,bj
459 INTEGER l
460 _RL resLoc, facLS
461 LOGICAL doLineSearch
462 C nVec :: size of the input vector(s)
463 C resTmp :: vector version of the residuals
464 INTEGER nVec
465 PARAMETER ( nVec = 2*sNx*sNy )
466 _RL resTmp (nVec,1,nSx,nSy)
467
468 CHARACTER*(MAX_LEN_MBUF) msgBuf
469 CEOP
470
471 C Initialise some local variables
472 l = 0
473 resLoc = JFNKresidual
474 facLS = 1. _d 0
475 doLineSearch = .TRUE.
476 DO WHILE ( doLineSearch )
477 C Create update
478 DO bj=myByLo(myThid),myByHi(myThid)
479 DO bi=myBxLo(myThid),myBxHi(myThid)
480 DO J=1-OLy,sNy+OLy
481 DO I=1-OLx,sNx+OLx
482 uIce(I,J,bi,bj) = uIce(I,J,bi,bj)+facLS*duIce(I,J,bi,bj)
483 vIce(I,J,bi,bj) = vIce(I,J,bi,bj)+facLS*dvIce(I,J,bi,bj)
484 ENDDO
485 ENDDO
486 ENDDO
487 ENDDO
488 C Compute current residual F(u), (includes re-computation of global
489 C variables DWATN, zeta, and eta, i.e. they are different after this)
490 CALL SEAICE_CALC_RESIDUAL(
491 I uIce, vIce,
492 O uIceRes, vIceRes,
493 I newtonIter, 0, myTime, myIter, myThid )
494 C Important: Compute the norm of the residual using the same scalar
495 C product that SEAICE_FGMRES does
496 CALL SEAICE_MAP2VEC(nVec,uIceRes,vIceRes,resTmp,.TRUE.,myThid)
497 CALL SEAICE_SCALPROD(nVec,1,1,1,resTmp,resTmp,resLoc,myThid)
498 resLoc = SQRT(resLoc)
499 C Determine, if we need more iterations
500 doLineSearch = resLoc .GE. JFNKresidual
501 C Limit the maximum number of iterations arbitrarily to four
502 doLineSearch = doLineSearch .AND. l .LT. 4
503 C For the first iteration du/vIce = 0 and there will be no
504 C improvement of the residual possible, so we do only the first
505 C iteration
506 IF ( newtonIter .EQ. 1 ) doLineSearch = .FALSE.
507 C Only start a linesearch after some Newton iterations
508 IF ( newtonIter .LE. SEAICE_JFNK_lsIter ) doLineSearch = .FALSE.
509 C Increment counter
510 l = l + 1
511 C some output diagnostics
512 IF ( debugLevel.GE.debLevA .AND. doLineSearch ) THEN
513 _BEGIN_MASTER( myThid )
514 WRITE(msgBuf,'(2A,2(1XI6),3E12.5)')
515 & ' S/R SEAICE_JFNK_UPDATE: Newton iter, LSiter, ',
516 & 'facLS, JFNKresidual, resLoc = ',
517 & newtonIter, l, facLS, JFNKresidual, resLoc
518 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
519 & SQUEEZE_RIGHT, myThid )
520 _END_MASTER( myThid )
521 ENDIF
522 C Get ready for the next iteration: after adding du/vIce in the first
523 C iteration, we substract 0.5*du/vIce from u/vIce in the next
524 C iterations, 0.25*du/vIce in the second, etc.
525 facLS = - 0.5 _d 0 * ABS(facLS)
526 ENDDO
527 C This is the new residual
528 JFNKresidual = resLoc
529
530 #endif /* SEAICE_ALLOW_JFNK */
531
532 RETURN
533 END

  ViewVC Help
Powered by ViewVC 1.1.22