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