/[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.9 - (show annotations) (download)
Mon Nov 12 09:46:38 2012 UTC (12 years, 8 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint64a
Changes since 1.8: +13 -3 lines
add calendar dump capability for the JFNK counters

1 C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_jfnk.F,v 1.8 2012/11/09 14:11:43 mlosch Exp $
2 C $Name: $
3
4 #include "SEAICE_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: SEAICE_JFNK
8 C !INTERFACE:
9 SUBROUTINE SEAICE_JFNK( myTime, myIter, myThid )
10
11 C !DESCRIPTION: \bv
12 C *==========================================================*
13 C | SUBROUTINE SEAICE_JFKF
14 C | o Ice dynamics using a Jacobian-free Newton-Krylov solver
15 C | following J.-F. Lemieux et al. Improving the numerical
16 C | convergence of viscous-plastic sea ice models with the
17 C | Jacobian-free Newton-Krylov method. J. Comp. Phys. 229,
18 C | 2840-2852 (2010).
19 C | o The logic follows JFs code.
20 C *==========================================================*
21 C | written by Martin Losch, Oct 2012
22 C *==========================================================*
23 C \ev
24
25 C !USES:
26 IMPLICIT NONE
27
28 C === Global variables ===
29 #include "SIZE.h"
30 #include "EEPARAMS.h"
31 #include "PARAMS.h"
32 #include "DYNVARS.h"
33 #include "GRID.h"
34 #include "SEAICE_SIZE.h"
35 #include "SEAICE_PARAMS.h"
36 #include "SEAICE.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
51 #if ( (defined SEAICE_CGRID) && \
52 (defined SEAICE_ALLOW_JFNK) && \
53 (defined SEAICE_ALLOW_DYNAMICS) )
54 C !FUNCTIONS:
55 LOGICAL DIFFERENT_MULTIPLE
56 EXTERNAL DIFFERENT_MULTIPLE
57
58 C i,j,bi,bj :: loop indices
59 INTEGER i,j,bi,bj
60 C loop indices
61 INTEGER newtonIter
62 INTEGER krylovIter, krylovFails
63 INTEGER totalKrylovItersLoc
64 C FGMRES flag that determines amount of output messages of fgmres
65 INTEGER iOutFGMRES
66 C FGMRES flag that indicates what fgmres wants us to do next
67 INTEGER iCode
68 _RL JFNKresidual, JFNKresidualTile(nSx,nSy)
69 _RL JFNKresidualKm1
70 C parameters to compute convergence criterion
71 _RL phi_e, alp_e, JFNKgamma_lin
72 _RL FGMRESeps
73 _RL JFNKtol
74 C
75 _RL recip_deltaT
76 LOGICAL JFNKconverged, krylovConverged
77 LOGICAL writeNow
78 CHARACTER*(MAX_LEN_MBUF) msgBuf
79 C
80 C u/vIceRes :: residual of sea-ice momentum equations
81 _RL uIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
82 _RL vIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
83 C du/vIce :: ice velocity increment to be added to u/vIce
84 _RL duIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
85 _RL dvIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
86 C precomputed (= constant per Newton iteration) versions of
87 C zeta, eta, and DWATN, press
88 _RL zetaPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
89 _RL etaPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
90 _RL etaZPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
91 _RL dwatPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
92 _RL pressPre(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
93 CEOP
94
95 C Initialise
96 newtonIter = 0
97 krylovFails = 0
98 totalKrylovItersLoc = 0
99 JFNKconverged = .FALSE.
100 JFNKtol = 0. _d 0
101 JFNKresidual = 0. _d 0
102 JFNKresidualKm1 = 0. _d 0
103 FGMRESeps = 0. _d 0
104 recip_deltaT = 1. _d 0 / SEAICE_deltaTdyn
105
106 iOutFGMRES=0
107 C iOutFgmres=1 gives a little bit of output
108 IF ( debugLevel.GE.debLevA .AND.
109 & DIFFERENT_MULTIPLE( SEAICE_monFreq, myTime, deltaTClock ) )
110 & iOutFGMRES=1
111
112 C
113 DO bj=myByLo(myThid),myByHi(myThid)
114 DO bi=myBxLo(myThid),myBxHi(myThid)
115 DO J=1-Oly,sNy+Oly
116 DO I=1-Olx,sNx+Olx
117 uIceRes(I,J,bi,bj) = 0. _d 0
118 vIceRes(I,J,bi,bj) = 0. _d 0
119 duIce (I,J,bi,bj) = 0. _d 0
120 dvIce (I,J,bi,bj) = 0. _d 0
121 uIceNm1(I,J,bi,bj) = uIce(I,J,bi,bj)
122 vIceNm1(I,J,bi,bj) = vIce(I,J,bi,bj)
123 ENDDO
124 ENDDO
125 C Compute things that do no change during the Newton iteration:
126 C sea-surface tilt and wind stress:
127 C FORCEX/Y0 - mass*(u/vIceNm1)/deltaT
128 DO J=1-Oly,sNy+Oly
129 DO I=1-Olx,sNx+Olx
130 FORCEX(I,J,bi,bj) = FORCEX0(I,J,bi,bj)
131 & + seaiceMassU(I,J,bi,bj)*uIceNm1(I,J,bi,bj)*recip_deltaT
132 FORCEY(I,J,bi,bj) = FORCEY0(I,J,bi,bj)
133 & + seaiceMassV(I,J,bi,bj)*vIceNm1(I,J,bi,bj)*recip_deltaT
134 ENDDO
135 ENDDO
136 ENDDO
137 ENDDO
138 C Start nonlinear Newton iteration: outer loop iteration
139 DO WHILE ( newtonIter.LT.SEAICEnewtonIterMax .AND.
140 & .NOT.JFNKconverged )
141 newtonIter = newtonIter + 1
142 C Compute initial residual F(u), (includes computation of global
143 C variables DWATN, zeta, and eta)
144 CALL SEAICE_CALC_RESIDUAL(
145 I uIce, vIce,
146 O uIceRes, vIceRes,
147 I newtonIter, 0, myTime, myIter, myThid )
148 CALL EXCH_UV_XY_RL( uIceRes, vIceRes,.TRUE.,myThid)
149 C local copies of precomputed coefficients that are to stay
150 C constant for the preconditioner
151 DO bj=myByLo(myThid),myByHi(myThid)
152 DO bi=myBxLo(myThid),myBxHi(myThid)
153 DO j=1-Oly,sNy+Oly
154 DO i=1-Olx,sNx+Olx
155 zetaPre(I,J,bi,bj) = zeta(I,J,bi,bj)
156 etaPre(I,J,bi,bj) = eta(I,J,bi,bj)
157 etaZPre(I,J,bi,bj) = etaZ(I,J,bi,bj)
158 dwatPre(I,J,bi,bj) = DWATN(I,J,bi,bj)
159 pressPre(I,J,bi,bj) = press(I,J,bi,bj)
160 ENDDO
161 ENDDO
162 ENDDO
163 ENDDO
164 C
165 DO bj=myByLo(myThid),myByHi(myThid)
166 DO bi=myBxLo(myThid),myBxHi(myThid)
167 JFNKresidualTile(bi,bj) = 0. _d 0
168 DO J=1,sNy
169 DO I=1,sNx
170 #ifdef CG2D_SINGLECPU_SUM
171 JFNKlocalBuf(I,J,bi,bj) =
172 #else
173 JFNKresidualTile(bi,bj) = JFNKresidualTile(bi,bj) +
174 #endif
175 & uIceRes(I,J,bi,bj)*uIceRes(I,J,bi,bj) +
176 & vIceRes(I,J,bi,bj)*vIceRes(I,J,bi,bj)
177 ENDDO
178 ENDDO
179 ENDDO
180 ENDDO
181 JFNKresidual = 0. _d 0
182 #ifdef CG2D_SINGLECPU_SUM
183 CALL GLOBAL_SUM_SINGLECPU_RL(
184 & JFNKlocalBuf,JFNKresidual, 0, 0, myThid)
185 #else
186 CALL GLOBAL_SUM_TILE_RL( JFNKresidualTile,JFNKresidual,myThid )
187 #endif
188 JFNKresidual = SQRT(JFNKresidual)
189 C compute convergence criterion for linear preconditioned FGMRES
190 JFNKgamma_lin = JFNKgamma_lin_max
191 IF ( newtonIter.GT.1.AND.newtonIter.LE.100
192 & .AND.JFNKresidual.LT.JFNKres_t ) THEN
193 C Eisenstat, 1996, equ.(2.6)
194 phi_e = 1. _d 0
195 alp_e = 1. _d 0
196 JFNKgamma_lin = phi_e*( JFNKresidual/JFNKresidualKm1 )**alp_e
197 JFNKgamma_lin = min(JFNKgamma_lin_max, JFNKgamma_lin)
198 JFNKgamma_lin = max(JFNKgamma_lin_min, JFNKgamma_lin)
199 ENDIF
200 C save the residual for the next iteration
201 JFNKresidualKm1 = JFNKresidual
202 C
203 C The Krylov iteration using FGMRES, the preconditioner is LSOR
204 C for now. The code is adapted from SEAICE_LSR, but heavily stripped
205 C down.
206 C krylovIter is mapped into "its" in seaice_fgmres and is incremented
207 C in that routine
208 krylovIter = 0
209 iCode = 0
210 IF ( debugLevel.GE.debLevA ) THEN
211 _BEGIN_MASTER( myThid )
212 WRITE(msgBuf,'(2A,2(1XI6),2E12.5)')
213 & ' S/R SEAICE_JFNK: newtonIter,',
214 & ' total newtonIter, JFNKgamma_lin, initial norm = ',
215 & newtonIter,SEAICEnewtonIterMax*(myIter-nIter0)+newtonIter,
216 & JFNKgamma_lin, JFNKresidual
217 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
218 & SQUEEZE_RIGHT, myThid )
219 _END_MASTER( myThid )
220 ENDIF
221 C
222 JFNKconverged = JFNKresidual.LT.JFNKtol
223 C
224 C do Krylov loop only if convergence is not reached
225 C
226 IF ( .NOT.JFNKconverged ) THEN
227 C
228 C start Krylov iteration (FGMRES)
229 C
230 krylovConverged = .FALSE.
231 FGMRESeps = JFNKgamma_lin * JFNKresidual
232 DO WHILE ( .NOT.krylovConverged )
233 C solution vector sol = du/vIce
234 C residual vector (rhs) Fu = u/vIceRes
235 C output work vectors wk1, -> input work vector wk2
236 C
237 CALL SEAICE_FGMRES_DRIVER(
238 I uIceRes, vIceRes,
239 U duIce, dvIce, iCode,
240 I FGMRESeps, iOutFGMRES,
241 I newtonIter, krylovIter, myTime, myIter, myThid )
242 C FGMRES returns iCode either asking for an new preconditioned vector
243 C or product of matrix (Jacobian) times vector. For iCode = 0, terminate
244 C iteration
245 IF (iCode.EQ.1) THEN
246 C Call preconditioner
247 IF ( SOLV_MAX_ITERS .GT. 0 )
248 & CALL SEAICE_PRECONDITIONER(
249 U duIce, dvIce,
250 I zetaPre, etaPre, etaZpre, dwatPre, pressPre,
251 I newtonIter, krylovIter, myTime, myIter, myThid )
252 ELSEIF (iCode.GE.2) THEN
253 C Compute Jacobian times vector
254 CALL SEAICE_JACVEC(
255 I uIce, vIce, uIceRes, vIceRes,
256 U duIce, dvIce,
257 I newtonIter, krylovIter, myTime, myIter, myThid )
258 ENDIF
259 krylovConverged = iCode.EQ.0
260 C End of Krylov iterate
261 ENDDO
262 totalKrylovItersLoc = totalKrylovItersLoc + krylovIter
263 C some output diagnostics
264 IF ( debugLevel.GE.debLevA ) THEN
265 _BEGIN_MASTER( myThid )
266 WRITE(msgBuf,'(3(A,I6))')
267 & ' S/R SEAICE_JFNK: Newton iterate / total = ', newtonIter,
268 & ' / ', SEAICEnewtonIterMax*(myIter-nIter0)+newtonIter,
269 & ', Nb. of FGMRES iterations = ', krylovIter
270 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
271 & SQUEEZE_RIGHT, myThid )
272 _END_MASTER( myThid )
273 ENDIF
274 IF ( krylovIter.EQ.SEAICEkrylovIterMax ) THEN
275 krylovFails = krylovFails + 1
276 ENDIF
277 C Update linear solution vector and return to Newton iteration
278 DO bj=myByLo(myThid),myByHi(myThid)
279 DO bi=myBxLo(myThid),myBxHi(myThid)
280 DO J=1-Oly,sNy+Oly
281 DO I=1-Olx,sNx+Olx
282 uIce(I,J,bi,bj) = uIce(I,J,bi,bj)+duIce(I,J,bi,bj)
283 vIce(I,J,bi,bj) = vIce(I,J,bi,bj)+dvIce(I,J,bi,bj)
284 C reset du/vIce here instead of setting sol = 0 in seaice_fgmres_driver
285 duIce(I,J,bi,bj)= 0. _d 0
286 dvIce(I,J,bi,bj)= 0. _d 0
287 ENDDO
288 ENDDO
289 ENDDO
290 ENDDO
291 C Set the stopping criterion for the Newton iteration
292 IF ( newtonIter .EQ. 1 ) JFNKtol=JFNKgamma_nonlin*JFNKresidual
293 ENDIF
294 C end of Newton iterate
295 ENDDO
296 C
297 C-- Output diagnostics
298 C
299 IF ( SEAICE_monFreq .GT. 0. _d 0 ) THEN
300 C Count iterations
301 totalJFNKtimeSteps = totalJFNKtimeSteps + 1
302 totalNewtonIters = totalNewtonIters + newtonIter
303 totalKrylovIters = totalKrylovIters + totalKrylovItersLoc
304 C Record failure
305 totalKrylovFails = totalKrylovFails + krylovFails
306 IF ( newtonIter .EQ. SEAICEnewtonIterMax ) THEN
307 totalNewtonFails = totalNewtonFails + 1
308 ENDIF
309 ENDIF
310 C Decide whether it is time to dump and reset the counter
311 writeNow = DIFFERENT_MULTIPLE(SEAICE_monFreq,
312 & myTime+deltaTClock, deltaTClock)
313 #ifdef ALLOW_CAL
314 IF ( useCAL ) THEN
315 CALL CAL_TIME2DUMP(
316 I zeroRL, SEAICE_monFreq, deltaTClock,
317 U writeNow,
318 I myTime+deltaTclock, myIter+1, myThid )
319 ENDIF
320 #endif
321 IF ( writeNow ) THEN
322 _BEGIN_MASTER( myThid )
323 WRITE(msgBuf,'(A)')
324 &' // ======================================================='
325 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
326 & SQUEEZE_RIGHT, myThid )
327 WRITE(msgBuf,'(A)') ' // Begin JFNK statistics'
328 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
329 & SQUEEZE_RIGHT, myThid )
330 WRITE(msgBuf,'(A)')
331 &' // ======================================================='
332 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
333 & SQUEEZE_RIGHT, myThid )
334 WRITE(msgBuf,'(A,I10)')
335 & ' %JFNK_MON: time step = ', myIter+1
336 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
337 & SQUEEZE_RIGHT, myThid )
338 WRITE(msgBuf,'(A,I10)')
339 & ' %JFNK_MON: Nb. of time steps = ', totalJFNKtimeSteps
340 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
341 & SQUEEZE_RIGHT, myThid )
342 WRITE(msgBuf,'(A,I10)')
343 & ' %JFNK_MON: Nb. of Newton steps = ', totalNewtonIters
344 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
345 & SQUEEZE_RIGHT, myThid )
346 WRITE(msgBuf,'(A,I10)')
347 & ' %JFNK_MON: Nb. of Krylov steps = ', totalKrylovIters
348 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
349 & SQUEEZE_RIGHT, myThid )
350 WRITE(msgBuf,'(A,I10)')
351 & ' %JFNK_MON: Nb. of Newton failures = ', totalNewtonFails
352 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
353 & SQUEEZE_RIGHT, myThid )
354 WRITE(msgBuf,'(A,I10)')
355 & ' %JFNK_MON: Nb. of Krylov failures = ', totalKrylovFails
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 WRITE(msgBuf,'(A)') ' // Begin JFNK statistics'
363 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
364 & SQUEEZE_RIGHT, myThid )
365 WRITE(msgBuf,'(A)')
366 &' // ======================================================='
367 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
368 & SQUEEZE_RIGHT, myThid )
369 _END_MASTER( myThid )
370 C reset and start again
371 totalJFNKtimeSteps = 0
372 totalNewtonIters = 0
373 totalKrylovIters = 0
374 totalKrylovFails = 0
375 totalNewtonFails = 0
376 ENDIF
377
378 C Print more debugging information
379 IF ( debugLevel.GE.debLevA ) THEN
380 IF ( newtonIter .EQ. SEAICEnewtonIterMax ) THEN
381 _BEGIN_MASTER( myThid )
382 WRITE(msgBuf,'(A,I10)')
383 & ' S/R SEAICE_JFNK: JFNK did not converge in timestep ',
384 & myIter+1
385 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
386 & SQUEEZE_RIGHT, myThid )
387 _END_MASTER( myThid )
388 ENDIF
389 IF ( krylovFails .GT. 0 ) THEN
390 _BEGIN_MASTER( myThid )
391 WRITE(msgBuf,'(A,I4,A,I10)')
392 & ' S/R SEAICE_JFNK: FGMRES did not converge ',
393 & krylovFails, ' times in timestep ', myIter+1
394 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
395 & SQUEEZE_RIGHT, myThid )
396 _END_MASTER( myThid )
397 ENDIF
398 _BEGIN_MASTER( myThid )
399 WRITE(msgBuf,'(A,I6,A,I10)')
400 & ' S/R SEAICE_JFNK: Total number FGMRES iterations = ',
401 & totalKrylovItersLoc, ' in timestep ', myIter+1
402 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
403 & SQUEEZE_RIGHT, myThid )
404 _END_MASTER( myThid )
405 ENDIF
406
407 #endif /* SEAICE_ALLOW_DYNAMICS and SEAICE_CGRID and SEAICE_ALLOW_JFNK */
408
409 RETURN
410 END

  ViewVC Help
Powered by ViewVC 1.1.22