/[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.22 - (hide annotations) (download)
Tue Apr 23 08:40:06 2013 UTC (12 years, 2 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint64g, checkpoint64h
Changes since 1.21: +5 -6 lines
turn some parameters for choosing the convergence criterion of the
inexact Newton method (JFNK) into runtime parameters for convenience

1 mlosch 1.22 C $Header: /u/gcmpack/MITgcm/pkg/seaice/seaice_jfnk.F,v 1.21 2013/04/04 07:02:51 mlosch Exp $
2 mlosch 1.1 C $Name: $
3    
4     #include "SEAICE_OPTIONS.h"
5    
6 mlosch 1.15 C-- File seaice_jfnk.F: seaice jfnk dynamical solver S/R:
7     C-- Contents
8     C-- o SEAICE_JFNK
9     C-- o SEAICE_JFNK_UPDATE
10    
11 mlosch 1.1 CBOP
12     C !ROUTINE: SEAICE_JFNK
13     C !INTERFACE:
14     SUBROUTINE SEAICE_JFNK( myTime, myIter, myThid )
15    
16     C !DESCRIPTION: \bv
17     C *==========================================================*
18 mlosch 1.15 C | SUBROUTINE SEAICE_JFNK
19 mlosch 1.1 C | o Ice dynamics using a Jacobian-free Newton-Krylov solver
20     C | following J.-F. Lemieux et al. Improving the numerical
21     C | convergence of viscous-plastic sea ice models with the
22     C | Jacobian-free Newton-Krylov method. J. Comp. Phys. 229,
23     C | 2840-2852 (2010).
24     C | o The logic follows JFs code.
25     C *==========================================================*
26     C | written by Martin Losch, Oct 2012
27     C *==========================================================*
28     C \ev
29    
30     C !USES:
31     IMPLICIT NONE
32    
33     C === Global variables ===
34     #include "SIZE.h"
35     #include "EEPARAMS.h"
36     #include "PARAMS.h"
37     #include "DYNVARS.h"
38     #include "GRID.h"
39     #include "SEAICE_SIZE.h"
40     #include "SEAICE_PARAMS.h"
41     #include "SEAICE.h"
42    
43     #ifdef ALLOW_AUTODIFF_TAMC
44     # include "tamc.h"
45     #endif
46    
47     C !INPUT/OUTPUT PARAMETERS:
48     C === Routine arguments ===
49     C myTime :: Simulation time
50     C myIter :: Simulation timestep number
51     C myThid :: my Thread Id. number
52     _RL myTime
53     INTEGER myIter
54     INTEGER myThid
55    
56 mlosch 1.21 #ifdef SEAICE_ALLOW_JFNK
57 mlosch 1.5 C !FUNCTIONS:
58     LOGICAL DIFFERENT_MULTIPLE
59     EXTERNAL DIFFERENT_MULTIPLE
60 mlosch 1.1
61 mlosch 1.16 C !LOCAL VARIABLES:
62     C === Local variables ===
63 mlosch 1.1 C i,j,bi,bj :: loop indices
64     INTEGER i,j,bi,bj
65     C loop indices
66 mlosch 1.5 INTEGER newtonIter
67     INTEGER krylovIter, krylovFails
68 mlosch 1.13 INTEGER totalKrylovItersLoc, totalNewtonItersLoc
69 mlosch 1.5 C FGMRES flag that determines amount of output messages of fgmres
70     INTEGER iOutFGMRES
71     C FGMRES flag that indicates what fgmres wants us to do next
72 mlosch 1.1 INTEGER iCode
73 mlosch 1.13 _RL JFNKresidual
74 mlosch 1.1 _RL JFNKresidualKm1
75     C parameters to compute convergence criterion
76 mlosch 1.22 _RL JFNKgamma_lin
77 mlosch 1.1 _RL FGMRESeps
78     _RL JFNKtol
79 jmc 1.20
80 mlosch 1.1 _RL recip_deltaT
81     LOGICAL JFNKconverged, krylovConverged
82 mlosch 1.9 LOGICAL writeNow
83 mlosch 1.1 CHARACTER*(MAX_LEN_MBUF) msgBuf
84 jmc 1.20
85 mlosch 1.1 C u/vIceRes :: residual of sea-ice momentum equations
86     _RL uIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
87     _RL vIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
88     C du/vIce :: ice velocity increment to be added to u/vIce
89     _RL duIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
90     _RL dvIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
91 jmc 1.20 C precomputed (= constant per Newton iteration) versions of
92 mlosch 1.2 C zeta, eta, and DWATN, press
93     _RL zetaPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
94     _RL etaPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
95 mlosch 1.8 _RL etaZPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
96 mlosch 1.2 _RL dwatPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
97 mlosch 1.1 CEOP
98    
99     C Initialise
100 mlosch 1.5 newtonIter = 0
101     krylovFails = 0
102     totalKrylovItersLoc = 0
103     JFNKconverged = .FALSE.
104     JFNKtol = 0. _d 0
105     JFNKresidual = 0. _d 0
106     JFNKresidualKm1 = 0. _d 0
107     FGMRESeps = 0. _d 0
108     recip_deltaT = 1. _d 0 / SEAICE_deltaTdyn
109    
110     iOutFGMRES=0
111 mlosch 1.12 C with iOutFgmres=1, seaice_fgmres prints the residual at each iteration
112     IF ( debugLevel.GE.debLevC .AND.
113 mlosch 1.5 & DIFFERENT_MULTIPLE( SEAICE_monFreq, myTime, deltaTClock ) )
114     & iOutFGMRES=1
115    
116 mlosch 1.1 DO bj=myByLo(myThid),myByHi(myThid)
117     DO bi=myBxLo(myThid),myBxHi(myThid)
118 jmc 1.20 DO J=1-OLy,sNy+OLy
119     DO I=1-OLx,sNx+OLx
120 mlosch 1.1 uIceRes(I,J,bi,bj) = 0. _d 0
121     vIceRes(I,J,bi,bj) = 0. _d 0
122     duIce (I,J,bi,bj) = 0. _d 0
123     dvIce (I,J,bi,bj) = 0. _d 0
124     uIceNm1(I,J,bi,bj) = uIce(I,J,bi,bj)
125     vIceNm1(I,J,bi,bj) = vIce(I,J,bi,bj)
126     ENDDO
127     ENDDO
128     C Compute things that do no change during the Newton iteration:
129 jmc 1.20 C sea-surface tilt and wind stress:
130 mlosch 1.1 C FORCEX/Y0 - mass*(u/vIceNm1)/deltaT
131 jmc 1.20 DO J=1-OLy,sNy+OLy
132     DO I=1-OLx,sNx+OLx
133 mlosch 1.1 FORCEX(I,J,bi,bj) = FORCEX0(I,J,bi,bj)
134 jmc 1.20 & + seaiceMassU(I,J,bi,bj)*uIceNm1(I,J,bi,bj)*recip_deltaT
135 mlosch 1.1 FORCEY(I,J,bi,bj) = FORCEY0(I,J,bi,bj)
136 jmc 1.20 & + seaiceMassV(I,J,bi,bj)*vIceNm1(I,J,bi,bj)*recip_deltaT
137 mlosch 1.1 ENDDO
138     ENDDO
139     ENDDO
140     ENDDO
141     C Start nonlinear Newton iteration: outer loop iteration
142     DO WHILE ( newtonIter.LT.SEAICEnewtonIterMax .AND.
143     & .NOT.JFNKconverged )
144     newtonIter = newtonIter + 1
145     C Compute initial residual F(u), (includes computation of global
146     C variables DWATN, zeta, and eta)
147 jmc 1.20 IF ( newtonIter .EQ. 1 ) CALL SEAICE_JFNK_UPDATE(
148     I duIce, dvIce,
149 mlosch 1.15 U uIce, vIce, JFNKresidual,
150     O uIceRes, vIceRes,
151     I newtonIter, myTime, myIter, myThid )
152 mlosch 1.1 C local copies of precomputed coefficients that are to stay
153     C constant for the preconditioner
154     DO bj=myByLo(myThid),myByHi(myThid)
155     DO bi=myBxLo(myThid),myBxHi(myThid)
156 jmc 1.20 DO j=1-OLy,sNy+OLy
157     DO i=1-OLx,sNx+OLx
158 mlosch 1.10 zetaPre(I,J,bi,bj) = zeta(I,J,bi,bj)
159     etaPre(I,J,bi,bj) = eta(I,J,bi,bj)
160     etaZPre(I,J,bi,bj) = etaZ(I,J,bi,bj)
161     dwatPre(I,J,bi,bj) = DWATN(I,J,bi,bj)
162 mlosch 1.1 ENDDO
163     ENDDO
164     ENDDO
165     ENDDO
166     C compute convergence criterion for linear preconditioned FGMRES
167     JFNKgamma_lin = JFNKgamma_lin_max
168 mlosch 1.18 IF ( newtonIter.GT.1.AND.newtonIter.LE.SEAICE_JFNK_tolIter
169 mlosch 1.1 & .AND.JFNKresidual.LT.JFNKres_t ) THEN
170 mlosch 1.22 C Eisenstat and Walker (1996), eq.(2.6)
171     JFNKgamma_lin = SEAICE_JFNKphi
172     & *( JFNKresidual/JFNKresidualKm1 )**SEAICE_JFNKalpha
173 mlosch 1.1 JFNKgamma_lin = min(JFNKgamma_lin_max, JFNKgamma_lin)
174     JFNKgamma_lin = max(JFNKgamma_lin_min, JFNKgamma_lin)
175     ENDIF
176     C save the residual for the next iteration
177     JFNKresidualKm1 = JFNKresidual
178 jmc 1.20
179 mlosch 1.1 C The Krylov iteration using FGMRES, the preconditioner is LSOR
180     C for now. The code is adapted from SEAICE_LSR, but heavily stripped
181     C down.
182     C krylovIter is mapped into "its" in seaice_fgmres and is incremented
183     C in that routine
184     krylovIter = 0
185     iCode = 0
186 jmc 1.20
187 mlosch 1.1 JFNKconverged = JFNKresidual.LT.JFNKtol
188 jmc 1.20
189 mlosch 1.1 C do Krylov loop only if convergence is not reached
190 jmc 1.20
191 mlosch 1.1 IF ( .NOT.JFNKconverged ) THEN
192 jmc 1.20
193 mlosch 1.1 C start Krylov iteration (FGMRES)
194 jmc 1.20
195 mlosch 1.1 krylovConverged = .FALSE.
196     FGMRESeps = JFNKgamma_lin * JFNKresidual
197 jmc 1.20 DO WHILE ( .NOT.krylovConverged )
198 mlosch 1.1 C solution vector sol = du/vIce
199     C residual vector (rhs) Fu = u/vIceRes
200 jmc 1.20 C output work vectors wk1, -> input work vector wk2
201    
202 mlosch 1.1 CALL SEAICE_FGMRES_DRIVER(
203 jmc 1.20 I uIceRes, vIceRes,
204 mlosch 1.1 U duIce, dvIce, iCode,
205 mlosch 1.5 I FGMRESeps, iOutFGMRES,
206 mlosch 1.1 I newtonIter, krylovIter, myTime, myIter, myThid )
207     C FGMRES returns iCode either asking for an new preconditioned vector
208     C or product of matrix (Jacobian) times vector. For iCode = 0, terminate
209     C iteration
210     IF (iCode.EQ.1) THEN
211 jmc 1.20 C Call preconditioner
212 mlosch 1.7 IF ( SOLV_MAX_ITERS .GT. 0 )
213 jmc 1.20 & CALL SEAICE_PRECONDITIONER(
214     U duIce, dvIce,
215     I zetaPre, etaPre, etaZpre, dwatPre,
216 mlosch 1.1 I newtonIter, krylovIter, myTime, myIter, myThid )
217     ELSEIF (iCode.GE.2) THEN
218     C Compute Jacobian times vector
219     CALL SEAICE_JACVEC(
220     I uIce, vIce, uIceRes, vIceRes,
221 jmc 1.20 U duIce, dvIce,
222 mlosch 1.1 I newtonIter, krylovIter, myTime, myIter, myThid )
223     ENDIF
224     krylovConverged = iCode.EQ.0
225     C End of Krylov iterate
226     ENDDO
227 mlosch 1.5 totalKrylovItersLoc = totalKrylovItersLoc + krylovIter
228 mlosch 1.1 C some output diagnostics
229     IF ( debugLevel.GE.debLevA ) THEN
230 mlosch 1.5 _BEGIN_MASTER( myThid )
231 jmc 1.20 totalNewtonItersLoc =
232 mlosch 1.13 & SEAICEnewtonIterMax*(myIter-nIter0)+newtonIter
233 jmc 1.20 WRITE(msgBuf,'(2A,2(1XI6),2E12.5)')
234 mlosch 1.13 & ' S/R SEAICE_JFNK: Newton iterate / total, ',
235     & 'JFNKgamma_lin, initial norm = ',
236     & newtonIter, totalNewtonItersLoc,
237     & JFNKgamma_lin,JFNKresidual
238     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
239     & SQUEEZE_RIGHT, myThid )
240 mlosch 1.1 WRITE(msgBuf,'(3(A,I6))')
241 jmc 1.20 & ' S/R SEAICE_JFNK: Newton iterate / total = ',newtonIter,
242 mlosch 1.13 & ' / ', totalNewtonItersLoc,
243 mlosch 1.1 & ', Nb. of FGMRES iterations = ', krylovIter
244     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
245     & SQUEEZE_RIGHT, myThid )
246 mlosch 1.5 _END_MASTER( myThid )
247 mlosch 1.1 ENDIF
248     IF ( krylovIter.EQ.SEAICEkrylovIterMax ) THEN
249 mlosch 1.5 krylovFails = krylovFails + 1
250 mlosch 1.1 ENDIF
251 mlosch 1.17 C Set the stopping criterion for the Newton iteration and the
252     C criterion for the transition from accurate to approximate FGMRES
253     IF ( newtonIter .EQ. 1 ) THEN
254     JFNKtol=JFNKgamma_nonlin*JFNKresidual
255     IF ( JFNKres_tFac .NE. UNSET_RL )
256     & JFNKres_t = JFNKresidual * JFNKres_tFac
257     ENDIF
258 mlosch 1.1 C Update linear solution vector and return to Newton iteration
259 mlosch 1.15 C Do a linesearch if necessary, and compute a new residual.
260     C Note that it should be possible to do the following operations
261     C at the beginning of the Newton iteration, thereby saving us from
262     C the extra call of seaice_jfnk_update, but unfortunately that
263     C changes the results, so we leave the stuff here for now.
264 jmc 1.20 CALL SEAICE_JFNK_UPDATE(
265     I duIce, dvIce,
266 mlosch 1.15 U uIce, vIce, JFNKresidual,
267     O uIceRes, vIceRes,
268     I newtonIter, myTime, myIter, myThid )
269     C reset du/vIce here instead of setting sol = 0 in seaice_fgmres_driver
270 mlosch 1.1 DO bj=myByLo(myThid),myByHi(myThid)
271     DO bi=myBxLo(myThid),myBxHi(myThid)
272 jmc 1.20 DO J=1-OLy,sNy+OLy
273     DO I=1-OLx,sNx+OLx
274 mlosch 1.4 duIce(I,J,bi,bj)= 0. _d 0
275     dvIce(I,J,bi,bj)= 0. _d 0
276 mlosch 1.1 ENDDO
277     ENDDO
278     ENDDO
279     ENDDO
280     ENDIF
281     C end of Newton iterate
282     ENDDO
283 jmc 1.20
284 mlosch 1.5 C-- Output diagnostics
285 jmc 1.20
286 mlosch 1.6 IF ( SEAICE_monFreq .GT. 0. _d 0 ) THEN
287 mlosch 1.5 C Count iterations
288 mlosch 1.6 totalJFNKtimeSteps = totalJFNKtimeSteps + 1
289     totalNewtonIters = totalNewtonIters + newtonIter
290     totalKrylovIters = totalKrylovIters + totalKrylovItersLoc
291 mlosch 1.5 C Record failure
292 mlosch 1.6 totalKrylovFails = totalKrylovFails + krylovFails
293     IF ( newtonIter .EQ. SEAICEnewtonIterMax ) THEN
294 jmc 1.20 totalNewtonFails = totalNewtonFails + 1
295 mlosch 1.6 ENDIF
296 mlosch 1.5 ENDIF
297     C Decide whether it is time to dump and reset the counter
298 mlosch 1.9 writeNow = DIFFERENT_MULTIPLE(SEAICE_monFreq,
299 jmc 1.20 & myTime+deltaTClock, deltaTClock)
300 mlosch 1.9 #ifdef ALLOW_CAL
301     IF ( useCAL ) THEN
302 jmc 1.20 CALL CAL_TIME2DUMP(
303 mlosch 1.9 I zeroRL, SEAICE_monFreq, deltaTClock,
304     U writeNow,
305     I myTime+deltaTclock, myIter+1, myThid )
306     ENDIF
307     #endif
308     IF ( writeNow ) THEN
309 mlosch 1.5 _BEGIN_MASTER( myThid )
310 jmc 1.20 WRITE(msgBuf,'(A)')
311 mlosch 1.5 &' // ======================================================='
312     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
313     & SQUEEZE_RIGHT, myThid )
314     WRITE(msgBuf,'(A)') ' // Begin JFNK statistics'
315     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
316     & SQUEEZE_RIGHT, myThid )
317 jmc 1.20 WRITE(msgBuf,'(A)')
318 mlosch 1.5 &' // ======================================================='
319     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
320     & SQUEEZE_RIGHT, myThid )
321 jmc 1.20 WRITE(msgBuf,'(A,I10)')
322 mlosch 1.5 & ' %JFNK_MON: time step = ', myIter+1
323     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
324     & SQUEEZE_RIGHT, myThid )
325 jmc 1.20 WRITE(msgBuf,'(A,I10)')
326 mlosch 1.5 & ' %JFNK_MON: Nb. of time steps = ', totalJFNKtimeSteps
327     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
328     & SQUEEZE_RIGHT, myThid )
329 jmc 1.20 WRITE(msgBuf,'(A,I10)')
330 mlosch 1.5 & ' %JFNK_MON: Nb. of Newton steps = ', totalNewtonIters
331     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
332     & SQUEEZE_RIGHT, myThid )
333 jmc 1.20 WRITE(msgBuf,'(A,I10)')
334 mlosch 1.5 & ' %JFNK_MON: Nb. of Krylov steps = ', totalKrylovIters
335     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
336     & SQUEEZE_RIGHT, myThid )
337 jmc 1.20 WRITE(msgBuf,'(A,I10)')
338 mlosch 1.5 & ' %JFNK_MON: Nb. of Newton failures = ', totalNewtonFails
339     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
340     & SQUEEZE_RIGHT, myThid )
341 jmc 1.20 WRITE(msgBuf,'(A,I10)')
342 mlosch 1.5 & ' %JFNK_MON: Nb. of Krylov failures = ', totalKrylovFails
343     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
344     & SQUEEZE_RIGHT, myThid )
345 jmc 1.20 WRITE(msgBuf,'(A)')
346 mlosch 1.5 &' // ======================================================='
347     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
348     & SQUEEZE_RIGHT, myThid )
349 mlosch 1.11 WRITE(msgBuf,'(A)') ' // End JFNK statistics'
350 mlosch 1.5 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
351     & SQUEEZE_RIGHT, myThid )
352 jmc 1.20 WRITE(msgBuf,'(A)')
353 mlosch 1.5 &' // ======================================================='
354     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
355     & SQUEEZE_RIGHT, myThid )
356     _END_MASTER( myThid )
357     C reset and start again
358     totalJFNKtimeSteps = 0
359     totalNewtonIters = 0
360     totalKrylovIters = 0
361     totalKrylovFails = 0
362     totalNewtonFails = 0
363     ENDIF
364    
365     C Print more debugging information
366 mlosch 1.1 IF ( debugLevel.GE.debLevA ) THEN
367     IF ( newtonIter .EQ. SEAICEnewtonIterMax ) THEN
368 mlosch 1.5 _BEGIN_MASTER( myThid )
369 jmc 1.20 WRITE(msgBuf,'(A,I10)')
370 mlosch 1.1 & ' S/R SEAICE_JFNK: JFNK did not converge in timestep ',
371 mlosch 1.5 & myIter+1
372 mlosch 1.1 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
373     & SQUEEZE_RIGHT, myThid )
374 mlosch 1.5 _END_MASTER( myThid )
375 mlosch 1.1 ENDIF
376 mlosch 1.5 IF ( krylovFails .GT. 0 ) THEN
377     _BEGIN_MASTER( myThid )
378 jmc 1.20 WRITE(msgBuf,'(A,I4,A,I10)')
379 mlosch 1.1 & ' S/R SEAICE_JFNK: FGMRES did not converge ',
380 mlosch 1.5 & krylovFails, ' times in timestep ', myIter+1
381 mlosch 1.1 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
382     & SQUEEZE_RIGHT, myThid )
383 mlosch 1.5 _END_MASTER( myThid )
384 mlosch 1.1 ENDIF
385 mlosch 1.5 _BEGIN_MASTER( myThid )
386 jmc 1.20 WRITE(msgBuf,'(A,I6,A,I10)')
387 mlosch 1.1 & ' S/R SEAICE_JFNK: Total number FGMRES iterations = ',
388 mlosch 1.5 & totalKrylovItersLoc, ' in timestep ', myIter+1
389     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
390     & SQUEEZE_RIGHT, myThid )
391     _END_MASTER( myThid )
392 mlosch 1.1 ENDIF
393    
394 mlosch 1.15 RETURN
395     END
396    
397 mlosch 1.16 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
398 mlosch 1.15 CBOP
399     C !ROUTINE: SEAICE_JFNK_UPDATE
400     C !INTERFACE:
401    
402 jmc 1.20 SUBROUTINE SEAICE_JFNK_UPDATE(
403     I duIce, dvIce,
404 mlosch 1.15 U uIce, vIce, JFNKresidual,
405     O uIceRes, vIceRes,
406     I newtonIter, myTime, myIter, myThid )
407    
408     C !DESCRIPTION: \bv
409     C *==========================================================*
410     C | SUBROUTINE SEAICE_JFNK_UPDATE
411     C | o Update velocities with incremental solutions of FGMRES
412     C | o compute residual of updated solutions and do
413     C | o linesearch:
414     C | reduce update until residual is smaller than previous
415     C | one (input)
416     C *==========================================================*
417     C | written by Martin Losch, Jan 2013
418     C *==========================================================*
419     C \ev
420    
421     C !USES:
422     IMPLICIT NONE
423    
424     C === Global variables ===
425     #include "SIZE.h"
426     #include "EEPARAMS.h"
427     #include "PARAMS.h"
428     #include "SEAICE_SIZE.h"
429     #include "SEAICE_PARAMS.h"
430    
431     C !INPUT/OUTPUT PARAMETERS:
432     C === Routine arguments ===
433     C myTime :: Simulation time
434     C myIter :: Simulation timestep number
435     C myThid :: my Thread Id. number
436     C newtonIter :: current iterate of Newton iteration
437     _RL myTime
438     INTEGER myIter
439     INTEGER myThid
440     INTEGER newtonIter
441     C JFNKresidual :: Residual at the beginning of the FGMRES iteration,
442     C changes with newtonIter (updated)
443     _RL JFNKresidual
444     C du/vIce :: ice velocity increment to be added to u/vIce (input)
445     _RL duIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
446     _RL dvIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
447     C u/vIce :: ice velocity increment to be added to u/vIce (updated)
448     _RL uIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
449     _RL vIce (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
450     C u/vIceRes :: residual of sea-ice momentum equations (output)
451     _RL uIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
452     _RL vIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
453    
454 mlosch 1.16 C !LOCAL VARIABLES:
455     C === Local variables ===
456 mlosch 1.15 C i,j,bi,bj :: loop indices
457     INTEGER i,j,bi,bj
458     INTEGER l
459     _RL resLoc, facLS
460     LOGICAL doLineSearch
461     C nVec :: size of the input vector(s)
462 jmc 1.20 C resTmp :: vector version of the residuals
463 mlosch 1.15 INTEGER nVec
464     PARAMETER ( nVec = 2*sNx*sNy )
465     _RL resTmp (nVec,1,nSx,nSy)
466 jmc 1.20
467 mlosch 1.15 CHARACTER*(MAX_LEN_MBUF) msgBuf
468     CEOP
469    
470     C Initialise some local variables
471     l = 0
472     resLoc = JFNKresidual
473     facLS = 1. _d 0
474     doLineSearch = .TRUE.
475     DO WHILE ( doLineSearch )
476     C Create update
477     DO bj=myByLo(myThid),myByHi(myThid)
478     DO bi=myBxLo(myThid),myBxHi(myThid)
479 jmc 1.20 DO J=1-OLy,sNy+OLy
480     DO I=1-OLx,sNx+OLx
481 mlosch 1.15 uIce(I,J,bi,bj) = uIce(I,J,bi,bj)+facLS*duIce(I,J,bi,bj)
482     vIce(I,J,bi,bj) = vIce(I,J,bi,bj)+facLS*dvIce(I,J,bi,bj)
483     ENDDO
484     ENDDO
485     ENDDO
486     ENDDO
487     C Compute current residual F(u), (includes re-computation of global
488     C variables DWATN, zeta, and eta, i.e. they are different after this)
489 jmc 1.20 CALL SEAICE_CALC_RESIDUAL(
490     I uIce, vIce,
491     O uIceRes, vIceRes,
492 mlosch 1.15 I newtonIter, 0, myTime, myIter, myThid )
493     C Important: Compute the norm of the residual using the same scalar
494     C product that SEAICE_FGMRES does
495     CALL SEAICE_MAP2VEC(nVec,uIceRes,vIceRes,resTmp,.TRUE.,myThid)
496     CALL SEAICE_SCALPROD(nVec,1,1,1,resTmp,resTmp,resLoc,myThid)
497     resLoc = SQRT(resLoc)
498 mlosch 1.19 C Determine, if we need more iterations
499 jmc 1.20 doLineSearch = resLoc .GE. JFNKresidual
500 mlosch 1.19 C Limit the maximum number of iterations arbitrarily to four
501 jmc 1.20 doLineSearch = doLineSearch .AND. l .LT. 4
502 mlosch 1.19 C For the first iteration du/vIce = 0 and there will be no
503     C improvement of the residual possible, so we do only the first
504     C iteration
505     IF ( newtonIter .EQ. 1 ) doLineSearch = .FALSE.
506     C Only start a linesearch after some Newton iterations
507     IF ( newtonIter .LE. SEAICE_JFNK_lsIter ) doLineSearch = .FALSE.
508     C Increment counter
509     l = l + 1
510 mlosch 1.15 C some output diagnostics
511     IF ( debugLevel.GE.debLevA .AND. doLineSearch ) THEN
512     _BEGIN_MASTER( myThid )
513 jmc 1.20 WRITE(msgBuf,'(2A,2(1XI6),3E12.5)')
514 mlosch 1.15 & ' S/R SEAICE_JFNK_UPDATE: Newton iter, LSiter, ',
515     & 'facLS, JFNKresidual, resLoc = ',
516     & newtonIter, l, facLS, JFNKresidual, resLoc
517     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
518     & SQUEEZE_RIGHT, myThid )
519     _END_MASTER( myThid )
520     ENDIF
521     C Get ready for the next iteration: after adding du/vIce in the first
522     C iteration, we substract 0.5*du/vIce from u/vIce in the next
523     C iterations, 0.25*du/vIce in the second, etc.
524     facLS = - 0.5 _d 0 * ABS(facLS)
525     ENDDO
526     C This is the new residual
527     JFNKresidual = resLoc
528    
529 mlosch 1.21 #endif /* SEAICE_ALLOW_JFNK */
530 mlosch 1.1
531     RETURN
532     END

  ViewVC Help
Powered by ViewVC 1.1.22