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

source: trunk/NEMO/OPA_SRC/restart.F90 @ 1191

Last change on this file since 1191 was 1146, checked in by rblod, 16 years ago

Add svn Id (first try), see ticket #210

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id
File size: 12.3 KB
RevLine 
[3]1MODULE restart
2   !!======================================================================
3   !!                     ***  MODULE  restart  ***
4   !! Ocean restart :  write the ocean restart file
[508]5   !!======================================================================
6   !! History :        !  99-11  (M. Imbard)  Original code
7   !!             8.5  !  02-08  (G. Madec)  F90: Free form
8   !!             9.0  !  05-11  (V. Garnier) Surface pressure gradient organization
9   !!             9.0  !  06-07  (S. Masson)  use IOM for restart
10   !!----------------------------------------------------------------------
[3]11
12   !!----------------------------------------------------------------------
[508]13   !!   rst_opn    : open the ocean restart file
14   !!   rst_write  : write the ocean restart file
15   !!   rst_read   : read the ocean restart file
[3]16   !!----------------------------------------------------------------------
17   USE dom_oce         ! ocean space and time domain
18   USE oce             ! ocean dynamics and tracers
19   USE phycst          ! physical constants
[467]20   USE cpl_oce, ONLY : lk_cpl              !
[508]21   USE in_out_manager  ! I/O manager
22   USE iom             ! I/O module
[900]23   USE c1d             ! re-initialization of u-v mask for the 1D configuration
[544]24   USE zpshde          ! partial step: hor. derivative (zps_hde routine)
25   USE eosbn2          ! equation of state            (eos bn2 routine)
[579]26   USE trdmld_oce      ! ocean active mixed layer tracers trends variables
[3]27
28   IMPLICIT NONE
29   PRIVATE
30
[508]31   PUBLIC   rst_opn    ! routine called by step module
32   PUBLIC   rst_write  ! routine called by step module
33   PUBLIC   rst_read   ! routine called by opa  module
[3]34
[521]35   LOGICAL, PUBLIC ::   lrst_oce                  !: logical to control the oce restart write
[579]36   INTEGER, PUBLIC ::   numror, numrow            !: logical unit for cean restart (read and write)
[508]37
38   !! * Substitutions
39#  include "vectopt_loop_substitute.h90"
[3]40   !!----------------------------------------------------------------------
[508]41   !!  OPA 9.0 , LOCEAN-IPSL (2006)
[888]42   !! $Id$
[508]43   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
[359]44   !!----------------------------------------------------------------------
[3]45
46CONTAINS
47
[508]48   SUBROUTINE rst_opn( kt )
49      !!---------------------------------------------------------------------
50      !!                   ***  ROUTINE rst_opn  ***
51      !!                     
52      !! ** Purpose : + initialization (should be read in the namelist) of nitrst
53      !!              + open the restart when we are one time step before nitrst
54      !!                   - restart header is defined when kt = nitrst-1
55      !!                   - restart data  are written when kt = nitrst
56      !!              + define lrst_oce to .TRUE. when we need to define or write the restart
57      !!----------------------------------------------------------------------
58      INTEGER, INTENT(in) ::   kt     ! ocean time-step
59      !!
60      CHARACTER(LEN=20)   ::   clkt     ! ocean time-step deine as a character
61      CHARACTER(LEN=50)   ::   clname   ! ice output restart file name
62      !!----------------------------------------------------------------------
63      !
[783]64      IF( kt == nit000 ) THEN   ! default definitions
65         lrst_oce = .FALSE.   
66         nitrst = nitend
[508]67      ENDIF
[783]68      IF( MOD( kt - 1, nstock ) == 0 ) THEN   
69         ! we use kt - 1 and not kt - nit000 to keep the same periodicity from the beginning of the experiment
70         nitrst = kt + nstock - 1                  ! define the next value of nitrst for restart writing
71         IF( nitrst > nitend )   nitrst = nitend   ! make sure we write a restart at the end of the run
72      ENDIF
[632]73
[783]74      ! to get better performances with NetCDF format:
75      ! we open and define the ocean restart file one time step before writing the data (-> at nitrst - 1)
76      ! except if we write ocean restart files every time step or if an ocean restart file was writen at nitend - 1
77      IF( kt == nitrst - 1 .OR. nstock == 1 .OR. ( kt == nitend .AND. .NOT. lrst_oce ) ) THEN
78         ! beware of the format used to write kt (default is i8.8, that should be large enough...)
79         IF( nitrst > 999999999 ) THEN   ;   WRITE(clkt, *       ) nitrst
80         ELSE                            ;   WRITE(clkt, '(i8.8)') nitrst
[508]81         ENDIF
82         ! create the file
83         clname = TRIM(cexper)//"_"//TRIM(ADJUSTL(clkt))//"_restart"
[611]84         IF(lwp) THEN
85            WRITE(numout,*)
[632]86            SELECT CASE ( jprstlib )
[783]87            CASE ( jprstdimg )   ;   WRITE(numout,*) '             open ocean restart binary file: '//clname
88            CASE DEFAULT         ;   WRITE(numout,*) '             open ocean restart NetCDF file: '//clname
[632]89            END SELECT
[1130]90            IF( kt == nitrst - 1 ) THEN   ;   WRITE(numout,*) '             kt = nitrst - 1 = ', kt
91            ELSE                          ;   WRITE(numout,*) '             kt = '             , kt
[611]92            ENDIF
93         ENDIF
94
[547]95         CALL iom_open( clname, numrow, ldwrt = .TRUE., kiolib = jprstlib )
[508]96         lrst_oce = .TRUE.
97      ENDIF
98      !
99   END SUBROUTINE rst_opn
100
101
[3]102   SUBROUTINE rst_write( kt )
103      !!---------------------------------------------------------------------
104      !!                   ***  ROUTINE rstwrite  ***
105      !!                     
[632]106      !! ** Purpose :   Write restart fields in the format corresponding to jprstlib
[3]107      !!
[508]108      !! ** Method  :   Write in numrow when kt == nitrst in NetCDF
[3]109      !!      file, save fields which are necessary for restart
110      !!----------------------------------------------------------------------
[508]111      INTEGER, INTENT(in) ::   kt   ! ocean time-step
[3]112      !!----------------------------------------------------------------------
[508]113      !                                                                     ! the begining of the run [s]
[544]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
[3]116
[508]117      ! prognostic variables
118      CALL iom_rstput( kt, nitrst, numrow, 'ub'     , ub      )   
119      CALL iom_rstput( kt, nitrst, numrow, 'vb'     , vb      )
120      CALL iom_rstput( kt, nitrst, numrow, 'tb'     , tb      )
121      CALL iom_rstput( kt, nitrst, numrow, 'sb'     , sb      )
122      CALL iom_rstput( kt, nitrst, numrow, 'rotb'   , rotb    )
123      CALL iom_rstput( kt, nitrst, numrow, 'hdivb'  , hdivb   )
124      CALL iom_rstput( kt, nitrst, numrow, 'un'     , un      )
125      CALL iom_rstput( kt, nitrst, numrow, 'vn'     , vn      )
[593]126      IF( lk_vvl ) CALL iom_rstput( kt, nitrst, numrow, 'wn'     , wn      )
[508]127      CALL iom_rstput( kt, nitrst, numrow, 'tn'     , tn      )
128      CALL iom_rstput( kt, nitrst, numrow, 'sn'     , sn      )
129      CALL iom_rstput( kt, nitrst, numrow, 'rotn'   , rotn    )
130      CALL iom_rstput( kt, nitrst, numrow, 'hdivn'  , hdivn   )
[632]131
[593]132      IF( nn_dynhpg_rst == 1 .OR. lk_vvl ) THEN
[544]133         CALL iom_rstput( kt, nitrst, numrow, 'rhd' , rhd  )
134         CALL iom_rstput( kt, nitrst, numrow, 'rhop', rhop )
135         IF( ln_zps ) THEN
136            CALL iom_rstput( kt, nitrst, numrow, 'gtu' , gtu )
137            CALL iom_rstput( kt, nitrst, numrow, 'gsu' , gsu )
138            CALL iom_rstput( kt, nitrst, numrow, 'gru' , gru )
139            CALL iom_rstput( kt, nitrst, numrow, 'gtv' , gtv )
140            CALL iom_rstput( kt, nitrst, numrow, 'gsv' , gsv )
141            CALL iom_rstput( kt, nitrst, numrow, 'grv' , grv )
142         ENDIF
143      ENDIF
144
[508]145      IF( kt == nitrst ) THEN
146         CALL iom_close( numrow )     ! close the restart file (only at last time step)
[579]147         IF( .NOT. lk_trdmld )   lrst_oce = .FALSE.
[3]148      ENDIF
[508]149      !
[3]150   END SUBROUTINE rst_write
151
152
153   SUBROUTINE rst_read
154      !!----------------------------------------------------------------------
155      !!                   ***  ROUTINE rst_read  ***
156      !!
[632]157      !! ** Purpose :   Read files for restart (format fixed by jprstlib )
[3]158      !!
[632]159      !! ** Method  :   Read the previous fields on the NetCDF/binary file
[3]160      !!      the first record indicates previous characterics
161      !!      after control with the present run, we read :
162      !!      - prognostic variables on the second record
163      !!      - elliptic solver arrays
[359]164      !!      - barotropic stream function arrays ("key_dynspg_rl" defined)
165      !!        or free surface arrays
[3]166      !!      - tke arrays (lk_zdftke=T)
167      !!      for this last three records,  the previous characteristics
168      !!      could be different with those used in the present run.
169      !!----------------------------------------------------------------------
[1130]170      REAL(wp) ::   zrdt, zrdttra1
[3]171      !!----------------------------------------------------------------------
172
[508]173      IF(lwp) THEN                                             ! Contol prints
174         WRITE(numout,*)
[632]175         SELECT CASE ( jprstlib )
[1130]176         CASE ( jpnf90    )   ;   WRITE(numout,*) 'rst_read : read oce NetCDF restart file'
177         CASE ( jprstdimg )   ;   WRITE(numout,*) 'rst_read : read oce binary restart file'
[632]178         END SELECT
[508]179         WRITE(numout,*) '~~~~~~~~'
[3]180      ENDIF
181
[547]182      CALL iom_open( 'restart', numror, kiolib = jprstlib )
[3]183
[544]184      ! Check dynamics and tracer time-step consistency and force Euler restart if changed
[746]185      IF( iom_varid( numror, 'rdt', ldstop = .FALSE. ) > 0 )   THEN
[544]186         CALL iom_get( numror, 'rdt', zrdt )
187         IF( zrdt /= rdt )   neuler = 0
188      ENDIF
[746]189      IF( iom_varid( numror, 'rdttra1', ldstop = .FALSE. ) > 0 )   THEN
[544]190         CALL iom_get( numror, 'rdttra1', zrdttra1 )
191         IF( zrdttra1 /= rdttra(1) )   neuler = 0
192      ENDIF
[508]193      !                                                       ! Read prognostic variables
[683]194      CALL iom_get( numror, jpdom_autoglo, 'ub'   , ub    )        ! before i-component velocity
195      CALL iom_get( numror, jpdom_autoglo, 'vb'   , vb    )        ! before j-component velocity
196      CALL iom_get( numror, jpdom_autoglo, 'tb'   , tb    )        ! before temperature
197      CALL iom_get( numror, jpdom_autoglo, 'sb'   , sb    )        ! before salinity
198      CALL iom_get( numror, jpdom_autoglo, 'rotb' , rotb  )        ! before curl
199      CALL iom_get( numror, jpdom_autoglo, 'hdivb', hdivb )        ! before horizontal divergence
200      CALL iom_get( numror, jpdom_autoglo, 'un'   , un    )        ! now    i-component velocity
201      CALL iom_get( numror, jpdom_autoglo, 'vn'   , vn    )        ! now    j-component velocity
202      IF( lk_vvl ) CALL iom_get( numror, jpdom_autoglo, 'wn'   , wn    )        ! now    k-component velocity
203      CALL iom_get( numror, jpdom_autoglo, 'tn'   , tn    )        ! now    temperature
204      CALL iom_get( numror, jpdom_autoglo, 'sn'   , sn    )        ! now    salinity
205      CALL iom_get( numror, jpdom_autoglo, 'rotn' , rotn  )        ! now    curl
206      CALL iom_get( numror, jpdom_autoglo, 'hdivn', hdivn )        ! now    horizontal divergence
[508]207
208
209      IF( neuler == 0 ) THEN                                  ! Euler restart (neuler=0)
210         tb   (:,:,:) = tn   (:,:,:)                             ! all before fields set to now field values
211         sb   (:,:,:) = sn   (:,:,:)
212         ub   (:,:,:) = un   (:,:,:)
213         vb   (:,:,:) = vn   (:,:,:)
214         rotb (:,:,:) = rotn (:,:,:)
215         hdivb(:,:,:) = hdivn(:,:,:)
[3]216      ENDIF
[508]217
[746]218      IF( iom_varid( numror, 'rhd', ldstop = .FALSE. ) > 0 ) THEN
[683]219         CALL iom_get( numror, jpdom_autoglo, 'rhd' , rhd  )
220         CALL iom_get( numror, jpdom_autoglo, 'rhop', rhop )
[544]221      ELSE
222         CALL eos( tb, sb, rhd, rhop )        ! before potential and in situ densities
223      ENDIF
[899]224      IF( ln_zps .AND. .NOT. lk_c1d ) THEN
[746]225         IF( iom_varid( numror, 'gtu', ldstop = .FALSE. ) > 0 ) THEN
[683]226            CALL iom_get( numror, jpdom_autoglo, 'gtu' , gtu )
227            CALL iom_get( numror, jpdom_autoglo, 'gsu' , gsu )
228            CALL iom_get( numror, jpdom_autoglo, 'gru' , gru )
229            CALL iom_get( numror, jpdom_autoglo, 'gtv' , gtv )
230            CALL iom_get( numror, jpdom_autoglo, 'gsv' , gsv )
231            CALL iom_get( numror, jpdom_autoglo, 'grv' , grv )
[544]232         ELSE
233            CALL zps_hde( nit000, tb , sb , rhd,   &  ! Partial steps: before Horizontal DErivative
234               &                  gtu, gsu, gru,   &  ! of t, s, rd at the bottom ocean level
235               &                  gtv, gsv, grv )
236         ENDIF
237      ENDIF
[508]238      !
239   END SUBROUTINE rst_read
[473]240
[3]241
242   !!=====================================================================
243END MODULE restart
Note: See TracBrowser for help on using the repository browser.