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

source: branches/2012/dev_r3337_NOCS10_ICB/NEMOGCM/NEMO/OPA_SRC/ICB/icbrst.F90 @ 3339

Last change on this file since 3339 was 3339, checked in by sga, 12 years ago

NEMO branch dev_r3337_NOCS10_ICB: add new iceberg sub-directory ICB

File size: 19.9 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   !!   icebergs_read_restart    : initialise
17   !!   icebergs_write_restart   : generate test icebergs
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 icb_oce        ! define iceberg arrays
25   USE icbutl         ! iceberg utility routines
26
27   IMPLICIT NONE
28   PRIVATE
29
30   REAL(wp), DIMENSION(:,:,:), PRIVATE, POINTER :: griddata => NULL()    ! need 2d array to read in with
31   INTEGER,  DIMENSION(3),     PRIVATE          :: strt3, lngth3
32   INTEGER,                    PRIVATE          :: lonid, latid, xid, yid, uvelid, vvelid
33   INTEGER,                    PRIVATE          :: massid, thicknessid, widthid, lengthid
34   INTEGER,                    PRIVATE          :: yearid, dayid
35   INTEGER,                    PRIVATE          :: scaling_id, mass_of_bits_id, heat_density_id, numberid
36   INTEGER,                    PRIVATE          :: siceid, sheatid, calvid, calvhid, kountid
37
38   PUBLIC   icebergs_read_restart       ! routine called in xxx.F90 module
39   PUBLIC   icebergs_write_restart      ! routine called in xxx.F90 module
40
41CONTAINS
42
43   SUBROUTINE icebergs_read_restart()
44      ! sga - read a restart file
45      !       for this version, we just read back in the restart for this processor
46      !       so we cannot change the processor layout currently with iceberg code
47
48      ! Local variables
49      INTEGER                               :: ndim, nvar, natt
50      INTEGER                               :: k, iret, ncid, unlim_dim, x_dim, y_dim, c_dim, nbergs_in_file
51      INTEGER                               :: mclass
52      INTEGER, DIMENSION(1)                 :: strt, lngth, idata
53      INTEGER, DIMENSION(2)                 :: strt2, lngth2
54      INTEGER, DIMENSION(nkounts)           :: idata2
55      REAL(wp), DIMENSION(1)                :: data                                         ! need 1d array to read in with
56                                                                                            ! start and count arrays
57      LOGICAL                               :: lres, found_restart, multiPErestart=.FALSE.
58      CHARACTER(len=80)                     :: filename
59      CHARACTER(len=NF90_MAX_NAME)          :: dname
60      TYPE(iceberg)                         :: localberg ! NOT a pointer but an actual local variable
61      TYPE(point)                           :: localpt   ! NOT a pointer but an actual local variable
62
63      IF( .NOT. ASSOCIATED(griddata) ) ALLOCATE( griddata(jpi,jpj,1) )
64
65      ! Find a restart file
66      multiPErestart=.FALSE.
67      DO
68         filename = ' '
69         filename = 'restart_icebergs.nc'
70         INQUIRE(file=TRIM(filename),exist=found_restart)
71         IF (found_restart) EXIT
72         filename = ' '
73         WRITE(filename,'("restart_icebergs_",I4.4,".nc")') narea-1
74         INQUIRE(file=TRIM(filename),exist=found_restart)
75         IF (found_restart) THEN
76            multiPErestart=.TRUE.
77            EXIT
78         ENDIF
79         IF (nn_verbose_level >= 0 .AND. lwp) &
80            WRITE(numout,'(a)') 'read_restart_bergs: no restart file found'
81         multiPErestart=.TRUE. ! force sanity checking in a multi-PE mode if no file found on this PE
82         EXIT
83      ENDDO
84
85      IF ( .NOT. found_restart) THEN ! only do the following if a file was found
86         CALL ctl_stop('icebergs: no restart file found')
87      ENDIF
88
89      IF (nn_verbose_level >= 0 .AND. lwp)  &
90         WRITE(numout,'(2a)') 'icebergs, read_restart_bergs: found restart file = ',TRIM(filename)
91
92      iret = NF90_OPEN(TRIM(filename), NF90_NOWRITE, ncid)
93      IF (iret .ne. NF90_NOERR) CALL ctl_stop('icebergs, read_restart_bergs: nf_open failed')
94
95      iret = nf90_inquire(ncid, ndim, nvar, natt, unlim_dim)
96      IF (iret .ne. NF90_NOERR) CALL ctl_stop('icebergs, read_restart_bergs: nf_inquire failed')
97
98      IF( unlim_dim .NE. -1) THEN
99
100         iret = nf90_inquire_dimension(ncid, unlim_dim, dname, nbergs_in_file)
101         IF (iret .ne. NF90_NOERR) CALL ctl_stop('icebergs, read_restart_bergs: nf_inq_dimlen failed')
102
103         iret = NF90_INQ_VARID(ncid, 'number', numberid)
104         iret = NF90_INQ_VARID(ncid, 'mass_scaling', scaling_id)
105         iret = NF90_INQ_VARID(ncid, 'xi', xid)
106         iret = NF90_INQ_VARID(ncid, 'yj', yid)
107         iret = NF90_INQ_VARID(ncid, 'lon', lonid)
108         iret = NF90_INQ_VARID(ncid, 'lat', latid)
109         iret = NF90_INQ_VARID(ncid, 'uvel', uvelid)
110         iret = NF90_INQ_VARID(ncid, 'vvel', vvelid)
111         iret = NF90_INQ_VARID(ncid, 'mass', massid)
112         iret = NF90_INQ_VARID(ncid, 'thickness', thicknessid)
113         iret = NF90_INQ_VARID(ncid, 'width', widthid)
114         iret = NF90_INQ_VARID(ncid, 'length', lengthid)
115         iret = NF90_INQ_VARID(ncid, 'year', yearid)
116         iret = NF90_INQ_VARID(ncid, 'day', dayid)
117         iret = NF90_INQ_VARID(ncid, 'mass_of_bits', mass_of_bits_id)
118         iret = NF90_INQ_VARID(ncid, 'heat_density', heat_density_id)
119
120         lngth(1) = 1
121         strt2(1) = 1
122         lngth2(1) = nkounts
123         lngth2(2) = 1
124         DO k=1, nbergs_in_file
125
126            strt(1) = k
127            strt2(2) = k
128
129            iret = NF90_GET_VAR(ncid, numberid, idata2, strt2, lngth2 )
130            localberg%number(:) = idata2(:)
131
132            iret = NF90_GET_VAR(ncid, scaling_id, data, strt, lngth )
133            localberg%mass_scaling = data(1)
134
135            iret = NF90_GET_VAR(ncid, lonid, data, strt, lngth)
136            localpt%lon = data(1)
137            iret = NF90_GET_VAR(ncid, latid, data, strt, lngth)
138            localpt%lat = data(1)
139            IF (nn_verbose_level >= 2 .AND. lwp) THEN
140               WRITE(numout,'(a,i5,a,2f10.4,a,i5)') 'icebergs, read_restart_bergs: berg ',k,' is at ', &
141                                              localpt%lon,localpt%lat,' on PE ',narea-1
142            ENDIF
143            iret = NF90_GET_VAR(ncid, xid, data, strt, lngth)
144            localpt%xi = data(1)
145            iret = NF90_GET_VAR(ncid, yid, data, strt, lngth)
146            localpt%yj = data(1)
147            iret = NF90_GET_VAR(ncid, uvelid, data, strt, lngth )
148            localpt%uvel = data(1)
149            iret = NF90_GET_VAR(ncid, vvelid, data, strt, lngth )
150            localpt%vvel = data(1)
151            iret = NF90_GET_VAR(ncid, massid, data, strt, lngth )
152            localpt%mass = data(1)
153            iret = NF90_GET_VAR(ncid, thicknessid, data, strt, lngth )
154            localpt%thickness = data(1)
155            iret = NF90_GET_VAR(ncid, widthid, data, strt, lngth )
156            localpt%width = data(1)
157            iret = NF90_GET_VAR(ncid, lengthid, data, strt, lngth )
158            localpt%length = data(1)
159            iret = NF90_GET_VAR(ncid, yearid, idata, strt, lngth )
160            localpt%year = idata(1)
161            iret = NF90_GET_VAR(ncid, dayid, data, strt, lngth )
162            localpt%day = data(1)
163            iret = NF90_GET_VAR(ncid, mass_of_bits_id, data, strt, lngth )
164            localpt%mass_of_bits = data(1)
165            iret = NF90_GET_VAR(ncid, heat_density_id, data, strt, lngth )
166            localpt%heat_density = data(1)
167            !
168            CALL add_new_berg_to_list( localberg, localpt )
169         ENDDO
170
171      ENDIF
172
173      iret = NF90_INQ_DIMID( ncid, 'c', c_dim )
174      IF (iret .ne. NF90_NOERR) CALL ctl_stop('icebergs, read_restart: nf_inq_dimid c failed')
175
176      iret = NF90_INQUIRE_DIMENSION( ncid, c_dim, dname, mclass )
177      IF (iret .ne. NF90_NOERR) CALL ctl_stop('icebergs, read_restart: nf_inquire_dimension failed')
178
179      iret = NF90_INQ_VARID(ncid, 'kount'       , kountid)
180      iret = NF90_INQ_VARID(ncid, 'calving'     , calvid)
181      iret = NF90_INQ_VARID(ncid, 'calving_hflx', calvhid)
182      iret = NF90_INQ_VARID(ncid, 'stored_ice'  , siceid)
183      iret = NF90_INQ_VARID(ncid, 'stored_heat' , sheatid)
184
185      strt3(1) = 1
186      strt3(2) = 1
187      lngth3(1) = jpi
188      lngth3(2) = jpj
189      lngth3(3) = 1
190
191      DO k=1,mclass
192
193         strt3(3) = k
194
195         iret = NF90_GET_VAR( ncid, siceid , griddata, strt3, lngth3 )
196         berg_grid%stored_ice(:,:,k) = griddata(:,:,1)
197
198      ENDDO
199
200      iret = NF90_GET_VAR( ncid, calvid , p_calving(:,:) )
201      iret = NF90_GET_VAR( ncid, calvhid, p_calving_hflx(:,:) )
202      iret = NF90_GET_VAR( ncid, sheatid, berg_grid%stored_heat(:,:) )
203      iret = NF90_GET_VAR( ncid, kountid, idata2(:) )
204      kount_bergs(:) = idata2(:)
205
206      ! Finish up
207      iret = NF90_CLOSE(ncid)
208      IF (iret .ne. NF90_NOERR) CALL ctl_stop('icebergs, read_restart: nf_close failed')
209
210      ! Sanity check
211      k = count_bergs()
212      IF (nn_verbose_level >= 0)   &
213         WRITE(numout,'(2(a,i5))') 'icebergs, read_restart_bergs: # bergs =',k,' on PE',narea-1
214      IF( lk_mpp ) THEN
215         IF (multiPErestart) CALL mpp_sum(nbergs_in_file) ! In case PE 0 didn't open a file
216         CALL mpp_sum(k)
217      ENDIF
218      IF ( lwp ) THEN
219         WRITE(numout,'(a,i5,a,i5,a)') 'icebergs, read_restart_bergs: there were',nbergs_in_file, &
220                                       ' bergs in the restart file and', k,' bergs have been read'
221      ENDIF
222
223      IF ( lwp .and. nn_verbose_level >= 0)  &
224         WRITE(numout,'(a)') 'icebergs, read_restart_bergs: completed'
225
226   END SUBROUTINE icebergs_read_restart
227
228! ##############################################################################
229
230   SUBROUTINE icebergs_write_restart( kt )
231
232      INTEGER, INTENT( in )                 :: kt
233      ! Local variables
234      INTEGER                               :: iret, ncid, k_dim, n_dim, x_dim, y_dim, c_dim, k
235      CHARACTER(len=80)                     :: filename
236      TYPE(iceberg), POINTER                :: this
237      TYPE(point)  , POINTER                :: pt
238
239      IF( .NOT. ASSOCIATED(griddata) ) ALLOCATE( griddata(jpi,jpj,1) )
240
241      WRITE(filename,'("icebergs_",I8.8,"_restart_",I4.4,".nc")') kt, narea-1
242      IF (nn_verbose_level >= 0) WRITE(numout,'(2a)') 'icebergs, write_restart: creating ',TRIM(filename)
243
244      iret = NF90_CREATE(TRIM(filename), NF90_CLOBBER, ncid)
245      IF (iret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_create failed')
246
247      ! Dimensions
248      iret = NF90_DEF_DIM(ncid, 'x', jpi, x_dim)
249      IF (iret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim x failed')
250
251      iret = NF90_DEF_DIM(ncid, 'y', jpj, y_dim)
252      IF (iret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim y failed')
253
254      iret = NF90_DEF_DIM(ncid, 'c', nclasses, c_dim)
255      IF (iret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim c failed')
256
257      iret = NF90_DEF_DIM(ncid, 'k', nkounts, k_dim)
258      IF (iret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim k failed')
259
260      IF (associated(first_berg)) then
261         iret = NF90_DEF_DIM(ncid, 'n', NF90_UNLIMITED, n_dim)
262         IF (iret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim n failed')
263      ENDIF
264
265      ! Variables
266      iret = NF90_DEF_VAR(ncid, 'kount'       , NF90_INT   , (/ k_dim /), kountid)
267      iret = NF90_DEF_VAR(ncid, 'calving'     , NF90_DOUBLE, (/ x_dim, y_dim /), calvid)
268      iret = NF90_DEF_VAR(ncid, 'calving_hflx', NF90_DOUBLE, (/ x_dim, y_dim /), calvhid)
269      iret = NF90_DEF_VAR(ncid, 'stored_ice'  , NF90_DOUBLE, (/ x_dim, y_dim, c_dim /), siceid)
270      iret = NF90_DEF_VAR(ncid, 'stored_heat' , NF90_DOUBLE, (/ x_dim, y_dim /), sheatid)
271
272      ! Attributes
273      iret = NF90_PUT_ATT(ncid, calvid , 'long_name', 'iceberg calving')
274      iret = NF90_PUT_ATT(ncid, calvid , 'units', 'some')
275      iret = NF90_PUT_ATT(ncid, calvhid, 'long_name', 'heat flux associated with iceberg calving')
276      iret = NF90_PUT_ATT(ncid, calvhid, 'units', 'some')
277      iret = NF90_PUT_ATT(ncid, siceid , 'long_name', 'stored ice used to calve icebergs')
278      iret = NF90_PUT_ATT(ncid, siceid , 'units', 'kg/s')
279      iret = NF90_PUT_ATT(ncid, sheatid, 'long_name', 'heat in stored ice used to calve icebergs')
280      iret = NF90_PUT_ATT(ncid, sheatid, 'units', 'J/kg/s')
281
282      IF ( ASSOCIATED(first_berg) ) THEN
283
284         ! Only add berg variables for this PE if we have anything to say
285
286         ! Variables
287         iret = NF90_DEF_VAR(ncid, 'lon', NF90_DOUBLE, n_dim, lonid)
288         iret = NF90_DEF_VAR(ncid, 'lat', NF90_DOUBLE, n_dim, latid)
289         iret = NF90_DEF_VAR(ncid, 'xi', NF90_DOUBLE, n_dim, xid)
290         iret = NF90_DEF_VAR(ncid, 'yj', NF90_DOUBLE, n_dim, yid)
291         iret = NF90_DEF_VAR(ncid, 'uvel', NF90_DOUBLE, n_dim, uvelid)
292         iret = NF90_DEF_VAR(ncid, 'vvel', NF90_DOUBLE, n_dim, vvelid)
293         iret = NF90_DEF_VAR(ncid, 'mass', NF90_DOUBLE, n_dim, massid)
294         iret = NF90_DEF_VAR(ncid, 'thickness', NF90_DOUBLE, n_dim, thicknessid)
295         iret = NF90_DEF_VAR(ncid, 'width', NF90_DOUBLE, n_dim, widthid)
296         iret = NF90_DEF_VAR(ncid, 'length', NF90_DOUBLE, n_dim, lengthid)
297         iret = NF90_DEF_VAR(ncid, 'number', NF90_INT, (/k_dim,n_dim/), numberid)
298         iret = NF90_DEF_VAR(ncid, 'year', NF90_INT, n_dim, yearid)
299         iret = NF90_DEF_VAR(ncid, 'day', NF90_DOUBLE, n_dim, dayid)
300         iret = NF90_DEF_VAR(ncid, 'mass_scaling', NF90_DOUBLE, n_dim, scaling_id)
301         iret = NF90_DEF_VAR(ncid, 'mass_of_bits', NF90_DOUBLE, n_dim, mass_of_bits_id)
302         iret = NF90_DEF_VAR(ncid, 'heat_density', NF90_DOUBLE, n_dim, heat_density_id)
303
304         ! Attributes
305         iret = NF90_PUT_ATT(ncid, lonid, 'long_name', 'longitude')
306         iret = NF90_PUT_ATT(ncid, lonid, 'units', 'degrees_E')
307         iret = NF90_PUT_ATT(ncid, latid, 'long_name', 'latitude')
308         iret = NF90_PUT_ATT(ncid, latid, 'units', 'degrees_N')
309         iret = NF90_PUT_ATT(ncid, xid, 'long_name', 'x grid box position')
310         iret = NF90_PUT_ATT(ncid, xid, 'units', 'fractional')
311         iret = NF90_PUT_ATT(ncid, yid, 'long_name', 'y grid box position')
312         iret = NF90_PUT_ATT(ncid, yid, 'units', 'fractional')
313         iret = NF90_PUT_ATT(ncid, uvelid, 'long_name', 'zonal velocity')
314         iret = NF90_PUT_ATT(ncid, uvelid, 'units', 'm/s')
315         iret = NF90_PUT_ATT(ncid, vvelid, 'long_name', 'meridional velocity')
316         iret = NF90_PUT_ATT(ncid, vvelid, 'units', 'm/s')
317         iret = NF90_PUT_ATT(ncid, massid, 'long_name', 'mass')
318         iret = NF90_PUT_ATT(ncid, massid, 'units', 'kg')
319         iret = NF90_PUT_ATT(ncid, thicknessid, 'long_name', 'thickness')
320         iret = NF90_PUT_ATT(ncid, thicknessid, 'units', 'm')
321         iret = NF90_PUT_ATT(ncid, widthid, 'long_name', 'width')
322         iret = NF90_PUT_ATT(ncid, widthid, 'units', 'm')
323         iret = NF90_PUT_ATT(ncid, lengthid, 'long_name', 'length')
324         iret = NF90_PUT_ATT(ncid, lengthid, 'units', 'm')
325         iret = NF90_PUT_ATT(ncid, numberid, 'long_name', 'iceberg number on this processor')
326         iret = NF90_PUT_ATT(ncid, numberid, 'units', 'count')
327         iret = NF90_PUT_ATT(ncid, yearid, 'long_name', 'calendar year of calving event')
328         iret = NF90_PUT_ATT(ncid, yearid, 'units', 'years')
329         iret = NF90_PUT_ATT(ncid, dayid, 'long_name', 'year day of calving event')
330         iret = NF90_PUT_ATT(ncid, dayid, 'units', 'days')
331         iret = NF90_PUT_ATT(ncid, scaling_id, 'long_name', 'scaling factor for mass of calving berg')
332         iret = NF90_PUT_ATT(ncid, scaling_id, 'units', 'none')
333         iret = NF90_PUT_ATT(ncid, mass_of_bits_id, 'long_name', 'mass of bergy bits')
334         iret = NF90_PUT_ATT(ncid, mass_of_bits_id, 'units', 'kg')
335         iret = NF90_PUT_ATT(ncid, heat_density_id, 'long_name', 'heat density')
336         iret = NF90_PUT_ATT(ncid, heat_density_id, 'units', 'J/kg')
337
338      ENDIF ! associated(first_berg)
339
340      ! End define mode
341      iret = NF90_ENDDEF(ncid)
342
343      ! --------------------------------
344      ! now write some data
345
346      strt3(1) = 1
347      strt3(2) = 1
348      lngth3(1) = jpi
349      lngth3(2) = jpj
350      lngth3(3) = 1
351
352      DO k=1,nclasses
353         griddata(:,:,1) = berg_grid%stored_ice(:,:,k)
354         strt3(3) = k
355         iret = NF90_PUT_VAR( ncid, siceid, griddata, strt3, lngth3 )
356         IF (iret .ne. NF90_NOERR) THEN
357            IF( lwp ) WRITE(numout,*) TRIM(NF90_STRERROR( iret ))
358            CALL ctl_stop('icebergs, write_restart: nf_put_var stored_ice failed')
359         ENDIF
360      ENDDO
361      IF( lwp ) WRITE(numout,*) 'file: ',TRIM(filename),' var: stored_ice  written'
362
363      iret = NF90_PUT_VAR( ncid, kountid, kount_bergs(:) )
364      IF (iret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var kount failed')
365
366      iret = NF90_PUT_VAR( ncid, sheatid, berg_grid%stored_heat(:,:) )
367      IF (iret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var stored_heat failed')
368      IF( lwp ) WRITE(numout,*) 'file: ',TRIM(filename),' var: stored_heat written'
369
370      iret = NF90_PUT_VAR( ncid, calvid , p_calving(:,:) )
371      IF (iret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var calving failed')
372      iret = NF90_PUT_VAR( ncid, calvhid, p_calving_hflx(:,:) )
373      IF (iret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var calving_hflx failed')
374      IF( lwp ) WRITE(numout,*) 'file: ',TRIM(filename),' var: calving written'
375
376      IF ( ASSOCIATED(first_berg) ) THEN
377
378         ! Write variables
379         ! sga - just write out the current point of the trajectory
380
381         this => first_berg
382         k = 0
383         DO WHILE (ASSOCIATED(this))
384            pt => this%current_point
385            k=k+1
386
387            iret = NF90_PUT_VAR(ncid, numberid, this%number, (/1,k/), (/nkounts,1/) )
388            iret = NF90_PUT_VAR(ncid, scaling_id, this%mass_scaling, (/ k /) )
389
390            iret = NF90_PUT_VAR(ncid, lonid, pt%lon, (/ k /) )
391            iret = NF90_PUT_VAR(ncid, latid, pt%lat, (/ k /) )
392            iret = NF90_PUT_VAR(ncid, xid, pt%xi, (/ k /) )
393            iret = NF90_PUT_VAR(ncid, yid, pt%yj, (/ k /) )
394            iret = NF90_PUT_VAR(ncid, uvelid, pt%uvel, (/ k /) )
395            iret = NF90_PUT_VAR(ncid, vvelid, pt%vvel, (/ k /) )
396            iret = NF90_PUT_VAR(ncid, massid, pt%mass, (/ k /) )
397            iret = NF90_PUT_VAR(ncid, thicknessid, pt%thickness, (/ k /) )
398            iret = NF90_PUT_VAR(ncid, widthid, pt%width, (/ k /) )
399            iret = NF90_PUT_VAR(ncid, lengthid, pt%length, (/ k /) )
400            iret = NF90_PUT_VAR(ncid, yearid, pt%year, (/ k /) )
401            iret = NF90_PUT_VAR(ncid, dayid, pt%day, (/ k /) )
402            iret = NF90_PUT_VAR(ncid, mass_of_bits_id, pt%mass_of_bits, (/ k /) )
403            iret = NF90_PUT_VAR(ncid, heat_density_id, pt%heat_density, (/ k /) )
404
405            this=>this%next
406         ENDDO
407
408      ENDIF ! associated(first_berg)
409
410      ! Finish up
411      iret = NF90_CLOSE(ncid)
412      IF (iret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_close failed')
413
414   END SUBROUTINE icebergs_write_restart
415
416END MODULE icbrst
Note: See TracBrowser for help on using the repository browser.