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.
icbrst.F90 in branches/UKMO/r8395_cpl_tauwav/NEMOGCM/NEMO/OPA_SRC/ICB – NEMO

source: branches/UKMO/r8395_cpl_tauwav/NEMOGCM/NEMO/OPA_SRC/ICB/icbrst.F90 @ 12286

Last change on this file since 12286 was 12286, checked in by jcastill, 4 years ago

Remove svn keywords

File size: 18.6 KB
Line 
1MODULE icbrst
2
3   !!======================================================================
4   !!                       ***  MODULE  icbrst  ***
5   !! Ocean physics:  read and write iceberg restart files
6   !!======================================================================
7   !! History : 3.3.1 !  2010-01  (Martin&Adcroft) Original code
8   !!            -    !  2011-03  (Madec)          Part conversion to NEMO form
9   !!            -    !                            Removal of mapping from another grid
10   !!            -    !  2011-04  (Alderson)       Split into separate modules
11   !!            -    !  2011-04  (Alderson)       Restore restart routine
12   !!            -    !                            Currently needs a fixed processor
13   !!            -    !                            layout between restarts
14   !!            -    !  2015-11  Dave Storkey     Convert icb_rst_read to use IOM so can
15   !!                                              read single restart files
16   !!----------------------------------------------------------------------
17   !!----------------------------------------------------------------------
18   !!   icb_rst_read    : read restart file
19   !!   icb_rst_write   : write restart file
20   !!----------------------------------------------------------------------
21   USE par_oce        ! NEMO parameters
22   USE dom_oce        ! NEMO domain
23   USE in_out_manager ! NEMO IO routines
24   USE lib_mpp        ! NEMO MPI library, lk_mpp in particular
25   USE netcdf         ! netcdf routines for IO
26   USE iom
27   USE icb_oce        ! define iceberg arrays
28   USE icbutl         ! iceberg utility routines
29
30   IMPLICIT NONE
31   PRIVATE
32
33   PUBLIC   icb_rst_read    ! routine called in icbini.F90 module
34   PUBLIC   icb_rst_write   ! routine called in icbstp.F90 module
35   
36   INTEGER ::   nlonid, nlatid, nxid, nyid, nuvelid, nvvelid
37   INTEGER ::   nmassid, nthicknessid, nwidthid, nlengthid
38   INTEGER ::   nyearid, ndayid
39   INTEGER ::   nscaling_id, nmass_of_bits_id, nheat_density_id, numberid
40   INTEGER ::   nsiceid, nsheatid, ncalvid, ncalvhid, nkountid
41   INTEGER ::   nret, ncid, nc_dim
42   
43   INTEGER,  DIMENSION(3)                  :: nstrt3, nlngth3
44
45   !!----------------------------------------------------------------------
46   !! NEMO/OPA 3.3 , NEMO Consortium (2011)
47   !! $Id$
48   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
49   !!----------------------------------------------------------------------
50CONTAINS
51
52   SUBROUTINE icb_rst_read()
53      !!----------------------------------------------------------------------
54      !!                 ***  SUBROUTINE icb_rst_read  ***
55      !!
56      !! ** Purpose :   read a iceberg restart file
57      !!      NB: for this version, we just read back in the restart for this processor
58      !!      so we cannot change the processor layout currently with iceberg code
59      !!----------------------------------------------------------------------
60      INTEGER                      ::   idim, ivar, iatt
61      INTEGER                      ::   jn, iunlim_dim, ibergs_in_file
62      INTEGER                      ::   ii,ij,iclass
63      REAL(wp), DIMENSION(nkounts) ::   zdata     
64      LOGICAL                      ::   ll_found_restart
65      CHARACTER(len=256)           ::   cl_path
66      CHARACTER(len=256)           ::   cl_filename
67      CHARACTER(len=NF90_MAX_NAME) ::   cl_dname
68      TYPE(iceberg)                ::   localberg ! NOT a pointer but an actual local variable
69      TYPE(point)                  ::   localpt   ! NOT a pointer but an actual local variable
70      !!----------------------------------------------------------------------
71
72      ! Find a restart file. Assume iceberg restarts in same directory as ocean restarts
73      ! and are called TRIM(cn_ocerst)//'_icebergs'
74      cl_path = TRIM(cn_ocerst_indir)
75      IF( cl_path(LEN_TRIM(cl_path):) /= '/' ) cl_path = TRIM(cl_path) // '/'
76      cl_filename = TRIM(cn_ocerst_in)//'_icebergs'
77      CALL iom_open( TRIM(cl_path)//cl_filename, ncid )
78
79      IF( iom_file(ncid)%iduld .GE. 0) THEN
80
81         ibergs_in_file = iom_file(ncid)%lenuld
82         DO jn = 1,ibergs_in_file
83
84            ! iom_get treats the unlimited dimension as time. Here the unlimited dimension
85            ! is the iceberg index, but we can still use the ktime keyword to get the iceberg we want.
86
87            CALL iom_get( ncid, 'xi'     ,localpt%xi  , ktime=jn )
88            CALL iom_get( ncid, 'yj'     ,localpt%yj  , ktime=jn )
89
90            ii = INT( localpt%xi + 0.5 )
91            ij = INT( localpt%yj + 0.5 )
92            ! Only proceed if this iceberg is on the local processor (excluding halos).
93            IF ( ii .GE. nldi+nimpp-1 .AND. ii .LE. nlei+nimpp-1 .AND. &
94           &     ij .GE. nldj+njmpp-1 .AND. ij .LE. nlej+njmpp-1 ) THEN           
95
96               CALL iom_get( ncid, jpdom_unknown, 'number'       , zdata(:) , ktime=jn, kstart=(/1/), kcount=(/nkounts/) )
97               localberg%number(:) = INT(zdata(:))
98               CALL iom_get( ncid, 'mass_scaling' , localberg%mass_scaling, ktime=jn )
99               CALL iom_get( ncid, 'lon'          , localpt%lon           , ktime=jn )
100               CALL iom_get( ncid, 'lat'          , localpt%lat           , ktime=jn )
101               CALL iom_get( ncid, 'uvel'         , localpt%uvel          , ktime=jn )
102               CALL iom_get( ncid, 'vvel'         , localpt%vvel          , ktime=jn )
103               CALL iom_get( ncid, 'mass'         , localpt%mass          , ktime=jn )
104               CALL iom_get( ncid, 'thickness'    , localpt%thickness     , ktime=jn )
105               CALL iom_get( ncid, 'width'        , localpt%width         , ktime=jn )
106               CALL iom_get( ncid, 'length'       , localpt%length        , ktime=jn )
107               CALL iom_get( ncid, 'year'         , zdata(1)              , ktime=jn )
108               localpt%year = INT(zdata(1))
109               CALL iom_get( ncid, 'day'          , localpt%day           , ktime=jn )
110               CALL iom_get( ncid, 'mass_of_bits' , localpt%mass_of_bits  , ktime=jn )
111               CALL iom_get( ncid, 'heat_density' , localpt%heat_density  , ktime=jn )
112
113               !
114               CALL icb_utl_add( localberg, localpt )
115
116            ENDIF
117
118         END DO
119
120      ENDIF 
121
122      ! Gridded variables
123      CALL iom_get( ncid, jpdom_autoglo,    'calving'     , src_calving  )
124      CALL iom_get( ncid, jpdom_autoglo,    'calving_hflx', src_calving_hflx  )
125      CALL iom_get( ncid, jpdom_autoglo,    'stored_heat' , berg_grid%stored_heat  )
126      CALL iom_get( ncid, jpdom_autoglo_xy, 'stored_ice'  , berg_grid%stored_ice, kstart=(/1,1,1/), kcount=(/1,1,nclasses/) )
127     
128      CALL iom_get( ncid, jpdom_unknown, 'kount' , zdata(:) )
129      num_bergs(:) = INT(zdata(:))
130
131      ! Sanity check
132      jn = icb_utl_count()
133      IF (nn_verbose_level >= 0)   &
134         WRITE(numout,'(2(a,i5))') 'icebergs, read_restart_bergs: # bergs =',jn,' on PE',narea-1
135      IF( lk_mpp ) THEN
136         ! Only mpp_sum ibergs_in_file if we are reading from multiple restart files.
137         IF( INDEX(iom_file(ncid)%name,'icebergs.nc' ) .EQ. 0 ) CALL mpp_sum(ibergs_in_file)
138         CALL mpp_sum(jn)
139      ENDIF
140      IF(lwp)   WRITE(numout,'(a,i5,a,i5,a)') 'icebergs, read_restart_bergs: there were',ibergs_in_file,   &
141         &                                    ' bergs in the restart file and', jn,' bergs have been read'
142      !
143      ! Finish up
144      CALL iom_close( ncid )
145      !
146      IF( lwp .and. nn_verbose_level >= 0)  WRITE(numout,'(a)') 'icebergs, read_restart_bergs: completed'
147      !
148   END SUBROUTINE icb_rst_read
149
150
151   SUBROUTINE icb_rst_write( kt )
152      !!----------------------------------------------------------------------
153      !!                 ***  SUBROUTINE icb_rst_write  ***
154      !!
155      !!----------------------------------------------------------------------
156      INTEGER, INTENT( in ) :: kt
157      !
158      INTEGER ::   jn   ! dummy loop index
159      INTEGER ::   ix_dim, iy_dim, ik_dim, in_dim
160      CHARACTER(len=256)     :: cl_path
161      CHARACTER(len=256)     :: cl_filename
162      TYPE(iceberg), POINTER :: this
163      TYPE(point)  , POINTER :: pt
164      !!----------------------------------------------------------------------
165
166      ! Assume we write iceberg restarts to same directory as ocean restarts.
167      cl_path = TRIM(cn_ocerst_outdir)
168      IF( cl_path(LEN_TRIM(cl_path):) /= '/' ) cl_path = TRIM(cl_path) // '/'
169      IF( lk_mpp ) THEN
170         WRITE(cl_filename,'(A,"_icebergs_",I8.8,"_restart_",I4.4,".nc")') TRIM(cexper), kt, narea-1
171      ELSE
172         WRITE(cl_filename,'(A,"_icebergs_",I8.8,"_restart.nc")') TRIM(cexper), kt
173      ENDIF
174      IF (nn_verbose_level >= 0) WRITE(numout,'(2a)') 'icebergs, write_restart: creating ',TRIM(cl_path)//TRIM(cl_filename)
175
176      nret = NF90_CREATE(TRIM(cl_path)//TRIM(cl_filename), NF90_CLOBBER, ncid)
177      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_create failed')
178
179      ! Dimensions
180      nret = NF90_DEF_DIM(ncid, 'x', jpi, ix_dim)
181      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim x failed')
182
183      nret = NF90_DEF_DIM(ncid, 'y', jpj, iy_dim)
184      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim y failed')
185
186      nret = NF90_DEF_DIM(ncid, 'c', nclasses, nc_dim)
187      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim c failed')
188
189      nret = NF90_DEF_DIM(ncid, 'k', nkounts, ik_dim)
190      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim k failed')
191
192      ! global attributes
193      IF( lk_mpp ) THEN
194         ! Set domain parameters (assume jpdom_local_full)
195         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_number_total'   , jpnij              )
196         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_number'         , narea-1            )
197         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_dimensions_ids' , (/1     , 2     /) )
198         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_size_global'    , (/jpiglo, jpjglo/) )
199         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_size_local'     , (/jpi   , jpj   /) )
200         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_position_first' , (/nimpp , njmpp /) )
201         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_position_last'  , (/nimpp + jpi - 1 , njmpp + jpj - 1  /) )
202         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_halo_size_start', (/nldi - 1        , nldj - 1         /) )
203         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_halo_size_end'  , (/jpi - nlei      , jpj - nlej       /) )
204         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_type'           , 'BOX'              )
205      ENDIF
206     
207      IF (associated(first_berg)) then
208         nret = NF90_DEF_DIM(ncid, 'n', NF90_UNLIMITED, in_dim)
209         IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim n failed')
210      ENDIF
211
212      ! Variables
213      nret = NF90_DEF_VAR(ncid, 'kount'       , NF90_INT   , (/ ik_dim /), nkountid)
214      nret = NF90_DEF_VAR(ncid, 'calving'     , NF90_DOUBLE, (/ ix_dim, iy_dim /), ncalvid)
215      nret = NF90_DEF_VAR(ncid, 'calving_hflx', NF90_DOUBLE, (/ ix_dim, iy_dim /), ncalvhid)
216      nret = NF90_DEF_VAR(ncid, 'stored_ice'  , NF90_DOUBLE, (/ ix_dim, iy_dim, nc_dim /), nsiceid)
217      nret = NF90_DEF_VAR(ncid, 'stored_heat' , NF90_DOUBLE, (/ ix_dim, iy_dim /), nsheatid)
218
219      ! Attributes
220      nret = NF90_PUT_ATT(ncid, ncalvid , 'long_name', 'iceberg calving')
221      nret = NF90_PUT_ATT(ncid, ncalvid , 'units', 'some')
222      nret = NF90_PUT_ATT(ncid, ncalvhid, 'long_name', 'heat flux associated with iceberg calving')
223      nret = NF90_PUT_ATT(ncid, ncalvhid, 'units', 'some')
224      nret = NF90_PUT_ATT(ncid, nsiceid , 'long_name', 'stored ice used to calve icebergs')
225      nret = NF90_PUT_ATT(ncid, nsiceid , 'units', 'kg/s')
226      nret = NF90_PUT_ATT(ncid, nsheatid, 'long_name', 'heat in stored ice used to calve icebergs')
227      nret = NF90_PUT_ATT(ncid, nsheatid, 'units', 'J/kg/s')
228
229      IF ( ASSOCIATED(first_berg) ) THEN
230
231         ! Only add berg variables for this PE if we have anything to say
232
233         ! Variables
234         nret = NF90_DEF_VAR(ncid, 'lon', NF90_DOUBLE, in_dim, nlonid)
235         nret = NF90_DEF_VAR(ncid, 'lat', NF90_DOUBLE, in_dim, nlatid)
236         nret = NF90_DEF_VAR(ncid, 'xi', NF90_DOUBLE, in_dim, nxid)
237         nret = NF90_DEF_VAR(ncid, 'yj', NF90_DOUBLE, in_dim, nyid)
238         nret = NF90_DEF_VAR(ncid, 'uvel', NF90_DOUBLE, in_dim, nuvelid)
239         nret = NF90_DEF_VAR(ncid, 'vvel', NF90_DOUBLE, in_dim, nvvelid)
240         nret = NF90_DEF_VAR(ncid, 'mass', NF90_DOUBLE, in_dim, nmassid)
241         nret = NF90_DEF_VAR(ncid, 'thickness', NF90_DOUBLE, in_dim, nthicknessid)
242         nret = NF90_DEF_VAR(ncid, 'width', NF90_DOUBLE, in_dim, nwidthid)
243         nret = NF90_DEF_VAR(ncid, 'length', NF90_DOUBLE, in_dim, nlengthid)
244         nret = NF90_DEF_VAR(ncid, 'number', NF90_INT, (/ik_dim,in_dim/), numberid)
245         nret = NF90_DEF_VAR(ncid, 'year', NF90_INT, in_dim, nyearid)
246         nret = NF90_DEF_VAR(ncid, 'day', NF90_DOUBLE, in_dim, ndayid)
247         nret = NF90_DEF_VAR(ncid, 'mass_scaling', NF90_DOUBLE, in_dim, nscaling_id)
248         nret = NF90_DEF_VAR(ncid, 'mass_of_bits', NF90_DOUBLE, in_dim, nmass_of_bits_id)
249         nret = NF90_DEF_VAR(ncid, 'heat_density', NF90_DOUBLE, in_dim, nheat_density_id)
250
251         ! Attributes
252         nret = NF90_PUT_ATT(ncid, nlonid, 'long_name', 'longitude')
253         nret = NF90_PUT_ATT(ncid, nlonid, 'units', 'degrees_E')
254         nret = NF90_PUT_ATT(ncid, nlatid, 'long_name', 'latitude')
255         nret = NF90_PUT_ATT(ncid, nlatid, 'units', 'degrees_N')
256         nret = NF90_PUT_ATT(ncid, nxid, 'long_name', 'x grid box position')
257         nret = NF90_PUT_ATT(ncid, nxid, 'units', 'fractional')
258         nret = NF90_PUT_ATT(ncid, nyid, 'long_name', 'y grid box position')
259         nret = NF90_PUT_ATT(ncid, nyid, 'units', 'fractional')
260         nret = NF90_PUT_ATT(ncid, nuvelid, 'long_name', 'zonal velocity')
261         nret = NF90_PUT_ATT(ncid, nuvelid, 'units', 'm/s')
262         nret = NF90_PUT_ATT(ncid, nvvelid, 'long_name', 'meridional velocity')
263         nret = NF90_PUT_ATT(ncid, nvvelid, 'units', 'm/s')
264         nret = NF90_PUT_ATT(ncid, nmassid, 'long_name', 'mass')
265         nret = NF90_PUT_ATT(ncid, nmassid, 'units', 'kg')
266         nret = NF90_PUT_ATT(ncid, nthicknessid, 'long_name', 'thickness')
267         nret = NF90_PUT_ATT(ncid, nthicknessid, 'units', 'm')
268         nret = NF90_PUT_ATT(ncid, nwidthid, 'long_name', 'width')
269         nret = NF90_PUT_ATT(ncid, nwidthid, 'units', 'm')
270         nret = NF90_PUT_ATT(ncid, nlengthid, 'long_name', 'length')
271         nret = NF90_PUT_ATT(ncid, nlengthid, 'units', 'm')
272         nret = NF90_PUT_ATT(ncid, numberid, 'long_name', 'iceberg number on this processor')
273         nret = NF90_PUT_ATT(ncid, numberid, 'units', 'count')
274         nret = NF90_PUT_ATT(ncid, nyearid, 'long_name', 'calendar year of calving event')
275         nret = NF90_PUT_ATT(ncid, nyearid, 'units', 'years')
276         nret = NF90_PUT_ATT(ncid, ndayid, 'long_name', 'year day of calving event')
277         nret = NF90_PUT_ATT(ncid, ndayid, 'units', 'days')
278         nret = NF90_PUT_ATT(ncid, nscaling_id, 'long_name', 'scaling factor for mass of calving berg')
279         nret = NF90_PUT_ATT(ncid, nscaling_id, 'units', 'none')
280         nret = NF90_PUT_ATT(ncid, nmass_of_bits_id, 'long_name', 'mass of bergy bits')
281         nret = NF90_PUT_ATT(ncid, nmass_of_bits_id, 'units', 'kg')
282         nret = NF90_PUT_ATT(ncid, nheat_density_id, 'long_name', 'heat density')
283         nret = NF90_PUT_ATT(ncid, nheat_density_id, 'units', 'J/kg')
284
285      ENDIF ! associated(first_berg)
286
287      ! End define mode
288      nret = NF90_ENDDEF(ncid)
289
290      ! --------------------------------
291      ! now write some data
292
293      nstrt3(1) = 1
294      nstrt3(2) = 1
295      nlngth3(1) = jpi
296      nlngth3(2) = jpj
297      nlngth3(3) = 1
298
299      DO jn=1,nclasses
300         griddata(:,:,1) = berg_grid%stored_ice(:,:,jn)
301         nstrt3(3) = jn
302         nret = NF90_PUT_VAR( ncid, nsiceid, griddata, nstrt3, nlngth3 )
303         IF (nret .ne. NF90_NOERR) THEN
304            IF( lwp ) WRITE(numout,*) TRIM(NF90_STRERROR( nret ))
305            CALL ctl_stop('icebergs, write_restart: nf_put_var stored_ice failed')
306         ENDIF
307      ENDDO
308      IF( lwp ) WRITE(numout,*) 'file: ',TRIM(cl_path)//TRIM(cl_filename),' var: stored_ice  written'
309
310      nret = NF90_PUT_VAR( ncid, nkountid, num_bergs(:) )
311      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var kount failed')
312
313      nret = NF90_PUT_VAR( ncid, nsheatid, berg_grid%stored_heat(:,:) )
314      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var stored_heat failed')
315      IF( lwp ) WRITE(numout,*) 'file: ',TRIM(cl_path)//TRIM(cl_filename),' var: stored_heat written'
316
317      nret = NF90_PUT_VAR( ncid, ncalvid , src_calving(:,:) )
318      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var calving failed')
319      nret = NF90_PUT_VAR( ncid, ncalvhid, src_calving_hflx(:,:) )
320      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var calving_hflx failed')
321      IF( lwp ) WRITE(numout,*) 'file: ',TRIM(cl_path)//TRIM(cl_filename),' var: calving written'
322
323      IF ( ASSOCIATED(first_berg) ) THEN
324
325         ! Write variables
326         ! just write out the current point of the trajectory
327
328         this => first_berg
329         jn = 0
330         DO WHILE (ASSOCIATED(this))
331            pt => this%current_point
332            jn=jn+1
333
334            nret = NF90_PUT_VAR(ncid, numberid, this%number, (/1,jn/), (/nkounts,1/) )
335            nret = NF90_PUT_VAR(ncid, nscaling_id, this%mass_scaling, (/ jn /) )
336
337            nret = NF90_PUT_VAR(ncid, nlonid, pt%lon, (/ jn /) )
338            nret = NF90_PUT_VAR(ncid, nlatid, pt%lat, (/ jn /) )
339            nret = NF90_PUT_VAR(ncid, nxid, pt%xi, (/ jn /) )
340            nret = NF90_PUT_VAR(ncid, nyid, pt%yj, (/ jn /) )
341            nret = NF90_PUT_VAR(ncid, nuvelid, pt%uvel, (/ jn /) )
342            nret = NF90_PUT_VAR(ncid, nvvelid, pt%vvel, (/ jn /) )
343            nret = NF90_PUT_VAR(ncid, nmassid, pt%mass, (/ jn /) )
344            nret = NF90_PUT_VAR(ncid, nthicknessid, pt%thickness, (/ jn /) )
345            nret = NF90_PUT_VAR(ncid, nwidthid, pt%width, (/ jn /) )
346            nret = NF90_PUT_VAR(ncid, nlengthid, pt%length, (/ jn /) )
347            nret = NF90_PUT_VAR(ncid, nyearid, pt%year, (/ jn /) )
348            nret = NF90_PUT_VAR(ncid, ndayid, pt%day, (/ jn /) )
349            nret = NF90_PUT_VAR(ncid, nmass_of_bits_id, pt%mass_of_bits, (/ jn /) )
350            nret = NF90_PUT_VAR(ncid, nheat_density_id, pt%heat_density, (/ jn /) )
351
352            this=>this%next
353         END DO
354         !
355      ENDIF ! associated(first_berg)
356
357      ! Finish up
358      nret = NF90_CLOSE(ncid)
359      IF (nret /= NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_close failed')
360      !
361   END SUBROUTINE icb_rst_write
362   !
363END MODULE icbrst
Note: See TracBrowser for help on using the repository browser.