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/dev_r5518_Surge_Modelling/NEMOGCM/NEMO/OPA_SRC/IOM – NEMO

source: branches/UKMO/dev_r5518_Surge_Modelling/NEMOGCM/NEMO/OPA_SRC/IOM/restart.F90 @ 12963

Last change on this file since 12963 was 5954, checked in by rfurner, 9 years ago

added surge code from 2014_Surge_Modelling branch

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