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

source: branches/UKMO/dev_r5518_AMM15_package/NEMOGCM/NEMO/OPA_SRC/ICB/icbrst.F90 @ 10253

Last change on this file since 10253 was 10253, checked in by kingr, 5 years ago

Merged AMM15_v3_6_STABLE_package_collate@10237

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