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 on Ticket #374 – Attachment – NEMO

Ticket #374: restart.F90

File restart.F90, 14.1 KB (added by ed.blockley, 15 years ago)

restart.F90 with modified iom_open call

Line 
1MODULE restart
2   !!======================================================================
3   !!                     ***  MODULE  restart  ***
4   !! Ocean restart :  write the ocean restart file
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   !!----------------------------------------------------------------------
11
12   !!----------------------------------------------------------------------
13   !!   rst_opn    : open the ocean restart file
14   !!   rst_write  : write the ocean restart file
15   !!   rst_read   : read the ocean restart file
16   !!----------------------------------------------------------------------
17   USE dom_oce         ! ocean space and time domain
18   USE oce             ! ocean dynamics and tracers
19   USE phycst          ! physical constants
20   USE in_out_manager  ! I/O manager
21   USE iom             ! I/O module
22   USE c1d             ! re-initialization of u-v mask for the 1D configuration
23   USE zpshde          ! partial step: hor. derivative (zps_hde routine)
24   USE eosbn2          ! equation of state            (eos bn2 routine)
25   USE trdmld_oce      ! ocean active mixed layer tracers trends variables
26#if defined key_zdftke2
27   USE zdf_oce
28#endif
29
30   IMPLICIT NONE
31   PRIVATE
32
33   PUBLIC   rst_opn    ! routine called by step module
34   PUBLIC   rst_write  ! routine called by step module
35   PUBLIC   rst_read   ! routine called by opa  module
36
37   LOGICAL, PUBLIC ::   lrst_oce                  !: logical to control the oce restart write
38   INTEGER, PUBLIC ::   numror, numrow            !: logical unit for cean restart (read and write)
39
40   !! * Substitutions
41#  include "vectopt_loop_substitute.h90"
42   !!----------------------------------------------------------------------
43   !!  OPA 9.0 , LOCEAN-IPSL (2006)
44   !! $Id: restart.F90 1239 2008-12-31 09:35:35Z ctlod $
45   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
46   !!----------------------------------------------------------------------
47
48CONTAINS
49
50   SUBROUTINE rst_opn( kt )
51      !!---------------------------------------------------------------------
52      !!                   ***  ROUTINE rst_opn  ***
53      !!                     
54      !! ** Purpose : + initialization (should be read in the namelist) of nitrst
55      !!              + open the restart when we are one time step before nitrst
56      !!                   - restart header is defined when kt = nitrst-1
57      !!                   - restart data  are written when kt = nitrst
58      !!              + define lrst_oce to .TRUE. when we need to define or write the restart
59      !!----------------------------------------------------------------------
60      INTEGER, INTENT(in) ::   kt     ! ocean time-step
61      !!
62      CHARACTER(LEN=20)   ::   clkt     ! ocean time-step deine as a character
63      CHARACTER(LEN=50)   ::   clname   ! ice output restart file name
64      !!----------------------------------------------------------------------
65      !
66      IF( kt == nit000 ) THEN   ! default definitions
67         lrst_oce = .FALSE.   
68         nitrst = nitend
69#if defined key_zdftke2
70         nitrst_tke2 = nitrst + 1
71#endif
72      ENDIF
73      IF( MOD( kt - 1, nstock ) == 0 ) THEN   
74         ! we use kt - 1 and not kt - nit000 to keep the same periodicity from the beginning of the experiment
75         nitrst = kt + nstock - 1                  ! define the next value of nitrst for restart writing
76         IF( nitrst > nitend )   nitrst = nitend   ! make sure we write a restart at the end of the run
77      ENDIF
78#if defined key_zdftke2
79      IF ( nitrst_tke2 .NE. kt ) nitrst_tke2 = nitrst + 1
80#endif
81      ! to get better performances with NetCDF format:
82      ! we open and define the ocean restart file one time step before writing the data (-> at nitrst - 1)
83      ! except if we write ocean restart files every time step or if an ocean restart file was writen at nitend - 1
84      IF( kt == nitrst - 1 .OR. nstock == 1 .OR. ( kt == nitend .AND. .NOT. lrst_oce ) ) THEN
85         ! beware of the format used to write kt (default is i8.8, that should be large enough...)
86         IF( nitrst > 999999999 ) THEN   ;   WRITE(clkt, *       ) nitrst
87         ELSE                            ;   WRITE(clkt, '(i8.8)') nitrst
88         ENDIF
89         ! create the file
90         clname = TRIM(cexper)//"_"//TRIM(ADJUSTL(clkt))//"_"//TRIM(cn_ocerst_out)
91         IF(lwp) THEN
92            WRITE(numout,*)
93            SELECT CASE ( jprstlib )
94            CASE ( jprstdimg )   ;   WRITE(numout,*) '             open ocean restart binary file: '//clname
95            CASE DEFAULT         ;   WRITE(numout,*) '             open ocean restart NetCDF file: '//clname
96            END SELECT
97            IF( kt == nitrst - 1 ) THEN   ;   WRITE(numout,*) '             kt = nitrst - 1 = ', kt
98            ELSE                          ;   WRITE(numout,*) '             kt = '             , kt
99            ENDIF
100         ENDIF
101
102         CALL iom_open( clname, numrow, ldwrt = .TRUE., kiolib = jprstlib, ldclobber=.TRUE. )
103         lrst_oce = .TRUE.
104      ENDIF
105      !
106   END SUBROUTINE rst_opn
107
108
109   SUBROUTINE rst_write( kt )
110      !!---------------------------------------------------------------------
111      !!                   ***  ROUTINE rstwrite  ***
112      !!                     
113      !! ** Purpose :   Write restart fields in the format corresponding to jprstlib
114      !!
115      !! ** Method  :   Write in numrow when kt == nitrst in NetCDF
116      !!      file, save fields which are necessary for restart
117      !!----------------------------------------------------------------------
118      INTEGER, INTENT(in) ::   kt   ! ocean time-step
119      !!----------------------------------------------------------------------
120
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
127      !                                                                     ! the begining of the run [s]
128      CALL iom_rstput( kt, nitrst, numrow, 'rdt'    , rdt               )   ! dynamics time step
129      CALL iom_rstput( kt, nitrst, numrow, 'rdttra1', rdttra(1)         )   ! surface tracer time step
130
131      ! prognostic variables
132      CALL iom_rstput( kt, nitrst, numrow, 'ub'     , ub      )   
133      CALL iom_rstput( kt, nitrst, numrow, 'vb'     , vb      )
134      CALL iom_rstput( kt, nitrst, numrow, 'tb'     , tb      )
135      CALL iom_rstput( kt, nitrst, numrow, 'sb'     , sb      )
136      CALL iom_rstput( kt, nitrst, numrow, 'rotb'   , rotb    )
137      CALL iom_rstput( kt, nitrst, numrow, 'hdivb'  , hdivb   )
138      CALL iom_rstput( kt, nitrst, numrow, 'un'     , un      )
139      CALL iom_rstput( kt, nitrst, numrow, 'vn'     , vn      )
140      IF( lk_vvl ) CALL iom_rstput( kt, nitrst, numrow, 'wn'     , wn      )
141      CALL iom_rstput( kt, nitrst, numrow, 'tn'     , tn      )
142      CALL iom_rstput( kt, nitrst, numrow, 'sn'     , sn      )
143      CALL iom_rstput( kt, nitrst, numrow, 'rotn'   , rotn    )
144      CALL iom_rstput( kt, nitrst, numrow, 'hdivn'  , hdivn   )
145
146#if defined key_zdftke2
147         CALL iom_rstput( kt, nitrst, numrow, 'rhop'   , rhop    )
148#endif
149      IF( nn_dynhpg_rst == 1 .OR. lk_vvl ) THEN
150         CALL iom_rstput( kt, nitrst, numrow, 'rhd' , rhd  )
151#if defined key_zdftke2
152         CALL iom_rstput( kt, nitrst, numrow, 'rn2' , rn2  )
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
168      IF( kt == nitrst ) THEN
169         CALL iom_close( numrow )     ! close the restart file (only at last time step)
170         IF( .NOT. lk_trdmld )   lrst_oce = .FALSE.
171      ENDIF
172#endif
173      !
174   END SUBROUTINE rst_write
175
176
177   SUBROUTINE rst_read
178      !!----------------------------------------------------------------------
179      !!                   ***  ROUTINE rst_read  ***
180      !!
181      !! ** Purpose :   Read files for restart (format fixed by jprstlib )
182      !!
183      !! ** Method  :   Read the previous fields on the NetCDF/binary file
184      !!      the first record indicates previous characterics
185      !!      after control with the present run, we read :
186      !!      - prognostic variables on the second record
187      !!      - elliptic solver arrays
188      !!      - barotropic stream function arrays ("key_dynspg_rl" defined)
189      !!        or free surface arrays
190      !!      - tke arrays (lk_zdftke=T .OR. lk_zdftke2=T)
191      !!      for this last three records,  the previous characteristics
192      !!      could be different with those used in the present run.
193      !!----------------------------------------------------------------------
194      REAL(wp) ::   zrdt, zrdttra1
195      !!----------------------------------------------------------------------
196
197      IF(lwp) THEN                                             ! Contol prints
198         WRITE(numout,*)
199         SELECT CASE ( jprstlib )
200         CASE ( jpnf90    )   ;   WRITE(numout,*) 'rst_read : read oce NetCDF restart file'
201         CASE ( jprstdimg )   ;   WRITE(numout,*) 'rst_read : read oce binary restart file'
202         END SELECT
203         WRITE(numout,*) '~~~~~~~~'
204      ENDIF
205
206      CALL iom_open( cn_ocerst_in, numror, kiolib = jprstlib )
207
208      ! Check dynamics and tracer time-step consistency and force Euler restart if changed
209      IF( iom_varid( numror, 'rdt', ldstop = .FALSE. ) > 0 )   THEN
210         CALL iom_get( numror, 'rdt', zrdt )
211         IF( zrdt /= rdt )   neuler = 0
212      ENDIF
213      IF( iom_varid( numror, 'rdttra1', ldstop = .FALSE. ) > 0 )   THEN
214         CALL iom_get( numror, 'rdttra1', zrdttra1 )
215         IF( zrdttra1 /= rdttra(1) )   neuler = 0
216      ENDIF
217      !                                                       ! Read prognostic variables
218      CALL iom_get( numror, jpdom_autoglo, 'ub'   , ub    )        ! before i-component velocity
219      CALL iom_get( numror, jpdom_autoglo, 'vb'   , vb    )        ! before j-component velocity
220      CALL iom_get( numror, jpdom_autoglo, 'tb'   , tb    )        ! before temperature
221      CALL iom_get( numror, jpdom_autoglo, 'sb'   , sb    )        ! before salinity
222      CALL iom_get( numror, jpdom_autoglo, 'rotb' , rotb  )        ! before curl
223      CALL iom_get( numror, jpdom_autoglo, 'hdivb', hdivb )        ! before horizontal divergence
224      CALL iom_get( numror, jpdom_autoglo, 'un'   , un    )        ! now    i-component velocity
225      CALL iom_get( numror, jpdom_autoglo, 'vn'   , vn    )        ! now    j-component velocity
226      IF( lk_vvl ) CALL iom_get( numror, jpdom_autoglo, 'wn'   , wn    )        ! now    k-component velocity
227      CALL iom_get( numror, jpdom_autoglo, 'tn'   , tn    )        ! now    temperature
228      CALL iom_get( numror, jpdom_autoglo, 'sn'   , sn    )        ! now    salinity
229      CALL iom_get( numror, jpdom_autoglo, 'rotn' , rotn  )        ! now    curl
230      CALL iom_get( numror, jpdom_autoglo, 'hdivn', hdivn )        ! now    horizontal divergence
231
232
233      IF( neuler == 0 ) THEN                                  ! Euler restart (neuler=0)
234         tb   (:,:,:) = tn   (:,:,:)                             ! all before fields set to now field values
235         sb   (:,:,:) = sn   (:,:,:)
236         ub   (:,:,:) = un   (:,:,:)
237         vb   (:,:,:) = vn   (:,:,:)
238         rotb (:,:,:) = rotn (:,:,:)
239         hdivb(:,:,:) = hdivn(:,:,:)
240      ENDIF
241
242#if defined key_zdftke2
243      CALL eos( tb, sb, rhd, rhop )        ! before potential and in situ densities
244      IF( iom_varid( numror, 'rhd', ldstop = .FALSE. ) > 0 ) THEN
245         CALL iom_get( numror, jpdom_autoglo, 'rhd' , rhd  )
246      ENDIF
247      IF( iom_varid( numror, 'rhop', ldstop = .FALSE. ) > 0 ) THEN
248         CALL iom_get( numror, jpdom_autoglo, 'rhop', rhop )
249      ENDIF
250#else
251      IF( iom_varid( numror, 'rhd', ldstop = .FALSE. ) > 0 ) THEN
252         CALL iom_get( numror, jpdom_autoglo, 'rhd' , rhd  )
253         CALL iom_get( numror, jpdom_autoglo, 'rhop', rhop )
254      ELSE
255         CALL eos( tb, sb, rhd, rhop )        ! before potential and in situ densities
256      ENDIF
257#endif
258#if defined key_zdftke2
259      CALL eos_init            ! usefull to get the equation state type neos parameter
260      IF( iom_varid( numror, 'rn2', ldstop = .FALSE. ) > 0 ) THEN
261         CALL iom_get( numror, jpdom_autoglo, 'rn2' , rn2  )
262      ELSE
263         IF ( ln_dynhpg_imp ) THEN
264            CALL bn2( tb, sb, rn2 )  ! before Brunt-Vaisala frequency   
265         ENDIF     
266      ENDIF
267      IF( iom_varid( numror, 'avt', ldstop = .FALSE. ) > 0 ) THEN
268         CALL iom_get( numror, jpdom_autoglo, 'avt' , avt  )
269      ELSE
270         IF ( ln_dynhpg_imp ) avt (:,:,:) = 1.2e-5 * tmask(:,:,:)
271      ENDIF
272#endif
273
274      IF( ln_zps .AND. .NOT. lk_c1d ) THEN
275         IF( iom_varid( numror, 'gtu', ldstop = .FALSE. ) > 0 ) THEN
276            CALL iom_get( numror, jpdom_autoglo, 'gtu' , gtu )
277            CALL iom_get( numror, jpdom_autoglo, 'gsu' , gsu )
278            CALL iom_get( numror, jpdom_autoglo, 'gru' , gru )
279            CALL iom_get( numror, jpdom_autoglo, 'gtv' , gtv )
280            CALL iom_get( numror, jpdom_autoglo, 'gsv' , gsv )
281            CALL iom_get( numror, jpdom_autoglo, 'grv' , grv )
282         ELSE
283            CALL zps_hde( nit000, tb , sb , rhd,   &  ! Partial steps: before Horizontal DErivative
284               &                  gtu, gsu, gru,   &  ! of t, s, rd at the bottom ocean level
285               &                  gtv, gsv, grv )
286         ENDIF
287      ENDIF
288      !
289   END SUBROUTINE rst_read
290
291
292   !!=====================================================================
293END MODULE restart