/[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.16 - (hide annotations) (download)
Thu Jan 17 08:51:15 2013 UTC (12 years, 6 months ago) by mlosch
Branch: MAIN
Changes since 1.15: +8 -29 lines
fix previous check-in:
- remove superfluous code and commented code
- add comments
- call new s/r seaice_jfnk_update only for newtonIter = 1 at the
beginning of the Newton loop

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

  ViewVC Help
Powered by ViewVC 1.1.22