MODULE restart !!====================================================================== !! *** MODULE restart *** !! Ocean restart : write the ocean restart file !!====================================================================== !! History : ! 99-11 (M. Imbard) Original code !! 8.5 ! 02-08 (G. Madec) F90: Free form !! 9.0 ! 05-11 (V. Garnier) Surface pressure gradient organization !! 9.0 ! 06-07 (S. Masson) use IOM for restart !!---------------------------------------------------------------------- !!---------------------------------------------------------------------- !! rst_opn : open the ocean restart file !! rst_write : write the ocean restart file !! rst_read : read the ocean restart file !!---------------------------------------------------------------------- USE dom_oce ! ocean space and time domain USE oce ! ocean dynamics and tracers USE phycst ! physical constants USE in_out_manager ! I/O manager USE iom ! I/O module USE c1d ! re-initialization of u-v mask for the 1D configuration USE zpshde ! partial step: hor. derivative (zps_hde routine) USE eosbn2 ! equation of state (eos bn2 routine) USE trdmld_oce ! ocean active mixed layer tracers trends variables #if defined key_zdftke2 USE zdf_oce #endif IMPLICIT NONE PRIVATE PUBLIC rst_opn ! routine called by step module PUBLIC rst_write ! routine called by step module PUBLIC rst_read ! routine called by opa module LOGICAL, PUBLIC :: lrst_oce !: logical to control the oce restart write INTEGER, PUBLIC :: numror, numrow !: logical unit for cean restart (read and write) !! * Substitutions # include "vectopt_loop_substitute.h90" !!---------------------------------------------------------------------- !! OPA 9.0 , LOCEAN-IPSL (2006) !! $Id: restart.F90 1239 2008-12-31 09:35:35Z ctlod $ !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt) !!---------------------------------------------------------------------- CONTAINS SUBROUTINE rst_opn( kt ) !!--------------------------------------------------------------------- !! *** ROUTINE rst_opn *** !! !! ** Purpose : + initialization (should be read in the namelist) of nitrst !! + open the restart when we are one time step before nitrst !! - restart header is defined when kt = nitrst-1 !! - restart data are written when kt = nitrst !! + define lrst_oce to .TRUE. when we need to define or write the restart !!---------------------------------------------------------------------- INTEGER, INTENT(in) :: kt ! ocean time-step !! CHARACTER(LEN=20) :: clkt ! ocean time-step deine as a character CHARACTER(LEN=50) :: clname ! ice output restart file name !!---------------------------------------------------------------------- ! IF( kt == nit000 ) THEN ! default definitions lrst_oce = .FALSE. nitrst = nitend #if defined key_zdftke2 nitrst_tke2 = nitrst + 1 #endif ENDIF IF( MOD( kt - 1, nstock ) == 0 ) THEN ! we use kt - 1 and not kt - nit000 to keep the same periodicity from the beginning of the experiment nitrst = kt + nstock - 1 ! define the next value of nitrst for restart writing IF( nitrst > nitend ) nitrst = nitend ! make sure we write a restart at the end of the run ENDIF #if defined key_zdftke2 IF ( nitrst_tke2 .NE. kt ) nitrst_tke2 = nitrst + 1 #endif ! to get better performances with NetCDF format: ! we open and define the ocean restart file one time step before writing the data (-> at nitrst - 1) ! except if we write ocean restart files every time step or if an ocean restart file was writen at nitend - 1 IF( kt == nitrst - 1 .OR. nstock == 1 .OR. ( kt == nitend .AND. .NOT. lrst_oce ) ) THEN ! beware of the format used to write kt (default is i8.8, that should be large enough...) IF( nitrst > 999999999 ) THEN ; WRITE(clkt, * ) nitrst ELSE ; WRITE(clkt, '(i8.8)') nitrst ENDIF ! create the file clname = TRIM(cexper)//"_"//TRIM(ADJUSTL(clkt))//"_"//TRIM(cn_ocerst_out) IF(lwp) THEN WRITE(numout,*) SELECT CASE ( jprstlib ) CASE ( jprstdimg ) ; WRITE(numout,*) ' open ocean restart binary file: '//clname CASE DEFAULT ; WRITE(numout,*) ' open ocean restart NetCDF file: '//clname END SELECT IF( kt == nitrst - 1 ) THEN ; WRITE(numout,*) ' kt = nitrst - 1 = ', kt ELSE ; WRITE(numout,*) ' kt = ' , kt ENDIF ENDIF CALL iom_open( clname, numrow, ldwrt = .TRUE., kiolib = jprstlib, ldclobber=.TRUE. ) lrst_oce = .TRUE. ENDIF ! END SUBROUTINE rst_opn SUBROUTINE rst_write( kt ) !!--------------------------------------------------------------------- !! *** ROUTINE rstwrite *** !! !! ** Purpose : Write restart fields in the format corresponding to jprstlib !! !! ** Method : Write in numrow when kt == nitrst in NetCDF !! file, save fields which are necessary for restart !!---------------------------------------------------------------------- INTEGER, INTENT(in) :: kt ! ocean time-step !!---------------------------------------------------------------------- #if defined key_zdftke2 IF( kt == nitrst_tke2 ) THEN CALL iom_close( numrow ) ! close the restart file (only at last time step) IF( .NOT. lk_trdmld ) lrst_oce = .FALSE. ELSE #endif ! ! the begining of the run [s] CALL iom_rstput( kt, nitrst, numrow, 'rdt' , rdt ) ! dynamics time step CALL iom_rstput( kt, nitrst, numrow, 'rdttra1', rdttra(1) ) ! surface tracer time step ! prognostic variables CALL iom_rstput( kt, nitrst, numrow, 'ub' , ub ) CALL iom_rstput( kt, nitrst, numrow, 'vb' , vb ) CALL iom_rstput( kt, nitrst, numrow, 'tb' , tb ) CALL iom_rstput( kt, nitrst, numrow, 'sb' , sb ) CALL iom_rstput( kt, nitrst, numrow, 'rotb' , rotb ) CALL iom_rstput( kt, nitrst, numrow, 'hdivb' , hdivb ) CALL iom_rstput( kt, nitrst, numrow, 'un' , un ) CALL iom_rstput( kt, nitrst, numrow, 'vn' , vn ) IF( lk_vvl ) CALL iom_rstput( kt, nitrst, numrow, 'wn' , wn ) CALL iom_rstput( kt, nitrst, numrow, 'tn' , tn ) CALL iom_rstput( kt, nitrst, numrow, 'sn' , sn ) CALL iom_rstput( kt, nitrst, numrow, 'rotn' , rotn ) CALL iom_rstput( kt, nitrst, numrow, 'hdivn' , hdivn ) #if defined key_zdftke2 CALL iom_rstput( kt, nitrst, numrow, 'rhop' , rhop ) #endif IF( nn_dynhpg_rst == 1 .OR. lk_vvl ) THEN CALL iom_rstput( kt, nitrst, numrow, 'rhd' , rhd ) #if defined key_zdftke2 CALL iom_rstput( kt, nitrst, numrow, 'rn2' , rn2 ) CALL iom_rstput( kt, nitrst, numrow, 'avt' , avt ) ENDIF #else CALL iom_rstput( kt, nitrst, numrow, 'rhop', rhop ) #endif IF( ln_zps ) THEN CALL iom_rstput( kt, nitrst, numrow, 'gtu' , gtu ) CALL iom_rstput( kt, nitrst, numrow, 'gsu' , gsu ) CALL iom_rstput( kt, nitrst, numrow, 'gru' , gru ) CALL iom_rstput( kt, nitrst, numrow, 'gtv' , gtv ) CALL iom_rstput( kt, nitrst, numrow, 'gsv' , gsv ) CALL iom_rstput( kt, nitrst, numrow, 'grv' , grv ) ENDIF ENDIF #if ! defined key_zdftke2 IF( kt == nitrst ) THEN CALL iom_close( numrow ) ! close the restart file (only at last time step) IF( .NOT. lk_trdmld ) lrst_oce = .FALSE. ENDIF #endif ! END SUBROUTINE rst_write SUBROUTINE rst_read !!---------------------------------------------------------------------- !! *** ROUTINE rst_read *** !! !! ** Purpose : Read files for restart (format fixed by jprstlib ) !! !! ** Method : Read the previous fields on the NetCDF/binary file !! the first record indicates previous characterics !! after control with the present run, we read : !! - prognostic variables on the second record !! - elliptic solver arrays !! - barotropic stream function arrays ("key_dynspg_rl" defined) !! or free surface arrays !! - tke arrays (lk_zdftke=T .OR. lk_zdftke2=T) !! for this last three records, the previous characteristics !! could be different with those used in the present run. !!---------------------------------------------------------------------- REAL(wp) :: zrdt, zrdttra1 !!---------------------------------------------------------------------- IF(lwp) THEN ! Contol prints WRITE(numout,*) SELECT CASE ( jprstlib ) CASE ( jpnf90 ) ; WRITE(numout,*) 'rst_read : read oce NetCDF restart file' CASE ( jprstdimg ) ; WRITE(numout,*) 'rst_read : read oce binary restart file' END SELECT WRITE(numout,*) '~~~~~~~~' ENDIF CALL iom_open( cn_ocerst_in, numror, kiolib = jprstlib ) ! Check dynamics and tracer time-step consistency and force Euler restart if changed IF( iom_varid( numror, 'rdt', ldstop = .FALSE. ) > 0 ) THEN CALL iom_get( numror, 'rdt', zrdt ) IF( zrdt /= rdt ) neuler = 0 ENDIF IF( iom_varid( numror, 'rdttra1', ldstop = .FALSE. ) > 0 ) THEN CALL iom_get( numror, 'rdttra1', zrdttra1 ) IF( zrdttra1 /= rdttra(1) ) neuler = 0 ENDIF ! ! Read prognostic variables CALL iom_get( numror, jpdom_autoglo, 'ub' , ub ) ! before i-component velocity CALL iom_get( numror, jpdom_autoglo, 'vb' , vb ) ! before j-component velocity CALL iom_get( numror, jpdom_autoglo, 'tb' , tb ) ! before temperature CALL iom_get( numror, jpdom_autoglo, 'sb' , sb ) ! before salinity CALL iom_get( numror, jpdom_autoglo, 'rotb' , rotb ) ! before curl CALL iom_get( numror, jpdom_autoglo, 'hdivb', hdivb ) ! before horizontal divergence CALL iom_get( numror, jpdom_autoglo, 'un' , un ) ! now i-component velocity CALL iom_get( numror, jpdom_autoglo, 'vn' , vn ) ! now j-component velocity IF( lk_vvl ) CALL iom_get( numror, jpdom_autoglo, 'wn' , wn ) ! now k-component velocity CALL iom_get( numror, jpdom_autoglo, 'tn' , tn ) ! now temperature CALL iom_get( numror, jpdom_autoglo, 'sn' , sn ) ! now salinity CALL iom_get( numror, jpdom_autoglo, 'rotn' , rotn ) ! now curl CALL iom_get( numror, jpdom_autoglo, 'hdivn', hdivn ) ! now horizontal divergence IF( neuler == 0 ) THEN ! Euler restart (neuler=0) tb (:,:,:) = tn (:,:,:) ! all before fields set to now field values sb (:,:,:) = sn (:,:,:) ub (:,:,:) = un (:,:,:) vb (:,:,:) = vn (:,:,:) rotb (:,:,:) = rotn (:,:,:) hdivb(:,:,:) = hdivn(:,:,:) ENDIF #if defined key_zdftke2 CALL eos( tb, sb, rhd, rhop ) ! before potential and in situ densities IF( iom_varid( numror, 'rhd', ldstop = .FALSE. ) > 0 ) THEN CALL iom_get( numror, jpdom_autoglo, 'rhd' , rhd ) ENDIF IF( iom_varid( numror, 'rhop', ldstop = .FALSE. ) > 0 ) THEN CALL iom_get( numror, jpdom_autoglo, 'rhop', rhop ) ENDIF #else IF( iom_varid( numror, 'rhd', ldstop = .FALSE. ) > 0 ) THEN CALL iom_get( numror, jpdom_autoglo, 'rhd' , rhd ) CALL iom_get( numror, jpdom_autoglo, 'rhop', rhop ) ELSE CALL eos( tb, sb, rhd, rhop ) ! before potential and in situ densities ENDIF #endif #if defined key_zdftke2 CALL eos_init ! usefull to get the equation state type neos parameter IF( iom_varid( numror, 'rn2', ldstop = .FALSE. ) > 0 ) THEN CALL iom_get( numror, jpdom_autoglo, 'rn2' , rn2 ) ELSE IF ( ln_dynhpg_imp ) THEN CALL bn2( tb, sb, rn2 ) ! before Brunt-Vaisala frequency ENDIF ENDIF IF( iom_varid( numror, 'avt', ldstop = .FALSE. ) > 0 ) THEN CALL iom_get( numror, jpdom_autoglo, 'avt' , avt ) ELSE IF ( ln_dynhpg_imp ) avt (:,:,:) = 1.2e-5 * tmask(:,:,:) ENDIF #endif IF( ln_zps .AND. .NOT. lk_c1d ) THEN IF( iom_varid( numror, 'gtu', ldstop = .FALSE. ) > 0 ) THEN CALL iom_get( numror, jpdom_autoglo, 'gtu' , gtu ) CALL iom_get( numror, jpdom_autoglo, 'gsu' , gsu ) CALL iom_get( numror, jpdom_autoglo, 'gru' , gru ) CALL iom_get( numror, jpdom_autoglo, 'gtv' , gtv ) CALL iom_get( numror, jpdom_autoglo, 'gsv' , gsv ) CALL iom_get( numror, jpdom_autoglo, 'grv' , grv ) ELSE CALL zps_hde( nit000, tb , sb , rhd, & ! Partial steps: before Horizontal DErivative & gtu, gsu, gru, & ! of t, s, rd at the bottom ocean level & gtv, gsv, grv ) ENDIF ENDIF ! END SUBROUTINE rst_read !!===================================================================== END MODULE restart