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

source: branches/UKMO/2014_Surge_Modelling/NEMOGCM/NEMO/OPA_SRC/ICB/icbrst.F90 @ 5746

Last change on this file since 5746 was 5746, checked in by rfurner, 9 years ago

changes to name restart files with date instead of time stamp

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