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

source: branches/UKMO/dev_r5107_kara_mld/NEMOGCM/NEMO/OPA_SRC/ICB/icbrst.F90 @ 5244

Last change on this file since 5244 was 5244, checked in by davestorkey, 9 years ago

UKMO Kara MLD branch: remove svn keyword property and clear keywords.

File size: 20.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   !!   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 icb_oce        ! define iceberg arrays
25   USE icbutl         ! iceberg utility routines
26
27   IMPLICIT NONE
28   PRIVATE
29
30   PUBLIC   icb_rst_read    ! routine called in icbini.F90 module
31   PUBLIC   icb_rst_write   ! routine called in icbstp.F90 module
32   
33   INTEGER ::   nlonid, nlatid, nxid, nyid, nuvelid, nvvelid
34   INTEGER ::   nmassid, nthicknessid, nwidthid, nlengthid
35   INTEGER ::   nyearid, ndayid
36   INTEGER ::   nscaling_id, nmass_of_bits_id, nheat_density_id, numberid
37   INTEGER ::   nsiceid, nsheatid, ncalvid, ncalvhid, nkountid
38   INTEGER ::   nret, ncid, nc_dim
39   
40   INTEGER,  DIMENSION(3)                  :: nstrt3, nlngth3
41
42   !!----------------------------------------------------------------------
43   !! NEMO/OPA 3.3 , NEMO Consortium (2011)
44   !! $Id$
45   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
46   !!----------------------------------------------------------------------
47CONTAINS
48
49   SUBROUTINE icb_rst_read()
50      !!----------------------------------------------------------------------
51      !!                 ***  SUBROUTINE icb_rst_read  ***
52      !!
53      !! ** Purpose :   read a iceberg restart file
54      !!      NB: for this version, we just read back in the restart for this processor
55      !!      so we cannot change the processor layout currently with iceberg code
56      !!----------------------------------------------------------------------
57      INTEGER                      ::   idim, ivar, iatt
58      INTEGER                      ::   jn, iunlim_dim, ibergs_in_file
59      INTEGER                      ::   iclass
60      INTEGER, DIMENSION(1)        ::   istrt, ilngth, idata
61      INTEGER, DIMENSION(2)        ::   istrt2, ilngth2
62      INTEGER, DIMENSION(nkounts)  ::   idata2
63      REAL(wp), DIMENSION(1)       ::   zdata                                         ! need 1d array to read in with
64                                                                                            ! start and count arrays
65      LOGICAL                      ::   ll_found_restart
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
73      cl_filename = ' '
74      IF ( lk_mpp ) THEN
75         cl_filename = ' '
76         WRITE( cl_filename, '("restart_icebergs_",I4.4,".nc")' ) narea-1
77         INQUIRE( file=TRIM(cl_filename), exist=ll_found_restart )
78      ELSE
79         cl_filename = 'restart_icebergs.nc'
80         INQUIRE( file=TRIM(cl_filename), exist=ll_found_restart )
81      ENDIF
82
83      IF ( .NOT. ll_found_restart) THEN                     ! only do the following if a file was found
84         CALL ctl_stop('icebergs: no restart file found')
85      ENDIF
86
87      IF (nn_verbose_level >= 0 .AND. lwp)  &
88         WRITE(numout,'(2a)') 'icebergs, read_restart_bergs: found restart file = ',TRIM(cl_filename)
89
90      nret = NF90_OPEN(TRIM(cl_filename), NF90_NOWRITE, ncid)
91      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, read_restart_bergs: nf_open failed')
92
93      nret = nf90_inquire(ncid, idim, ivar, iatt, iunlim_dim)
94      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, read_restart_bergs: nf_inquire failed')
95
96      IF( iunlim_dim .NE. -1) THEN
97
98         nret = nf90_inquire_dimension(ncid, iunlim_dim, cl_dname, ibergs_in_file)
99         IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, read_restart_bergs: nf_inq_dimlen failed')
100
101         nret = NF90_INQ_VARID(ncid, 'number', numberid)
102         nret = NF90_INQ_VARID(ncid, 'mass_scaling', nscaling_id)
103         nret = NF90_INQ_VARID(ncid, 'xi', nxid)
104         nret = NF90_INQ_VARID(ncid, 'yj', nyid)
105         nret = NF90_INQ_VARID(ncid, 'lon', nlonid)
106         nret = NF90_INQ_VARID(ncid, 'lat', nlatid)
107         nret = NF90_INQ_VARID(ncid, 'uvel', nuvelid)
108         nret = NF90_INQ_VARID(ncid, 'vvel', nvvelid)
109         nret = NF90_INQ_VARID(ncid, 'mass', nmassid)
110         nret = NF90_INQ_VARID(ncid, 'thickness', nthicknessid)
111         nret = NF90_INQ_VARID(ncid, 'width', nwidthid)
112         nret = NF90_INQ_VARID(ncid, 'length', nlengthid)
113         nret = NF90_INQ_VARID(ncid, 'year', nyearid)
114         nret = NF90_INQ_VARID(ncid, 'day', ndayid)
115         nret = NF90_INQ_VARID(ncid, 'mass_of_bits', nmass_of_bits_id)
116         nret = NF90_INQ_VARID(ncid, 'heat_density', nheat_density_id)
117
118         ilngth(1) = 1
119         istrt2(1) = 1
120         ilngth2(1) = nkounts
121         ilngth2(2) = 1
122         DO jn=1, ibergs_in_file
123
124            istrt(1) = jn
125            istrt2(2) = jn
126
127            nret = NF90_GET_VAR(ncid, numberid, idata2, istrt2, ilngth2 )
128            localberg%number(:) = idata2(:)
129
130            nret = NF90_GET_VAR(ncid, nscaling_id, zdata, istrt, ilngth )
131            localberg%mass_scaling = zdata(1)
132
133            nret = NF90_GET_VAR(ncid, nlonid, zdata, istrt, ilngth)
134            localpt%lon = zdata(1)
135            nret = NF90_GET_VAR(ncid, nlatid, zdata, istrt, ilngth)
136            localpt%lat = zdata(1)
137            IF (nn_verbose_level >= 2 .AND. lwp) THEN
138               WRITE(numout,'(a,i5,a,2f10.4,a,i5)') 'icebergs, read_restart_bergs: berg ',jn,' is at ', &
139                                              localpt%lon,localpt%lat,' on PE ',narea-1
140            ENDIF
141            nret = NF90_GET_VAR(ncid, nxid, zdata, istrt, ilngth)
142            localpt%xi = zdata(1)
143            nret = NF90_GET_VAR(ncid, nyid, zdata, istrt, ilngth)
144            localpt%yj = zdata(1)
145            nret = NF90_GET_VAR(ncid, nuvelid, zdata, istrt, ilngth )
146            localpt%uvel = zdata(1)
147            nret = NF90_GET_VAR(ncid, nvvelid, zdata, istrt, ilngth )
148            localpt%vvel = zdata(1)
149            nret = NF90_GET_VAR(ncid, nmassid, zdata, istrt, ilngth )
150            localpt%mass = zdata(1)
151            nret = NF90_GET_VAR(ncid, nthicknessid, zdata, istrt, ilngth )
152            localpt%thickness = zdata(1)
153            nret = NF90_GET_VAR(ncid, nwidthid, zdata, istrt, ilngth )
154            localpt%width = zdata(1)
155            nret = NF90_GET_VAR(ncid, nlengthid, zdata, istrt, ilngth )
156            localpt%length = zdata(1)
157            nret = NF90_GET_VAR(ncid, nyearid, idata, istrt, ilngth )
158            localpt%year = idata(1)
159            nret = NF90_GET_VAR(ncid, ndayid, zdata, istrt, ilngth )
160            localpt%day = zdata(1)
161            nret = NF90_GET_VAR(ncid, nmass_of_bits_id, zdata, istrt, ilngth )
162            localpt%mass_of_bits = zdata(1)
163            nret = NF90_GET_VAR(ncid, nheat_density_id, zdata, istrt, ilngth )
164            localpt%heat_density = zdata(1)
165            !
166            CALL icb_utl_add( localberg, localpt )
167         END DO
168         !
169      ENDIF
170
171      nret = NF90_INQ_DIMID( ncid, 'c', nc_dim )
172      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, read_restart: nf_inq_dimid c failed')
173
174      nret = NF90_INQUIRE_DIMENSION( ncid, nc_dim, cl_dname, iclass )
175      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, read_restart: nf_inquire_dimension failed')
176
177      nret = NF90_INQ_VARID(ncid, 'kount'       , nkountid)
178      nret = NF90_INQ_VARID(ncid, 'calving'     , ncalvid)
179      nret = NF90_INQ_VARID(ncid, 'calving_hflx', ncalvhid)
180      nret = NF90_INQ_VARID(ncid, 'stored_ice'  , nsiceid)
181      nret = NF90_INQ_VARID(ncid, 'stored_heat' , nsheatid)
182
183      nstrt3(1) = 1
184      nstrt3(2) = 1
185      nlngth3(1) = jpi
186      nlngth3(2) = jpj
187      nlngth3(3) = 1
188
189      DO jn = 1, iclass
190         nstrt3(3) = jn
191         nret      = NF90_GET_VAR( ncid, nsiceid , griddata, nstrt3, nlngth3 )
192         berg_grid%stored_ice(:,:,jn) = griddata(:,:,1)
193      END DO
194
195      nret = NF90_GET_VAR( ncid, ncalvid , src_calving          (:,:) )
196      nret = NF90_GET_VAR( ncid, ncalvhid, src_calving_hflx     (:,:) )
197      nret = NF90_GET_VAR( ncid, nsheatid, berg_grid%stored_heat(:,:) )
198      nret = NF90_GET_VAR( ncid, nkountid, idata2(:) )
199      num_bergs(:) = idata2(:)
200
201      ! Finish up
202      nret = NF90_CLOSE(ncid)
203      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, read_restart: nf_close failed')
204
205      ! Sanity check
206      jn = icb_utl_count()
207      IF (nn_verbose_level >= 0)   &
208         WRITE(numout,'(2(a,i5))') 'icebergs, read_restart_bergs: # bergs =',jn,' on PE',narea-1
209      IF( lk_mpp ) THEN
210         CALL mpp_sum(ibergs_in_file)
211         CALL mpp_sum(jn)
212      ENDIF
213      IF(lwp)   WRITE(numout,'(a,i5,a,i5,a)') 'icebergs, read_restart_bergs: there were',ibergs_in_file,   &
214         &                                    ' bergs in the restart file and', jn,' bergs have been read'
215      !
216      IF( lwp .and. nn_verbose_level >= 0)  WRITE(numout,'(a)') 'icebergs, read_restart_bergs: completed'
217      !
218   END SUBROUTINE icb_rst_read
219
220
221   SUBROUTINE icb_rst_write( kt )
222      !!----------------------------------------------------------------------
223      !!                 ***  SUBROUTINE icb_rst_write  ***
224      !!
225      !!----------------------------------------------------------------------
226      INTEGER, INTENT( in ) :: kt
227      !
228      INTEGER ::   jn   ! dummy loop index
229      INTEGER ::   ix_dim, iy_dim, ik_dim, in_dim
230      CHARACTER(len=256)     :: cl_filename
231      TYPE(iceberg), POINTER :: this
232      TYPE(point)  , POINTER :: pt
233      !!----------------------------------------------------------------------
234
235      IF( lk_mpp ) THEN
236         WRITE(cl_filename,'("icebergs_",I8.8,"_restart_",I4.4,".nc")') kt, narea-1
237      ELSE
238         WRITE(cl_filename,'("icebergs_",I8.8,"_restart.nc")') kt
239      ENDIF
240      IF (nn_verbose_level >= 0) WRITE(numout,'(2a)') 'icebergs, write_restart: creating ',TRIM(cl_filename)
241
242      nret = NF90_CREATE(TRIM(cl_filename), NF90_CLOBBER, ncid)
243      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_create failed')
244
245      ! Dimensions
246      nret = NF90_DEF_DIM(ncid, 'x', jpi, ix_dim)
247      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim x failed')
248
249      nret = NF90_DEF_DIM(ncid, 'y', jpj, iy_dim)
250      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim y failed')
251
252      nret = NF90_DEF_DIM(ncid, 'c', nclasses, nc_dim)
253      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim c failed')
254
255      nret = NF90_DEF_DIM(ncid, 'k', nkounts, ik_dim)
256      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim k failed')
257
258      ! global attributes
259      IF( lk_mpp ) THEN
260         ! Set domain parameters (assume jpdom_local_full)
261         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_number_total'   , jpnij              )
262         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_number'         , narea-1            )
263         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_dimensions_ids' , (/1     , 2     /) )
264         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_size_global'    , (/jpiglo, jpjglo/) )
265         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_size_local'     , (/jpi   , jpj   /) )
266         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_position_first' , (/nimpp , njmpp /) )
267         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_position_last'  , (/nimpp + jpi - 1 , njmpp + jpj - 1  /) )
268         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_halo_size_start', (/nldi - 1        , nldj - 1         /) )
269         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_halo_size_end'  , (/jpi - nlei      , jpj - nlej       /) )
270         nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_type'           , 'BOX'              )
271      ENDIF
272     
273      IF (associated(first_berg)) then
274         nret = NF90_DEF_DIM(ncid, 'n', NF90_UNLIMITED, in_dim)
275         IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim n failed')
276      ENDIF
277
278      ! Variables
279      nret = NF90_DEF_VAR(ncid, 'kount'       , NF90_INT   , (/ ik_dim /), nkountid)
280      nret = NF90_DEF_VAR(ncid, 'calving'     , NF90_DOUBLE, (/ ix_dim, iy_dim /), ncalvid)
281      nret = NF90_DEF_VAR(ncid, 'calving_hflx', NF90_DOUBLE, (/ ix_dim, iy_dim /), ncalvhid)
282      nret = NF90_DEF_VAR(ncid, 'stored_ice'  , NF90_DOUBLE, (/ ix_dim, iy_dim, nc_dim /), nsiceid)
283      nret = NF90_DEF_VAR(ncid, 'stored_heat' , NF90_DOUBLE, (/ ix_dim, iy_dim /), nsheatid)
284
285      ! Attributes
286      nret = NF90_PUT_ATT(ncid, ncalvid , 'long_name', 'iceberg calving')
287      nret = NF90_PUT_ATT(ncid, ncalvid , 'units', 'some')
288      nret = NF90_PUT_ATT(ncid, ncalvhid, 'long_name', 'heat flux associated with iceberg calving')
289      nret = NF90_PUT_ATT(ncid, ncalvhid, 'units', 'some')
290      nret = NF90_PUT_ATT(ncid, nsiceid , 'long_name', 'stored ice used to calve icebergs')
291      nret = NF90_PUT_ATT(ncid, nsiceid , 'units', 'kg/s')
292      nret = NF90_PUT_ATT(ncid, nsheatid, 'long_name', 'heat in stored ice used to calve icebergs')
293      nret = NF90_PUT_ATT(ncid, nsheatid, 'units', 'J/kg/s')
294
295      IF ( ASSOCIATED(first_berg) ) THEN
296
297         ! Only add berg variables for this PE if we have anything to say
298
299         ! Variables
300         nret = NF90_DEF_VAR(ncid, 'lon', NF90_DOUBLE, in_dim, nlonid)
301         nret = NF90_DEF_VAR(ncid, 'lat', NF90_DOUBLE, in_dim, nlatid)
302         nret = NF90_DEF_VAR(ncid, 'xi', NF90_DOUBLE, in_dim, nxid)
303         nret = NF90_DEF_VAR(ncid, 'yj', NF90_DOUBLE, in_dim, nyid)
304         nret = NF90_DEF_VAR(ncid, 'uvel', NF90_DOUBLE, in_dim, nuvelid)
305         nret = NF90_DEF_VAR(ncid, 'vvel', NF90_DOUBLE, in_dim, nvvelid)
306         nret = NF90_DEF_VAR(ncid, 'mass', NF90_DOUBLE, in_dim, nmassid)
307         nret = NF90_DEF_VAR(ncid, 'thickness', NF90_DOUBLE, in_dim, nthicknessid)
308         nret = NF90_DEF_VAR(ncid, 'width', NF90_DOUBLE, in_dim, nwidthid)
309         nret = NF90_DEF_VAR(ncid, 'length', NF90_DOUBLE, in_dim, nlengthid)
310         nret = NF90_DEF_VAR(ncid, 'number', NF90_INT, (/ik_dim,in_dim/), numberid)
311         nret = NF90_DEF_VAR(ncid, 'year', NF90_INT, in_dim, nyearid)
312         nret = NF90_DEF_VAR(ncid, 'day', NF90_DOUBLE, in_dim, ndayid)
313         nret = NF90_DEF_VAR(ncid, 'mass_scaling', NF90_DOUBLE, in_dim, nscaling_id)
314         nret = NF90_DEF_VAR(ncid, 'mass_of_bits', NF90_DOUBLE, in_dim, nmass_of_bits_id)
315         nret = NF90_DEF_VAR(ncid, 'heat_density', NF90_DOUBLE, in_dim, nheat_density_id)
316
317         ! Attributes
318         nret = NF90_PUT_ATT(ncid, nlonid, 'long_name', 'longitude')
319         nret = NF90_PUT_ATT(ncid, nlonid, 'units', 'degrees_E')
320         nret = NF90_PUT_ATT(ncid, nlatid, 'long_name', 'latitude')
321         nret = NF90_PUT_ATT(ncid, nlatid, 'units', 'degrees_N')
322         nret = NF90_PUT_ATT(ncid, nxid, 'long_name', 'x grid box position')
323         nret = NF90_PUT_ATT(ncid, nxid, 'units', 'fractional')
324         nret = NF90_PUT_ATT(ncid, nyid, 'long_name', 'y grid box position')
325         nret = NF90_PUT_ATT(ncid, nyid, 'units', 'fractional')
326         nret = NF90_PUT_ATT(ncid, nuvelid, 'long_name', 'zonal velocity')
327         nret = NF90_PUT_ATT(ncid, nuvelid, 'units', 'm/s')
328         nret = NF90_PUT_ATT(ncid, nvvelid, 'long_name', 'meridional velocity')
329         nret = NF90_PUT_ATT(ncid, nvvelid, 'units', 'm/s')
330         nret = NF90_PUT_ATT(ncid, nmassid, 'long_name', 'mass')
331         nret = NF90_PUT_ATT(ncid, nmassid, 'units', 'kg')
332         nret = NF90_PUT_ATT(ncid, nthicknessid, 'long_name', 'thickness')
333         nret = NF90_PUT_ATT(ncid, nthicknessid, 'units', 'm')
334         nret = NF90_PUT_ATT(ncid, nwidthid, 'long_name', 'width')
335         nret = NF90_PUT_ATT(ncid, nwidthid, 'units', 'm')
336         nret = NF90_PUT_ATT(ncid, nlengthid, 'long_name', 'length')
337         nret = NF90_PUT_ATT(ncid, nlengthid, 'units', 'm')
338         nret = NF90_PUT_ATT(ncid, numberid, 'long_name', 'iceberg number on this processor')
339         nret = NF90_PUT_ATT(ncid, numberid, 'units', 'count')
340         nret = NF90_PUT_ATT(ncid, nyearid, 'long_name', 'calendar year of calving event')
341         nret = NF90_PUT_ATT(ncid, nyearid, 'units', 'years')
342         nret = NF90_PUT_ATT(ncid, ndayid, 'long_name', 'year day of calving event')
343         nret = NF90_PUT_ATT(ncid, ndayid, 'units', 'days')
344         nret = NF90_PUT_ATT(ncid, nscaling_id, 'long_name', 'scaling factor for mass of calving berg')
345         nret = NF90_PUT_ATT(ncid, nscaling_id, 'units', 'none')
346         nret = NF90_PUT_ATT(ncid, nmass_of_bits_id, 'long_name', 'mass of bergy bits')
347         nret = NF90_PUT_ATT(ncid, nmass_of_bits_id, 'units', 'kg')
348         nret = NF90_PUT_ATT(ncid, nheat_density_id, 'long_name', 'heat density')
349         nret = NF90_PUT_ATT(ncid, nheat_density_id, 'units', 'J/kg')
350
351      ENDIF ! associated(first_berg)
352
353      ! End define mode
354      nret = NF90_ENDDEF(ncid)
355
356      ! --------------------------------
357      ! now write some data
358
359      nstrt3(1) = 1
360      nstrt3(2) = 1
361      nlngth3(1) = jpi
362      nlngth3(2) = jpj
363      nlngth3(3) = 1
364
365      DO jn=1,nclasses
366         griddata(:,:,1) = berg_grid%stored_ice(:,:,jn)
367         nstrt3(3) = jn
368         nret = NF90_PUT_VAR( ncid, nsiceid, griddata, nstrt3, nlngth3 )
369         IF (nret .ne. NF90_NOERR) THEN
370            IF( lwp ) WRITE(numout,*) TRIM(NF90_STRERROR( nret ))
371            CALL ctl_stop('icebergs, write_restart: nf_put_var stored_ice failed')
372         ENDIF
373      ENDDO
374      IF( lwp ) WRITE(numout,*) 'file: ',TRIM(cl_filename),' var: stored_ice  written'
375
376      nret = NF90_PUT_VAR( ncid, nkountid, num_bergs(:) )
377      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var kount failed')
378
379      nret = NF90_PUT_VAR( ncid, nsheatid, berg_grid%stored_heat(:,:) )
380      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var stored_heat failed')
381      IF( lwp ) WRITE(numout,*) 'file: ',TRIM(cl_filename),' var: stored_heat written'
382
383      nret = NF90_PUT_VAR( ncid, ncalvid , src_calving(:,:) )
384      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var calving failed')
385      nret = NF90_PUT_VAR( ncid, ncalvhid, src_calving_hflx(:,:) )
386      IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var calving_hflx failed')
387      IF( lwp ) WRITE(numout,*) 'file: ',TRIM(cl_filename),' var: calving written'
388
389      IF ( ASSOCIATED(first_berg) ) THEN
390
391         ! Write variables
392         ! just write out the current point of the trajectory
393
394         this => first_berg
395         jn = 0
396         DO WHILE (ASSOCIATED(this))
397            pt => this%current_point
398            jn=jn+1
399
400            nret = NF90_PUT_VAR(ncid, numberid, this%number, (/1,jn/), (/nkounts,1/) )
401            nret = NF90_PUT_VAR(ncid, nscaling_id, this%mass_scaling, (/ jn /) )
402
403            nret = NF90_PUT_VAR(ncid, nlonid, pt%lon, (/ jn /) )
404            nret = NF90_PUT_VAR(ncid, nlatid, pt%lat, (/ jn /) )
405            nret = NF90_PUT_VAR(ncid, nxid, pt%xi, (/ jn /) )
406            nret = NF90_PUT_VAR(ncid, nyid, pt%yj, (/ jn /) )
407            nret = NF90_PUT_VAR(ncid, nuvelid, pt%uvel, (/ jn /) )
408            nret = NF90_PUT_VAR(ncid, nvvelid, pt%vvel, (/ jn /) )
409            nret = NF90_PUT_VAR(ncid, nmassid, pt%mass, (/ jn /) )
410            nret = NF90_PUT_VAR(ncid, nthicknessid, pt%thickness, (/ jn /) )
411            nret = NF90_PUT_VAR(ncid, nwidthid, pt%width, (/ jn /) )
412            nret = NF90_PUT_VAR(ncid, nlengthid, pt%length, (/ jn /) )
413            nret = NF90_PUT_VAR(ncid, nyearid, pt%year, (/ jn /) )
414            nret = NF90_PUT_VAR(ncid, ndayid, pt%day, (/ jn /) )
415            nret = NF90_PUT_VAR(ncid, nmass_of_bits_id, pt%mass_of_bits, (/ jn /) )
416            nret = NF90_PUT_VAR(ncid, nheat_density_id, pt%heat_density, (/ jn /) )
417
418            this=>this%next
419         END DO
420         !
421      ENDIF ! associated(first_berg)
422
423      ! Finish up
424      nret = NF90_CLOSE(ncid)
425      IF (nret /= NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_close failed')
426      !
427   END SUBROUTINE icb_rst_write
428   !
429END MODULE icbrst
Note: See TracBrowser for help on using the repository browser.