/[MITgcm]/MITgcm/pkg/seaice/seaice_jfnk.F
ViewVC logotype

Annotation of /MITgcm/pkg/seaice/seaice_jfnk.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.31 - (hide annotations) (download)
Thu Jun 8 15:33:50 2017 UTC (8 years, 1 month ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, HEAD
Changes since 1.30: +2 -1 lines
prevent JFNK and KRYLOV solvers from trying to reduce a residual of 0

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

  ViewVC Help
Powered by ViewVC 1.1.22