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 @ 5865

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

IOM-inification of reading of iceberg restarts.

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