/[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.13 - (hide annotations) (download)
Mon Dec 17 15:06:02 2012 UTC (12 years, 6 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint64b
Changes since 1.12: +21 -40 lines
  - add a metric based on grid cell area to SEAICE_SCALPROD

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

  ViewVC Help
Powered by ViewVC 1.1.22