/[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.4 - (hide annotations) (download)
Tue Nov 6 13:18:14 2012 UTC (12 years, 8 months ago) by mlosch
Branch: MAIN
Changes since 1.3: +4 -1 lines
move resetting initial guess for fgmres outside of seaice_fgmres_driver
in order to make restarts with sol .ne. zero (fgmres with restarts)

1 mlosch 1.4 C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_jfnk.F,v 1.3 2012/11/06 12:53:14 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    
55     C i,j,bi,bj :: loop indices
56     INTEGER i,j,bi,bj
57     C loop indices
58     INTEGER newtonIter, newtonIterFail
59     INTEGER krylovIter, krylovIterFail
60     INTEGER totalKrylovIter
61     C FGMRES flag that indicates what to do next
62     INTEGER iCode
63     _RL JFNKresidual, JFNKresidualTile(nSx,nSy)
64     _RL JFNKresidualKm1
65     C parameters to compute convergence criterion
66     _RL phi_e, alp_e, JFNKgamma_lin
67     _RL FGMRESeps
68     _RL JFNKtol
69     C
70     _RL recip_deltaT
71     LOGICAL JFNKconverged, krylovConverged
72     CHARACTER*(MAX_LEN_MBUF) msgBuf
73     C
74     C u/vIceRes :: residual of sea-ice momentum equations
75     _RL uIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
76     _RL vIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
77     C du/vIce :: ice velocity increment to be added to u/vIce
78     _RL duIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
79     _RL dvIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
80     C precomputed (= constant per Newton iteration) versions of
81 mlosch 1.2 C zeta, eta, and DWATN, press
82     _RL zetaPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
83     _RL etaPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
84     _RL dwatPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
85     _RL pressPre(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
86 mlosch 1.1 CEOP
87    
88     C Initialise
89     newtonIter = 0
90     newtonIterFail = 0
91     krylovIterFail = 0
92     totalKrylovIter = 0
93     JFNKconverged = .FALSE.
94     JFNKtol = 0. _d 0
95     JFNKresidual = 0. _d 0
96     JFNKresidualKm1 = 0. _d 0
97     FGMRESeps = 0. _d 0
98     recip_deltaT = 1. _d 0 / SEAICE_deltaTdyn
99     C
100     DO bj=myByLo(myThid),myByHi(myThid)
101     DO bi=myBxLo(myThid),myBxHi(myThid)
102     DO J=1-Oly,sNy+Oly
103     DO I=1-Olx,sNx+Olx
104     uIceRes(I,J,bi,bj) = 0. _d 0
105     vIceRes(I,J,bi,bj) = 0. _d 0
106     duIce (I,J,bi,bj) = 0. _d 0
107     dvIce (I,J,bi,bj) = 0. _d 0
108     uIceNm1(I,J,bi,bj) = uIce(I,J,bi,bj)
109     vIceNm1(I,J,bi,bj) = vIce(I,J,bi,bj)
110     ENDDO
111     ENDDO
112     C Compute things that do no change during the Newton iteration:
113     C sea-surface tilt and wind stress:
114     C FORCEX/Y0 - mass*(u/vIceNm1)/deltaT
115     DO J=1-Oly,sNy+Oly
116     DO I=1-Olx,sNx+Olx
117     FORCEX(I,J,bi,bj) = FORCEX0(I,J,bi,bj)
118     & + seaiceMassU(I,J,bi,bj)*uIceNm1(I,J,bi,bj)*recip_deltaT
119     FORCEY(I,J,bi,bj) = FORCEY0(I,J,bi,bj)
120     & + seaiceMassV(I,J,bi,bj)*vIceNm1(I,J,bi,bj)*recip_deltaT
121     ENDDO
122     ENDDO
123     ENDDO
124     ENDDO
125     C Start nonlinear Newton iteration: outer loop iteration
126     DO WHILE ( newtonIter.LT.SEAICEnewtonIterMax .AND.
127     & .NOT.JFNKconverged )
128     newtonIter = newtonIter + 1
129     C Compute initial residual F(u), (includes computation of global
130     C variables DWATN, zeta, and eta)
131     CALL SEAICE_CALC_RESIDUAL(
132     I uIce, vIce,
133     O uIceRes, vIceRes,
134     I newtonIter, 0, myTime, myIter, myThid )
135 mlosch 1.3 CALL EXCH_UV_XY_RL( uIceRes, vIceRes,.TRUE.,myThid)
136 mlosch 1.1 C local copies of precomputed coefficients that are to stay
137     C constant for the preconditioner
138     DO bj=myByLo(myThid),myByHi(myThid)
139     DO bi=myBxLo(myThid),myBxHi(myThid)
140     DO j=1-Oly,sNy+Oly
141     DO i=1-Olx,sNx+Olx
142 mlosch 1.2 zetaPre(I,J,bi,bj) = zeta(I,J,bi,bj)
143     etaPre(I,J,bi,bj) = eta(I,J,bi,bj)
144     dwatPre(I,J,bi,bj) = DWATN(I,J,bi,bj)
145     pressPre(I,J,bi,bj) = press(I,J,bi,bj)
146 mlosch 1.1 ENDDO
147     ENDDO
148     ENDDO
149     ENDDO
150     C
151     DO bj=myByLo(myThid),myByHi(myThid)
152     DO bi=myBxLo(myThid),myBxHi(myThid)
153     JFNKresidualTile(bi,bj) = 0. _d 0
154     DO J=1,sNy
155     DO I=1,sNx
156     #ifdef CG2D_SINGLECPU_SUM
157     JFNKlocalBuf(I,J,bi,bj) =
158     #else
159     JFNKresidualTile(bi,bj) = JFNKresidualTile(bi,bj) +
160     #endif
161     & uIceRes(I,J,bi,bj)*uIceRes(I,J,bi,bj) +
162     & vIceRes(I,J,bi,bj)*vIceRes(I,J,bi,bj)
163     ENDDO
164     ENDDO
165     ENDDO
166     ENDDO
167     JFNKresidual = 0. _d 0
168     #ifdef CG2D_SINGLECPU_SUM
169     CALL GLOBAL_SUM_SINGLECPU_RL(
170     & JFNKlocalBuf,JFNKresidual, 0, 0, myThid)
171     #else
172     CALL GLOBAL_SUM_TILE_RL( JFNKresidualTile,JFNKresidual,myThid )
173     #endif
174     JFNKresidual = SQRT(JFNKresidual)
175     C compute convergence criterion for linear preconditioned FGMRES
176     JFNKgamma_lin = JFNKgamma_lin_max
177     IF ( newtonIter.GT.1.AND.newtonIter.LE.100
178     & .AND.JFNKresidual.LT.JFNKres_t ) THEN
179     C Eisenstat, 1996, equ.(2.6)
180     phi_e = 1. _d 0
181     alp_e = 1. _d 0
182     JFNKgamma_lin = phi_e*( JFNKresidual/JFNKresidualKm1 )**alp_e
183     JFNKgamma_lin = min(JFNKgamma_lin_max, JFNKgamma_lin)
184     JFNKgamma_lin = max(JFNKgamma_lin_min, JFNKgamma_lin)
185     ENDIF
186     C save the residual for the next iteration
187     JFNKresidualKm1 = JFNKresidual
188     C
189     C The Krylov iteration using FGMRES, the preconditioner is LSOR
190     C for now. The code is adapted from SEAICE_LSR, but heavily stripped
191     C down.
192     C krylovIter is mapped into "its" in seaice_fgmres and is incremented
193     C in that routine
194     krylovIter = 0
195     iCode = 0
196     IF ( debugLevel.GE.debLevA ) THEN
197     WRITE(msgBuf,'(2A,2(1XI6),2E12.5)')
198     & ' S/R SEAICE_JFNK: newtonIter,',
199     & ' total newtonIter, JFNKgamma_lin, initial norm = ',
200 mlosch 1.3 & newtonIter,SEAICEnewtonIterMax*(myIter-nIter0)+newtonIter,
201 mlosch 1.1 & JFNKgamma_lin, JFNKresidual
202     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
203     & SQUEEZE_RIGHT, myThid )
204     ENDIF
205     C
206     JFNKconverged = JFNKresidual.LT.JFNKtol
207     C
208     C do Krylov loop only if convergence is not reached
209     C
210     IF ( .NOT.JFNKconverged ) THEN
211     C
212     C start Krylov iteration (FGMRES)
213     C
214     krylovConverged = .FALSE.
215     FGMRESeps = JFNKgamma_lin * JFNKresidual
216     DO WHILE ( .NOT.krylovConverged )
217     C solution vector sol = du/vIce
218     C residual vector (rhs) Fu = u/vIceRes
219     C output work vectors wk1, -> input work vector wk2
220     C
221     CALL SEAICE_FGMRES_DRIVER(
222     I uIceRes, vIceRes,
223     U duIce, dvIce, iCode,
224     I FGMRESeps,
225     I newtonIter, krylovIter, myTime, myIter, myThid )
226     C FGMRES returns iCode either asking for an new preconditioned vector
227     C or product of matrix (Jacobian) times vector. For iCode = 0, terminate
228     C iteration
229     IF (iCode.EQ.1) THEN
230     C Call preconditioner
231     CALL SEAICE_PRECONDITIONER(
232     U duIce, dvIce,
233 mlosch 1.2 I zetaPre, etaPre, dwatPre, pressPre,
234 mlosch 1.1 I newtonIter, krylovIter, myTime, myIter, myThid )
235     ELSEIF (iCode.GE.2) THEN
236     C Compute Jacobian times vector
237     CALL SEAICE_JACVEC(
238     I uIce, vIce, uIceRes, vIceRes,
239     U duIce, dvIce,
240     I newtonIter, krylovIter, myTime, myIter, myThid )
241     ENDIF
242     krylovConverged = iCode.EQ.0
243     C End of Krylov iterate
244     ENDDO
245     totalKrylovIter = totalKrylovIter + krylovIter
246     C some output diagnostics
247     IF ( debugLevel.GE.debLevA ) THEN
248     WRITE(msgBuf,'(3(A,I6))')
249     & ' S/R SEAICE_JFNK: Newton iterate / total = ', newtonIter,
250 mlosch 1.3 & ' / ', SEAICEnewtonIterMax*(myIter-nIter0)+newtonIter,
251 mlosch 1.1 & ', Nb. of FGMRES iterations = ', krylovIter
252     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
253     & SQUEEZE_RIGHT, myThid )
254     ENDIF
255     IF ( krylovIter.EQ.SEAICEkrylovIterMax ) THEN
256     krylovIterFail = krylovIterFail + 1
257     ENDIF
258     C Update linear solution vector and return to Newton iteration
259     DO bj=myByLo(myThid),myByHi(myThid)
260     DO bi=myBxLo(myThid),myBxHi(myThid)
261     DO J=1-Oly,sNy+Oly
262     DO I=1-Olx,sNx+Olx
263     uIce(I,J,bi,bj) = uIce(I,J,bi,bj)+duIce(I,J,bi,bj)
264     vIce(I,J,bi,bj) = vIce(I,J,bi,bj)+dvIce(I,J,bi,bj)
265 mlosch 1.4 C reset du/vIce here instead of setting sol = 0 in seaice_fgmres_driver
266     duIce(I,J,bi,bj)= 0. _d 0
267     dvIce(I,J,bi,bj)= 0. _d 0
268 mlosch 1.1 ENDDO
269     ENDDO
270     ENDDO
271     ENDDO
272     C Set the stopping criterion for the Newton iteration
273     IF ( newtonIter .EQ. 1 ) JFNKtol=JFNKgamma_nonlin*JFNKresidual
274     ENDIF
275     C end of Newton iterate
276     ENDDO
277     C some output diagnostics
278     IF ( debugLevel.GE.debLevA ) THEN
279     C Record failure
280     IF ( newtonIter .EQ. SEAICEnewtonIterMax ) THEN
281     newtonIterFail = newtonIterFail + 1
282     WRITE(msgBuf,'(A,I10)')
283     & ' S/R SEAICE_JFNK: JFNK did not converge in timestep ',
284     & myIter
285     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
286     & SQUEEZE_RIGHT, myThid )
287     ENDIF
288     IF ( krylovIterFail .GT. 0 ) THEN
289     WRITE(msgBuf,'(A,I4,A,I10)')
290     & ' S/R SEAICE_JFNK: FGMRES did not converge ',
291     & krylovIterFail, ' times in timestep ', myIter
292     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
293     & SQUEEZE_RIGHT, myThid )
294     ENDIF
295     WRITE(msgBuf,'(A,I6)')
296     & ' S/R SEAICE_JFNK: Total number FGMRES iterations = ',
297     & totalKrylovIter
298     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
299     & SQUEEZE_RIGHT, myThid )
300    
301     ENDIF
302    
303     #endif /* SEAICE_ALLOW_DYNAMICS and SEAICE_CGRID and SEAICE_ALLOW_JFNK */
304    
305     RETURN
306     END

  ViewVC Help
Powered by ViewVC 1.1.22