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/icebergs_restart_single_file/NEMOGCM/NEMO/OPA_SRC/ICB – NEMO

source: branches/UKMO/icebergs_restart_single_file/NEMOGCM/NEMO/OPA_SRC/ICB/icbrst.F90 @ 5911

Last change on this file since 5911 was 5911, checked in by davestorkey, 8 years ago

UKMO/icebergs_restart_single_file branch: Add rebuild script and comment.

File size: 19.1 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                      ::   iclass
63      INTEGER, DIMENSION(1)        ::   istrt, ilngth, idata
64      INTEGER, DIMENSION(2)        ::   istrt2, ilngth2
65      INTEGER, DIMENSION(nkounts)  ::   idata2
66      REAL(wp), DIMENSION(1)       ::   zdata                                         ! need 1d array to read in with
67                                                                                            ! start and count arrays
68      LOGICAL                      ::   ll_found_restart
69      CHARACTER(len=256)           ::   cl_path
70      CHARACTER(len=256)           ::   cl_filename
71      CHARACTER(len=NF90_MAX_NAME) ::   cl_dname
72      TYPE(iceberg)                ::   localberg ! NOT a pointer but an actual local variable
73      TYPE(point)                  ::   localpt   ! NOT a pointer but an actual local variable
74      REAL(wp), DIMENSION(1,nkounts) :: zdata_in
75      INTEGER                      ::   ii,ij
76      !!----------------------------------------------------------------------
77
78      ! Find a restart file. Assume iceberg restarts in same directory as ocean restarts
79      ! and are called TRIM(cn_ocerst)//'_icebergs'
80      cl_path = TRIM(cn_ocerst_indir)
81      IF( cl_path(LEN_TRIM(cl_path):) /= '/' ) cl_path = TRIM(cl_path) // '/'
82      cl_filename = TRIM(cn_ocerst_in)//'_icebergs'
83      CALL iom_open( TRIM(cl_path)//cl_filename, ncid )
84
85      IF( iom_file(ncid)%iduld .GE. 0) THEN
86
87         ibergs_in_file = iom_file(ncid)%lenuld
88         DO jn = 1,ibergs_in_file
89
90            ! iom_get treats the unlimited dimension as time. Here the unlimited dimension
91            ! is the iceberg index, but we can still use the ktime keyword to get the iceberg we want.
92
93            CALL iom_get( ncid, 'xi'     ,localpt%xi  , ktime=jn )
94            CALL iom_get( ncid, 'yj'     ,localpt%yj  , ktime=jn )
95
96            ii = INT( localpt%xi + 0.5 )
97            ij = INT( localpt%yj + 0.5 )
98            ! Only proceed if this iceberg is on the local processor (excluding halos).
99            IF ( ii .GE. nldi+nimpp-1 .AND. ii .LE. nlei+nimpp-1 .AND. &
100           &     ij .GE. nldj+njmpp-1 .AND. ij .LE. nlej+njmpp-1 ) THEN           
101
102               CALL iom_get( ncid, jpdom_unknown, 'number'       , (/zdata_in(1,:)/) , ktime=jn, kstart=(/1/), kcount=(/nkounts/) )
103               localberg%number(:) = INT(zdata_in(1,:))
104               CALL iom_get( ncid, 'mass_scaling' , localberg%mass_scaling, ktime=jn )
105               CALL iom_get( ncid, 'lon'          , localpt%lon           , ktime=jn )
106               CALL iom_get( ncid, 'lat'          , localpt%lat           , ktime=jn )
107               CALL iom_get( ncid, 'uvel'         , localpt%uvel          , ktime=jn )
108               CALL iom_get( ncid, 'vvel'         , localpt%vvel          , ktime=jn )
109               CALL iom_get( ncid, 'mass'         , localpt%mass          , ktime=jn )
110               CALL iom_get( ncid, 'thickness'    , localpt%thickness     , ktime=jn )
111               CALL iom_get( ncid, 'width'        , localpt%width         , ktime=jn )
112               CALL iom_get( ncid, 'length'       , localpt%length        , ktime=jn )
113               CALL iom_get( ncid, 'year'         , zdata_in(1,1)         , ktime=jn )
114               localpt%year = INT(zdata_in(1,1))
115               CALL iom_get( ncid, 'day'          , localpt%day           , ktime=jn )
116               CALL iom_get( ncid, 'mass_of_bits' , localpt%mass_of_bits  , ktime=jn )
117               CALL iom_get( ncid, 'heat_density' , localpt%heat_density  , ktime=jn )
118
119               !
120               CALL icb_utl_add( localberg, localpt )
121
122            ENDIF
123
124         END DO
125
126      ENDIF 
127
128      ! Gridded variables
129      CALL iom_get( ncid, jpdom_autoglo,    'calving'     , src_calving  )
130      CALL iom_get( ncid, jpdom_autoglo,    'calving_hflx', src_calving_hflx  )
131      CALL iom_get( ncid, jpdom_autoglo,    'stored_heat' , berg_grid%stored_heat  )
132      CALL iom_get( ncid, jpdom_autoglo_xy, 'stored_ice'  , berg_grid%stored_ice, kstart=(/1,1,1/), kcount=(/1,1,nclasses/) )
133     
134      CALL iom_get( ncid, jpdom_unknown, 'kount' , (/zdata_in(1,:)/) )
135      num_bergs(:) = INT(zdata_in(1,:))
136
137      ! Sanity check
138      jn = icb_utl_count()
139      IF (nn_verbose_level >= 0)   &
140         WRITE(numout,'(2(a,i5))') 'icebergs, read_restart_bergs: # bergs =',jn,' on PE',narea-1
141      IF( lk_mpp ) THEN
142         ! Only mpp_sum ibergs_in_file if we are reading from multiple restart files.
143         IF( INDEX(iom_file(ncid)%name,'icebergs.nc' ) .EQ. 0 ) CALL mpp_sum(ibergs_in_file)
144         CALL mpp_sum(jn)
145      ENDIF
146      IF(lwp)   WRITE(numout,'(a,i5,a,i5,a)') 'icebergs, read_restart_bergs: there were',ibergs_in_file,   &
147         &                                    ' bergs in the restart file and', jn,' bergs have been read'
148      !
149      ! Finish up
150      CALL iom_close( ncid )
151      !
152      IF( lwp .and. nn_verbose_level >= 0)  WRITE(numout,'(a)') 'icebergs, read_restart_bergs: completed'
153      !
154   END SUBROUTINE icb_rst_read
155
156
157   SUBROUTINE icb_rst_write( kt )
158      !!----------------------------------------------------------------------
159      !!                 ***  SUBROUTINE icb_rst_write  ***
160      !!
161      !!----------------------------------------------------------------------
162      INTEGER, INTENT( in ) :: kt
163      !
164      INTEGER ::   jn   ! dummy loop index
165      INTEGER ::   ix_dim, iy_dim, ik_dim, in_dim
166      CHARACTER(len=256)     :: cl_path
167      CHARACTER(len=256)     :: cl_filename
168      TYPE(iceberg), POINTER :: this
169      TYPE(point)  , POINTER :: pt
170      !!----------------------------------------------------------------------
171
172      ! Assume we write iceberg restarts to same directory as ocean restarts.
173      cl_path = TRIM(cn_ocerst_outdir)
174      IF( cl_path(LEN_TRIM(cl_path):) /= '/' ) cl_path = TRIM(cl_path) // '/'
175      IF( lk_mpp ) THEN
176         WRITE(cl_filename,'(A,"_icebergs_",I8.8,"_restart_",I4.4,".nc")') TRIM(cexper), kt, narea-1
177      ELSE
178         WRITE(cl_filename,'(A,"_icebergs_",I8.8,"_restart.nc")') TRIM(cexper), kt
179      ENDIF
180      IF (nn_verbose_level >= 0) WRITE(numout,'(2a)') 'icebergs, write_restart: creating ',TRIM(cl_path)//TRIM(cl_filename)
181
182      nret = NF90_CREATE(TRIM(cl_path)//TRIM(cl_filename), NF90_CLOBBER, ncid)
183      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_create failed')
184
185      ! Dimensions
186      nret = NF90_DEF_DIM(ncid, 'x', jpi, ix_dim)
187      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim x failed')
188
189      nret = NF90_DEF_DIM(ncid, 'y', jpj, iy_dim)
190      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim y failed')
191
192      nret = NF90_DEF_DIM(ncid, 'c', nclasses, nc_dim)
193      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim c failed')
194
195      nret = NF90_DEF_DIM(ncid, 'k', nkounts, ik_dim)
196      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim k failed')
197
198      ! global attributes
199      IF( lk_mpp ) THEN
200         ! Set domain parameters (assume jpdom_local_full)
201         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_number_total'   , jpnij              )
202         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_number'         , narea-1            )
203         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_dimensions_ids' , (/1     , 2     /) )
204         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_size_global'    , (/jpiglo, jpjglo/) )
205         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_size_local'     , (/jpi   , jpj   /) )
206         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_position_first' , (/nimpp , njmpp /) )
207         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_position_last'  , (/nimpp + jpi - 1 , njmpp + jpj - 1  /) )
208         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_halo_size_start', (/nldi - 1        , nldj - 1         /) )
209         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_halo_size_end'  , (/jpi - nlei      , jpj - nlej       /) )
210         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_type'           , 'BOX'              )
211      ENDIF
212     
213      IF (associated(first_berg)) then
214         nret = NF90_DEF_DIM(ncid, 'n', NF90_UNLIMITED, in_dim)
215         IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim n failed')
216      ENDIF
217
218      ! Variables
219      nret = NF90_DEF_VAR(ncid, 'kount'       , NF90_INT   , (/ ik_dim /), nkountid)
220      nret = NF90_DEF_VAR(ncid, 'calving'     , NF90_DOUBLE, (/ ix_dim, iy_dim /), ncalvid)
221      nret = NF90_DEF_VAR(ncid, 'calving_hflx', NF90_DOUBLE, (/ ix_dim, iy_dim /), ncalvhid)
222      nret = NF90_DEF_VAR(ncid, 'stored_ice'  , NF90_DOUBLE, (/ ix_dim, iy_dim, nc_dim /), nsiceid)
223      nret = NF90_DEF_VAR(ncid, 'stored_heat' , NF90_DOUBLE, (/ ix_dim, iy_dim /), nsheatid)
224
225      ! Attributes
226      nret = NF90_PUT_ATT(ncid, ncalvid , 'long_name', 'iceberg calving')
227      nret = NF90_PUT_ATT(ncid, ncalvid , 'units', 'some')
228      nret = NF90_PUT_ATT(ncid, ncalvhid, 'long_name', 'heat flux associated with iceberg calving')
229      nret = NF90_PUT_ATT(ncid, ncalvhid, 'units', 'some')
230      nret = NF90_PUT_ATT(ncid, nsiceid , 'long_name', 'stored ice used to calve icebergs')
231      nret = NF90_PUT_ATT(ncid, nsiceid , 'units', 'kg/s')
232      nret = NF90_PUT_ATT(ncid, nsheatid, 'long_name', 'heat in stored ice used to calve icebergs')
233      nret = NF90_PUT_ATT(ncid, nsheatid, 'units', 'J/kg/s')
234
235      IF ( ASSOCIATED(first_berg) ) THEN
236
237         ! Only add berg variables for this PE if we have anything to say
238
239         ! Variables
240         nret = NF90_DEF_VAR(ncid, 'lon', NF90_DOUBLE, in_dim, nlonid)
241         nret = NF90_DEF_VAR(ncid, 'lat', NF90_DOUBLE, in_dim, nlatid)
242         nret = NF90_DEF_VAR(ncid, 'xi', NF90_DOUBLE, in_dim, nxid)
243         nret = NF90_DEF_VAR(ncid, 'yj', NF90_DOUBLE, in_dim, nyid)
244         nret = NF90_DEF_VAR(ncid, 'uvel', NF90_DOUBLE, in_dim, nuvelid)
245         nret = NF90_DEF_VAR(ncid, 'vvel', NF90_DOUBLE, in_dim, nvvelid)
246         nret = NF90_DEF_VAR(ncid, 'mass', NF90_DOUBLE, in_dim, nmassid)
247         nret = NF90_DEF_VAR(ncid, 'thickness', NF90_DOUBLE, in_dim, nthicknessid)
248         nret = NF90_DEF_VAR(ncid, 'width', NF90_DOUBLE, in_dim, nwidthid)
249         nret = NF90_DEF_VAR(ncid, 'length', NF90_DOUBLE, in_dim, nlengthid)
250         nret = NF90_DEF_VAR(ncid, 'number', NF90_INT, (/ik_dim,in_dim/), numberid)
251         nret = NF90_DEF_VAR(ncid, 'year', NF90_INT, in_dim, nyearid)
252         nret = NF90_DEF_VAR(ncid, 'day', NF90_DOUBLE, in_dim, ndayid)
253         nret = NF90_DEF_VAR(ncid, 'mass_scaling', NF90_DOUBLE, in_dim, nscaling_id)
254         nret = NF90_DEF_VAR(ncid, 'mass_of_bits', NF90_DOUBLE, in_dim, nmass_of_bits_id)
255         nret = NF90_DEF_VAR(ncid, 'heat_density', NF90_DOUBLE, in_dim, nheat_density_id)
256
257         ! Attributes
258         nret = NF90_PUT_ATT(ncid, nlonid, 'long_name', 'longitude')
259         nret = NF90_PUT_ATT(ncid, nlonid, 'units', 'degrees_E')
260         nret = NF90_PUT_ATT(ncid, nlatid, 'long_name', 'latitude')
261         nret = NF90_PUT_ATT(ncid, nlatid, 'units', 'degrees_N')
262         nret = NF90_PUT_ATT(ncid, nxid, 'long_name', 'x grid box position')
263         nret = NF90_PUT_ATT(ncid, nxid, 'units', 'fractional')
264         nret = NF90_PUT_ATT(ncid, nyid, 'long_name', 'y grid box position')
265         nret = NF90_PUT_ATT(ncid, nyid, 'units', 'fractional')
266         nret = NF90_PUT_ATT(ncid, nuvelid, 'long_name', 'zonal velocity')
267         nret = NF90_PUT_ATT(ncid, nuvelid, 'units', 'm/s')
268         nret = NF90_PUT_ATT(ncid, nvvelid, 'long_name', 'meridional velocity')
269         nret = NF90_PUT_ATT(ncid, nvvelid, 'units', 'm/s')
270         nret = NF90_PUT_ATT(ncid, nmassid, 'long_name', 'mass')
271         nret = NF90_PUT_ATT(ncid, nmassid, 'units', 'kg')
272         nret = NF90_PUT_ATT(ncid, nthicknessid, 'long_name', 'thickness')
273         nret = NF90_PUT_ATT(ncid, nthicknessid, 'units', 'm')
274         nret = NF90_PUT_ATT(ncid, nwidthid, 'long_name', 'width')
275         nret = NF90_PUT_ATT(ncid, nwidthid, 'units', 'm')
276         nret = NF90_PUT_ATT(ncid, nlengthid, 'long_name', 'length')
277         nret = NF90_PUT_ATT(ncid, nlengthid, 'units', 'm')
278         nret = NF90_PUT_ATT(ncid, numberid, 'long_name', 'iceberg number on this processor')
279         nret = NF90_PUT_ATT(ncid, numberid, 'units', 'count')
280         nret = NF90_PUT_ATT(ncid, nyearid, 'long_name', 'calendar year of calving event')
281         nret = NF90_PUT_ATT(ncid, nyearid, 'units', 'years')
282         nret = NF90_PUT_ATT(ncid, ndayid, 'long_name', 'year day of calving event')
283         nret = NF90_PUT_ATT(ncid, ndayid, 'units', 'days')
284         nret = NF90_PUT_ATT(ncid, nscaling_id, 'long_name', 'scaling factor for mass of calving berg')
285         nret = NF90_PUT_ATT(ncid, nscaling_id, 'units', 'none')
286         nret = NF90_PUT_ATT(ncid, nmass_of_bits_id, 'long_name', 'mass of bergy bits')
287         nret = NF90_PUT_ATT(ncid, nmass_of_bits_id, 'units', 'kg')
288         nret = NF90_PUT_ATT(ncid, nheat_density_id, 'long_name', 'heat density')
289         nret = NF90_PUT_ATT(ncid, nheat_density_id, 'units', 'J/kg')
290
291      ENDIF ! associated(first_berg)
292
293      ! End define mode
294      nret = NF90_ENDDEF(ncid)
295
296      ! --------------------------------
297      ! now write some data
298
299      nstrt3(1) = 1
300      nstrt3(2) = 1
301      nlngth3(1) = jpi
302      nlngth3(2) = jpj
303      nlngth3(3) = 1
304
305      DO jn=1,nclasses
306         griddata(:,:,1) = berg_grid%stored_ice(:,:,jn)
307         nstrt3(3) = jn
308         nret = NF90_PUT_VAR( ncid, nsiceid, griddata, nstrt3, nlngth3 )
309         IF (nret .ne. NF90_NOERR) THEN
310            IF( lwp ) WRITE(numout,*) TRIM(NF90_STRERROR( nret ))
311            CALL ctl_stop('icebergs, write_restart: nf_put_var stored_ice failed')
312         ENDIF
313      ENDDO
314      IF( lwp ) WRITE(numout,*) 'file: ',TRIM(cl_path)//TRIM(cl_filename),' var: stored_ice  written'
315
316      nret = NF90_PUT_VAR( ncid, nkountid, num_bergs(:) )
317      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var kount failed')
318
319      nret = NF90_PUT_VAR( ncid, nsheatid, berg_grid%stored_heat(:,:) )
320      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var stored_heat failed')
321      IF( lwp ) WRITE(numout,*) 'file: ',TRIM(cl_path)//TRIM(cl_filename),' var: stored_heat written'
322
323      nret = NF90_PUT_VAR( ncid, ncalvid , src_calving(:,:) )
324      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var calving failed')
325      nret = NF90_PUT_VAR( ncid, ncalvhid, src_calving_hflx(:,:) )
326      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var calving_hflx failed')
327      IF( lwp ) WRITE(numout,*) 'file: ',TRIM(cl_path)//TRIM(cl_filename),' var: calving written'
328
329      IF ( ASSOCIATED(first_berg) ) THEN
330
331         ! Write variables
332         ! just write out the current point of the trajectory
333
334         this => first_berg
335         jn = 0
336         DO WHILE (ASSOCIATED(this))
337            pt => this%current_point
338            jn=jn+1
339
340            nret = NF90_PUT_VAR(ncid, numberid, this%number, (/1,jn/), (/nkounts,1/) )
341            nret = NF90_PUT_VAR(ncid, nscaling_id, this%mass_scaling, (/ jn /) )
342
343            nret = NF90_PUT_VAR(ncid, nlonid, pt%lon, (/ jn /) )
344            nret = NF90_PUT_VAR(ncid, nlatid, pt%lat, (/ jn /) )
345            nret = NF90_PUT_VAR(ncid, nxid, pt%xi, (/ jn /) )
346            nret = NF90_PUT_VAR(ncid, nyid, pt%yj, (/ jn /) )
347            nret = NF90_PUT_VAR(ncid, nuvelid, pt%uvel, (/ jn /) )
348            nret = NF90_PUT_VAR(ncid, nvvelid, pt%vvel, (/ jn /) )
349            nret = NF90_PUT_VAR(ncid, nmassid, pt%mass, (/ jn /) )
350            nret = NF90_PUT_VAR(ncid, nthicknessid, pt%thickness, (/ jn /) )
351            nret = NF90_PUT_VAR(ncid, nwidthid, pt%width, (/ jn /) )
352            nret = NF90_PUT_VAR(ncid, nlengthid, pt%length, (/ jn /) )
353            nret = NF90_PUT_VAR(ncid, nyearid, pt%year, (/ jn /) )
354            nret = NF90_PUT_VAR(ncid, ndayid, pt%day, (/ jn /) )
355            nret = NF90_PUT_VAR(ncid, nmass_of_bits_id, pt%mass_of_bits, (/ jn /) )
356            nret = NF90_PUT_VAR(ncid, nheat_density_id, pt%heat_density, (/ jn /) )
357
358            this=>this%next
359         END DO
360         !
361      ENDIF ! associated(first_berg)
362
363      ! Finish up
364      nret = NF90_CLOSE(ncid)
365      IF (nret /= NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_close failed')
366      !
367   END SUBROUTINE icb_rst_write
368   !
369END MODULE icbrst
Note: See TracBrowser for help on using the repository browser.