/[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.2 - (hide annotations) (download)
Wed Oct 17 14:53:51 2012 UTC (12 years, 9 months ago) by mlosch
Branch: MAIN
Changes since 1.1: +11 -9 lines
save also precomputed pressure and pass it to the preconditioner

1 mlosch 1.2 C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_jfnk.F,v 1.1 2012/10/16 07:00:21 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     C local copies of precomputed coefficients that are to stay
136     C constant for the preconditioner
137     DO bj=myByLo(myThid),myByHi(myThid)
138     DO bi=myBxLo(myThid),myBxHi(myThid)
139     DO j=1-Oly,sNy+Oly
140     DO i=1-Olx,sNx+Olx
141 mlosch 1.2 zetaPre(I,J,bi,bj) = zeta(I,J,bi,bj)
142     etaPre(I,J,bi,bj) = eta(I,J,bi,bj)
143     dwatPre(I,J,bi,bj) = DWATN(I,J,bi,bj)
144     pressPre(I,J,bi,bj) = press(I,J,bi,bj)
145 mlosch 1.1 ENDDO
146     ENDDO
147     ENDDO
148     ENDDO
149     C
150     DO bj=myByLo(myThid),myByHi(myThid)
151     DO bi=myBxLo(myThid),myBxHi(myThid)
152     JFNKresidualTile(bi,bj) = 0. _d 0
153     DO J=1,sNy
154     DO I=1,sNx
155     #ifdef CG2D_SINGLECPU_SUM
156     JFNKlocalBuf(I,J,bi,bj) =
157     #else
158     JFNKresidualTile(bi,bj) = JFNKresidualTile(bi,bj) +
159     #endif
160     & uIceRes(I,J,bi,bj)*uIceRes(I,J,bi,bj) +
161     & vIceRes(I,J,bi,bj)*vIceRes(I,J,bi,bj)
162     ENDDO
163     ENDDO
164     ENDDO
165     ENDDO
166     JFNKresidual = 0. _d 0
167     #ifdef CG2D_SINGLECPU_SUM
168     CALL GLOBAL_SUM_SINGLECPU_RL(
169     & JFNKlocalBuf,JFNKresidual, 0, 0, myThid)
170     #else
171     CALL GLOBAL_SUM_TILE_RL( JFNKresidualTile,JFNKresidual,myThid )
172     #endif
173     JFNKresidual = SQRT(JFNKresidual)
174     C compute convergence criterion for linear preconditioned FGMRES
175     JFNKgamma_lin = JFNKgamma_lin_max
176     IF ( newtonIter.GT.1.AND.newtonIter.LE.100
177     & .AND.JFNKresidual.LT.JFNKres_t ) THEN
178     C Eisenstat, 1996, equ.(2.6)
179     phi_e = 1. _d 0
180     alp_e = 1. _d 0
181     JFNKgamma_lin = phi_e*( JFNKresidual/JFNKresidualKm1 )**alp_e
182     JFNKgamma_lin = min(JFNKgamma_lin_max, JFNKgamma_lin)
183     JFNKgamma_lin = max(JFNKgamma_lin_min, JFNKgamma_lin)
184     ENDIF
185     C save the residual for the next iteration
186     JFNKresidualKm1 = JFNKresidual
187     C
188     C The Krylov iteration using FGMRES, the preconditioner is LSOR
189     C for now. The code is adapted from SEAICE_LSR, but heavily stripped
190     C down.
191     C krylovIter is mapped into "its" in seaice_fgmres and is incremented
192     C in that routine
193     krylovIter = 0
194     iCode = 0
195     IF ( debugLevel.GE.debLevA ) THEN
196     WRITE(msgBuf,'(2A,2(1XI6),2E12.5)')
197     & ' S/R SEAICE_JFNK: newtonIter,',
198     & ' total newtonIter, JFNKgamma_lin, initial norm = ',
199     & newtonIter, SEAICEnewtonIterMax*myIter+newtonIter,
200     & JFNKgamma_lin, JFNKresidual
201     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
202     & SQUEEZE_RIGHT, myThid )
203     ENDIF
204     C
205     JFNKconverged = JFNKresidual.LT.JFNKtol
206     C
207     C do Krylov loop only if convergence is not reached
208     C
209     IF ( .NOT.JFNKconverged ) THEN
210     C
211     C start Krylov iteration (FGMRES)
212     C
213     krylovConverged = .FALSE.
214     FGMRESeps = JFNKgamma_lin * JFNKresidual
215     DO WHILE ( .NOT.krylovConverged )
216     C solution vector sol = du/vIce
217     C residual vector (rhs) Fu = u/vIceRes
218     C output work vectors wk1, -> input work vector wk2
219     C
220     CALL EXCH_UV_XY_RL( uIceRes, vIceRes,.TRUE.,myThid)
221     CALL EXCH_UV_XY_RL( duIce, dvIce,.TRUE.,myThid)
222     CALL SEAICE_FGMRES_DRIVER(
223     I uIceRes, vIceRes,
224     U duIce, dvIce, iCode,
225     I FGMRESeps,
226     I newtonIter, krylovIter, myTime, myIter, myThid )
227     C FGMRES returns iCode either asking for an new preconditioned vector
228     C or product of matrix (Jacobian) times vector. For iCode = 0, terminate
229     C iteration
230     IF (iCode.EQ.1) THEN
231     C Call preconditioner
232     CALL SEAICE_PRECONDITIONER(
233     U duIce, dvIce,
234 mlosch 1.2 I zetaPre, etaPre, dwatPre, pressPre,
235 mlosch 1.1 I newtonIter, krylovIter, myTime, myIter, myThid )
236     ELSEIF (iCode.GE.2) THEN
237     C Compute Jacobian times vector
238     CALL SEAICE_JACVEC(
239     I uIce, vIce, uIceRes, vIceRes,
240     U duIce, dvIce,
241     I newtonIter, krylovIter, myTime, myIter, myThid )
242     ENDIF
243     krylovConverged = iCode.EQ.0
244     C End of Krylov iterate
245     ENDDO
246     totalKrylovIter = totalKrylovIter + krylovIter
247     C some output diagnostics
248     IF ( debugLevel.GE.debLevA ) THEN
249     WRITE(msgBuf,'(3(A,I6))')
250     & ' S/R SEAICE_JFNK: Newton iterate / total = ', newtonIter,
251     & ' / ', SEAICEnewtonIterMax*myIter+newtonIter,
252     & ', Nb. of FGMRES iterations = ', krylovIter
253     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
254     & SQUEEZE_RIGHT, myThid )
255     ENDIF
256     IF ( krylovIter.EQ.SEAICEkrylovIterMax ) THEN
257     krylovIterFail = krylovIterFail + 1
258     ENDIF
259     C Update linear solution vector and return to Newton iteration
260     DO bj=myByLo(myThid),myByHi(myThid)
261     DO bi=myBxLo(myThid),myBxHi(myThid)
262     DO J=1-Oly,sNy+Oly
263     DO I=1-Olx,sNx+Olx
264     uIce(I,J,bi,bj) = uIce(I,J,bi,bj)+duIce(I,J,bi,bj)
265     vIce(I,J,bi,bj) = vIce(I,J,bi,bj)+dvIce(I,J,bi,bj)
266     ENDDO
267     ENDDO
268     ENDDO
269     ENDDO
270     C Set the stopping criterion for the Newton iteration
271     IF ( newtonIter .EQ. 1 ) JFNKtol=JFNKgamma_nonlin*JFNKresidual
272     ENDIF
273     C end of Newton iterate
274     ENDDO
275     C some output diagnostics
276     IF ( debugLevel.GE.debLevA ) THEN
277     C Record failure
278     IF ( newtonIter .EQ. SEAICEnewtonIterMax ) THEN
279     newtonIterFail = newtonIterFail + 1
280     WRITE(msgBuf,'(A,I10)')
281     & ' S/R SEAICE_JFNK: JFNK did not converge in timestep ',
282     & myIter
283     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
284     & SQUEEZE_RIGHT, myThid )
285     ENDIF
286     IF ( krylovIterFail .GT. 0 ) THEN
287     WRITE(msgBuf,'(A,I4,A,I10)')
288     & ' S/R SEAICE_JFNK: FGMRES did not converge ',
289     & krylovIterFail, ' times in timestep ', myIter
290     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
291     & SQUEEZE_RIGHT, myThid )
292     ENDIF
293     WRITE(msgBuf,'(A,I6)')
294     & ' S/R SEAICE_JFNK: Total number FGMRES iterations = ',
295     & totalKrylovIter
296     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
297     & SQUEEZE_RIGHT, myThid )
298    
299     ENDIF
300    
301     #endif /* SEAICE_ALLOW_DYNAMICS and SEAICE_CGRID and SEAICE_ALLOW_JFNK */
302    
303     RETURN
304     END

  ViewVC Help
Powered by ViewVC 1.1.22