1 | MODULE icbrst |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE icbrst *** |
---|
4 | !! Ocean physics: read and write iceberg restart files |
---|
5 | !!====================================================================== |
---|
6 | !! History : 3.3.1 ! 2010-01 (Martin&Adcroft) Original code |
---|
7 | !! - ! 2011-03 (Madec) Part conversion to NEMO form |
---|
8 | !! - ! Removal of mapping from another grid |
---|
9 | !! - ! 2011-04 (Alderson) Split into separate modules |
---|
10 | !! - ! 2011-04 (Alderson) Restore restart routine |
---|
11 | !! - ! Currently needs a fixed processor |
---|
12 | !! - ! layout between restarts |
---|
13 | !! - ! 2015-11 Dave Storkey Convert icb_rst_read to use IOM so can |
---|
14 | !! read single restart files |
---|
15 | !!---------------------------------------------------------------------- |
---|
16 | |
---|
17 | !!---------------------------------------------------------------------- |
---|
18 | !! icb_rst_read : read restart file |
---|
19 | !! icb_rst_write : write restart file |
---|
20 | !!---------------------------------------------------------------------- |
---|
21 | USE par_oce ! NEMO parameters |
---|
22 | USE dom_oce ! NEMO domain |
---|
23 | USE in_out_manager ! NEMO IO routines |
---|
24 | USE lib_mpp ! NEMO MPI library, lk_mpp in particular |
---|
25 | USE netcdf ! netcdf routines for IO |
---|
26 | USE iom |
---|
27 | USE ioipsl, ONLY : ju2ymds ! for calendar |
---|
28 | USE icb_oce ! define iceberg arrays |
---|
29 | USE icbutl ! iceberg utility routines |
---|
30 | |
---|
31 | IMPLICIT NONE |
---|
32 | PRIVATE |
---|
33 | |
---|
34 | PUBLIC icb_rst_read ! routine called in icbini.F90 module |
---|
35 | PUBLIC icb_rst_write ! routine called in icbstp.F90 module |
---|
36 | |
---|
37 | INTEGER :: nlonid, nlatid, nxid, nyid, nuvelid, nvvelid |
---|
38 | INTEGER :: nmassid, nthicknessid, nwidthid, nlengthid |
---|
39 | INTEGER :: nyearid, ndayid |
---|
40 | INTEGER :: nscaling_id, nmass_of_bits_id, nheat_density_id, numberid |
---|
41 | INTEGER :: nsiceid, nsheatid, ncalvid, ncalvhid, nkountid |
---|
42 | INTEGER :: nret, ncid, nc_dim |
---|
43 | |
---|
44 | INTEGER, DIMENSION(3) :: nstrt3, nlngth3 |
---|
45 | |
---|
46 | !!---------------------------------------------------------------------- |
---|
47 | !! NEMO/OCE 4.0 , NEMO Consortium (2018) |
---|
48 | !! $Id$ |
---|
49 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
50 | !!---------------------------------------------------------------------- |
---|
51 | CONTAINS |
---|
52 | |
---|
53 | SUBROUTINE icb_rst_read() |
---|
54 | !!---------------------------------------------------------------------- |
---|
55 | !! *** SUBROUTINE icb_rst_read *** |
---|
56 | !! |
---|
57 | !! ** Purpose : read a iceberg restart file |
---|
58 | !! NB: for this version, we just read back in the restart for this processor |
---|
59 | !! so we cannot change the processor layout currently with iceberg code |
---|
60 | !!---------------------------------------------------------------------- |
---|
61 | INTEGER :: idim, ivar, iatt |
---|
62 | INTEGER :: jn, iunlim_dim, ibergs_in_file |
---|
63 | INTEGER :: ii, ij, iclass, ibase_err, imax_icb |
---|
64 | REAL(wp), DIMENSION(nkounts) :: zdata |
---|
65 | LOGICAL :: ll_found_restart |
---|
66 | CHARACTER(len=256) :: cl_path |
---|
67 | CHARACTER(len=256) :: 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 | ! Find a restart file. Assume iceberg restarts in same directory as ocean restarts |
---|
74 | ! and are called TRIM(cn_ocerst)//'_icebergs' |
---|
75 | cl_path = TRIM(cn_ocerst_indir) |
---|
76 | IF( cl_path(LEN_TRIM(cl_path):) /= '/' ) cl_path = TRIM(cl_path) // '/' |
---|
77 | cl_filename = TRIM(cn_ocerst_in)//'_icebergs' |
---|
78 | CALL iom_open( TRIM(cl_path)//cl_filename, ncid ) |
---|
79 | |
---|
80 | imax_icb = 0 |
---|
81 | IF( iom_file(ncid)%iduld .GE. 0) THEN |
---|
82 | |
---|
83 | ibergs_in_file = iom_file(ncid)%lenuld |
---|
84 | DO jn = 1,ibergs_in_file |
---|
85 | |
---|
86 | ! iom_get treats the unlimited dimension as time. Here the unlimited dimension |
---|
87 | ! is the iceberg index, but we can still use the ktime keyword to get the iceberg we want. |
---|
88 | |
---|
89 | CALL iom_get( ncid, 'xi' ,localpt%xi , ktime=jn ) |
---|
90 | CALL iom_get( ncid, 'yj' ,localpt%yj , ktime=jn ) |
---|
91 | |
---|
92 | ii = INT( localpt%xi + 0.5 ) |
---|
93 | ij = INT( localpt%yj + 0.5 ) |
---|
94 | ! Only proceed if this iceberg is on the local processor (excluding halos). |
---|
95 | IF ( ii .GE. nldi+nimpp-1 .AND. ii .LE. nlei+nimpp-1 .AND. & |
---|
96 | & ij .GE. nldj+njmpp-1 .AND. ij .LE. nlej+njmpp-1 ) THEN |
---|
97 | |
---|
98 | CALL iom_get( ncid, jpdom_unknown, 'number' , zdata(:) , ktime=jn, kstart=(/1/), kcount=(/nkounts/) ) |
---|
99 | localberg%number(:) = INT(zdata(:)) |
---|
100 | imax_icb = MAX( imax_icb, INT(zdata(1)) ) |
---|
101 | CALL iom_get( ncid, 'mass_scaling' , localberg%mass_scaling, ktime=jn ) |
---|
102 | CALL iom_get( ncid, 'lon' , localpt%lon , ktime=jn ) |
---|
103 | CALL iom_get( ncid, 'lat' , localpt%lat , ktime=jn ) |
---|
104 | CALL iom_get( ncid, 'uvel' , localpt%uvel , ktime=jn ) |
---|
105 | CALL iom_get( ncid, 'vvel' , localpt%vvel , ktime=jn ) |
---|
106 | CALL iom_get( ncid, 'mass' , localpt%mass , ktime=jn ) |
---|
107 | CALL iom_get( ncid, 'thickness' , localpt%thickness , ktime=jn ) |
---|
108 | CALL iom_get( ncid, 'width' , localpt%width , ktime=jn ) |
---|
109 | CALL iom_get( ncid, 'length' , localpt%length , ktime=jn ) |
---|
110 | CALL iom_get( ncid, 'year' , zdata(1) , ktime=jn ) |
---|
111 | localpt%year = INT(zdata(1)) |
---|
112 | CALL iom_get( ncid, 'day' , localpt%day , ktime=jn ) |
---|
113 | CALL iom_get( ncid, 'mass_of_bits' , localpt%mass_of_bits , ktime=jn ) |
---|
114 | CALL iom_get( ncid, 'heat_density' , localpt%heat_density , ktime=jn ) |
---|
115 | ! |
---|
116 | CALL icb_utl_add( localberg, localpt ) |
---|
117 | ! |
---|
118 | ENDIF |
---|
119 | ! |
---|
120 | END DO |
---|
121 | ! |
---|
122 | ELSE |
---|
123 | ibergs_in_file = 0 |
---|
124 | ENDIF |
---|
125 | |
---|
126 | ! Gridded variables |
---|
127 | CALL iom_get( ncid, jpdom_autoglo, 'calving' , src_calving ) |
---|
128 | CALL iom_get( ncid, jpdom_autoglo, 'calving_hflx', src_calving_hflx ) |
---|
129 | CALL iom_get( ncid, jpdom_autoglo, 'stored_heat' , berg_grid%stored_heat ) |
---|
130 | CALL iom_get( ncid, jpdom_autoglo_xy, 'stored_ice' , berg_grid%stored_ice, kstart=(/1,1,1/), kcount=(/1,1,nclasses/) ) |
---|
131 | |
---|
132 | CALL iom_get( ncid, jpdom_unknown, 'kount' , zdata(:) ) |
---|
133 | num_bergs(:) = INT(zdata(:)) |
---|
134 | ! Close file |
---|
135 | CALL iom_close( ncid ) |
---|
136 | ! |
---|
137 | |
---|
138 | ! Sanity checks |
---|
139 | jn = icb_utl_count() |
---|
140 | IF ( lwp .AND. nn_verbose_level >= 0 ) & |
---|
141 | WRITE(numout,'(2(a,i5))') 'icebergs, read_restart_bergs: # bergs =',jn,' on PE',narea-1 |
---|
142 | IF( lk_mpp ) THEN |
---|
143 | ! Only mpp_sum ibergs_in_file if we are reading from multiple restart files. |
---|
144 | IF( INDEX(iom_file(ncid)%name,'icebergs.nc' ) .EQ. 0 ) CALL mpp_sum('icbrst', ibergs_in_file) |
---|
145 | CALL mpp_sum('icbrst', jn) |
---|
146 | ENDIF |
---|
147 | IF( lwp ) WRITE(numout,'(a,i5,a,i5,a)') 'icebergs, icb_rst_read: there were',ibergs_in_file, & |
---|
148 | & ' bergs in the restart file and', jn,' bergs have been read' |
---|
149 | ! |
---|
150 | ! Confirm that all areas have a suitable base for assigning new iceberg |
---|
151 | ! numbers. This will not be the case if restarting from a collated dataset |
---|
152 | ! (even if using the same processor decomposition) |
---|
153 | ! |
---|
154 | ibase_err = 0 |
---|
155 | IF( num_bergs(1) < 0 .AND. num_bergs(1) /= narea - jpnij ) THEN |
---|
156 | ! If this area has never calved a new berg then the base should be |
---|
157 | ! set to narea - jpnij. If it is negative but something else then |
---|
158 | ! a new base will be needed to guarantee unique, future iceberg numbers |
---|
159 | ibase_err = 1 |
---|
160 | ELSEIF( MOD( num_bergs(1) - narea , jpnij ) /= 0 ) THEN |
---|
161 | ! If this area has a base which is not in the set {narea + N*jpnij} |
---|
162 | ! for positive integers N then a new base will be needed to guarantee |
---|
163 | ! unique, future iceberg numbers |
---|
164 | ibase_err = 1 |
---|
165 | ENDIF |
---|
166 | IF( lk_mpp ) THEN |
---|
167 | CALL mpp_sum('icbrst', ibase_err) |
---|
168 | ENDIF |
---|
169 | IF( ibase_err > 0 ) THEN |
---|
170 | ! |
---|
171 | ! A new base is needed. The only secure solution is to set bases such that |
---|
172 | ! all future icebergs numbers will be greater than the current global maximum |
---|
173 | IF( lk_mpp ) THEN |
---|
174 | CALL mpp_max('icbrst', imax_icb) |
---|
175 | ENDIF |
---|
176 | num_bergs(1) = imax_icb - jpnij + narea |
---|
177 | ENDIF |
---|
178 | ! |
---|
179 | IF( lwp .AND. nn_verbose_level >= 0 ) WRITE(numout,'(a)') 'icebergs, icb_rst_read: completed' |
---|
180 | ! |
---|
181 | END SUBROUTINE icb_rst_read |
---|
182 | |
---|
183 | |
---|
184 | SUBROUTINE icb_rst_write( kt ) |
---|
185 | !!---------------------------------------------------------------------- |
---|
186 | !! *** SUBROUTINE icb_rst_write *** |
---|
187 | !! |
---|
188 | !!---------------------------------------------------------------------- |
---|
189 | INTEGER, INTENT( in ) :: kt |
---|
190 | ! |
---|
191 | INTEGER :: jn ! dummy loop index |
---|
192 | INTEGER :: ix_dim, iy_dim, ik_dim, in_dim |
---|
193 | INTEGER :: iyear, imonth, iday |
---|
194 | REAL (wp) :: zsec |
---|
195 | REAL (wp) :: zfjulday |
---|
196 | CHARACTER(len=256) :: cl_path |
---|
197 | CHARACTER(len=256) :: cl_filename |
---|
198 | CHARACTER(LEN=20) :: clkt ! ocean time-step deine as a character |
---|
199 | TYPE(iceberg), POINTER :: this |
---|
200 | TYPE(point) , POINTER :: pt |
---|
201 | !!---------------------------------------------------------------------- |
---|
202 | RETURN |
---|
203 | ! Following the normal restart procedure, this routine will be called |
---|
204 | ! the timestep before a restart stage as well as the restart timestep. |
---|
205 | ! This is a performance step enabling the file to be opened and contents |
---|
206 | ! defined in advance of the write. This is not possible with icebergs |
---|
207 | ! since the number of bergs to be written could change between timesteps |
---|
208 | IF( kt == nitrst ) THEN |
---|
209 | ! Only operate on the restart timestep itself. |
---|
210 | ! Assume we write iceberg restarts to same directory as ocean restarts. |
---|
211 | cl_path = TRIM(cn_ocerst_outdir) |
---|
212 | IF( cl_path(LEN_TRIM(cl_path):) /= '/' ) cl_path = TRIM(cl_path) // '/' |
---|
213 | IF ( ln_rstdate ) THEN |
---|
214 | zfjulday = fjulday + rdt / rday |
---|
215 | IF( ABS(zfjulday - REAL(NINT(zfjulday),wp)) < 0.1 / rday ) zfjulday = REAL(NINT(zfjulday),wp) ! avoid truncation error |
---|
216 | CALL ju2ymds( zfjulday, iyear, imonth, iday, zsec ) |
---|
217 | WRITE(clkt, '(i4.4,2i2.2)') iyear, imonth, iday |
---|
218 | ELSE |
---|
219 | IF( kt > 999999999 ) THEN ; WRITE(clkt, * ) kt |
---|
220 | ELSE ; WRITE(clkt, '(i8.8)') kt |
---|
221 | ENDIF |
---|
222 | ENDIF |
---|
223 | IF( lk_mpp ) THEN |
---|
224 | WRITE(cl_filename,'(A,"_icebergs_",A,"_restart_",I6.6,".nc")') TRIM(cexper), TRIM(ADJUSTL(clkt)), narea-1 |
---|
225 | ELSE |
---|
226 | WRITE(cl_filename,'(A,"_icebergs_",A,"_restart.nc")') TRIM(cexper), TRIM(ADJUSTL(clkt)) |
---|
227 | ENDIF |
---|
228 | IF ( lwp .AND. nn_verbose_level >= 0) WRITE(numout,'(2a)') 'icebergs, write_restart: creating ', & |
---|
229 | & TRIM(cl_path)//TRIM(cl_filename) |
---|
230 | |
---|
231 | nret = NF90_CREATE(TRIM(cl_path)//TRIM(cl_filename), NF90_CLOBBER, ncid) |
---|
232 | IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_create failed') |
---|
233 | |
---|
234 | ! Dimensions |
---|
235 | nret = NF90_DEF_DIM(ncid, 'x', jpi, ix_dim) |
---|
236 | IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim x failed') |
---|
237 | |
---|
238 | nret = NF90_DEF_DIM(ncid, 'y', jpj, iy_dim) |
---|
239 | IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim y failed') |
---|
240 | |
---|
241 | nret = NF90_DEF_DIM(ncid, 'c', nclasses, nc_dim) |
---|
242 | IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim c failed') |
---|
243 | |
---|
244 | nret = NF90_DEF_DIM(ncid, 'k', nkounts, ik_dim) |
---|
245 | IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim k failed') |
---|
246 | |
---|
247 | ! global attributes |
---|
248 | IF( lk_mpp ) THEN |
---|
249 | ! Set domain parameters (assume jpdom_local_full) |
---|
250 | nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_number_total' , jpnij ) |
---|
251 | nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_number' , narea-1 ) |
---|
252 | nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_dimensions_ids' , (/1 , 2 /) ) |
---|
253 | nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_size_global' , (/jpiglo, jpjglo/) ) |
---|
254 | nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_size_local' , (/jpi , jpj /) ) |
---|
255 | nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_position_first' , (/nimpp , njmpp /) ) |
---|
256 | nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_position_last' , (/nimpp + jpi - 1 , njmpp + jpj - 1 /) ) |
---|
257 | nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_halo_size_start', (/nldi - 1 , nldj - 1 /) ) |
---|
258 | nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_halo_size_end' , (/jpi - nlei , jpj - nlej /) ) |
---|
259 | nret = NF90_PUT_ATT( ncid, NF90_GLOBAL, 'DOMAIN_type' , 'BOX' ) |
---|
260 | ENDIF |
---|
261 | |
---|
262 | IF (associated(first_berg)) then |
---|
263 | nret = NF90_DEF_DIM(ncid, 'n', NF90_UNLIMITED, in_dim) |
---|
264 | IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_def_dim n failed') |
---|
265 | ENDIF |
---|
266 | |
---|
267 | ! Variables |
---|
268 | nret = NF90_DEF_VAR(ncid, 'kount' , NF90_INT , (/ ik_dim /), nkountid) |
---|
269 | nret = NF90_DEF_VAR(ncid, 'calving' , NF90_DOUBLE, (/ ix_dim, iy_dim /), ncalvid) |
---|
270 | nret = NF90_DEF_VAR(ncid, 'calving_hflx', NF90_DOUBLE, (/ ix_dim, iy_dim /), ncalvhid) |
---|
271 | nret = NF90_DEF_VAR(ncid, 'stored_ice' , NF90_DOUBLE, (/ ix_dim, iy_dim, nc_dim /), nsiceid) |
---|
272 | nret = NF90_DEF_VAR(ncid, 'stored_heat' , NF90_DOUBLE, (/ ix_dim, iy_dim /), nsheatid) |
---|
273 | |
---|
274 | ! Attributes |
---|
275 | nret = NF90_PUT_ATT(ncid, ncalvid , 'long_name', 'iceberg calving') |
---|
276 | nret = NF90_PUT_ATT(ncid, ncalvid , 'units', 'some') |
---|
277 | nret = NF90_PUT_ATT(ncid, ncalvhid, 'long_name', 'heat flux associated with iceberg calving') |
---|
278 | nret = NF90_PUT_ATT(ncid, ncalvhid, 'units', 'some') |
---|
279 | nret = NF90_PUT_ATT(ncid, nsiceid , 'long_name', 'stored ice used to calve icebergs') |
---|
280 | nret = NF90_PUT_ATT(ncid, nsiceid , 'units', 'kg/s') |
---|
281 | nret = NF90_PUT_ATT(ncid, nsheatid, 'long_name', 'heat in stored ice used to calve icebergs') |
---|
282 | nret = NF90_PUT_ATT(ncid, nsheatid, 'units', 'J/kg/s') |
---|
283 | |
---|
284 | IF ( ASSOCIATED(first_berg) ) THEN |
---|
285 | |
---|
286 | ! Only add berg variables for this PE if we have anything to say |
---|
287 | |
---|
288 | ! Variables |
---|
289 | nret = NF90_DEF_VAR(ncid, 'lon', NF90_DOUBLE, in_dim, nlonid) |
---|
290 | nret = NF90_DEF_VAR(ncid, 'lat', NF90_DOUBLE, in_dim, nlatid) |
---|
291 | nret = NF90_DEF_VAR(ncid, 'xi', NF90_DOUBLE, in_dim, nxid) |
---|
292 | nret = NF90_DEF_VAR(ncid, 'yj', NF90_DOUBLE, in_dim, nyid) |
---|
293 | nret = NF90_DEF_VAR(ncid, 'uvel', NF90_DOUBLE, in_dim, nuvelid) |
---|
294 | nret = NF90_DEF_VAR(ncid, 'vvel', NF90_DOUBLE, in_dim, nvvelid) |
---|
295 | nret = NF90_DEF_VAR(ncid, 'mass', NF90_DOUBLE, in_dim, nmassid) |
---|
296 | nret = NF90_DEF_VAR(ncid, 'thickness', NF90_DOUBLE, in_dim, nthicknessid) |
---|
297 | nret = NF90_DEF_VAR(ncid, 'width', NF90_DOUBLE, in_dim, nwidthid) |
---|
298 | nret = NF90_DEF_VAR(ncid, 'length', NF90_DOUBLE, in_dim, nlengthid) |
---|
299 | nret = NF90_DEF_VAR(ncid, 'number', NF90_INT, (/ik_dim,in_dim/), numberid) |
---|
300 | nret = NF90_DEF_VAR(ncid, 'year', NF90_INT, in_dim, nyearid) |
---|
301 | nret = NF90_DEF_VAR(ncid, 'day', NF90_DOUBLE, in_dim, ndayid) |
---|
302 | nret = NF90_DEF_VAR(ncid, 'mass_scaling', NF90_DOUBLE, in_dim, nscaling_id) |
---|
303 | nret = NF90_DEF_VAR(ncid, 'mass_of_bits', NF90_DOUBLE, in_dim, nmass_of_bits_id) |
---|
304 | nret = NF90_DEF_VAR(ncid, 'heat_density', NF90_DOUBLE, in_dim, nheat_density_id) |
---|
305 | |
---|
306 | ! Attributes |
---|
307 | nret = NF90_PUT_ATT(ncid, nlonid, 'long_name', 'longitude') |
---|
308 | nret = NF90_PUT_ATT(ncid, nlonid, 'units', 'degrees_E') |
---|
309 | nret = NF90_PUT_ATT(ncid, nlatid, 'long_name', 'latitude') |
---|
310 | nret = NF90_PUT_ATT(ncid, nlatid, 'units', 'degrees_N') |
---|
311 | nret = NF90_PUT_ATT(ncid, nxid, 'long_name', 'x grid box position') |
---|
312 | nret = NF90_PUT_ATT(ncid, nxid, 'units', 'fractional') |
---|
313 | nret = NF90_PUT_ATT(ncid, nyid, 'long_name', 'y grid box position') |
---|
314 | nret = NF90_PUT_ATT(ncid, nyid, 'units', 'fractional') |
---|
315 | nret = NF90_PUT_ATT(ncid, nuvelid, 'long_name', 'zonal velocity') |
---|
316 | nret = NF90_PUT_ATT(ncid, nuvelid, 'units', 'm/s') |
---|
317 | nret = NF90_PUT_ATT(ncid, nvvelid, 'long_name', 'meridional velocity') |
---|
318 | nret = NF90_PUT_ATT(ncid, nvvelid, 'units', 'm/s') |
---|
319 | nret = NF90_PUT_ATT(ncid, nmassid, 'long_name', 'mass') |
---|
320 | nret = NF90_PUT_ATT(ncid, nmassid, 'units', 'kg') |
---|
321 | nret = NF90_PUT_ATT(ncid, nthicknessid, 'long_name', 'thickness') |
---|
322 | nret = NF90_PUT_ATT(ncid, nthicknessid, 'units', 'm') |
---|
323 | nret = NF90_PUT_ATT(ncid, nwidthid, 'long_name', 'width') |
---|
324 | nret = NF90_PUT_ATT(ncid, nwidthid, 'units', 'm') |
---|
325 | nret = NF90_PUT_ATT(ncid, nlengthid, 'long_name', 'length') |
---|
326 | nret = NF90_PUT_ATT(ncid, nlengthid, 'units', 'm') |
---|
327 | nret = NF90_PUT_ATT(ncid, numberid, 'long_name', 'iceberg number on this processor') |
---|
328 | nret = NF90_PUT_ATT(ncid, numberid, 'units', 'count') |
---|
329 | nret = NF90_PUT_ATT(ncid, nyearid, 'long_name', 'calendar year of calving event') |
---|
330 | nret = NF90_PUT_ATT(ncid, nyearid, 'units', 'years') |
---|
331 | nret = NF90_PUT_ATT(ncid, ndayid, 'long_name', 'year day of calving event') |
---|
332 | nret = NF90_PUT_ATT(ncid, ndayid, 'units', 'days') |
---|
333 | nret = NF90_PUT_ATT(ncid, nscaling_id, 'long_name', 'scaling factor for mass of calving berg') |
---|
334 | nret = NF90_PUT_ATT(ncid, nscaling_id, 'units', 'none') |
---|
335 | nret = NF90_PUT_ATT(ncid, nmass_of_bits_id, 'long_name', 'mass of bergy bits') |
---|
336 | nret = NF90_PUT_ATT(ncid, nmass_of_bits_id, 'units', 'kg') |
---|
337 | nret = NF90_PUT_ATT(ncid, nheat_density_id, 'long_name', 'heat density') |
---|
338 | nret = NF90_PUT_ATT(ncid, nheat_density_id, 'units', 'J/kg') |
---|
339 | |
---|
340 | ENDIF ! associated(first_berg) |
---|
341 | |
---|
342 | ! End define mode |
---|
343 | nret = NF90_ENDDEF(ncid) |
---|
344 | |
---|
345 | ! -------------------------------- |
---|
346 | ! now write some data |
---|
347 | |
---|
348 | nstrt3(1) = 1 |
---|
349 | nstrt3(2) = 1 |
---|
350 | nlngth3(1) = jpi |
---|
351 | nlngth3(2) = jpj |
---|
352 | nlngth3(3) = 1 |
---|
353 | |
---|
354 | DO jn=1,nclasses |
---|
355 | griddata(:,:,1) = berg_grid%stored_ice(:,:,jn) |
---|
356 | nstrt3(3) = jn |
---|
357 | nret = NF90_PUT_VAR( ncid, nsiceid, griddata, nstrt3, nlngth3 ) |
---|
358 | IF (nret .ne. NF90_NOERR) THEN |
---|
359 | IF( lwp ) WRITE(numout,*) TRIM(NF90_STRERROR( nret )) |
---|
360 | CALL ctl_stop('icebergs, write_restart: nf_put_var stored_ice failed') |
---|
361 | ENDIF |
---|
362 | ENDDO |
---|
363 | IF( lwp ) WRITE(numout,*) 'file: ',TRIM(cl_path)//TRIM(cl_filename),' var: stored_ice written' |
---|
364 | |
---|
365 | nret = NF90_PUT_VAR( ncid, nkountid, num_bergs(:) ) |
---|
366 | IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var kount failed') |
---|
367 | |
---|
368 | nret = NF90_PUT_VAR( ncid, nsheatid, berg_grid%stored_heat(:,:) ) |
---|
369 | IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var stored_heat failed') |
---|
370 | IF( lwp ) WRITE(numout,*) 'file: ',TRIM(cl_path)//TRIM(cl_filename),' var: stored_heat written' |
---|
371 | |
---|
372 | nret = NF90_PUT_VAR( ncid, ncalvid , src_calving(:,:) ) |
---|
373 | IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var calving failed') |
---|
374 | nret = NF90_PUT_VAR( ncid, ncalvhid, src_calving_hflx(:,:) ) |
---|
375 | IF (nret .ne. NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_put_var calving_hflx failed') |
---|
376 | IF( lwp ) WRITE(numout,*) 'file: ',TRIM(cl_path)//TRIM(cl_filename),' var: calving written' |
---|
377 | |
---|
378 | IF ( ASSOCIATED(first_berg) ) THEN |
---|
379 | |
---|
380 | ! Write variables |
---|
381 | ! just write out the current point of the trajectory |
---|
382 | |
---|
383 | this => first_berg |
---|
384 | jn = 0 |
---|
385 | DO WHILE (ASSOCIATED(this)) |
---|
386 | pt => this%current_point |
---|
387 | jn=jn+1 |
---|
388 | |
---|
389 | nret = NF90_PUT_VAR(ncid, numberid, this%number, (/1,jn/), (/nkounts,1/) ) |
---|
390 | nret = NF90_PUT_VAR(ncid, nscaling_id, this%mass_scaling, (/ jn /) ) |
---|
391 | |
---|
392 | nret = NF90_PUT_VAR(ncid, nlonid, pt%lon, (/ jn /) ) |
---|
393 | nret = NF90_PUT_VAR(ncid, nlatid, pt%lat, (/ jn /) ) |
---|
394 | nret = NF90_PUT_VAR(ncid, nxid, pt%xi, (/ jn /) ) |
---|
395 | nret = NF90_PUT_VAR(ncid, nyid, pt%yj, (/ jn /) ) |
---|
396 | nret = NF90_PUT_VAR(ncid, nuvelid, pt%uvel, (/ jn /) ) |
---|
397 | nret = NF90_PUT_VAR(ncid, nvvelid, pt%vvel, (/ jn /) ) |
---|
398 | nret = NF90_PUT_VAR(ncid, nmassid, pt%mass, (/ jn /) ) |
---|
399 | nret = NF90_PUT_VAR(ncid, nthicknessid, pt%thickness, (/ jn /) ) |
---|
400 | nret = NF90_PUT_VAR(ncid, nwidthid, pt%width, (/ jn /) ) |
---|
401 | nret = NF90_PUT_VAR(ncid, nlengthid, pt%length, (/ jn /) ) |
---|
402 | nret = NF90_PUT_VAR(ncid, nyearid, pt%year, (/ jn /) ) |
---|
403 | nret = NF90_PUT_VAR(ncid, ndayid, pt%day, (/ jn /) ) |
---|
404 | nret = NF90_PUT_VAR(ncid, nmass_of_bits_id, pt%mass_of_bits, (/ jn /) ) |
---|
405 | nret = NF90_PUT_VAR(ncid, nheat_density_id, pt%heat_density, (/ jn /) ) |
---|
406 | |
---|
407 | this=>this%next |
---|
408 | END DO |
---|
409 | ! |
---|
410 | ENDIF ! associated(first_berg) |
---|
411 | |
---|
412 | ! Finish up |
---|
413 | nret = NF90_CLOSE(ncid) |
---|
414 | IF (nret /= NF90_NOERR) CALL ctl_stop('icebergs, write_restart: nf_close failed') |
---|
415 | |
---|
416 | ! Sanity check |
---|
417 | jn = icb_utl_count() |
---|
418 | IF ( lwp .AND. nn_verbose_level >= 0) & |
---|
419 | WRITE(numout,'(2(a,i5))') 'icebergs, icb_rst_write: # bergs =',jn,' on PE',narea-1 |
---|
420 | IF( lk_mpp ) THEN |
---|
421 | CALL mpp_sum('icbrst', jn) |
---|
422 | ENDIF |
---|
423 | IF(lwp) WRITE(numout,'(a,i5,a,i5,a)') 'icebergs, icb_rst_write: ', jn, & |
---|
424 | & ' bergs in total have been written at timestep ', kt |
---|
425 | ! |
---|
426 | ! Finish up |
---|
427 | ! |
---|
428 | ENDIF |
---|
429 | END SUBROUTINE icb_rst_write |
---|
430 | ! |
---|
431 | !!====================================================================== |
---|
432 | END MODULE icbrst |
---|