/[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.25 - (show annotations) (download)
Fri Feb 7 14:27:21 2014 UTC (11 years, 5 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint64u
Changes since 1.24: +2 -2 lines
fix a comment

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

  ViewVC Help
Powered by ViewVC 1.1.22