/[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.19 - (show annotations) (download)
Mon Feb 25 10:44:10 2013 UTC (12 years, 4 months ago) by mlosch
Branch: MAIN
Changes since 1.18: +13 -13 lines
fix line search in s/r seaice_jfnk_update: move the logic to after the
computation of the new residual (and now the line search actually
improves the convergence)

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

  ViewVC Help
Powered by ViewVC 1.1.22