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 branches/UKMO/r6232_collate_bgc_diagnostics/NEMOGCM/NEMO/OPA_SRC/IOM – NEMO

source: branches/UKMO/r6232_collate_bgc_diagnostics/NEMOGCM/NEMO/OPA_SRC/IOM/restart.F90 @ 11134

Last change on this file since 11134 was 11134, checked in by jcastill, 5 years ago

Full set of changes as in the original branch

File size: 14.4 KB
RevLine 
[3]1MODULE restart
2   !!======================================================================
3   !!                     ***  MODULE  restart  ***
4   !! Ocean restart :  write the ocean restart file
[508]5   !!======================================================================
[2528]6   !! History :  OPA  !  1999-11  (M. Imbard)  Original code
7   !!   NEMO     1.0  !  2002-08  (G. Madec)  F90: Free form
8   !!            2.0  !  2006-07  (S. Masson)  use IOM for restart
9   !!            3.3  !  2010-04  (M. Leclair, G. Madec)  modified LF-RA
10   !!            - -  !  2010-10  (C. Ethe, G. Madec) TRC-TRA merge (T-S in 4D)
[508]11   !!----------------------------------------------------------------------
[3]12
13   !!----------------------------------------------------------------------
[508]14   !!   rst_opn    : open the ocean restart file
15   !!   rst_write  : write the ocean restart file
16   !!   rst_read   : read the ocean restart file
[3]17   !!----------------------------------------------------------------------
[2528]18   USE oce             ! ocean dynamics and tracers
[3]19   USE dom_oce         ! ocean space and time domain
20   USE phycst          ! physical constants
[508]21   USE in_out_manager  ! I/O manager
22   USE iom             ! I/O module
[11134]23   USE ioipsl, ONLY : ju2ymds    ! for calendar
[544]24   USE eosbn2          ! equation of state            (eos bn2 routine)
[4990]25   USE trdmxl_oce      ! ocean active mixed layer tracers trends variables
[3680]26   USE divcur          ! hor. divergence and curl      (div & cur routines)
[3]27
28   IMPLICIT NONE
29   PRIVATE
30
[4292]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 istate module
34   PUBLIC   rst_read_open   ! routine called in rst_read and (possibly) in dom_vvl_init
[3]35
[508]36   !! * Substitutions
[2528]37#  include "domzgr_substitute.h90"
[508]38#  include "vectopt_loop_substitute.h90"
[3]39   !!----------------------------------------------------------------------
[2528]40   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
[888]41   !! $Id$
[2528]42   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[359]43   !!----------------------------------------------------------------------
[3]44CONTAINS
45
[508]46   SUBROUTINE rst_opn( kt )
47      !!---------------------------------------------------------------------
48      !!                   ***  ROUTINE rst_opn  ***
49      !!                     
50      !! ** Purpose : + initialization (should be read in the namelist) of nitrst
51      !!              + open the restart when we are one time step before nitrst
52      !!                   - restart header is defined when kt = nitrst-1
53      !!                   - restart data  are written when kt = nitrst
54      !!              + define lrst_oce to .TRUE. when we need to define or write the restart
55      !!----------------------------------------------------------------------
56      INTEGER, INTENT(in) ::   kt     ! ocean time-step
[11134]57      INTEGER             ::   iyear, imonth, iday
58      REAL (wp)           ::   zsec
59      REAL (wp)           ::   zfjulday
[508]60      !!
61      CHARACTER(LEN=20)   ::   clkt     ! ocean time-step deine as a character
[5341]62      CHARACTER(LEN=50)   ::   clname   ! ocean output restart file name
[11134]63      CHARACTER(LEN=150)  ::   clpath   ! full path to ocean output restart file
[508]64      !!----------------------------------------------------------------------
65      !
[783]66      IF( kt == nit000 ) THEN   ! default definitions
67         lrst_oce = .FALSE.   
[5341]68         IF( ln_rst_list ) THEN
69            nrst_lst = 1
70            nitrst = nstocklist( nrst_lst )
71         ELSE
72            nitrst = nitend
73         ENDIF
[508]74      ENDIF
[5341]75
76      ! frequency-based restart dumping (nn_stock)
77      IF( .NOT. ln_rst_list .AND. MOD( kt - 1, nstock ) == 0 ) THEN   
[783]78         ! we use kt - 1 and not kt - nit000 to keep the same periodicity from the beginning of the experiment
79         nitrst = kt + nstock - 1                  ! define the next value of nitrst for restart writing
80         IF( nitrst > nitend )   nitrst = nitend   ! make sure we write a restart at the end of the run
81      ENDIF
82      ! to get better performances with NetCDF format:
83      ! we open and define the ocean restart file one time step before writing the data (-> at nitrst - 1)
84      ! except if we write ocean restart files every time step or if an ocean restart file was writen at nitend - 1
85      IF( kt == nitrst - 1 .OR. nstock == 1 .OR. ( kt == nitend .AND. .NOT. lrst_oce ) ) THEN
[5341]86         IF( nitrst <= nitend .AND. nitrst > 0 ) THEN
[11134]87            IF ( ln_rstdate ) THEN
88               zfjulday = fjulday + rdttra(1) / rday
89               IF( ABS(zfjulday - REAL(NINT(zfjulday),wp)) < 0.1 / rday )   zfjulday = REAL(NINT(zfjulday),wp)   ! avoid truncation error
90               CALL ju2ymds( zfjulday, iyear, imonth, iday, zsec )           
91               WRITE(clkt, '(i4.4,2i2.2)') iyear, imonth, iday
92            ELSE
93               ! beware of the format used to write kt (default is i8.8, that should be large enough...)
94               IF( nitrst > 999999999 ) THEN   ;   WRITE(clkt, *       ) nitrst
95               ELSE                            ;   WRITE(clkt, '(i8.8)') nitrst
96               ENDIF
[611]97            ENDIF
[5341]98            ! create the file
99            clname = TRIM(cexper)//"_"//TRIM(ADJUSTL(clkt))//"_"//TRIM(cn_ocerst_out)
100            clpath = TRIM(cn_ocerst_outdir)
101            IF( clpath(LEN_TRIM(clpath):) /= '/' ) clpath = TRIM(clpath) // '/'
102            IF(lwp) THEN
103               WRITE(numout,*)
104               SELECT CASE ( jprstlib )
105               CASE ( jprstdimg )   ;   WRITE(numout,*)                            &
106                   '             open ocean restart binary file: ',TRIM(clpath)//clname
107               CASE DEFAULT         ;   WRITE(numout,*)                            &
108                   '             open ocean restart NetCDF file: ',TRIM(clpath)//clname
109               END SELECT
110               IF ( snc4set%luse )      WRITE(numout,*) '             opened for NetCDF4 chunking and compression'
111               IF( kt == nitrst - 1 ) THEN   ;   WRITE(numout,*) '             kt = nitrst - 1 = ', kt
112               ELSE                          ;   WRITE(numout,*) '             kt = '             , kt
113               ENDIF
114            ENDIF
115            !
116            CALL iom_open( TRIM(clpath)//TRIM(clname), numrow, ldwrt = .TRUE., kiolib = jprstlib )
117            lrst_oce = .TRUE.
[611]118         ENDIF
[508]119      ENDIF
120      !
121   END SUBROUTINE rst_opn
122
123
[3]124   SUBROUTINE rst_write( kt )
125      !!---------------------------------------------------------------------
126      !!                   ***  ROUTINE rstwrite  ***
127      !!                     
[632]128      !! ** Purpose :   Write restart fields in the format corresponding to jprstlib
[3]129      !!
[508]130      !! ** Method  :   Write in numrow when kt == nitrst in NetCDF
[2528]131      !!              file, save fields which are necessary for restart
[3]132      !!----------------------------------------------------------------------
[508]133      INTEGER, INTENT(in) ::   kt   ! ocean time-step
[3]134      !!----------------------------------------------------------------------
[1239]135
[2528]136                     CALL iom_rstput( kt, nitrst, numrow, 'rdt'    , rdt       )   ! dynamics time step
137                     CALL iom_rstput( kt, nitrst, numrow, 'rdttra1', rdttra(1) )   ! surface tracer time step
[3]138
[2528]139                     CALL iom_rstput( kt, nitrst, numrow, 'ub'     , ub        )     ! before fields
140                     CALL iom_rstput( kt, nitrst, numrow, 'vb'     , vb        )
[3294]141                     CALL iom_rstput( kt, nitrst, numrow, 'tb'     , tsb(:,:,:,jp_tem) )
142                     CALL iom_rstput( kt, nitrst, numrow, 'sb'     , tsb(:,:,:,jp_sal) )
[2528]143                     CALL iom_rstput( kt, nitrst, numrow, 'rotb'   , rotb      )
144                     CALL iom_rstput( kt, nitrst, numrow, 'hdivb'  , hdivb     )
145                     CALL iom_rstput( kt, nitrst, numrow, 'sshb'   , sshb      )
[4990]146                     !
[2528]147                     CALL iom_rstput( kt, nitrst, numrow, 'un'     , un        )     ! now fields
148                     CALL iom_rstput( kt, nitrst, numrow, 'vn'     , vn        )
[3294]149                     CALL iom_rstput( kt, nitrst, numrow, 'tn'     , tsn(:,:,:,jp_tem) )
150                     CALL iom_rstput( kt, nitrst, numrow, 'sn'     , tsn(:,:,:,jp_sal) )
[2528]151                     CALL iom_rstput( kt, nitrst, numrow, 'rotn'   , rotn      )
152                     CALL iom_rstput( kt, nitrst, numrow, 'hdivn'  , hdivn     )
153                     CALL iom_rstput( kt, nitrst, numrow, 'sshn'   , sshn      )
154                     CALL iom_rstput( kt, nitrst, numrow, 'rhop'   , rhop      )
[1545]155#if defined key_zdfkpp
[2528]156                     CALL iom_rstput( kt, nitrst, numrow, 'rhd'    , rhd       )
[1483]157#endif
[508]158      IF( kt == nitrst ) THEN
159         CALL iom_close( numrow )     ! close the restart file (only at last time step)
[4990]160!!gm         IF( .NOT. lk_trdmld )   lrst_oce = .FALSE.
161!!gm  not sure what to do here   ===>>>  ask to Sebastian
162         lrst_oce = .FALSE.
[5341]163            IF( ln_rst_list ) THEN
164               nrst_lst = MIN(nrst_lst + 1, SIZE(nstocklist,1))
165               nitrst = nstocklist( nrst_lst )
166            ENDIF
167            lrst_oce = .FALSE.
[3]168      ENDIF
[508]169      !
[3]170   END SUBROUTINE rst_write
171
[4990]172
[4292]173   SUBROUTINE rst_read_open
174      !!----------------------------------------------------------------------
175      !!                   ***  ROUTINE rst_read_open  ***
176      !!
177      !! ** Purpose :   Open read files for restart (format fixed by jprstlib )
178      !!
179      !! ** Method  :   Use a non-zero, positive value of numror to assess whether or not
180      !!                the file has already been opened
181      !!----------------------------------------------------------------------
[5341]182      INTEGER        ::   jlibalt = jprstlib
183      LOGICAL        ::   llok
184      CHARACTER(lc)  ::   clpath   ! full path to ocean output restart file
[4292]185      !!----------------------------------------------------------------------
[4990]186      !
187      IF( numror <= 0 ) THEN
[4292]188         IF(lwp) THEN                                             ! Contol prints
189            WRITE(numout,*)
190            SELECT CASE ( jprstlib )
191            CASE ( jpnf90    )   ;   WRITE(numout,*) 'rst_read : read oce NetCDF restart file'
192            CASE ( jprstdimg )   ;   WRITE(numout,*) 'rst_read : read oce binary restart file'
193            END SELECT
194            IF ( snc4set%luse )      WRITE(numout,*) 'rst_read : configured with NetCDF4 support'
195            WRITE(numout,*) '~~~~~~~~'
196         ENDIF
197
[5341]198         clpath = TRIM(cn_ocerst_indir)
199         IF( clpath(LEN_TRIM(clpath):) /= '/' ) clpath = TRIM(clpath) // '/'
[4292]200         IF ( jprstlib == jprstdimg ) THEN
201           ! eventually read netcdf file (monobloc)  for restarting on different number of processors
202           ! if {cn_ocerst_in}.nc exists, then set jlibalt to jpnf90
[5341]203           INQUIRE( FILE = TRIM(cn_ocerst_indir)//'/'//TRIM(cn_ocerst_in)//'.nc', EXIST = llok )
[4292]204           IF ( llok ) THEN ; jlibalt = jpnf90  ; ELSE ; jlibalt = jprstlib ; ENDIF
205         ENDIF
[5341]206         CALL iom_open( TRIM(clpath)//cn_ocerst_in, numror, kiolib = jlibalt )
[4292]207      ENDIF
208   END SUBROUTINE rst_read_open
209
[3]210   SUBROUTINE rst_read
211      !!----------------------------------------------------------------------
212      !!                   ***  ROUTINE rst_read  ***
213      !!
[632]214      !! ** Purpose :   Read files for restart (format fixed by jprstlib )
[3]215      !!
[1531]216      !! ** Method  :   Read in restart.nc file fields which are necessary for restart
[3]217      !!----------------------------------------------------------------------
[1130]218      REAL(wp) ::   zrdt, zrdttra1
[4292]219      INTEGER  ::   jk
[1473]220      LOGICAL  ::   llok
[3]221      !!----------------------------------------------------------------------
222
[4292]223      CALL rst_read_open           ! open restart for reading (if not already opened)
[3]224
[544]225      ! Check dynamics and tracer time-step consistency and force Euler restart if changed
[746]226      IF( iom_varid( numror, 'rdt', ldstop = .FALSE. ) > 0 )   THEN
[544]227         CALL iom_get( numror, 'rdt', zrdt )
228         IF( zrdt /= rdt )   neuler = 0
229      ENDIF
[746]230      IF( iom_varid( numror, 'rdttra1', ldstop = .FALSE. ) > 0 )   THEN
[544]231         CALL iom_get( numror, 'rdttra1', zrdttra1 )
232         IF( zrdttra1 /= rdttra(1) )   neuler = 0
233      ENDIF
[1607]234      !
[3680]235      IF( iom_varid( numror, 'ub', ldstop = .FALSE. ) > 0 ) THEN
236         CALL iom_get( numror, jpdom_autoglo, 'ub'     , ub      )   ! before fields
237         CALL iom_get( numror, jpdom_autoglo, 'vb'     , vb      )
238         CALL iom_get( numror, jpdom_autoglo, 'tb'     , tsb(:,:,:,jp_tem) )
239         CALL iom_get( numror, jpdom_autoglo, 'sb'     , tsb(:,:,:,jp_sal) )
240         CALL iom_get( numror, jpdom_autoglo, 'rotb'   , rotb    )
241         CALL iom_get( numror, jpdom_autoglo, 'hdivb'  , hdivb   )
242         CALL iom_get( numror, jpdom_autoglo, 'sshb'   , sshb    )
243      ELSE
244         neuler = 0
245      ENDIF
246      !
247      CALL iom_get( numror, jpdom_autoglo, 'un'     , un      )   ! now    fields
248      CALL iom_get( numror, jpdom_autoglo, 'vn'     , vn      )
249      CALL iom_get( numror, jpdom_autoglo, 'tn'     , tsn(:,:,:,jp_tem) )
250      CALL iom_get( numror, jpdom_autoglo, 'sn'     , tsn(:,:,:,jp_sal) )
251      CALL iom_get( numror, jpdom_autoglo, 'sshn'   , sshn    )
252      IF( iom_varid( numror, 'rotn', ldstop = .FALSE. ) > 0 ) THEN
253         CALL iom_get( numror, jpdom_autoglo, 'rotn'   , rotn    )
254         CALL iom_get( numror, jpdom_autoglo, 'hdivn'  , hdivn   )
255      ELSE
256         CALL div_cur( 0 )                              ! Horizontal divergence & Relative vorticity
257      ENDIF
258      IF( iom_varid( numror, 'rhop', ldstop = .FALSE. ) > 0 ) THEN
259         CALL iom_get( numror, jpdom_autoglo, 'rhop'   , rhop    )   ! now    potential density
260      ELSE
[4313]261         CALL eos    ( tsn, rhd, rhop, fsdept_n(:,:,:) )   
[3680]262      ENDIF
[1545]263#if defined key_zdfkpp
264      IF( iom_varid( numror, 'rhd', ldstop = .FALSE. ) > 0 ) THEN
[3680]265         CALL iom_get( numror, jpdom_autoglo, 'rhd'    , rhd     )   ! now    in situ density anomaly
[1545]266      ELSE
[4313]267         CALL eos( tsn, rhd, fsdept_n(:,:,:) )   ! compute rhd
[1545]268      ENDIF
[1483]269#endif
[2528]270      !
[508]271      IF( neuler == 0 ) THEN                                  ! Euler restart (neuler=0)
[3294]272         tsb  (:,:,:,:) = tsn  (:,:,:,:)                             ! all before fields set to now values
273         ub   (:,:,:)   = un   (:,:,:)
274         vb   (:,:,:)   = vn   (:,:,:)
275         rotb (:,:,:)   = rotn (:,:,:)
276         hdivb(:,:,:)   = hdivn(:,:,:)
277         sshb (:,:)     = sshn (:,:)
[4990]278
279         IF( lk_vvl ) THEN
[4689]280            DO jk = 1, jpk
281               fse3t_b(:,:,jk) = fse3t_n(:,:,jk)
282            END DO
283         ENDIF
[4990]284
[3]285      ENDIF
[508]286      !
287   END SUBROUTINE rst_read
[473]288
[3]289   !!=====================================================================
290END MODULE restart
Note: See TracBrowser for help on using the repository browser.