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 1481 for trunk/NEMO/OPA_SRC/restart.F90 – NEMO

Ignore:
Timestamp:
2009-07-03T17:07:08+02:00 (15 years ago)
Author:
ctlod
Message:

ensure restartability of TKE2 in coupled mode + few minor bugs corrections, see ticket: #466

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMO/OPA_SRC/restart.F90

    r1473 r1481  
    2323   USE zpshde          ! partial step: hor. derivative (zps_hde routine) 
    2424   USE eosbn2          ! equation of state            (eos bn2 routine) 
     25   USE zdfddm          ! double diffusion mixing  
     26   USE zdfmxl          ! mixed layer depth 
    2527   USE trdmld_oce      ! ocean active mixed layer tracers trends variables 
    26 #if defined key_zdftke2 
    27    USE zdf_oce 
    28 #endif 
    2928 
    3029   IMPLICIT NONE 
     
    6766         lrst_oce = .FALSE.    
    6867         nitrst = nitend 
    69 #if defined key_zdftke2 
    70          nitrst_tke2 = nitrst + 1 
    71 #endif 
    7268      ENDIF 
    7369      IF( MOD( kt - 1, nstock ) == 0 ) THEN    
     
    7672         IF( nitrst > nitend )   nitrst = nitend   ! make sure we write a restart at the end of the run 
    7773      ENDIF 
    78 #if defined key_zdftke2 
    79       IF ( nitrst_tke2 .NE. kt ) nitrst_tke2 = nitrst + 1 
    80 #endif 
    8174      ! to get better performances with NetCDF format: 
    8275      ! we open and define the ocean restart file one time step before writing the data (-> at nitrst - 1) 
     
    119112      !!---------------------------------------------------------------------- 
    120113 
    121 #if defined key_zdftke2 
    122       IF( kt == nitrst_tke2 ) THEN 
    123          CALL iom_close( numrow )     ! close the restart file (only at last time step) 
    124          IF( .NOT. lk_trdmld )   lrst_oce = .FALSE. 
    125       ELSE 
    126 #endif 
    127114      !                                                                     ! the begining of the run [s] 
    128115      CALL iom_rstput( kt, nitrst, numrow, 'rdt'    , rdt               )   ! dynamics time step 
     
    138125      CALL iom_rstput( kt, nitrst, numrow, 'un'     , un      ) 
    139126      CALL iom_rstput( kt, nitrst, numrow, 'vn'     , vn      ) 
    140       IF( lk_vvl ) CALL iom_rstput( kt, nitrst, numrow, 'wn'     , wn      ) 
    141127      CALL iom_rstput( kt, nitrst, numrow, 'tn'     , tn      ) 
    142128      CALL iom_rstput( kt, nitrst, numrow, 'sn'     , sn      ) 
    143129      CALL iom_rstput( kt, nitrst, numrow, 'rotn'   , rotn    ) 
    144130      CALL iom_rstput( kt, nitrst, numrow, 'hdivn'  , hdivn   ) 
    145       CALL iom_rstput( kt, nitrst, numrow, 'rn2' , rn2  ) 
    146  
    147 #if defined key_zdftke2 
    148          CALL iom_rstput( kt, nitrst, numrow, 'rhop'   , rhop    ) 
    149 #endif 
    150       IF( nn_dynhpg_rst == 1 .OR. lk_vvl ) THEN 
    151          CALL iom_rstput( kt, nitrst, numrow, 'rhd' , rhd  ) 
    152 #if defined key_zdftke2 
    153          CALL iom_rstput( kt, nitrst, numrow, 'avt' , avt  ) 
    154       ENDIF 
    155 #else 
    156          CALL iom_rstput( kt, nitrst, numrow, 'rhop', rhop ) 
    157 #endif 
    158          IF( ln_zps ) THEN 
    159             CALL iom_rstput( kt, nitrst, numrow, 'gtu' , gtu ) 
    160             CALL iom_rstput( kt, nitrst, numrow, 'gsu' , gsu ) 
    161             CALL iom_rstput( kt, nitrst, numrow, 'gru' , gru ) 
    162             CALL iom_rstput( kt, nitrst, numrow, 'gtv' , gtv ) 
    163             CALL iom_rstput( kt, nitrst, numrow, 'gsv' , gsv ) 
    164             CALL iom_rstput( kt, nitrst, numrow, 'grv' , grv ) 
    165          ENDIF 
    166       ENDIF 
    167 #if ! defined key_zdftke2 
     131 
     132      CALL iom_rstput( kt, nitrst, numrow, 'rhop'  , rhop   ) 
     133 
     134      IF( lk_zdfddm )   CALL iom_rstput( kt, nitrst, numrow, 'rrau' , rrau  ) 
     135 
    168136      IF( kt == nitrst ) THEN 
    169137         CALL iom_close( numrow )     ! close the restart file (only at last time step) 
    170138         IF( .NOT. lk_trdmld )   lrst_oce = .FALSE. 
    171139      ENDIF 
    172 #endif 
    173140      ! 
    174141   END SUBROUTINE rst_write 
     
    232199      CALL iom_get( numror, jpdom_autoglo, 'un'   , un    )        ! now    i-component velocity 
    233200      CALL iom_get( numror, jpdom_autoglo, 'vn'   , vn    )        ! now    j-component velocity 
    234       IF( lk_vvl ) CALL iom_get( numror, jpdom_autoglo, 'wn'   , wn    )        ! now    k-component velocity 
    235201      CALL iom_get( numror, jpdom_autoglo, 'tn'   , tn    )        ! now    temperature 
    236202      CALL iom_get( numror, jpdom_autoglo, 'sn'   , sn    )        ! now    salinity 
     
    238204      CALL iom_get( numror, jpdom_autoglo, 'hdivn', hdivn )        ! now    horizontal divergence 
    239205 
     206      CALL iom_get( numror, jpdom_autoglo, 'rhop', rhop )          ! now    mixed layer depth 
     207 
     208      IF( lk_zdfddm ) THEN 
     209         IF( iom_varid( numror, 'rrau', ldstop = .FALSE. ) > 0 ) THEN 
     210            CALL iom_get( numror, jpdom_autoglo, 'rrau' , rrau  ) 
     211         ELSE 
     212            CALL eos_init                   ! read equation state type neos parameter 
     213            CALL eos( tb, sb, rhd, rhop )   ! compute rrau 
     214         ENDIF 
     215      ENDIF  
    240216 
    241217      IF( neuler == 0 ) THEN                                  ! Euler restart (neuler=0) 
     
    247223         hdivb(:,:,:) = hdivn(:,:,:) 
    248224      ENDIF 
    249       CALL eos_init            ! usefull to get the equation state type neos parameter 
    250       IF( iom_varid( numror, 'rn2', ldstop = .FALSE. ) > 0 ) THEN 
    251          CALL iom_get( numror, jpdom_autoglo, 'rn2' , rn2  ) 
    252       ELSE 
    253          CALL bn2( tb, sb, rn2 )  ! before Brunt-Vaisala frequency    
    254       ENDIF 
    255  
    256 #if defined key_zdftke2 
    257       CALL eos( tb, sb, rhd, rhop )        ! before potential and in situ densities 
    258       IF( iom_varid( numror, 'rhd', ldstop = .FALSE. ) > 0 ) THEN 
    259          CALL iom_get( numror, jpdom_autoglo, 'rhd' , rhd  ) 
    260       ENDIF 
    261       IF( iom_varid( numror, 'rhop', ldstop = .FALSE. ) > 0 ) THEN 
    262          CALL iom_get( numror, jpdom_autoglo, 'rhop', rhop ) 
    263       ENDIF 
    264 #else 
    265       IF( iom_varid( numror, 'rhd', ldstop = .FALSE. ) > 0 ) THEN 
    266          CALL iom_get( numror, jpdom_autoglo, 'rhd' , rhd  ) 
    267          CALL iom_get( numror, jpdom_autoglo, 'rhop', rhop ) 
    268       ELSE 
    269          CALL eos( tb, sb, rhd, rhop )        ! before potential and in situ densities 
    270       ENDIF 
    271 #endif 
    272 #if defined key_zdftke2 
    273       IF( iom_varid( numror, 'avt', ldstop = .FALSE. ) > 0 ) THEN 
    274          CALL iom_get( numror, jpdom_autoglo, 'avt' , avt  ) 
    275       ELSE 
    276          IF ( ln_dynhpg_imp ) avt (:,:,:) = 1.2e-5 * tmask(:,:,:) 
    277       ENDIF 
    278 #endif 
    279  
    280       IF( ln_zps .AND. .NOT. lk_c1d ) THEN 
    281          IF( iom_varid( numror, 'gtu', ldstop = .FALSE. ) > 0 ) THEN 
    282             CALL iom_get( numror, jpdom_autoglo, 'gtu' , gtu ) 
    283             CALL iom_get( numror, jpdom_autoglo, 'gsu' , gsu ) 
    284             CALL iom_get( numror, jpdom_autoglo, 'gru' , gru ) 
    285             CALL iom_get( numror, jpdom_autoglo, 'gtv' , gtv ) 
    286             CALL iom_get( numror, jpdom_autoglo, 'gsv' , gsv ) 
    287             CALL iom_get( numror, jpdom_autoglo, 'grv' , grv ) 
    288          ELSE 
    289             CALL zps_hde( nit000, tb , sb , rhd,   &  ! Partial steps: before Horizontal DErivative 
    290                &                  gtu, gsu, gru,   &  ! of t, s, rd at the bottom ocean level 
    291                &                  gtv, gsv, grv ) 
    292          ENDIF 
    293       ENDIF 
    294225      ! 
    295226   END SUBROUTINE rst_read 
Note: See TracChangeset for help on using the changeset viewer.