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 1607 – NEMO

Changeset 1607


Ignore:
Timestamp:
2009-08-12T10:01:24+02:00 (15 years ago)
Author:
ctlod
Message:

recover the sea-ice LIM3.0 restartability, see ticket: #529

Location:
trunk/NEMO/OPA_SRC
Files:
2 edited

Legend:

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

    r1565 r1607  
    1010   !!   ssh_wzv        : after ssh & now vertical velocity 
    1111   !!   ssh_nxt        : filter ans swap the ssh arrays 
    12    !!   ssh_rst        : read/write ssh restart fields in the ocean restart file 
    1312   !!---------------------------------------------------------------------- 
    1413   USE oce             ! ocean dynamics and tracers variables 
     
    7877         IF(lwp) WRITE(numout,*) 'ssh_wzv : after sea surface height and now vertical velocity ' 
    7978         IF(lwp) WRITE(numout,*) '~~~~~~~ ' 
    80          ! 
    81          CALL ssh_rst( nit000, 'READ' )       ! read or initialize sshb and sshn 
    8279         ! 
    8380         wn(:,:,jpk) = 0.e0                   ! bottom boundary condition: w=0 (set once for all) 
     
    111108 
    112109      !                                           !------------------------------! 
    113       !                                           !  Update Now Vertical coord.  ! 
    114       !                                           !------------------------------! 
    115       IF( lk_vvl ) THEN                           ! only in vvl case) 
    116          !                                             ! now local depth and scale factors (stored in fse3. arrays) 
     110      IF( lk_vvl ) THEN                           !  Update Now Vertical coord.  !   (only in vvl case) 
     111      !                                           !------------------------------! 
    117112         DO jk = 1, jpkm1 
    118             fsdept(:,:,jk) = fsdept_n(:,:,jk)               ! depths 
     113            fsdept(:,:,jk) = fsdept_n(:,:,jk)          ! now local depths stored in fsdep. arrays 
    119114            fsdepw(:,:,jk) = fsdepw_n(:,:,jk) 
    120115            fsde3w(:,:,jk) = fsde3w_n(:,:,jk) 
    121116            ! 
    122             fse3t (:,:,jk) = fse3t_n (:,:,jk)               ! vertical scale factors 
     117            fse3t (:,:,jk) = fse3t_n (:,:,jk)          ! vertical scale factors stored in fse3. arrays 
    123118            fse3u (:,:,jk) = fse3u_n (:,:,jk) 
    124119            fse3v (:,:,jk) = fse3v_n (:,:,jk) 
     
    128123            fse3vw(:,:,jk) = fse3vw_n(:,:,jk) 
    129124         END DO 
    130          !                                             ! ocean depth (at u- and v-points) 
     125         !                                             ! now ocean depth (at u- and v-points) 
    131126         hu(:,:) = hu_0(:,:) + sshu_n(:,:) 
    132127         hv(:,:) = hv_0(:,:) + sshv_n(:,:) 
    133          !                                             ! masked inverse of the ocean depth (at u- and v-points) 
     128         !                                             ! now masked inverse of the ocean depth (at u- and v-points) 
    134129         hur(:,:) = umask(:,:,1) / ( hu(:,:) + 1.e0 - umask(:,:,1) ) 
    135130         hvr(:,:) = vmask(:,:,1) / ( hv(:,:) + 1.e0 - vmask(:,:,1) ) 
     
    142137      ! set time step size (Euler/Leapfrog) 
    143138      z2dt = 2. * rdt  
    144       IF( neuler == 0 .AND. kt == nit000 ) z2dt =rdt 
     139      IF( neuler == 0 .AND. kt == nit000 )   z2dt =rdt 
    145140 
    146141      zraur = 1. / rauw 
     
    178173                  &                                  * ( e1t(ji,jj  ) * e2t(ji,jj  ) * ssha(ji,jj  )     & 
    179174                  &                                    + e1t(ji,jj+1) * e2t(ji,jj+1) * ssha(ji,jj+1) ) 
    180                sshf_a(ji,jj) = 0.25 * umask(ji,jj,1) * umask (ji,jj+1,1)   &   ! Caution : fmask not used 
     175               sshf_a(ji,jj) = 0.25 * umask(ji,jj,1) * umask (ji,jj+1,1)                                 &  
    181176                  &                                  * ( ssha(ji  ,jj) + ssha(ji  ,jj+1)                 & 
    182177                  &                                    + ssha(ji+1,jj) + ssha(ji+1,jj+1) ) 
     
    198193      END DO 
    199194      ! 
    200       CALL iom_put( "woce", wn   )                     ! vert. current 
     195      CALL iom_put( "woce", wn   )                     ! vertical velocity 
    201196      CALL iom_put( "ssh" , sshn )                     ! sea surface height 
    202197      ! 
     
    270265         ! 
    271266      ENDIF 
    272  
    273       !                      ! write filtered free surface arrays in restart file 
    274       IF( lrst_oce )   CALL ssh_rst( kt, 'WRITE' ) 
    275267      ! 
    276268      IF(ln_ctl)   CALL prt_ctl(tab2d_1=sshb    , clinfo1=' sshb  - : ', mask1=tmask, ovlap=1 ) 
    277269      ! 
    278270   END SUBROUTINE ssh_nxt 
    279  
    280  
    281    SUBROUTINE ssh_rst( kt, cdrw ) 
    282       !!--------------------------------------------------------------------- 
    283       !!                   ***  ROUTINE ssh_rst  *** 
    284       !! 
    285       !! ** Purpose :   Read or write Sea Surface Height arrays in restart file 
    286       !! 
    287       !! ** action  : - cdrw = READ  :  sshb, sshn  read    in ocean restart file 
    288       !!              - cdrw = WRITE :  sshb, sshn  written in ocean restart file 
    289       !!---------------------------------------------------------------------- 
    290       INTEGER         , INTENT(in) ::   kt         ! ocean time-step 
    291       CHARACTER(len=*), INTENT(in) ::   cdrw       ! "READ"/"WRITE" flag 
    292       !!---------------------------------------------------------------------- 
    293       ! 
    294       IF( TRIM(cdrw) == 'READ' ) THEN 
    295          IF( iom_varid( numror, 'sshn', ldstop = .FALSE. ) > 0 ) THEN 
    296             CALL iom_get( numror, jpdom_autoglo, 'sshb'  , sshb(:,:)   ) 
    297             CALL iom_get( numror, jpdom_autoglo, 'sshn'  , sshn(:,:)   ) 
    298             IF( neuler == 0 )   sshb(:,:) = sshn(:,:) 
    299          ELSE 
    300             IF( nn_rstssh == 1 ) THEN 
    301                sshb(:,:) = 0.e0 
    302                sshn(:,:) = 0.e0 
    303             ENDIF 
    304          ENDIF 
    305       ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN 
    306          CALL iom_rstput( kt, nitrst, numrow, 'sshb'  , sshb(:,:) ) 
    307          CALL iom_rstput( kt, nitrst, numrow, 'sshn'  , sshn(:,:) ) 
    308       ENDIF 
    309       ! 
    310    END SUBROUTINE ssh_rst 
    311271 
    312272   !!====================================================================== 
  • trunk/NEMO/OPA_SRC/restart.F90

    r1545 r1607  
    112112      !!---------------------------------------------------------------------- 
    113113 
    114       !                                                                     ! the begining of the run [s] 
    115       CALL iom_rstput( kt, nitrst, numrow, 'rdt'    , rdt               )   ! dynamics time step 
    116       CALL iom_rstput( kt, nitrst, numrow, 'rdttra1', rdttra(1)         )   ! surface tracer time step 
    117  
    118       ! prognostic variables 
    119       CALL iom_rstput( kt, nitrst, numrow, 'ub'     , ub      )    
     114      CALL iom_rstput( kt, nitrst, numrow, 'rdt'    , rdt       )   ! dynamics time step 
     115      CALL iom_rstput( kt, nitrst, numrow, 'rdttra1', rdttra(1) )   ! surface tracer time step 
     116 
     117      CALL iom_rstput( kt, nitrst, numrow, 'ub'     , ub      )     ! before fields 
    120118      CALL iom_rstput( kt, nitrst, numrow, 'vb'     , vb      ) 
    121119      CALL iom_rstput( kt, nitrst, numrow, 'tb'     , tb      ) 
     
    123121      CALL iom_rstput( kt, nitrst, numrow, 'rotb'   , rotb    ) 
    124122      CALL iom_rstput( kt, nitrst, numrow, 'hdivb'  , hdivb   ) 
    125       CALL iom_rstput( kt, nitrst, numrow, 'un'     , un      ) 
     123      CALL iom_rstput( kt, nitrst, numrow, 'sshb'   , sshb    ) 
     124      ! 
     125      CALL iom_rstput( kt, nitrst, numrow, 'un'     , un      )     ! now fields 
    126126      CALL iom_rstput( kt, nitrst, numrow, 'vn'     , vn      ) 
    127127      CALL iom_rstput( kt, nitrst, numrow, 'tn'     , tn      ) 
     
    129129      CALL iom_rstput( kt, nitrst, numrow, 'rotn'   , rotn    ) 
    130130      CALL iom_rstput( kt, nitrst, numrow, 'hdivn'  , hdivn   ) 
     131      CALL iom_rstput( kt, nitrst, numrow, 'sshn'   , sshn    ) 
    131132 
    132133      CALL iom_rstput( kt, nitrst, numrow, 'rhop'   , rhop    ) 
     
    182183         IF( zrdttra1 /= rdttra(1) )   neuler = 0 
    183184      ENDIF 
    184       !                                                       ! Read prognostic variables 
    185       CALL iom_get( numror, jpdom_autoglo, 'ub'   , ub    )        ! before i-component velocity 
    186       CALL iom_get( numror, jpdom_autoglo, 'vb'   , vb    )        ! before j-component velocity 
    187       CALL iom_get( numror, jpdom_autoglo, 'tb'   , tb    )        ! before temperature 
    188       CALL iom_get( numror, jpdom_autoglo, 'sb'   , sb    )        ! before salinity 
    189       CALL iom_get( numror, jpdom_autoglo, 'rotb' , rotb  )        ! before curl 
    190       CALL iom_get( numror, jpdom_autoglo, 'hdivb', hdivb )        ! before horizontal divergence 
    191       CALL iom_get( numror, jpdom_autoglo, 'un'   , un    )        ! now    i-component velocity 
    192       CALL iom_get( numror, jpdom_autoglo, 'vn'   , vn    )        ! now    j-component velocity 
    193       CALL iom_get( numror, jpdom_autoglo, 'tn'   , tn    )        ! now    temperature 
    194       CALL iom_get( numror, jpdom_autoglo, 'sn'   , sn    )        ! now    salinity 
    195       CALL iom_get( numror, jpdom_autoglo, 'rotn' , rotn  )        ! now    curl 
    196       CALL iom_get( numror, jpdom_autoglo, 'hdivn', hdivn )        ! now    horizontal divergence 
     185      !  
     186      CALL iom_get( numror, jpdom_autoglo, 'ub'   , ub    )        ! before fields 
     187      CALL iom_get( numror, jpdom_autoglo, 'vb'   , vb    ) 
     188      CALL iom_get( numror, jpdom_autoglo, 'tb'   , tb    ) 
     189      CALL iom_get( numror, jpdom_autoglo, 'sb'   , sb    ) 
     190      CALL iom_get( numror, jpdom_autoglo, 'rotb' , rotb  ) 
     191      CALL iom_get( numror, jpdom_autoglo, 'hdivb', hdivb ) 
     192      CALL iom_get( numror, jpdom_autoglo, 'sshb' , sshb  ) 
     193      ! 
     194      CALL iom_get( numror, jpdom_autoglo, 'un'   , un    )        ! now    fields 
     195      CALL iom_get( numror, jpdom_autoglo, 'vn'   , vn    ) 
     196      CALL iom_get( numror, jpdom_autoglo, 'tn'   , tn    ) 
     197      CALL iom_get( numror, jpdom_autoglo, 'sn'   , sn    ) 
     198      CALL iom_get( numror, jpdom_autoglo, 'rotn' , rotn  ) 
     199      CALL iom_get( numror, jpdom_autoglo, 'hdivn', hdivn ) 
     200      CALL iom_get( numror, jpdom_autoglo, 'sshn' , sshn  ) 
    197201 
    198202      CALL iom_get( numror, jpdom_autoglo, 'rhop' , rhop  )        ! now    potential density 
     
    206210 
    207211      IF( neuler == 0 ) THEN                                  ! Euler restart (neuler=0) 
    208          tb   (:,:,:) = tn   (:,:,:)                             ! all before fields set to now field values 
     212         tb   (:,:,:) = tn   (:,:,:)                             ! all before fields set to now values 
    209213         sb   (:,:,:) = sn   (:,:,:) 
    210214         ub   (:,:,:) = un   (:,:,:) 
     
    212216         rotb (:,:,:) = rotn (:,:,:) 
    213217         hdivb(:,:,:) = hdivn(:,:,:) 
     218         sshb (:,:)   = sshn (:,:) 
    214219      ENDIF 
    215220      ! 
    216221   END SUBROUTINE rst_read 
    217  
    218222 
    219223   !!===================================================================== 
Note: See TracChangeset for help on using the changeset viewer.