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

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

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

revision 1.13 by mlosch, Mon Dec 17 15:06:02 2012 UTC revision 1.17 by mlosch, Thu Jan 17 10:42:43 2013 UTC
# Line 3  C $Name$ Line 3  C $Name$
3    
4  #include "SEAICE_OPTIONS.h"  #include "SEAICE_OPTIONS.h"
5    
6    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  CBOP  CBOP
12  C     !ROUTINE: SEAICE_JFNK  C     !ROUTINE: SEAICE_JFNK
13  C     !INTERFACE:  C     !INTERFACE:
# Line 10  C     !INTERFACE: Line 15  C     !INTERFACE:
15    
16  C     !DESCRIPTION: \bv  C     !DESCRIPTION: \bv
17  C     *==========================================================*  C     *==========================================================*
18  C     | SUBROUTINE SEAICE_JFKF  C     | SUBROUTINE SEAICE_JFNK
19  C     | o Ice dynamics using a Jacobian-free Newton-Krylov solver  C     | o Ice dynamics using a Jacobian-free Newton-Krylov solver
20  C     |   following J.-F. Lemieux et al. Improving the numerical  C     |   following J.-F. Lemieux et al. Improving the numerical
21  C     |   convergence of viscous-plastic sea ice models with the  C     |   convergence of viscous-plastic sea ice models with the
# Line 55  C     !FUNCTIONS: Line 60  C     !FUNCTIONS:
60        LOGICAL  DIFFERENT_MULTIPLE        LOGICAL  DIFFERENT_MULTIPLE
61        EXTERNAL DIFFERENT_MULTIPLE        EXTERNAL DIFFERENT_MULTIPLE
62    
63    C     !LOCAL VARIABLES:
64    C     === Local variables ===
65  C     i,j,bi,bj :: loop indices  C     i,j,bi,bj :: loop indices
66        INTEGER i,j,bi,bj        INTEGER i,j,bi,bj
67  C     loop indices  C     loop indices
# Line 80  C Line 87  C
87  C     u/vIceRes :: residual of sea-ice momentum equations  C     u/vIceRes :: residual of sea-ice momentum equations
88        _RL uIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _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)        _RL vIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
90    C     vector version of the residuals
91          _RL resTmp (nVec,1,nSx,nSy)
92  C     du/vIce   :: ice velocity increment to be added to u/vIce  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)        _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)        _RL dvIce  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
# Line 90  C     zeta, eta, and DWATN, press Line 99  C     zeta, eta, and DWATN, press
99        _RL etaZPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL etaZPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
100        _RL dwatPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL dwatPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
101  CEOP  CEOP
       INTEGER n  
       _RL resTmp (nVec,1,nSx,nSy)  
102    
103  C     Initialise  C     Initialise
104        newtonIter          = 0        newtonIter          = 0
# Line 142  C     Start nonlinear Newton iteration: Line 149  C     Start nonlinear Newton iteration:
149         newtonIter = newtonIter + 1         newtonIter = newtonIter + 1
150  C     Compute initial residual F(u), (includes computation of global  C     Compute initial residual F(u), (includes computation of global
151  C     variables DWATN, zeta, and eta)  C     variables DWATN, zeta, and eta)
152         CALL SEAICE_CALC_RESIDUAL(         IF ( newtonIter .EQ. 1 ) CALL SEAICE_JFNK_UPDATE(
153       I      uIce, vIce,       I      duIce, dvIce,
154       O      uIceRes, vIceRes,       U      uIce, vIce, JFNKresidual,
155       I      newtonIter, 0, myTime, myIter, myThid )       O      uIceRes, vIceRes,
156  C     probably not necessary, will be removed later:       I      newtonIter, myTime, myIter, myThid )
        CALL EXCH_UV_XY_RL( uIceRes, vIceRes,.TRUE.,myThid)  
157  C     local copies of precomputed coefficients that are to stay  C     local copies of precomputed coefficients that are to stay
158  C     constant for the preconditioner  C     constant for the preconditioner
159         DO bj=myByLo(myThid),myByHi(myThid)         DO bj=myByLo(myThid),myByHi(myThid)
# Line 162  C     constant for the preconditioner Line 168  C     constant for the preconditioner
168           ENDDO           ENDDO
169          ENDDO          ENDDO
170         ENDDO         ENDDO
 C     Important: Compute the norm of the residual using the same scalar  
 C     product that SEAICE_FGMRES does  
        CALL SEAICE_MAP2VEC(nVec,uIceRes,vIceRes,resTmp,.TRUE.,myThid)  
        CALL SEAICE_SCALPROD(  
      &      nVec,1,1,1,resTmp,resTmp,JFNKresidual,myThid)  
        JFNKresidual = SQRT(JFNKresidual)  
171  C     compute convergence criterion for linear preconditioned FGMRES  C     compute convergence criterion for linear preconditioned FGMRES
172         JFNKgamma_lin = JFNKgamma_lin_max         JFNKgamma_lin = JFNKgamma_lin_max
173         IF ( newtonIter.GT.1.AND.newtonIter.LE.100         IF ( newtonIter.GT.1.AND.newtonIter.LE.100
# Line 254  C     some output diagnostics Line 254  C     some output diagnostics
254          IF ( krylovIter.EQ.SEAICEkrylovIterMax ) THEN          IF ( krylovIter.EQ.SEAICEkrylovIterMax ) THEN
255           krylovFails = krylovFails + 1           krylovFails = krylovFails + 1
256          ENDIF          ENDIF
257    C     Set the stopping criterion for the Newton iteration and the
258    C     criterion for the transition from accurate to approximate FGMRES
259            IF ( newtonIter .EQ. 1 ) THEN
260             JFNKtol=JFNKgamma_nonlin*JFNKresidual
261             IF ( JFNKres_tFac .NE. UNSET_RL )
262         &        JFNKres_t = JFNKresidual * JFNKres_tFac
263            ENDIF
264  C     Update linear solution vector and return to Newton iteration  C     Update linear solution vector and return to Newton iteration
265    C     Do a linesearch if necessary, and compute a new residual.
266    C     Note that it should be possible to do the following operations
267    C     at the beginning of the Newton iteration, thereby saving us from
268    C     the extra call of seaice_jfnk_update, but unfortunately that
269    C     changes the results, so we leave the stuff here for now.
270            CALL SEAICE_JFNK_UPDATE(
271         I       duIce, dvIce,
272         U       uIce, vIce, JFNKresidual,
273         O       uIceRes, vIceRes,
274         I       newtonIter, myTime, myIter, myThid )
275    C     reset du/vIce here instead of setting sol = 0 in seaice_fgmres_driver
276          DO bj=myByLo(myThid),myByHi(myThid)          DO bj=myByLo(myThid),myByHi(myThid)
277           DO bi=myBxLo(myThid),myBxHi(myThid)           DO bi=myBxLo(myThid),myBxHi(myThid)
278            DO J=1-Oly,sNy+Oly            DO J=1-Oly,sNy+Oly
279             DO I=1-Olx,sNx+Olx             DO I=1-Olx,sNx+Olx
             uIce(I,J,bi,bj) = uIce(I,J,bi,bj)+duIce(I,J,bi,bj)  
             vIce(I,J,bi,bj) = vIce(I,J,bi,bj)+dvIce(I,J,bi,bj)  
 C     reset du/vIce here instead of setting sol = 0 in seaice_fgmres_driver  
280              duIce(I,J,bi,bj)= 0. _d 0              duIce(I,J,bi,bj)= 0. _d 0
281              dvIce(I,J,bi,bj)= 0. _d 0              dvIce(I,J,bi,bj)= 0. _d 0
282             ENDDO             ENDDO
283            ENDDO            ENDDO
284           ENDDO           ENDDO
285          ENDDO          ENDDO
 C     Set the stopping criterion for the Newton iteration  
         IF ( newtonIter .EQ. 1 ) JFNKtol=JFNKgamma_nonlin*JFNKresidual  
286         ENDIF         ENDIF
287  C     end of Newton iterate  C     end of Newton iterate
288        ENDDO        ENDDO
# Line 384  C     Print more debugging information Line 397  C     Print more debugging information
397         _END_MASTER( myThid )         _END_MASTER( myThid )
398        ENDIF        ENDIF
399    
400          RETURN
401          END
402    
403    C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
404    CBOP
405    C     !ROUTINE: SEAICE_JFNK_UPDATE
406    C     !INTERFACE:
407    
408          SUBROUTINE SEAICE_JFNK_UPDATE(
409         I     duIce, dvIce,
410         U     uIce, vIce, JFNKresidual,
411         O     uIceRes, vIceRes,
412         I     newtonIter, myTime, myIter, myThid )
413    
414    C     !DESCRIPTION: \bv
415    C     *==========================================================*
416    C     | SUBROUTINE SEAICE_JFNK_UPDATE
417    C     | o Update velocities with incremental solutions of FGMRES
418    C     | o compute residual of updated solutions and do
419    C     | o linesearch:
420    C     |   reduce update until residual is smaller than previous
421    C     |   one (input)
422    C     *==========================================================*
423    C     | written by Martin Losch, Jan 2013
424    C     *==========================================================*
425    C     \ev
426    
427    C     !USES:
428          IMPLICIT NONE
429    
430    C     === Global variables ===
431    #include "SIZE.h"
432    #include "EEPARAMS.h"
433    #include "PARAMS.h"
434    #include "SEAICE_SIZE.h"
435    #include "SEAICE_PARAMS.h"
436    
437    C     !INPUT/OUTPUT PARAMETERS:
438    C     === Routine arguments ===
439    C     myTime :: Simulation time
440    C     myIter :: Simulation timestep number
441    C     myThid :: my Thread Id. number
442    C     newtonIter :: current iterate of Newton iteration
443          _RL     myTime
444          INTEGER myIter
445          INTEGER myThid
446          INTEGER newtonIter
447    C     JFNKresidual :: Residual at the beginning of the FGMRES iteration,
448    C                     changes with newtonIter (updated)
449          _RL     JFNKresidual
450    C     du/vIce   :: ice velocity increment to be added to u/vIce (input)
451          _RL duIce  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
452          _RL dvIce  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
453    C     u/vIce    :: ice velocity increment to be added to u/vIce (updated)
454          _RL uIce   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
455          _RL vIce   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
456    C     u/vIceRes :: residual of sea-ice momentum equations (output)
457          _RL uIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
458          _RL vIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
459    
460    C     !LOCAL VARIABLES:
461    C     === Local variables ===
462    C     i,j,bi,bj :: loop indices
463          INTEGER i,j,bi,bj
464          INTEGER l
465          _RL     resLoc, facLS
466          LOGICAL doLineSearch
467    C     nVec    :: size of the input vector(s)
468    C     vector version of the residuals
469          INTEGER nVec
470          PARAMETER ( nVec  = 2*sNx*sNy )
471          _RL resTmp (nVec,1,nSx,nSy)
472    C
473          CHARACTER*(MAX_LEN_MBUF) msgBuf
474    CEOP
475    
476    C     Initialise some local variables
477          l = 0
478          resLoc = JFNKresidual
479          facLS = 1. _d 0
480          doLineSearch = .TRUE.
481          DO WHILE ( doLineSearch )
482    C     Determine, if we need more iterations
483           doLineSearch = resLoc .GE. JFNKresidual
484    C     Limit the maximum number of iterations arbitrarily to four
485           doLineSearch = doLineSearch .AND. l .LE. 4
486    C     For the first iteration du/vIce = 0 and there will be no
487    C     improvement of the residual possible, so we do only the first
488    C     iteration
489           IF ( newtonIter .EQ. 1 ) doLineSearch = .FALSE.
490    C     Only start a linesearch after some Newton iterations
491           IF ( newtonIter .LE. SEAICE_JFNK_lsIter ) doLineSearch = .FALSE.
492    C     Increment counter
493           l = l + 1
494    C     Create update
495           DO bj=myByLo(myThid),myByHi(myThid)
496            DO bi=myBxLo(myThid),myBxHi(myThid)
497             DO J=1-Oly,sNy+Oly
498              DO I=1-Olx,sNx+Olx
499               uIce(I,J,bi,bj) = uIce(I,J,bi,bj)+facLS*duIce(I,J,bi,bj)
500               vIce(I,J,bi,bj) = vIce(I,J,bi,bj)+facLS*dvIce(I,J,bi,bj)
501              ENDDO
502             ENDDO
503            ENDDO
504           ENDDO
505    C     Compute current residual F(u), (includes re-computation of global
506    C     variables DWATN, zeta, and eta, i.e. they are different after this)
507           CALL SEAICE_CALC_RESIDUAL(
508         I      uIce, vIce,
509         O      uIceRes, vIceRes,
510         I      newtonIter, 0, myTime, myIter, myThid )
511    C     Important: Compute the norm of the residual using the same scalar
512    C     product that SEAICE_FGMRES does
513           CALL SEAICE_MAP2VEC(nVec,uIceRes,vIceRes,resTmp,.TRUE.,myThid)
514           CALL SEAICE_SCALPROD(nVec,1,1,1,resTmp,resTmp,resLoc,myThid)
515           resLoc = SQRT(resLoc)
516    C     some output diagnostics
517           IF ( debugLevel.GE.debLevA .AND. doLineSearch ) THEN
518            _BEGIN_MASTER( myThid )
519            WRITE(msgBuf,'(2A,2(1XI6),3E12.5)')
520         &       ' S/R SEAICE_JFNK_UPDATE: Newton iter, LSiter, ',
521         &       'facLS, JFNKresidual, resLoc = ',
522         &        newtonIter, l, facLS, JFNKresidual, resLoc
523            CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
524         &       SQUEEZE_RIGHT, myThid )
525            _END_MASTER( myThid )
526           ENDIF
527    C     Get ready for the next iteration: after adding du/vIce in the first
528    C     iteration, we substract 0.5*du/vIce from u/vIce in the next
529    C     iterations, 0.25*du/vIce in the second, etc.
530           facLS = - 0.5 _d 0 * ABS(facLS)
531          ENDDO
532    C     This is the new residual
533          JFNKresidual = resLoc
534    
535  #endif /* SEAICE_ALLOW_DYNAMICS and SEAICE_CGRID and SEAICE_ALLOW_JFNK */  #endif /* SEAICE_ALLOW_DYNAMICS and SEAICE_CGRID and SEAICE_ALLOW_JFNK */
536    
537        RETURN        RETURN

Legend:
Removed from v.1.13  
changed lines
  Added in v.1.17

  ViewVC Help
Powered by ViewVC 1.1.22