/[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.26 - (hide annotations) (download)
Thu Mar 20 09:24:49 2014 UTC (11 years, 3 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64w, checkpoint64v, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65d, checkpoint65e, checkpoint65
Changes since 1.25: +4 -3 lines
comment out if-statement related to IMEX to avoid suprises

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

  ViewVC Help
Powered by ViewVC 1.1.22