/[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.20 - (hide annotations) (download)
Sat Mar 2 04:35:05 2013 UTC (12 years, 4 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64e, checkpoint64f
Changes since 1.19: +68 -71 lines
remove unused variable

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

  ViewVC Help
Powered by ViewVC 1.1.22