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

source: branches/UKMO/dev_r5518_Surge_Modelling/NEMOGCM/NEMO/OPA_SRC/ICB/icbrst.F90 @ 6389

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

added surge code from 2014_Surge_Modelling branch

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