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

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

NEMO branch dev_r3337_NOCS10_ICB: change all routine names and create more Gurvanistic havoc

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