Line data Source code
1 : #ifdef ncdf
2 : #define USE_NETCDF
3 : #endif
4 : !=======================================================================
5 :
6 : ! Read and write ice model restart files using netCDF or binary
7 : ! interfaces.
8 : ! authors David A Bailey, NCAR
9 :
10 : module ice_restart
11 :
12 : use ice_broadcast
13 : use ice_communicate, only: my_task, master_task
14 : use ice_kinds_mod
15 : #ifdef USE_NETCDF
16 : use netcdf
17 : #endif
18 : use ice_restart_shared, only: &
19 : restart_ext, restart_dir, restart_file, pointer_file, & ! LCOV_EXCL_LINE
20 : runid, use_restart_time, lcdf64, lenstr, restart_coszen
21 : use ice_fileunits, only: nu_diag, nu_rst_pointer
22 : use ice_exit, only: abort_ice
23 : use icepack_intfc, only: icepack_query_parameters
24 : use icepack_intfc, only: icepack_query_tracer_sizes
25 : use icepack_intfc, only: icepack_query_tracer_flags
26 : use icepack_intfc, only: icepack_warnings_flush, icepack_warnings_aborted
27 :
28 : implicit none
29 : private
30 : public :: init_restart_write, init_restart_read, &
31 : read_restart_field, write_restart_field, final_restart, & ! LCOV_EXCL_LINE
32 : query_field
33 :
34 : integer (kind=int_kind) :: ncid
35 :
36 : !=======================================================================
37 :
38 : contains
39 :
40 : !=======================================================================
41 :
42 : ! Sets up restart file for reading.
43 : ! author David A Bailey, NCAR
44 :
45 21 : subroutine init_restart_read(ice_ic)
46 :
47 : use ice_calendar, only: msec, mmonth, mday, myear, &
48 : istep0, istep1, npt
49 :
50 : character(len=char_len_long), intent(in), optional :: ice_ic
51 :
52 : ! local variables
53 :
54 : character(len=char_len_long) :: &
55 : filename, filename0
56 :
57 : integer (kind=int_kind) :: status, status1
58 :
59 : character(len=*), parameter :: subname = '(init_restart_read)'
60 :
61 : #ifdef USE_NETCDF
62 29 : if (present(ice_ic)) then
63 17 : filename = trim(ice_ic)
64 : else
65 12 : if (my_task == master_task) then
66 2 : open(nu_rst_pointer,file=pointer_file)
67 2 : read(nu_rst_pointer,'(a)') filename0
68 2 : filename = trim(filename0)
69 2 : close(nu_rst_pointer)
70 2 : write(nu_diag,*) 'Read ',pointer_file(1:lenstr(pointer_file))
71 : endif
72 12 : call broadcast_scalar(filename, master_task)
73 : endif
74 :
75 29 : if (my_task == master_task) then
76 6 : write(nu_diag,*) 'Using restart dump=', trim(filename)
77 :
78 6 : status = nf90_open(trim(filename), nf90_nowrite, ncid)
79 6 : if (status /= nf90_noerr) call abort_ice(subname// &
80 0 : 'ERROR: reading restart ncfile '//trim(filename))
81 :
82 6 : if (use_restart_time) then
83 2 : status1 = nf90_noerr
84 2 : status = nf90_get_att(ncid, nf90_global, 'istep1', istep0)
85 2 : if (status /= nf90_noerr) status1 = status
86 : ! status = nf90_get_att(ncid, nf90_global, 'time', time)
87 : ! status = nf90_get_att(ncid, nf90_global, 'time_forc', time_forc)
88 2 : status = nf90_get_att(ncid, nf90_global, 'myear', myear)
89 2 : if (status /= nf90_noerr) status = nf90_get_att(ncid, nf90_global, 'nyr', myear)
90 2 : if (status /= nf90_noerr) status1 = status
91 2 : status = nf90_get_att(ncid, nf90_global, 'mmonth', mmonth)
92 2 : if (status /= nf90_noerr) status = nf90_get_att(ncid, nf90_global, 'month', mmonth)
93 2 : if (status /= nf90_noerr) status1 = status
94 2 : status = nf90_get_att(ncid, nf90_global, 'mday', mday)
95 2 : if (status /= nf90_noerr) status1 = status
96 2 : status = nf90_get_att(ncid, nf90_global, 'msec', msec)
97 2 : if (status /= nf90_noerr) status = nf90_get_att(ncid, nf90_global, 'sec', msec)
98 2 : if (status /= nf90_noerr) status1 = status
99 2 : if (status1 /= nf90_noerr) call abort_ice(subname// &
100 0 : 'ERROR: reading restart time '//trim(filename))
101 : endif ! use namelist values if use_restart_time = F
102 :
103 : endif
104 :
105 29 : call broadcast_scalar(istep0,master_task)
106 : ! call broadcast_scalar(time,master_task)
107 29 : call broadcast_scalar(myear,master_task)
108 29 : call broadcast_scalar(mmonth,master_task)
109 29 : call broadcast_scalar(mday,master_task)
110 29 : call broadcast_scalar(msec,master_task)
111 : ! call broadcast_scalar(time_forc,master_task)
112 :
113 29 : istep1 = istep0
114 :
115 : ! if runid is bering then need to correct npt for istep0
116 29 : if (trim(runid) == 'bering') then
117 0 : npt = npt - istep0
118 : endif
119 : #else
120 : call abort_ice(subname//'ERROR: USE_NETCDF cpp not defined for '//trim(ice_ic), &
121 : file=__FILE__, line=__LINE__)
122 : #endif
123 :
124 29 : end subroutine init_restart_read
125 :
126 : !=======================================================================
127 :
128 : ! Sets up restart file for writing.
129 : ! author David A Bailey, NCAR
130 :
131 37 : subroutine init_restart_write(filename_spec)
132 :
133 8 : use ice_blocks, only: nghost
134 : use ice_calendar, only: msec, mmonth, mday, myear, istep1
135 : use ice_domain_size, only: nx_global, ny_global, ncat, nilyr, nslyr, &
136 : n_iso, n_aero, nblyr, n_zaero, n_algae, n_doc, & ! LCOV_EXCL_LINE
137 : n_dic, n_don, n_fed, n_fep, nfsd
138 : use ice_arrays_column, only: oceanmixed_ice
139 : use ice_dyn_shared, only: kdyn
140 : use ice_grid, only: grid_ice
141 :
142 : character(len=char_len_long), intent(in), optional :: filename_spec
143 :
144 : ! local variables
145 :
146 : logical (kind=log_kind) :: &
147 : skl_bgc, z_tracers, tr_fsd, & ! LCOV_EXCL_LINE
148 : tr_iage, tr_FY, tr_lvl, tr_iso, tr_aero, & ! LCOV_EXCL_LINE
149 : tr_pond_topo, tr_pond_lvl, tr_brine, tr_snow, & ! LCOV_EXCL_LINE
150 : tr_bgc_N, tr_bgc_C, tr_bgc_Nit, & ! LCOV_EXCL_LINE
151 : tr_bgc_Sil, tr_bgc_DMS, & ! LCOV_EXCL_LINE
152 : tr_bgc_chl, tr_bgc_Am, & ! LCOV_EXCL_LINE
153 : tr_bgc_PON, tr_bgc_DON, & ! LCOV_EXCL_LINE
154 : tr_zaero, tr_bgc_Fe, & ! LCOV_EXCL_LINE
155 : tr_bgc_hum
156 :
157 : integer (kind=int_kind) :: &
158 : k, n, & ! index ! LCOV_EXCL_LINE
159 : nx, ny, & ! global array size ! LCOV_EXCL_LINE
160 : nbtrcr ! number of bgc tracers
161 :
162 : character(len=char_len_long) :: filename
163 :
164 49 : integer (kind=int_kind), allocatable :: dims(:)
165 :
166 : integer (kind=int_kind) :: &
167 : dimid_ni, & ! netCDF identifiers ! LCOV_EXCL_LINE
168 : dimid_nj, & ! ! LCOV_EXCL_LINE
169 : dimid_ncat, & ! ! LCOV_EXCL_LINE
170 : iflag, & ! netCDF creation flag ! LCOV_EXCL_LINE
171 : status ! status variable from netCDF routine
172 :
173 : character (len=3) :: nchar, ncharb
174 :
175 : character(len=*), parameter :: subname = '(init_restart_write)'
176 :
177 : #ifdef USE_NETCDF
178 : call icepack_query_parameters( &
179 49 : skl_bgc_out=skl_bgc, z_tracers_out=z_tracers)
180 : call icepack_query_tracer_sizes( &
181 49 : nbtrcr_out=nbtrcr)
182 : call icepack_query_tracer_flags( &
183 : tr_iage_out=tr_iage, tr_FY_out=tr_FY, tr_lvl_out=tr_lvl, tr_fsd_out=tr_fsd, & ! LCOV_EXCL_LINE
184 : tr_iso_out=tr_iso, tr_aero_out=tr_aero, & ! LCOV_EXCL_LINE
185 : tr_pond_topo_out=tr_pond_topo, tr_pond_lvl_out=tr_pond_lvl, & ! LCOV_EXCL_LINE
186 : tr_snow_out=tr_snow, tr_brine_out=tr_brine, & ! LCOV_EXCL_LINE
187 : tr_bgc_N_out=tr_bgc_N, tr_bgc_C_out=tr_bgc_C, tr_bgc_Nit_out=tr_bgc_Nit, & ! LCOV_EXCL_LINE
188 : tr_bgc_Sil_out=tr_bgc_Sil, tr_bgc_DMS_out=tr_bgc_DMS, & ! LCOV_EXCL_LINE
189 : tr_bgc_chl_out=tr_bgc_chl, tr_bgc_Am_out=tr_bgc_Am, & ! LCOV_EXCL_LINE
190 : tr_bgc_PON_out=tr_bgc_PON, tr_bgc_DON_out=tr_bgc_DON, & ! LCOV_EXCL_LINE
191 : tr_zaero_out=tr_zaero, tr_bgc_Fe_out=tr_bgc_Fe, & ! LCOV_EXCL_LINE
192 49 : tr_bgc_hum_out=tr_bgc_hum)
193 49 : call icepack_warnings_flush(nu_diag)
194 49 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
195 0 : file=__FILE__, line=__LINE__)
196 :
197 : ! construct path/file
198 49 : if (present(filename_spec)) then
199 0 : filename = trim(filename_spec)
200 : else
201 12 : write(filename,'(a,a,a,i4.4,a,i2.2,a,i2.2,a,i5.5)') &
202 : restart_dir(1:lenstr(restart_dir)), & ! LCOV_EXCL_LINE
203 : restart_file(1:lenstr(restart_file)),'.', & ! LCOV_EXCL_LINE
204 98 : myear,'-',mmonth,'-',mday,'-',msec
205 : end if
206 :
207 : ! write pointer (path/file)
208 49 : if (my_task == master_task) then
209 9 : filename = trim(filename) // '.nc'
210 9 : open(nu_rst_pointer,file=pointer_file)
211 9 : write(nu_rst_pointer,'(a)') filename
212 9 : close(nu_rst_pointer)
213 :
214 9 : iflag = 0
215 9 : if (lcdf64) iflag = nf90_64bit_offset
216 9 : status = nf90_create(trim(filename), iflag, ncid)
217 9 : if (status /= nf90_noerr) call abort_ice(subname// &
218 0 : 'ERROR: creating restart ncfile '//trim(filename))
219 :
220 9 : status = nf90_put_att(ncid,nf90_global,'istep1',istep1)
221 : ! status = nf90_put_att(ncid,nf90_global,'time',time)
222 : ! status = nf90_put_att(ncid,nf90_global,'time_forc',time_forc)
223 9 : status = nf90_put_att(ncid,nf90_global,'myear',myear)
224 9 : status = nf90_put_att(ncid,nf90_global,'mmonth',mmonth)
225 9 : status = nf90_put_att(ncid,nf90_global,'mday',mday)
226 9 : status = nf90_put_att(ncid,nf90_global,'msec',msec)
227 :
228 9 : nx = nx_global
229 9 : ny = ny_global
230 9 : if (restart_ext) then
231 0 : nx = nx_global + 2*nghost
232 0 : ny = ny_global + 2*nghost
233 : endif
234 9 : status = nf90_def_dim(ncid,'ni',nx,dimid_ni)
235 9 : status = nf90_def_dim(ncid,'nj',ny,dimid_nj)
236 :
237 9 : status = nf90_def_dim(ncid,'ncat',ncat,dimid_ncat)
238 :
239 : !-----------------------------------------------------------------
240 : ! 2D restart fields
241 : !-----------------------------------------------------------------
242 :
243 9 : allocate(dims(2))
244 :
245 9 : dims(1) = dimid_ni
246 9 : dims(2) = dimid_nj
247 :
248 9 : call define_rest_field(ncid,'uvel',dims)
249 9 : call define_rest_field(ncid,'vvel',dims)
250 :
251 9 : if (grid_ice == 'CD') then
252 0 : call define_rest_field(ncid,'uvelE',dims)
253 0 : call define_rest_field(ncid,'vvelE',dims)
254 0 : call define_rest_field(ncid,'uvelN',dims)
255 0 : call define_rest_field(ncid,'vvelN',dims)
256 : endif
257 :
258 9 : if (grid_ice == 'C') then
259 0 : call define_rest_field(ncid,'uvelE',dims)
260 0 : call define_rest_field(ncid,'vvelN',dims)
261 : endif
262 :
263 9 : if (restart_coszen) call define_rest_field(ncid,'coszen',dims)
264 :
265 9 : call define_rest_field(ncid,'scale_factor',dims)
266 9 : call define_rest_field(ncid,'swvdr',dims)
267 9 : call define_rest_field(ncid,'swvdf',dims)
268 9 : call define_rest_field(ncid,'swidr',dims)
269 9 : call define_rest_field(ncid,'swidf',dims)
270 :
271 9 : call define_rest_field(ncid,'strocnxT',dims)
272 9 : call define_rest_field(ncid,'strocnyT',dims)
273 :
274 9 : call define_rest_field(ncid,'stressp_1',dims)
275 9 : call define_rest_field(ncid,'stressp_2',dims)
276 9 : call define_rest_field(ncid,'stressp_3',dims)
277 9 : call define_rest_field(ncid,'stressp_4',dims)
278 :
279 9 : call define_rest_field(ncid,'stressm_1',dims)
280 9 : call define_rest_field(ncid,'stressm_2',dims)
281 9 : call define_rest_field(ncid,'stressm_3',dims)
282 9 : call define_rest_field(ncid,'stressm_4',dims)
283 :
284 9 : call define_rest_field(ncid,'stress12_1',dims)
285 9 : call define_rest_field(ncid,'stress12_2',dims)
286 9 : call define_rest_field(ncid,'stress12_3',dims)
287 9 : call define_rest_field(ncid,'stress12_4',dims)
288 :
289 9 : call define_rest_field(ncid,'iceumask',dims)
290 :
291 9 : if (grid_ice == 'CD' .or. grid_ice == 'C') then
292 0 : call define_rest_field(ncid,'stresspT' ,dims)
293 0 : call define_rest_field(ncid,'stressmT' ,dims)
294 0 : call define_rest_field(ncid,'stress12T',dims)
295 0 : call define_rest_field(ncid,'stresspU' ,dims)
296 0 : call define_rest_field(ncid,'stressmU' ,dims)
297 0 : call define_rest_field(ncid,'stress12U',dims)
298 0 : call define_rest_field(ncid,'icenmask',dims)
299 0 : call define_rest_field(ncid,'iceemask',dims)
300 : endif
301 :
302 :
303 9 : if (oceanmixed_ice) then
304 9 : call define_rest_field(ncid,'sst',dims)
305 9 : call define_rest_field(ncid,'frzmlt',dims)
306 : endif
307 :
308 9 : if (tr_FY) then
309 9 : call define_rest_field(ncid,'frz_onset',dims)
310 : endif
311 :
312 9 : if (kdyn == 2) then
313 0 : call define_rest_field(ncid,'a11_1',dims)
314 0 : call define_rest_field(ncid,'a11_2',dims)
315 0 : call define_rest_field(ncid,'a11_3',dims)
316 0 : call define_rest_field(ncid,'a11_4',dims)
317 0 : call define_rest_field(ncid,'a12_1',dims)
318 0 : call define_rest_field(ncid,'a12_2',dims)
319 0 : call define_rest_field(ncid,'a12_3',dims)
320 0 : call define_rest_field(ncid,'a12_4',dims)
321 : endif
322 :
323 9 : if (tr_pond_lvl) then
324 9 : call define_rest_field(ncid,'fsnow',dims)
325 : endif
326 :
327 9 : if (nbtrcr > 0) then
328 0 : if (tr_bgc_N) then
329 0 : do k=1,n_algae
330 0 : write(nchar,'(i3.3)') k
331 0 : call define_rest_field(ncid,'algalN'//trim(nchar),dims)
332 : enddo
333 : endif
334 0 : if (tr_bgc_C) then
335 0 : do k=1,n_doc
336 0 : write(nchar,'(i3.3)') k
337 0 : call define_rest_field(ncid,'doc'//trim(nchar),dims)
338 : enddo
339 0 : do k=1,n_dic
340 0 : write(nchar,'(i3.3)') k
341 0 : call define_rest_field(ncid,'dic'//trim(nchar),dims)
342 : enddo
343 : endif
344 0 : call define_rest_field(ncid,'nit' ,dims)
345 0 : if (tr_bgc_Am) &
346 0 : call define_rest_field(ncid,'amm' ,dims)
347 0 : if (tr_bgc_Sil) &
348 0 : call define_rest_field(ncid,'sil' ,dims)
349 0 : if (tr_bgc_hum) &
350 0 : call define_rest_field(ncid,'hum' ,dims)
351 0 : if (tr_bgc_DMS) then
352 0 : call define_rest_field(ncid,'dmsp' ,dims)
353 0 : call define_rest_field(ncid,'dms' ,dims)
354 : endif
355 0 : if (tr_bgc_DON) then
356 0 : do k=1,n_don
357 0 : write(nchar,'(i3.3)') k
358 0 : call define_rest_field(ncid,'don'//trim(nchar),dims)
359 : enddo
360 : endif
361 0 : if (tr_bgc_Fe ) then
362 0 : do k=1,n_fed
363 0 : write(nchar,'(i3.3)') k
364 0 : call define_rest_field(ncid,'fed'//trim(nchar),dims)
365 : enddo
366 0 : do k=1,n_fep
367 0 : write(nchar,'(i3.3)') k
368 0 : call define_rest_field(ncid,'fep'//trim(nchar),dims)
369 : enddo
370 : endif
371 0 : if (tr_zaero) then
372 0 : do k=1,n_zaero
373 0 : write(nchar,'(i3.3)') k
374 0 : call define_rest_field(ncid,'zaeros'//trim(nchar),dims)
375 : enddo
376 : endif
377 : endif !nbtrcr
378 :
379 9 : deallocate(dims)
380 :
381 : !-----------------------------------------------------------------
382 : ! 3D restart fields (ncat)
383 : !-----------------------------------------------------------------
384 :
385 9 : allocate(dims(3))
386 :
387 9 : dims(1) = dimid_ni
388 9 : dims(2) = dimid_nj
389 9 : dims(3) = dimid_ncat
390 :
391 9 : call define_rest_field(ncid,'aicen',dims)
392 9 : call define_rest_field(ncid,'vicen',dims)
393 9 : call define_rest_field(ncid,'vsnon',dims)
394 9 : call define_rest_field(ncid,'Tsfcn',dims)
395 :
396 9 : if (tr_iage) then
397 9 : call define_rest_field(ncid,'iage',dims)
398 : end if
399 :
400 9 : if (tr_FY) then
401 9 : call define_rest_field(ncid,'FY',dims)
402 : end if
403 :
404 9 : if (tr_lvl) then
405 9 : call define_rest_field(ncid,'alvl',dims)
406 9 : call define_rest_field(ncid,'vlvl',dims)
407 : end if
408 :
409 9 : if (tr_pond_topo) then
410 0 : call define_rest_field(ncid,'apnd',dims)
411 0 : call define_rest_field(ncid,'hpnd',dims)
412 0 : call define_rest_field(ncid,'ipnd',dims)
413 : end if
414 :
415 9 : if (tr_pond_lvl) then
416 9 : call define_rest_field(ncid,'apnd',dims)
417 9 : call define_rest_field(ncid,'hpnd',dims)
418 9 : call define_rest_field(ncid,'ipnd',dims)
419 9 : call define_rest_field(ncid,'dhs',dims)
420 9 : call define_rest_field(ncid,'ffrac',dims)
421 : end if
422 :
423 9 : if (tr_brine) then
424 0 : call define_rest_field(ncid,'fbrn',dims)
425 0 : call define_rest_field(ncid,'first_ice',dims)
426 : endif
427 :
428 9 : if (skl_bgc) then
429 0 : do k = 1, n_algae
430 0 : write(nchar,'(i3.3)') k
431 0 : call define_rest_field(ncid,'bgc_N'//trim(nchar) ,dims)
432 : enddo
433 0 : if (tr_bgc_C) then
434 : ! do k = 1, n_algae
435 : ! write(nchar,'(i3.3)') k
436 : ! call define_rest_field(ncid,'bgc_C'//trim(nchar) ,dims)
437 : ! enddo
438 0 : do k = 1, n_doc
439 0 : write(nchar,'(i3.3)') k
440 0 : call define_rest_field(ncid,'bgc_DOC'//trim(nchar) ,dims)
441 : enddo
442 0 : do k = 1, n_dic
443 0 : write(nchar,'(i3.3)') k
444 0 : call define_rest_field(ncid,'bgc_DIC'//trim(nchar) ,dims)
445 : enddo
446 : endif
447 0 : if (tr_bgc_chl) then
448 0 : do k = 1, n_algae
449 0 : write(nchar,'(i3.3)') k
450 0 : call define_rest_field(ncid,'bgc_chl'//trim(nchar) ,dims)
451 : enddo
452 : endif
453 0 : call define_rest_field(ncid,'bgc_Nit' ,dims)
454 0 : if (tr_bgc_Am) &
455 0 : call define_rest_field(ncid,'bgc_Am' ,dims)
456 0 : if (tr_bgc_Sil) &
457 0 : call define_rest_field(ncid,'bgc_Sil' ,dims)
458 0 : if (tr_bgc_hum) &
459 0 : call define_rest_field(ncid,'bgc_hum' ,dims)
460 0 : if (tr_bgc_DMS) then
461 0 : call define_rest_field(ncid,'bgc_DMSPp',dims)
462 0 : call define_rest_field(ncid,'bgc_DMSPd',dims)
463 0 : call define_rest_field(ncid,'bgc_DMS' ,dims)
464 : endif
465 0 : if (tr_bgc_PON) &
466 0 : call define_rest_field(ncid,'bgc_PON' ,dims)
467 0 : if (tr_bgc_DON) then
468 0 : do k = 1, n_don
469 0 : write(nchar,'(i3.3)') k
470 0 : call define_rest_field(ncid,'bgc_DON'//trim(nchar) ,dims)
471 : enddo
472 : endif
473 0 : if (tr_bgc_Fe ) then
474 0 : do k = 1, n_fed
475 0 : write(nchar,'(i3.3)') k
476 0 : call define_rest_field(ncid,'bgc_Fed'//trim(nchar) ,dims)
477 : enddo
478 0 : do k = 1, n_fep
479 0 : write(nchar,'(i3.3)') k
480 0 : call define_rest_field(ncid,'bgc_Fep'//trim(nchar) ,dims)
481 : enddo
482 : endif
483 : endif !skl_bgc
484 :
485 : !-----------------------------------------------------------------
486 : ! 4D restart fields, written as layers of 3D
487 : !-----------------------------------------------------------------
488 :
489 72 : do k=1,nilyr
490 63 : write(nchar,'(i3.3)') k
491 63 : call define_rest_field(ncid,'sice'//trim(nchar),dims)
492 72 : call define_rest_field(ncid,'qice'//trim(nchar),dims)
493 : enddo
494 :
495 18 : do k=1,nslyr
496 9 : write(nchar,'(i3.3)') k
497 18 : call define_rest_field(ncid,'qsno'//trim(nchar),dims)
498 : enddo
499 :
500 9 : if (tr_snow) then
501 0 : do k=1,nslyr
502 0 : write(nchar,'(i3.3)') k
503 0 : call define_rest_field(ncid,'smice'//trim(nchar),dims)
504 0 : call define_rest_field(ncid,'smliq'//trim(nchar),dims)
505 0 : call define_rest_field(ncid, 'rhos'//trim(nchar),dims)
506 0 : call define_rest_field(ncid, 'rsnw'//trim(nchar),dims)
507 : enddo
508 : endif
509 :
510 9 : if (tr_fsd) then
511 0 : do k=1,nfsd
512 0 : write(nchar,'(i3.3)') k
513 0 : call define_rest_field(ncid,'fsd'//trim(nchar),dims)
514 : enddo
515 : endif
516 :
517 9 : if (tr_iso) then
518 0 : do k=1,n_iso
519 0 : write(nchar,'(i3.3)') k
520 0 : call define_rest_field(ncid,'isosno'//trim(nchar),dims)
521 0 : call define_rest_field(ncid,'isoice'//trim(nchar),dims)
522 : enddo
523 : endif
524 :
525 9 : if (tr_aero) then
526 0 : do k=1,n_aero
527 0 : write(nchar,'(i3.3)') k
528 0 : call define_rest_field(ncid,'aerosnossl'//trim(nchar),dims)
529 0 : call define_rest_field(ncid,'aerosnoint'//trim(nchar),dims)
530 0 : call define_rest_field(ncid,'aeroicessl'//trim(nchar),dims)
531 0 : call define_rest_field(ncid,'aeroiceint'//trim(nchar),dims)
532 : enddo
533 : endif
534 :
535 9 : if (z_tracers) then
536 0 : if (tr_zaero) then
537 0 : do n = 1, n_zaero
538 0 : write(ncharb,'(i3.3)') n
539 0 : do k = 1, nblyr+3
540 0 : write(nchar,'(i3.3)') k
541 0 : call define_rest_field(ncid,'zaero'//trim(ncharb)//trim(nchar),dims)
542 : enddo !k
543 : enddo !n
544 : endif !tr_zaero
545 0 : if (tr_bgc_Nit) then
546 0 : do k = 1, nblyr+3
547 0 : write(nchar,'(i3.3)') k
548 0 : call define_rest_field(ncid,'bgc_Nit'//trim(nchar),dims)
549 : enddo
550 : endif
551 0 : if (tr_bgc_N) then
552 0 : do n = 1, n_algae
553 0 : write(ncharb,'(i3.3)') n
554 0 : do k = 1, nblyr+3
555 0 : write(nchar,'(i3.3)') k
556 0 : call define_rest_field(ncid,'bgc_N'//trim(ncharb)//trim(nchar),dims)
557 : enddo
558 : enddo
559 : endif
560 0 : if (tr_bgc_C) then
561 : ! do n = 1, n_algae
562 : ! write(ncharb,'(i3.3)') n
563 : ! do k = 1, nblyr+3
564 : ! write(nchar,'(i3.3)') k
565 : ! call define_rest_field(ncid,'bgc_C'//trim(ncharb)//trim(nchar),dims)
566 : ! enddo
567 : ! enddo
568 0 : do n = 1, n_doc
569 0 : write(ncharb,'(i3.3)') n
570 0 : do k = 1, nblyr+3
571 0 : write(nchar,'(i3.3)') k
572 0 : call define_rest_field(ncid,'bgc_DOC'//trim(ncharb)//trim(nchar),dims)
573 : enddo
574 : enddo
575 0 : do n = 1, n_dic
576 0 : write(ncharb,'(i3.3)') n
577 0 : do k = 1, nblyr+3
578 0 : write(nchar,'(i3.3)') k
579 0 : call define_rest_field(ncid,'bgc_DIC'//trim(ncharb)//trim(nchar),dims)
580 : enddo
581 : enddo
582 : endif
583 0 : if (tr_bgc_chl) then
584 0 : do n = 1, n_algae
585 0 : write(ncharb,'(i3.3)') n
586 0 : do k = 1, nblyr+3
587 0 : write(nchar,'(i3.3)') k
588 0 : call define_rest_field(ncid,'bgc_chl'//trim(ncharb)//trim(nchar),dims)
589 : enddo
590 : enddo
591 : endif
592 0 : if (tr_bgc_Am) then
593 0 : do k = 1, nblyr+3
594 0 : write(nchar,'(i3.3)') k
595 0 : call define_rest_field(ncid,'bgc_Am'//trim(nchar),dims)
596 : enddo
597 : endif
598 0 : if (tr_bgc_Sil) then
599 0 : do k = 1, nblyr+3
600 0 : write(nchar,'(i3.3)') k
601 0 : call define_rest_field(ncid,'bgc_Sil'//trim(nchar),dims)
602 : enddo
603 : endif
604 0 : if (tr_bgc_hum) then
605 0 : do k = 1, nblyr+3
606 0 : write(nchar,'(i3.3)') k
607 0 : call define_rest_field(ncid,'bgc_hum'//trim(nchar),dims)
608 : enddo
609 : endif
610 0 : if (tr_bgc_DMS) then
611 0 : do k = 1, nblyr+3
612 0 : write(nchar,'(i3.3)') k
613 0 : call define_rest_field(ncid,'bgc_DMSPp'//trim(nchar),dims)
614 0 : call define_rest_field(ncid,'bgc_DMSPd'//trim(nchar),dims)
615 0 : call define_rest_field(ncid,'bgc_DMS'//trim(nchar),dims)
616 : enddo
617 : endif
618 0 : if (tr_bgc_PON) then
619 0 : do k = 1, nblyr+3
620 0 : write(nchar,'(i3.3)') k
621 0 : call define_rest_field(ncid,'bgc_PON'//trim(nchar),dims)
622 : enddo
623 : endif
624 0 : if (tr_bgc_DON) then
625 0 : do n = 1, n_don
626 0 : write(ncharb,'(i3.3)') n
627 0 : do k = 1, nblyr+3
628 0 : write(nchar,'(i3.3)') k
629 0 : call define_rest_field(ncid,'bgc_DON'//trim(ncharb)//trim(nchar),dims)
630 : enddo
631 : enddo
632 : endif
633 0 : if (tr_bgc_Fe ) then
634 0 : do n = 1, n_fed
635 0 : write(ncharb,'(i3.3)') n
636 0 : do k = 1, nblyr+3
637 0 : write(nchar,'(i3.3)') k
638 0 : call define_rest_field(ncid,'bgc_Fed'//trim(ncharb)//trim(nchar),dims)
639 : enddo
640 : enddo
641 0 : do n = 1, n_fep
642 0 : write(ncharb,'(i3.3)') n
643 0 : do k = 1, nblyr+3
644 0 : write(nchar,'(i3.3)') k
645 0 : call define_rest_field(ncid,'bgc_Fep'//trim(ncharb)//trim(nchar),dims)
646 : enddo
647 : enddo
648 : endif
649 0 : do k = 1, nbtrcr
650 0 : write(nchar,'(i3.3)') k
651 0 : call define_rest_field(ncid,'zbgc_frac'//trim(nchar),dims)
652 : enddo
653 : endif !z_tracers
654 :
655 9 : deallocate(dims)
656 9 : status = nf90_enddef(ncid)
657 :
658 9 : write(nu_diag,*) 'Writing ',filename(1:lenstr(filename))
659 : endif ! master_task
660 :
661 : #else
662 : call abort_ice(subname//'ERROR: USE_NETCDF cpp not defined for '//trim(filename_spec), &
663 : file=__FILE__, line=__LINE__)
664 : #endif
665 :
666 49 : end subroutine init_restart_write
667 :
668 : !=======================================================================
669 :
670 : ! Reads a single restart field
671 : ! author David A Bailey, NCAR
672 :
673 991 : subroutine read_restart_field(nu,nrec,work,atype,vname,ndim3, &
674 : diag, field_loc, field_type)
675 :
676 12 : use ice_blocks, only: nx_block, ny_block
677 : use ice_domain_size, only: max_blocks, ncat
678 : use ice_read_write, only: ice_read_nc
679 :
680 : integer (kind=int_kind), intent(in) :: &
681 : nu , & ! unit number (not used for netcdf) ! LCOV_EXCL_LINE
682 : ndim3 , & ! third dimension ! LCOV_EXCL_LINE
683 : nrec ! record number (0 for sequential access)
684 :
685 : real (kind=dbl_kind), dimension(nx_block,ny_block,ndim3,max_blocks), intent(inout) :: &
686 : work ! input array (real, 8-byte)
687 :
688 : character (len=4), intent(in) :: &
689 : atype ! format for output array
690 : ! (real/integer, 4-byte/8-byte)
691 :
692 : logical (kind=log_kind), intent(in) :: &
693 : diag ! if true, write diagnostic output
694 :
695 : character (len=*), intent(in) :: vname
696 :
697 : integer (kind=int_kind), optional, intent(in) :: &
698 : field_loc, & ! location of field on staggered grid ! LCOV_EXCL_LINE
699 : field_type ! type of field (scalar, vector, angle)
700 :
701 : ! local variables
702 :
703 : real (kind=dbl_kind), dimension(nx_block,ny_block,max_blocks) :: &
704 1351058 : work2 ! input array (real, 8-byte)
705 :
706 : character(len=*), parameter :: subname = '(read_restart_field)'
707 :
708 : #ifdef USE_NETCDF
709 1379 : if (present(field_loc)) then
710 1379 : if (ndim3 == ncat) then
711 659 : if (restart_ext) then
712 : call ice_read_nc(ncid,1,vname,work,diag, &
713 0 : field_loc=field_loc,field_type=field_type,restart_ext=restart_ext)
714 : else
715 659 : call ice_read_nc(ncid,1,vname,work,diag,field_loc,field_type)
716 : endif
717 720 : elseif (ndim3 == 1) then
718 720 : if (restart_ext) then
719 : call ice_read_nc(ncid,1,vname,work2,diag, &
720 0 : field_loc=field_loc,field_type=field_type,restart_ext=restart_ext)
721 : else
722 720 : call ice_read_nc(ncid,1,vname,work2,diag,field_loc,field_type)
723 : endif
724 2294872 : work(:,:,1,:) = work2(:,:,:)
725 : else
726 0 : write(nu_diag,*) 'ndim3 not supported ',ndim3
727 : endif
728 : else
729 0 : if (ndim3 == ncat) then
730 0 : if (restart_ext) then
731 0 : call ice_read_nc(ncid, 1, vname, work, diag, restart_ext=restart_ext)
732 : else
733 0 : call ice_read_nc(ncid, 1, vname, work, diag)
734 : endif
735 0 : elseif (ndim3 == 1) then
736 0 : if (restart_ext) then
737 0 : call ice_read_nc(ncid, 1, vname, work2, diag, restart_ext=restart_ext)
738 : else
739 0 : call ice_read_nc(ncid, 1, vname, work2, diag)
740 : endif
741 0 : work(:,:,1,:) = work2(:,:,:)
742 : else
743 0 : write(nu_diag,*) 'ndim3 not supported ',ndim3
744 : endif
745 : endif
746 :
747 : #else
748 : call abort_ice(subname//'ERROR: USE_NETCDF cpp not defined', &
749 : file=__FILE__, line=__LINE__)
750 : #endif
751 :
752 1379 : end subroutine read_restart_field
753 :
754 : !=======================================================================
755 :
756 : ! Writes a single restart field.
757 : ! author David A Bailey, NCAR
758 :
759 1998 : subroutine write_restart_field(nu,nrec,work,atype,vname,ndim3,diag)
760 :
761 : use ice_blocks, only: nx_block, ny_block
762 : use ice_domain_size, only: max_blocks, ncat
763 : use ice_read_write, only: ice_write_nc
764 :
765 : integer (kind=int_kind), intent(in) :: &
766 : nu , & ! unit number ! LCOV_EXCL_LINE
767 : ndim3 , & ! third dimension ! LCOV_EXCL_LINE
768 : nrec ! record number (0 for sequential access)
769 :
770 : real (kind=dbl_kind), dimension(nx_block,ny_block,ndim3,max_blocks), intent(in) :: &
771 : work ! input array (real, 8-byte)
772 :
773 : character (len=4), intent(in) :: &
774 : atype ! format for output array
775 : ! (real/integer, 4-byte/8-byte)
776 :
777 : logical (kind=log_kind), intent(in) :: &
778 : diag ! if true, write diagnostic output
779 :
780 : character (len=*), intent(in) :: vname
781 :
782 : ! local variables
783 :
784 : integer (kind=int_kind) :: &
785 : varid, & ! variable id ! LCOV_EXCL_LINE
786 : status ! status variable from netCDF routine
787 :
788 : real (kind=dbl_kind), dimension(nx_block,ny_block,max_blocks) :: &
789 2257092 : work2 ! input array (real, 8-byte)
790 :
791 : character(len=*), parameter :: subname = '(write_restart_field)'
792 :
793 : #ifdef USE_NETCDF
794 2646 : status = nf90_inq_varid(ncid,trim(vname),varid)
795 2646 : if (ndim3 == ncat) then
796 1372 : if (restart_ext) then
797 0 : call ice_write_nc(ncid, 1, varid, work, diag, restart_ext, varname=trim(vname))
798 : else
799 1372 : call ice_write_nc(ncid, 1, varid, work, diag, varname=trim(vname))
800 : endif
801 1274 : elseif (ndim3 == 1) then
802 3767608 : work2(:,:,:) = work(:,:,1,:)
803 1274 : if (restart_ext) then
804 0 : call ice_write_nc(ncid, 1, varid, work2, diag, restart_ext, varname=trim(vname))
805 : else
806 1274 : call ice_write_nc(ncid, 1, varid, work2, diag, varname=trim(vname))
807 : endif
808 : else
809 0 : write(nu_diag,*) 'ndim3 not supported',ndim3
810 : endif
811 :
812 : #else
813 : call abort_ice(subname//'ERROR: USE_NETCDF cpp not defined', &
814 : file=__FILE__, line=__LINE__)
815 : #endif
816 :
817 2646 : end subroutine write_restart_field
818 :
819 : !=======================================================================
820 :
821 : ! Finalize the restart file.
822 : ! author David A Bailey, NCAR
823 :
824 49 : subroutine final_restart()
825 :
826 : use ice_calendar, only: istep1, myear, mmonth, mday, msec
827 :
828 : integer (kind=int_kind) :: status
829 :
830 : character(len=*), parameter :: subname = '(final_restart)'
831 :
832 : #ifdef USE_NETCDF
833 49 : status = nf90_close(ncid)
834 :
835 49 : if (my_task == master_task) &
836 9 : write(nu_diag,'(a,i8,4x,i4.4,a,i2.2,a,i2.2,a,i5.5)') 'Restart read/written ',istep1,myear,'-',mmonth,'-',mday,'-',msec
837 :
838 : #else
839 : call abort_ice(subname//'ERROR: USE_NETCDF cpp not defined', &
840 : file=__FILE__, line=__LINE__)
841 : #endif
842 :
843 49 : end subroutine final_restart
844 :
845 : !=======================================================================
846 :
847 : ! Defines a restart field
848 : ! author David A Bailey, NCAR
849 :
850 486 : subroutine define_rest_field(ncid, vname, dims)
851 :
852 : character (len=*) , intent(in) :: vname
853 : integer (kind=int_kind), intent(in) :: dims(:)
854 : integer (kind=int_kind), intent(in) :: ncid
855 :
856 : integer (kind=int_kind) :: varid
857 :
858 : integer (kind=int_kind) :: &
859 : status ! status variable from netCDF routine
860 :
861 : character(len=*), parameter :: subname = '(define_rest_field)'
862 :
863 : #ifdef USE_NETCDF
864 486 : status = nf90_def_var(ncid,trim(vname),nf90_double,dims,varid)
865 : #else
866 : call abort_ice(subname//'ERROR: USE_NETCDF cpp not defined', &
867 : file=__FILE__, line=__LINE__)
868 : #endif
869 :
870 486 : end subroutine define_rest_field
871 :
872 : !=======================================================================
873 :
874 : ! Inquire field existance
875 : ! author T. Craig
876 :
877 29 : logical function query_field(nu,vname)
878 :
879 : integer (kind=int_kind), intent(in) :: nu ! unit number
880 : character (len=*) , intent(in) :: vname ! variable name
881 :
882 : ! local variables
883 :
884 : integer (kind=int_kind) :: status, varid
885 : character(len=*), parameter :: subname = '(query_field)'
886 :
887 29 : query_field = .false.
888 : #ifdef USE_NETCDF
889 29 : if (my_task == master_task) then
890 6 : status = nf90_inq_varid(ncid,trim(vname),varid)
891 6 : if (status == nf90_noerr) query_field = .true.
892 : endif
893 29 : call broadcast_scalar(query_field,master_task)
894 : #else
895 : call abort_ice(subname//'ERROR: USE_NETCDF cpp not defined for '//trim(ice_ic), &
896 : file=__FILE__, line=__LINE__)
897 : #endif
898 :
899 58 : end function query_field
900 :
901 : !=======================================================================
902 :
903 : end module ice_restart
904 :
905 : !=======================================================================
|