/[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.14 by mlosch, Fri Jan 4 15:48:09 2013 UTC revision 1.15 by mlosch, Wed Jan 16 21:20:28 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 142  C     Start nonlinear Newton iteration: Line 147  C     Start nonlinear Newton iteration:
147         newtonIter = newtonIter + 1         newtonIter = newtonIter + 1
148  C     Compute initial residual F(u), (includes computation of global  C     Compute initial residual F(u), (includes computation of global
149  C     variables DWATN, zeta, and eta)  C     variables DWATN, zeta, and eta)
150         CALL SEAICE_CALC_RESIDUAL(  C     Update linear solution vector and return to Newton iteration
151       I      uIce, vIce,  C     Do the linesearch
152       O      uIceRes, vIceRes,         CALL SEAICE_JFNK_UPDATE(
153       I      newtonIter, 0, myTime, myIter, myThid )       I      duIce, dvIce,
154  C     probably not necessary, will be removed later:       U      uIce, vIce, JFNKresidual,
155         CALL EXCH_UV_XY_RL( uIceRes, vIceRes,.TRUE.,myThid)       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  C     local copies of precomputed coefficients that are to stay  C     local copies of precomputed coefficients that are to stay
182  C     constant for the preconditioner  C     constant for the preconditioner
183         DO bj=myByLo(myThid),myByHi(myThid)         DO bj=myByLo(myThid),myByHi(myThid)
# Line 162  C     constant for the preconditioner Line 192  C     constant for the preconditioner
192           ENDDO           ENDDO
193          ENDDO          ENDDO
194         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)  
195  C     compute convergence criterion for linear preconditioned FGMRES  C     compute convergence criterion for linear preconditioned FGMRES
196         JFNKgamma_lin = JFNKgamma_lin_max         JFNKgamma_lin = JFNKgamma_lin_max
197         IF ( newtonIter.GT.1.AND.newtonIter.LE.100         IF ( newtonIter.GT.1.AND.newtonIter.LE.100
# Line 254  C     some output diagnostics Line 278  C     some output diagnostics
278          IF ( krylovIter.EQ.SEAICEkrylovIterMax ) THEN          IF ( krylovIter.EQ.SEAICEkrylovIterMax ) THEN
279           krylovFails = krylovFails + 1           krylovFails = krylovFails + 1
280          ENDIF          ENDIF
281    C     Set the stopping criterion for the Newton iteration
282            IF ( newtonIter .EQ. 1 ) JFNKtol=JFNKgamma_nonlin*JFNKresidual
283  C     Update linear solution vector and return to Newton iteration  C     Update linear solution vector and return to Newton iteration
284    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          DO bj=myByLo(myThid),myByHi(myThid)          DO bj=myByLo(myThid),myByHi(myThid)
296           DO bi=myBxLo(myThid),myBxHi(myThid)           DO bi=myBxLo(myThid),myBxHi(myThid)
297            DO J=1-Oly,sNy+Oly            DO J=1-Oly,sNy+Oly
298             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  
299              duIce(I,J,bi,bj)= 0. _d 0              duIce(I,J,bi,bj)= 0. _d 0
300              dvIce(I,J,bi,bj)= 0. _d 0              dvIce(I,J,bi,bj)= 0. _d 0
301             ENDDO             ENDDO
302            ENDDO            ENDDO
303           ENDDO           ENDDO
304          ENDDO          ENDDO
 C     Set the stopping criterion for the Newton iteration  
         IF ( newtonIter .EQ. 1 ) JFNKtol=JFNKgamma_nonlin*JFNKresidual  
305         ENDIF         ENDIF
306  C     end of Newton iterate  C     end of Newton iterate
307        ENDDO        ENDDO
# Line 384  C     Print more debugging information Line 416  C     Print more debugging information
416         _END_MASTER( myThid )         _END_MASTER( myThid )
417        ENDIF        ENDIF
418    
419          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  #endif /* SEAICE_ALLOW_DYNAMICS and SEAICE_CGRID and SEAICE_ALLOW_JFNK */  #endif /* SEAICE_ALLOW_DYNAMICS and SEAICE_CGRID and SEAICE_ALLOW_JFNK */
552    
553        RETURN        RETURN

Legend:
Removed from v.1.14  
changed lines
  Added in v.1.15

  ViewVC Help
Powered by ViewVC 1.1.22