/[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.11 - (show annotations) (download)
Mon Dec 3 15:49:17 2012 UTC (12 years, 7 months ago) by mlosch
Branch: MAIN
Changes since 1.10: +2 -2 lines
fix output diagnostic's text

1 C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_jfnk.F,v 1.10 2012/11/26 08:04:50 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 CEOP
93
94 C Initialise
95 newtonIter = 0
96 krylovFails = 0
97 totalKrylovItersLoc = 0
98 JFNKconverged = .FALSE.
99 JFNKtol = 0. _d 0
100 JFNKresidual = 0. _d 0
101 JFNKresidualKm1 = 0. _d 0
102 FGMRESeps = 0. _d 0
103 recip_deltaT = 1. _d 0 / SEAICE_deltaTdyn
104
105 iOutFGMRES=0
106 C iOutFgmres=1 gives a little bit of output
107 IF ( debugLevel.GE.debLevA .AND.
108 & DIFFERENT_MULTIPLE( SEAICE_monFreq, myTime, deltaTClock ) )
109 & iOutFGMRES=1
110
111 C
112 DO bj=myByLo(myThid),myByHi(myThid)
113 DO bi=myBxLo(myThid),myBxHi(myThid)
114 DO J=1-Oly,sNy+Oly
115 DO I=1-Olx,sNx+Olx
116 uIceRes(I,J,bi,bj) = 0. _d 0
117 vIceRes(I,J,bi,bj) = 0. _d 0
118 duIce (I,J,bi,bj) = 0. _d 0
119 dvIce (I,J,bi,bj) = 0. _d 0
120 uIceNm1(I,J,bi,bj) = uIce(I,J,bi,bj)
121 vIceNm1(I,J,bi,bj) = vIce(I,J,bi,bj)
122 ENDDO
123 ENDDO
124 C Compute things that do no change during the Newton iteration:
125 C sea-surface tilt and wind stress:
126 C FORCEX/Y0 - mass*(u/vIceNm1)/deltaT
127 DO J=1-Oly,sNy+Oly
128 DO I=1-Olx,sNx+Olx
129 FORCEX(I,J,bi,bj) = FORCEX0(I,J,bi,bj)
130 & + seaiceMassU(I,J,bi,bj)*uIceNm1(I,J,bi,bj)*recip_deltaT
131 FORCEY(I,J,bi,bj) = FORCEY0(I,J,bi,bj)
132 & + seaiceMassV(I,J,bi,bj)*vIceNm1(I,J,bi,bj)*recip_deltaT
133 ENDDO
134 ENDDO
135 ENDDO
136 ENDDO
137 C Start nonlinear Newton iteration: outer loop iteration
138 DO WHILE ( newtonIter.LT.SEAICEnewtonIterMax .AND.
139 & .NOT.JFNKconverged )
140 newtonIter = newtonIter + 1
141 C Compute initial residual F(u), (includes computation of global
142 C variables DWATN, zeta, and eta)
143 CALL SEAICE_CALC_RESIDUAL(
144 I uIce, vIce,
145 O uIceRes, vIceRes,
146 I newtonIter, 0, myTime, myIter, myThid )
147 CALL EXCH_UV_XY_RL( uIceRes, vIceRes,.TRUE.,myThid)
148 C local copies of precomputed coefficients that are to stay
149 C constant for the preconditioner
150 DO bj=myByLo(myThid),myByHi(myThid)
151 DO bi=myBxLo(myThid),myBxHi(myThid)
152 DO j=1-Oly,sNy+Oly
153 DO i=1-Olx,sNx+Olx
154 zetaPre(I,J,bi,bj) = zeta(I,J,bi,bj)
155 etaPre(I,J,bi,bj) = eta(I,J,bi,bj)
156 etaZPre(I,J,bi,bj) = etaZ(I,J,bi,bj)
157 dwatPre(I,J,bi,bj) = DWATN(I,J,bi,bj)
158 ENDDO
159 ENDDO
160 ENDDO
161 ENDDO
162 C
163 DO bj=myByLo(myThid),myByHi(myThid)
164 DO bi=myBxLo(myThid),myBxHi(myThid)
165 JFNKresidualTile(bi,bj) = 0. _d 0
166 DO J=1,sNy
167 DO I=1,sNx
168 #ifdef CG2D_SINGLECPU_SUM
169 JFNKlocalBuf(I,J,bi,bj) =
170 #else
171 JFNKresidualTile(bi,bj) = JFNKresidualTile(bi,bj) +
172 #endif
173 & uIceRes(I,J,bi,bj)*uIceRes(I,J,bi,bj) +
174 & vIceRes(I,J,bi,bj)*vIceRes(I,J,bi,bj)
175 ENDDO
176 ENDDO
177 ENDDO
178 ENDDO
179 JFNKresidual = 0. _d 0
180 #ifdef CG2D_SINGLECPU_SUM
181 CALL GLOBAL_SUM_SINGLECPU_RL(
182 & JFNKlocalBuf,JFNKresidual, 0, 0, myThid)
183 #else
184 CALL GLOBAL_SUM_TILE_RL( JFNKresidualTile,JFNKresidual,myThid )
185 #endif
186 JFNKresidual = SQRT(JFNKresidual)
187 C compute convergence criterion for linear preconditioned FGMRES
188 JFNKgamma_lin = JFNKgamma_lin_max
189 IF ( newtonIter.GT.1.AND.newtonIter.LE.100
190 & .AND.JFNKresidual.LT.JFNKres_t ) THEN
191 C Eisenstat, 1996, equ.(2.6)
192 phi_e = 1. _d 0
193 alp_e = 1. _d 0
194 JFNKgamma_lin = phi_e*( JFNKresidual/JFNKresidualKm1 )**alp_e
195 JFNKgamma_lin = min(JFNKgamma_lin_max, JFNKgamma_lin)
196 JFNKgamma_lin = max(JFNKgamma_lin_min, JFNKgamma_lin)
197 ENDIF
198 C save the residual for the next iteration
199 JFNKresidualKm1 = JFNKresidual
200 C
201 C The Krylov iteration using FGMRES, the preconditioner is LSOR
202 C for now. The code is adapted from SEAICE_LSR, but heavily stripped
203 C down.
204 C krylovIter is mapped into "its" in seaice_fgmres and is incremented
205 C in that routine
206 krylovIter = 0
207 iCode = 0
208 IF ( debugLevel.GE.debLevA ) THEN
209 _BEGIN_MASTER( myThid )
210 WRITE(msgBuf,'(2A,2(1XI6),2E12.5)')
211 & ' S/R SEAICE_JFNK: newtonIter,',
212 & ' total newtonIter, JFNKgamma_lin, initial norm = ',
213 & newtonIter,SEAICEnewtonIterMax*(myIter-nIter0)+newtonIter,
214 & JFNKgamma_lin, JFNKresidual
215 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
216 & SQUEEZE_RIGHT, myThid )
217 _END_MASTER( myThid )
218 ENDIF
219 C
220 JFNKconverged = JFNKresidual.LT.JFNKtol
221 C
222 C do Krylov loop only if convergence is not reached
223 C
224 IF ( .NOT.JFNKconverged ) THEN
225 C
226 C start Krylov iteration (FGMRES)
227 C
228 krylovConverged = .FALSE.
229 FGMRESeps = JFNKgamma_lin * JFNKresidual
230 DO WHILE ( .NOT.krylovConverged )
231 C solution vector sol = du/vIce
232 C residual vector (rhs) Fu = u/vIceRes
233 C output work vectors wk1, -> input work vector wk2
234 C
235 CALL SEAICE_FGMRES_DRIVER(
236 I uIceRes, vIceRes,
237 U duIce, dvIce, iCode,
238 I FGMRESeps, iOutFGMRES,
239 I newtonIter, krylovIter, myTime, myIter, myThid )
240 C FGMRES returns iCode either asking for an new preconditioned vector
241 C or product of matrix (Jacobian) times vector. For iCode = 0, terminate
242 C iteration
243 IF (iCode.EQ.1) THEN
244 C Call preconditioner
245 IF ( SOLV_MAX_ITERS .GT. 0 )
246 & CALL SEAICE_PRECONDITIONER(
247 U duIce, dvIce,
248 I zetaPre, etaPre, etaZpre, dwatPre,
249 I newtonIter, krylovIter, myTime, myIter, myThid )
250 ELSEIF (iCode.GE.2) THEN
251 C Compute Jacobian times vector
252 CALL SEAICE_JACVEC(
253 I uIce, vIce, uIceRes, vIceRes,
254 U duIce, dvIce,
255 I newtonIter, krylovIter, myTime, myIter, myThid )
256 ENDIF
257 krylovConverged = iCode.EQ.0
258 C End of Krylov iterate
259 ENDDO
260 totalKrylovItersLoc = totalKrylovItersLoc + krylovIter
261 C some output diagnostics
262 IF ( debugLevel.GE.debLevA ) THEN
263 _BEGIN_MASTER( myThid )
264 WRITE(msgBuf,'(3(A,I6))')
265 & ' S/R SEAICE_JFNK: Newton iterate / total = ', newtonIter,
266 & ' / ', SEAICEnewtonIterMax*(myIter-nIter0)+newtonIter,
267 & ', Nb. of FGMRES iterations = ', krylovIter
268 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
269 & SQUEEZE_RIGHT, myThid )
270 _END_MASTER( myThid )
271 ENDIF
272 IF ( krylovIter.EQ.SEAICEkrylovIterMax ) THEN
273 krylovFails = krylovFails + 1
274 ENDIF
275 C Update linear solution vector and return to Newton iteration
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 uIce(I,J,bi,bj) = uIce(I,J,bi,bj)+duIce(I,J,bi,bj)
281 vIce(I,J,bi,bj) = vIce(I,J,bi,bj)+dvIce(I,J,bi,bj)
282 C reset du/vIce here instead of setting sol = 0 in seaice_fgmres_driver
283 duIce(I,J,bi,bj)= 0. _d 0
284 dvIce(I,J,bi,bj)= 0. _d 0
285 ENDDO
286 ENDDO
287 ENDDO
288 ENDDO
289 C Set the stopping criterion for the Newton iteration
290 IF ( newtonIter .EQ. 1 ) JFNKtol=JFNKgamma_nonlin*JFNKresidual
291 ENDIF
292 C end of Newton iterate
293 ENDDO
294 C
295 C-- Output diagnostics
296 C
297 IF ( SEAICE_monFreq .GT. 0. _d 0 ) THEN
298 C Count iterations
299 totalJFNKtimeSteps = totalJFNKtimeSteps + 1
300 totalNewtonIters = totalNewtonIters + newtonIter
301 totalKrylovIters = totalKrylovIters + totalKrylovItersLoc
302 C Record failure
303 totalKrylovFails = totalKrylovFails + krylovFails
304 IF ( newtonIter .EQ. SEAICEnewtonIterMax ) THEN
305 totalNewtonFails = totalNewtonFails + 1
306 ENDIF
307 ENDIF
308 C Decide whether it is time to dump and reset the counter
309 writeNow = DIFFERENT_MULTIPLE(SEAICE_monFreq,
310 & myTime+deltaTClock, deltaTClock)
311 #ifdef ALLOW_CAL
312 IF ( useCAL ) THEN
313 CALL CAL_TIME2DUMP(
314 I zeroRL, SEAICE_monFreq, deltaTClock,
315 U writeNow,
316 I myTime+deltaTclock, myIter+1, myThid )
317 ENDIF
318 #endif
319 IF ( writeNow ) THEN
320 _BEGIN_MASTER( myThid )
321 WRITE(msgBuf,'(A)')
322 &' // ======================================================='
323 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
324 & SQUEEZE_RIGHT, myThid )
325 WRITE(msgBuf,'(A)') ' // Begin JFNK statistics'
326 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
327 & SQUEEZE_RIGHT, myThid )
328 WRITE(msgBuf,'(A)')
329 &' // ======================================================='
330 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
331 & SQUEEZE_RIGHT, myThid )
332 WRITE(msgBuf,'(A,I10)')
333 & ' %JFNK_MON: time step = ', myIter+1
334 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
335 & SQUEEZE_RIGHT, myThid )
336 WRITE(msgBuf,'(A,I10)')
337 & ' %JFNK_MON: Nb. of time steps = ', totalJFNKtimeSteps
338 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
339 & SQUEEZE_RIGHT, myThid )
340 WRITE(msgBuf,'(A,I10)')
341 & ' %JFNK_MON: Nb. of Newton steps = ', totalNewtonIters
342 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
343 & SQUEEZE_RIGHT, myThid )
344 WRITE(msgBuf,'(A,I10)')
345 & ' %JFNK_MON: Nb. of Krylov steps = ', totalKrylovIters
346 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
347 & SQUEEZE_RIGHT, myThid )
348 WRITE(msgBuf,'(A,I10)')
349 & ' %JFNK_MON: Nb. of Newton failures = ', totalNewtonFails
350 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
351 & SQUEEZE_RIGHT, myThid )
352 WRITE(msgBuf,'(A,I10)')
353 & ' %JFNK_MON: Nb. of Krylov failures = ', totalKrylovFails
354 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
355 & SQUEEZE_RIGHT, myThid )
356 WRITE(msgBuf,'(A)')
357 &' // ======================================================='
358 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
359 & SQUEEZE_RIGHT, myThid )
360 WRITE(msgBuf,'(A)') ' // End JFNK statistics'
361 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
362 & SQUEEZE_RIGHT, myThid )
363 WRITE(msgBuf,'(A)')
364 &' // ======================================================='
365 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
366 & SQUEEZE_RIGHT, myThid )
367 _END_MASTER( myThid )
368 C reset and start again
369 totalJFNKtimeSteps = 0
370 totalNewtonIters = 0
371 totalKrylovIters = 0
372 totalKrylovFails = 0
373 totalNewtonFails = 0
374 ENDIF
375
376 C Print more debugging information
377 IF ( debugLevel.GE.debLevA ) THEN
378 IF ( newtonIter .EQ. SEAICEnewtonIterMax ) THEN
379 _BEGIN_MASTER( myThid )
380 WRITE(msgBuf,'(A,I10)')
381 & ' S/R SEAICE_JFNK: JFNK did not converge in timestep ',
382 & myIter+1
383 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
384 & SQUEEZE_RIGHT, myThid )
385 _END_MASTER( myThid )
386 ENDIF
387 IF ( krylovFails .GT. 0 ) THEN
388 _BEGIN_MASTER( myThid )
389 WRITE(msgBuf,'(A,I4,A,I10)')
390 & ' S/R SEAICE_JFNK: FGMRES did not converge ',
391 & krylovFails, ' times in timestep ', myIter+1
392 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
393 & SQUEEZE_RIGHT, myThid )
394 _END_MASTER( myThid )
395 ENDIF
396 _BEGIN_MASTER( myThid )
397 WRITE(msgBuf,'(A,I6,A,I10)')
398 & ' S/R SEAICE_JFNK: Total number FGMRES iterations = ',
399 & totalKrylovItersLoc, ' in timestep ', myIter+1
400 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
401 & SQUEEZE_RIGHT, myThid )
402 _END_MASTER( myThid )
403 ENDIF
404
405 #endif /* SEAICE_ALLOW_DYNAMICS and SEAICE_CGRID and SEAICE_ALLOW_JFNK */
406
407 RETURN
408 END

  ViewVC Help
Powered by ViewVC 1.1.22