/[MITgcm]/MITgcm/pkg/seaice/seaice_jfnk.F
ViewVC logotype

Contents of /MITgcm/pkg/seaice/seaice_jfnk.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.1 - (show annotations) (download)
Tue Oct 16 07:00:21 2012 UTC (12 years, 9 months ago) by mlosch
Branch: MAIN
add JFNK-solver routines, mostly parallel and mult-threaded (except
for FGMRES)

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

  ViewVC Help
Powered by ViewVC 1.1.22