/[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.8 - (hide annotations) (download)
Fri Nov 9 14:11:43 2012 UTC (12 years, 8 months ago) by mlosch
Branch: MAIN
Changes since 1.7: +4 -2 lines
compute etaZpre and modify call to preconditioner accordingly

1 mlosch 1.8 C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_jfnk.F,v 1.7 2012/11/09 12:56:00 mlosch Exp $
2 mlosch 1.1 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 mlosch 1.5 C !FUNCTIONS:
55     LOGICAL DIFFERENT_MULTIPLE
56     EXTERNAL DIFFERENT_MULTIPLE
57 mlosch 1.1
58     C i,j,bi,bj :: loop indices
59     INTEGER i,j,bi,bj
60     C loop indices
61 mlosch 1.5 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 mlosch 1.1 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 mlosch 1.2 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 mlosch 1.8 _RL etaZPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
90 mlosch 1.2 _RL dwatPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
91     _RL pressPre(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
92 mlosch 1.1 CEOP
93    
94     C Initialise
95 mlosch 1.5 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 mlosch 1.1 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 mlosch 1.3 CALL EXCH_UV_XY_RL( uIceRes, vIceRes,.TRUE.,myThid)
148 mlosch 1.1 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 mlosch 1.2 zetaPre(I,J,bi,bj) = zeta(I,J,bi,bj)
155     etaPre(I,J,bi,bj) = eta(I,J,bi,bj)
156 mlosch 1.8 etaZPre(I,J,bi,bj) = etaZ(I,J,bi,bj)
157 mlosch 1.2 dwatPre(I,J,bi,bj) = DWATN(I,J,bi,bj)
158     pressPre(I,J,bi,bj) = press(I,J,bi,bj)
159 mlosch 1.1 ENDDO
160     ENDDO
161     ENDDO
162     ENDDO
163     C
164     DO bj=myByLo(myThid),myByHi(myThid)
165     DO bi=myBxLo(myThid),myBxHi(myThid)
166     JFNKresidualTile(bi,bj) = 0. _d 0
167     DO J=1,sNy
168     DO I=1,sNx
169     #ifdef CG2D_SINGLECPU_SUM
170     JFNKlocalBuf(I,J,bi,bj) =
171     #else
172     JFNKresidualTile(bi,bj) = JFNKresidualTile(bi,bj) +
173     #endif
174     & uIceRes(I,J,bi,bj)*uIceRes(I,J,bi,bj) +
175     & vIceRes(I,J,bi,bj)*vIceRes(I,J,bi,bj)
176     ENDDO
177     ENDDO
178     ENDDO
179     ENDDO
180     JFNKresidual = 0. _d 0
181     #ifdef CG2D_SINGLECPU_SUM
182     CALL GLOBAL_SUM_SINGLECPU_RL(
183     & JFNKlocalBuf,JFNKresidual, 0, 0, myThid)
184     #else
185     CALL GLOBAL_SUM_TILE_RL( JFNKresidualTile,JFNKresidual,myThid )
186     #endif
187     JFNKresidual = SQRT(JFNKresidual)
188     C compute convergence criterion for linear preconditioned FGMRES
189     JFNKgamma_lin = JFNKgamma_lin_max
190     IF ( newtonIter.GT.1.AND.newtonIter.LE.100
191     & .AND.JFNKresidual.LT.JFNKres_t ) THEN
192     C Eisenstat, 1996, equ.(2.6)
193     phi_e = 1. _d 0
194     alp_e = 1. _d 0
195     JFNKgamma_lin = phi_e*( JFNKresidual/JFNKresidualKm1 )**alp_e
196     JFNKgamma_lin = min(JFNKgamma_lin_max, JFNKgamma_lin)
197     JFNKgamma_lin = max(JFNKgamma_lin_min, JFNKgamma_lin)
198     ENDIF
199     C save the residual for the next iteration
200     JFNKresidualKm1 = JFNKresidual
201     C
202     C The Krylov iteration using FGMRES, the preconditioner is LSOR
203     C for now. The code is adapted from SEAICE_LSR, but heavily stripped
204     C down.
205     C krylovIter is mapped into "its" in seaice_fgmres and is incremented
206     C in that routine
207     krylovIter = 0
208     iCode = 0
209     IF ( debugLevel.GE.debLevA ) THEN
210 mlosch 1.5 _BEGIN_MASTER( myThid )
211 mlosch 1.1 WRITE(msgBuf,'(2A,2(1XI6),2E12.5)')
212     & ' S/R SEAICE_JFNK: newtonIter,',
213     & ' total newtonIter, JFNKgamma_lin, initial norm = ',
214 mlosch 1.3 & newtonIter,SEAICEnewtonIterMax*(myIter-nIter0)+newtonIter,
215 mlosch 1.1 & JFNKgamma_lin, JFNKresidual
216     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
217     & SQUEEZE_RIGHT, myThid )
218 mlosch 1.5 _END_MASTER( myThid )
219 mlosch 1.1 ENDIF
220     C
221     JFNKconverged = JFNKresidual.LT.JFNKtol
222     C
223     C do Krylov loop only if convergence is not reached
224     C
225     IF ( .NOT.JFNKconverged ) THEN
226     C
227     C start Krylov iteration (FGMRES)
228     C
229     krylovConverged = .FALSE.
230     FGMRESeps = JFNKgamma_lin * JFNKresidual
231     DO WHILE ( .NOT.krylovConverged )
232     C solution vector sol = du/vIce
233     C residual vector (rhs) Fu = u/vIceRes
234     C output work vectors wk1, -> input work vector wk2
235     C
236     CALL SEAICE_FGMRES_DRIVER(
237     I uIceRes, vIceRes,
238     U duIce, dvIce, iCode,
239 mlosch 1.5 I FGMRESeps, iOutFGMRES,
240 mlosch 1.1 I newtonIter, krylovIter, myTime, myIter, myThid )
241     C FGMRES returns iCode either asking for an new preconditioned vector
242     C or product of matrix (Jacobian) times vector. For iCode = 0, terminate
243     C iteration
244     IF (iCode.EQ.1) THEN
245 mlosch 1.7 C Call preconditioner
246     IF ( SOLV_MAX_ITERS .GT. 0 )
247     & CALL SEAICE_PRECONDITIONER(
248 mlosch 1.1 U duIce, dvIce,
249 mlosch 1.8 I zetaPre, etaPre, etaZpre, dwatPre, pressPre,
250 mlosch 1.1 I newtonIter, krylovIter, myTime, myIter, myThid )
251     ELSEIF (iCode.GE.2) THEN
252     C Compute Jacobian times vector
253     CALL SEAICE_JACVEC(
254     I uIce, vIce, uIceRes, vIceRes,
255     U duIce, dvIce,
256     I newtonIter, krylovIter, myTime, myIter, myThid )
257     ENDIF
258     krylovConverged = iCode.EQ.0
259     C End of Krylov iterate
260     ENDDO
261 mlosch 1.5 totalKrylovItersLoc = totalKrylovItersLoc + krylovIter
262 mlosch 1.1 C some output diagnostics
263     IF ( debugLevel.GE.debLevA ) THEN
264 mlosch 1.5 _BEGIN_MASTER( myThid )
265 mlosch 1.1 WRITE(msgBuf,'(3(A,I6))')
266     & ' S/R SEAICE_JFNK: Newton iterate / total = ', newtonIter,
267 mlosch 1.3 & ' / ', SEAICEnewtonIterMax*(myIter-nIter0)+newtonIter,
268 mlosch 1.1 & ', Nb. of FGMRES iterations = ', krylovIter
269     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
270     & SQUEEZE_RIGHT, myThid )
271 mlosch 1.5 _END_MASTER( myThid )
272 mlosch 1.1 ENDIF
273     IF ( krylovIter.EQ.SEAICEkrylovIterMax ) THEN
274 mlosch 1.5 krylovFails = krylovFails + 1
275 mlosch 1.1 ENDIF
276     C Update linear solution vector and return to Newton iteration
277     DO bj=myByLo(myThid),myByHi(myThid)
278     DO bi=myBxLo(myThid),myBxHi(myThid)
279     DO J=1-Oly,sNy+Oly
280     DO I=1-Olx,sNx+Olx
281     uIce(I,J,bi,bj) = uIce(I,J,bi,bj)+duIce(I,J,bi,bj)
282     vIce(I,J,bi,bj) = vIce(I,J,bi,bj)+dvIce(I,J,bi,bj)
283 mlosch 1.4 C reset du/vIce here instead of setting sol = 0 in seaice_fgmres_driver
284     duIce(I,J,bi,bj)= 0. _d 0
285     dvIce(I,J,bi,bj)= 0. _d 0
286 mlosch 1.1 ENDDO
287     ENDDO
288     ENDDO
289     ENDDO
290     C Set the stopping criterion for the Newton iteration
291     IF ( newtonIter .EQ. 1 ) JFNKtol=JFNKgamma_nonlin*JFNKresidual
292     ENDIF
293     C end of Newton iterate
294     ENDDO
295 mlosch 1.5 C
296     C-- Output diagnostics
297     C
298 mlosch 1.6 IF ( SEAICE_monFreq .GT. 0. _d 0 ) THEN
299 mlosch 1.5 C Count iterations
300 mlosch 1.6 totalJFNKtimeSteps = totalJFNKtimeSteps + 1
301     totalNewtonIters = totalNewtonIters + newtonIter
302     totalKrylovIters = totalKrylovIters + totalKrylovItersLoc
303 mlosch 1.5 C Record failure
304 mlosch 1.6 totalKrylovFails = totalKrylovFails + krylovFails
305     IF ( newtonIter .EQ. SEAICEnewtonIterMax ) THEN
306     totalNewtonFails = totalNewtonFails + 1
307     ENDIF
308 mlosch 1.5 ENDIF
309     C Decide whether it is time to dump and reset the counter
310     IF ( DIFFERENT_MULTIPLE(SEAICE_monFreq,myTime+deltaTClock,
311     & deltaTClock) ) THEN
312     _BEGIN_MASTER( myThid )
313     WRITE(msgBuf,'(A)')
314     &' // ======================================================='
315     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
316     & SQUEEZE_RIGHT, myThid )
317     WRITE(msgBuf,'(A)') ' // Begin JFNK statistics'
318     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
319     & SQUEEZE_RIGHT, myThid )
320     WRITE(msgBuf,'(A)')
321     &' // ======================================================='
322     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
323     & SQUEEZE_RIGHT, myThid )
324     WRITE(msgBuf,'(A,I10)')
325     & ' %JFNK_MON: time step = ', myIter+1
326     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
327     & SQUEEZE_RIGHT, myThid )
328     WRITE(msgBuf,'(A,I10)')
329     & ' %JFNK_MON: Nb. of time steps = ', totalJFNKtimeSteps
330     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
331     & SQUEEZE_RIGHT, myThid )
332     WRITE(msgBuf,'(A,I10)')
333     & ' %JFNK_MON: Nb. of Newton steps = ', totalNewtonIters
334     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
335     & SQUEEZE_RIGHT, myThid )
336     WRITE(msgBuf,'(A,I10)')
337     & ' %JFNK_MON: Nb. of Krylov steps = ', totalKrylovIters
338     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
339     & SQUEEZE_RIGHT, myThid )
340     WRITE(msgBuf,'(A,I10)')
341     & ' %JFNK_MON: Nb. of Newton failures = ', totalNewtonFails
342     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
343     & SQUEEZE_RIGHT, myThid )
344     WRITE(msgBuf,'(A,I10)')
345     & ' %JFNK_MON: Nb. of Krylov failures = ', totalKrylovFails
346     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
347     & SQUEEZE_RIGHT, myThid )
348     WRITE(msgBuf,'(A)')
349     &' // ======================================================='
350     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
351     & SQUEEZE_RIGHT, myThid )
352     WRITE(msgBuf,'(A)') ' // Begin JFNK statistics'
353     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
354     & SQUEEZE_RIGHT, myThid )
355     WRITE(msgBuf,'(A)')
356     &' // ======================================================='
357     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
358     & SQUEEZE_RIGHT, myThid )
359     _END_MASTER( myThid )
360     C reset and start again
361     totalJFNKtimeSteps = 0
362     totalNewtonIters = 0
363     totalKrylovIters = 0
364     totalKrylovFails = 0
365     totalNewtonFails = 0
366     ENDIF
367    
368     C Print more debugging information
369 mlosch 1.1 IF ( debugLevel.GE.debLevA ) THEN
370     IF ( newtonIter .EQ. SEAICEnewtonIterMax ) THEN
371 mlosch 1.5 _BEGIN_MASTER( myThid )
372 mlosch 1.1 WRITE(msgBuf,'(A,I10)')
373     & ' S/R SEAICE_JFNK: JFNK did not converge in timestep ',
374 mlosch 1.5 & myIter+1
375 mlosch 1.1 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
376     & SQUEEZE_RIGHT, myThid )
377 mlosch 1.5 _END_MASTER( myThid )
378 mlosch 1.1 ENDIF
379 mlosch 1.5 IF ( krylovFails .GT. 0 ) THEN
380     _BEGIN_MASTER( myThid )
381 mlosch 1.1 WRITE(msgBuf,'(A,I4,A,I10)')
382     & ' S/R SEAICE_JFNK: FGMRES did not converge ',
383 mlosch 1.5 & krylovFails, ' times in timestep ', myIter+1
384 mlosch 1.1 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
385     & SQUEEZE_RIGHT, myThid )
386 mlosch 1.5 _END_MASTER( myThid )
387 mlosch 1.1 ENDIF
388 mlosch 1.5 _BEGIN_MASTER( myThid )
389     WRITE(msgBuf,'(A,I6,A,I10)')
390 mlosch 1.1 & ' S/R SEAICE_JFNK: Total number FGMRES iterations = ',
391 mlosch 1.5 & totalKrylovItersLoc, ' in timestep ', myIter+1
392     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
393     & SQUEEZE_RIGHT, myThid )
394     _END_MASTER( myThid )
395 mlosch 1.1 ENDIF
396    
397     #endif /* SEAICE_ALLOW_DYNAMICS and SEAICE_CGRID and SEAICE_ALLOW_JFNK */
398    
399     RETURN
400     END

  ViewVC Help
Powered by ViewVC 1.1.22