/[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.29 - (show annotations) (download)
Wed Jan 27 14:03:34 2016 UTC (9 years, 5 months ago) by mlosch
Branch: MAIN
Changes since 1.28: +38 -8 lines
put content of S/R SEAICE_FGMRES_DRIVER into S/R SEAICE_JFKN

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

  ViewVC Help
Powered by ViewVC 1.1.22