/[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.15 - (hide annotations) (download)
Wed Jan 16 21:20:28 2013 UTC (12 years, 6 months ago) by mlosch
Branch: MAIN
Changes since 1.14: +183 -19 lines
add a line search option to the JFNK solver

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

  ViewVC Help
Powered by ViewVC 1.1.22