New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
Changeset 508 for trunk/NEMO/OPA_SRC/DYN/dynspg_rl.F90 – NEMO

Ignore:
Timestamp:
2006-10-03T17:58:55+02:00 (18 years ago)
Author:
opalod
Message:

nemo_v1_update_071:RB: add iom for restart and reorganization of restart

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMO/OPA_SRC/DYN/dynspg_rl.F90

    r474 r508  
    44   !! Ocean dynamics:  surface pressure gradient trend 
    55   !!====================================================================== 
     6   !! History :  7.0  !  96-05  (G. Madec, M. Imbard, M. Guyon)  rewitting in 1 
     7   !!                           routine, without pointers, and with the same matrix 
     8   !!                           for sor and pcg, mpp exchange, and symmetric conditions 
     9   !!            " "  !  96-07  (A. Weaver)  Euler forward step 
     10   !!            " "  !  96-11  (A. Weaver)  correction to preconditioning: 
     11   !!            8.0  !  98-02  (M. Guyon)  FETI method 
     12   !!            " "  !  97-09  (J.-M. Molines)  Open boundaries 
     13   !!            8.5  !  02-08  (G. Madec)  F90: Free form and module 
     14   !!                 !  02-11  (C. Talandier, A-M Treguier) Open boundaries 
     15   !!            9.0  !  04-08  (C. Talandier)  New trends organization 
     16   !!            " "  !  05-11  (V. Garnier) Surface pressure gradient organization 
     17   !!            " "  !  06-07  (S. Masson)  distributed restart using iom 
     18   !!--------------------------------------------------------------------- 
    619#if   defined key_dynspg_rl   ||   defined key_esopa 
    720   !!---------------------------------------------------------------------- 
     
    1023   !!   dyn_spg_rl   : update the momentum trend with the surface pressure 
    1124   !!                  for the rigid-lid case. 
    12    !!---------------------------------------------------------------------- 
    13    !! * Modules used 
     25   !!   rl_rst       : read/write the rigid-lid restart fields in the ocean restart file 
     26   !!---------------------------------------------------------------------- 
    1427   USE oce             ! ocean dynamics and tracers 
    1528   USE dom_oce         ! ocean space and time domain 
     
    1932   USE zdf_oce         ! ocean vertical physics 
    2033   USE sol_oce         ! ocean elliptic solver 
     34   USE solver          ! solver initialization 
    2135   USE solpcg          ! preconditionned conjugate gradient solver 
    2236   USE solsor          ! Successive Over-relaxation solver 
     
    2842   USE lbclnk          ! ocean lateral boundary conditions (or mpp link) 
    2943   USE in_out_manager  ! I/O manager 
     44   USE iom 
     45   USE restart         ! only for lrst_oce 
    3046 
    3147   IMPLICIT NONE 
    3248   PRIVATE 
    3349 
    34    !! * Accessibility 
    3550   PUBLIC dyn_spg_rl   ! called by step.F90 
    3651 
     
    4257   !!   OPA 9.0 , LOCEAN-IPSL (2005)  
    4358   !! $Header$  
    44    !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt  
     59   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)  
    4560   !!---------------------------------------------------------------------- 
    4661 
     
    7691      !! ** Action : - Update (ua,va) with the surf. pressure gradient trend 
    7792      !! 
    78       !! References : 
    79       !!      Madec et al. 1988, ocean modelling, issue 78, 1-6. 
    80       !! 
    81       !! History : 
    82       !!        !  96-05  (G. Madec, M. Imbard, M. Guyon)  rewitting in 1 
    83       !!                  routine, without pointers, and with the same matrix 
    84       !!                  for sor and pcg, mpp exchange, and symmetric conditions 
    85       !!        !  96-07  (A. Weaver)  Euler forward step 
    86       !!        !  96-11  (A. Weaver)  correction to preconditioning: 
    87       !!        !  98-02  (M. Guyon)  FETI method 
    88       !!        !  98-05  (G. Roullet)  free surface 
    89       !!        !  97-09  (J.-M. Molines)  Open boundaries 
    90       !!   8.5  !  02-08  (G. Madec)  F90: Free form and module 
    91       !!        !  02-11  (C. Talandier, A-M Treguier) Open boundaries 
    92       !!   9.0  !  04-08  (C. Talandier)  New trends organization 
    93       !!    "   !  05-11  (V. Garnier) Surface pressure gradient organization 
     93      !! References : Madec et al. 1988, ocean modelling, issue 78, 1-6. 
    9494      !!--------------------------------------------------------------------- 
    95       !! * Arguments 
    9695      INTEGER, INTENT( in  ) ::   kt       ! ocean time-step index 
    97       INTEGER, INTENT( out ) ::   kindic   ! solver flag, take a negative value 
    98       !                                    ! when the solver doesnot converge 
    99       !! * Local declarations 
     96      INTEGER, INTENT( out ) ::   kindic   ! solver flag (<0 when the solver does not converge) 
    10097      INTEGER ::   ji, jj, jk    ! dummy loop indices 
    10198      REAL(wp) ::  zbsfa, zgcx, z2dt 
     
    114111 
    115112         ! set to zero rigid-lid specific arrays 
    116          spgu(:,:) = 0.e0      ! surface pressure gradient (i-direction)  
    117          spgv(:,:) = 0.e0      ! surface pressure gradient (j-direction) 
    118       ENDIF 
    119  
    120       ! 0. Initializations: 
    121       ! ------------------- 
    122 # if defined key_obc 
    123       ! space index on boundary arrays 
    124       ib = 1 
    125       ibm = 2 
    126       ibm2 = 3 
    127       ! time index on boundary arrays 
    128       it = 1 
    129       itm = 2 
    130       itm2 = 3 
    131 # endif 
    132  
     113         spgu(:,:) = 0.e0                   ! surface pressure gradient (i-direction)  
     114         spgv(:,:) = 0.e0                   ! surface pressure gradient (j-direction) 
     115 
     116         CALL solver_init( nit000 )         ! Elliptic solver initialisation 
     117 
     118         ! read rigid lid arrays in restart file 
     119         CALL rl_rst( nit000, 'READ' )      ! read or initialize the following fields: 
     120         !                                  ! gcx, gcxb, bsfb, bsfn, bsfd 
     121      ENDIF 
     122 
     123      !  Vertically averaged momentum trend 
     124      ! ------------------------------------ 
    133125      !                                                ! =============== 
    134126      DO jj = 2, jpjm1                                 !  Vertical slab 
    135127         !                                             ! =============== 
    136  
    137          ! 1. Vertically averaged momentum trend 
    138          ! ------------------------------------- 
    139          ! initialization to zero 
    140          spgu(:,jj) = 0. 
     128          
     129         spgu(:,jj) = 0.                          ! initialization to zero 
    141130         spgv(:,jj) = 0. 
    142          ! vertical sum 
    143          DO jk = 1, jpk 
     131         DO jk = 1, jpk                           ! vertical sum 
    144132            DO ji = 2, jpim1 
    145133               spgu(ji,jj) = spgu(ji,jj) + ua(ji,jj,jk) * fse3u(ji,jj,jk) * umask(ji,jj,jk) 
     
    147135            END DO  
    148136         END DO  
    149          ! divide by the depth 
    150          spgu(:,jj) = spgu(:,jj) * hur(:,jj) 
     137         spgu(:,jj) = spgu(:,jj) * hur(:,jj)      ! divide by the depth  
    151138         spgv(:,jj) = spgv(:,jj) * hvr(:,jj) 
    152139 
     
    155142      !                                                ! =============== 
    156143 
    157       !,,,,,,,,,,,,,,,,,,,,,,,,,,,,,synchro,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,, 
    158  
    159144      ! Boundary conditions on (spgu,spgv) 
    160145      CALL lbc_lnk( spgu, 'U', -1. ) 
    161146      CALL lbc_lnk( spgv, 'V', -1. ) 
    162147 
    163       !,,,,,,,,,,,,,,,,,,,,,,,,,,,,,synchro,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,, 
    164  
    165       ! 2. Barotropic streamfunction trend (bsfd) 
     148      !  Barotropic streamfunction trend (bsfd) 
    166149      ! ---------------------------------- 
    167  
    168       ! Curl of the vertically averaged velocity 
    169       DO jj = 2, jpjm1 
     150      DO jj = 2, jpjm1                            ! Curl of the vertically averaged velocity  
    170151         DO ji = 2, jpim1 
    171152            gcb(ji,jj) = -gcdprc(ji,jj)   & 
    172                        *(  ( e2v(ji+1,jj  )*spgv(ji+1,jj  ) - e2v(ji,jj)*spgv(ji,jj) )   & 
    173                           -( e1u(ji  ,jj+1)*spgu(ji  ,jj+1) - e1u(ji,jj)*spgu(ji,jj) )   )  
     153               &       *(  ( e2v(ji+1,jj  )*spgv(ji+1,jj  ) - e2v(ji,jj)*spgv(ji,jj) )   & 
     154               &          -( e1u(ji  ,jj+1)*spgu(ji  ,jj+1) - e1u(ji,jj)*spgu(ji,jj) )   )  
    174155         END DO 
    175156      END DO 
    176157 
    177158# if defined key_obc 
    178       ! Open boundary contribution 
    179       DO jj = 2, jpjm1 
     159      DO jj = 2, jpjm1                            ! Open boundary contribution  
    180160         DO ji = 2, jpim1 
    181161           gcb(ji,jj) = gcb(ji,jj) - gcdprc(ji,jj) * gcbob(ji,jj) 
     
    198178      ! applied the lateral boundary conditions 
    199179      IF( nsolv == 4)   CALL lbc_lnk_e( gcb, c_solver_pt, 1. )    
    200  
    201       !,,,,,,,,,,,,,,,,,,,,,,,,,,,,,synchro,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,, 
    202180 
    203181      ! Relative precision (computation on one processor) 
     
    234212      ENDIF 
    235213 
    236  
    237214      ! bsf trend update 
    238215      ! ---------------- 
    239  
    240216      bsfd(1:nlci,1:nlcj) = gcx(1:nlci,1:nlcj) 
    241  
    242217       
    243218      ! update bsf trend with islands trend 
    244219      ! ----------------------------------- 
    245  
    246220      IF( lk_isl )   CALL isl_dyn_spg         ! update bsfd 
    247  
    248221 
    249222# if defined key_obc 
     
    337310      ! 4. Barotropic stream function and array swap 
    338311      ! -------------------------------------------- 
    339  
    340312      ! Leap-frog time scheme, time filter and array swap 
    341313      IF( neuler == 0 .AND. kt == nit000 ) THEN 
     
    362334 
    363335# if defined key_obc 
     336      ib   = 1      ! space index on boundary arrays 
     337      ibm  = 2 
     338      ibm2 = 3 
     339      it   = 1      ! time index on boundary arrays 
     340      itm  = 2 
     341      itm2 = 3 
     342 
    364343      ! Swap of boundary arrays 
    365344      IF( lp_obc_east ) THEN 
     
    499478      ENDIF 
    500479# endif 
    501       ! 
    502       !,,,,,,,,,,,,,,,,,,,,,,,,,,,,,synchro,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,,, 
    503480       
    504481      !  add the surface pressure trend to the general trend 
    505482      ! ----------------------------------------------------- 
    506        
    507483      DO jj = 2, jpjm1 
    508  
    509484         ! update the surface pressure gradient with the barotropic trend 
    510485         DO ji = 2, jpim1 
     
    519494            END DO 
    520495         END DO 
    521  
    522       END DO 
     496      END DO 
     497 
     498      ! write rigid lid arrays in restart file 
     499      ! -------------------------------------- 
     500      IF( lrst_oce ) CALL rl_rst( kt, 'WRITE' ) 
    523501 
    524502   END SUBROUTINE dyn_spg_rl 
     503 
     504 
     505   SUBROUTINE rl_rst( kt, cdrw ) 
     506     !!--------------------------------------------------------------------- 
     507     !!                   ***  ROUTINE rl_rst  *** 
     508     !! 
     509     !! ** Purpose : Read or write rigid-lid arrays in ocean restart file 
     510     !!---------------------------------------------------------------------- 
     511     INTEGER         , INTENT(in) ::   kt         ! ocean time-step 
     512     CHARACTER(len=*), INTENT(in) ::   cdrw       ! "READ"/"WRITE" flag 
     513     !!---------------------------------------------------------------------- 
     514     ! 
     515     IF( TRIM(cdrw) == 'READ' ) THEN 
     516        IF( iom_varid( numror, 'gcx' ) > 0 ) THEN 
     517     ! Caution : extra-hallow 
     518     ! gcx and gcxb are defined as: DIMENSION(1-jpr2di:jpi+jpr2di,1-jpr2dj:jpj+jpr2dj) 
     519           CALL iom_get( numror, jpdom_local, 'gcx' , gcx (1:jpi,1:jpj) ) 
     520           CALL iom_get( numror, jpdom_local, 'gcxb', gcxb(1:jpi,1:jpj) ) 
     521           CALL iom_get( numror, jpdom_local, 'bsfb', bsfb(:,:)         ) 
     522           CALL iom_get( numror, jpdom_local, 'bsfn', bsfn(:,:)         ) 
     523           CALL iom_get( numror, jpdom_local, 'bsfd', bsfd(:,:)         ) 
     524           IF( neuler == 0 ) THEN 
     525              gcxb(:,:) = gcx (:,:) 
     526              bsfb(:,:) = bsfn(:,:) 
     527           ENDIF 
     528        ELSE 
     529           gcx (:,:) = 0.e0 
     530           gcxb(:,:) = 0.e0 
     531           bsfb(:,:) = 0.e0 
     532           bsfn(:,:) = 0.e0 
     533           bsfd(:,:) = 0.e0 
     534        ENDIF 
     535     ELSEIF(  TRIM(cdrw) == 'WRITE' ) THEN 
     536        ! Caution : extra-hallow, gcx and gcxb are defined as: DIMENSION(1-jpr2di:jpi+jpr2di,1-jpr2dj:jpj+jpr2dj) 
     537        CALL iom_rstput( kt, nitrst, numrow, 'gcx' , gcx (1:jpi,1:jpj) ) 
     538        CALL iom_rstput( kt, nitrst, numrow, 'gcxb', gcxb(1:jpi,1:jpj) ) 
     539        CALL iom_rstput( kt, nitrst, numrow, 'bsfb', bsfb(:,:)         ) 
     540        CALL iom_rstput( kt, nitrst, numrow, 'bsfn', bsfn(:,:)         ) 
     541        CALL iom_rstput( kt, nitrst, numrow, 'bsfd', bsfd(:,:)         ) 
     542     ENDIF 
     543     ! 
     544   END SUBROUTINE rl_rst 
    525545 
    526546#else 
Note: See TracChangeset for help on using the changeset viewer.