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

source: branches/UKMO/dev_r5518_GO6_package/NEMOGCM/NEMO/OPA_SRC/ICB/icbrst.F90 @ 6755

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

UKMO/dev_r5518_GO6_package branch: Merge back changes from UKMO/dev_r5518_GO6_package_MEDUSA_temporary
Custom merge into /branches/UKMO/dev_r5518_GO6_package/NEMOGCM: r6749 cf. r6618 of /branches/UKMO/dev_r5518_GO6_package_MEDUSA_temporary/NEMOGCM@6754

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