Line data Source code
1 : !=========================================================================
2 : !
3 : ! Restart routines for the column package.
4 : !
5 : ! author: Elizabeth C. Hunke, LANL
6 : !
7 : ! 2014: Moved subroutines from column package modules
8 :
9 : module ice_restart_column
10 :
11 : use ice_kinds_mod
12 : use ice_communicate, only: my_task, master_task
13 : use ice_constants, only: c0, c1, p5
14 : use ice_constants, only: field_loc_center, field_type_scalar
15 : use ice_domain_size, only: ncat, nslyr, nfsd, nblyr
16 : use ice_restart,only: read_restart_field, write_restart_field
17 : use ice_exit, only: abort_ice
18 : use ice_fileunits, only: nu_diag
19 : use icepack_intfc, only: icepack_warnings_flush, icepack_warnings_aborted
20 : use icepack_intfc, only: icepack_max_algae, icepack_max_doc, &
21 : icepack_max_don, icepack_max_dic, icepack_max_fe, icepack_max_aero
22 : use icepack_intfc, only: icepack_query_parameters, &
23 : icepack_query_tracer_sizes, icepack_query_tracer_flags, & ! LCOV_EXCL_LINE
24 : icepack_query_tracer_indices
25 :
26 : implicit none
27 :
28 : private
29 : public :: write_restart_age, read_restart_age, &
30 : write_restart_FY, read_restart_FY, & ! LCOV_EXCL_LINE
31 : write_restart_lvl, read_restart_lvl, & ! LCOV_EXCL_LINE
32 : write_restart_pond_lvl, read_restart_pond_lvl, & ! LCOV_EXCL_LINE
33 : write_restart_pond_topo, read_restart_pond_topo, & ! LCOV_EXCL_LINE
34 : write_restart_snow, read_restart_snow, & ! LCOV_EXCL_LINE
35 : write_restart_fsd, read_restart_fsd, & ! LCOV_EXCL_LINE
36 : write_restart_iso, read_restart_iso, & ! LCOV_EXCL_LINE
37 : write_restart_aero, read_restart_aero, & ! LCOV_EXCL_LINE
38 : write_restart_bgc, read_restart_bgc, & ! LCOV_EXCL_LINE
39 : write_restart_hbrine, read_restart_hbrine
40 :
41 : logical (kind=log_kind), public :: &
42 : restart_age , & ! if .true., read age tracer restart file ! LCOV_EXCL_LINE
43 : restart_FY , & ! if .true., read FY tracer restart file ! LCOV_EXCL_LINE
44 : restart_lvl , & ! if .true., read lvl tracer restart file ! LCOV_EXCL_LINE
45 : restart_pond_lvl , & ! if .true., read meltponds restart file ! LCOV_EXCL_LINE
46 : restart_pond_topo, & ! if .true., read meltponds restart file ! LCOV_EXCL_LINE
47 : restart_snow , & ! if .true., read snow tracer restart file ! LCOV_EXCL_LINE
48 : restart_fsd , & ! if .true., read floe size restart file ! LCOV_EXCL_LINE
49 : restart_iso , & ! if .true., read isotope tracer restart file ! LCOV_EXCL_LINE
50 : restart_aero , & ! if .true., read aerosol tracer restart file ! LCOV_EXCL_LINE
51 : restart_hbrine , & ! if .true., read hbrine from restart file ! LCOV_EXCL_LINE
52 : restart_bgc ! if .true., read bgc restart file
53 :
54 : !=======================================================================
55 :
56 : contains
57 :
58 : !=======================================================================
59 :
60 : ! Dumps all values needed for restarting
61 : ! author Elizabeth C. Hunke, LANL
62 :
63 49 : subroutine write_restart_age()
64 :
65 : use ice_fileunits, only: nu_dump_age
66 : use ice_state, only: trcrn
67 :
68 : ! local variables
69 :
70 : logical (kind=log_kind) :: diag
71 : integer (kind=int_kind) :: nt_iage
72 : character(len=*),parameter :: subname='(write_restart_age)'
73 :
74 49 : call icepack_query_tracer_indices(nt_iage_out=nt_iage)
75 49 : call icepack_warnings_flush(nu_diag)
76 49 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
77 0 : file=__FILE__, line=__LINE__)
78 :
79 49 : diag = .true.
80 :
81 : !-----------------------------------------------------------------
82 :
83 0 : call write_restart_field(nu_dump_age,0,trcrn(:,:,nt_iage,:,:),'ruf8', &
84 49 : 'iage',ncat,diag)
85 :
86 49 : end subroutine write_restart_age
87 :
88 : !=======================================================================
89 :
90 : ! Reads all values needed for an ice age restart
91 : ! author Elizabeth C. Hunke, LANL
92 :
93 12 : subroutine read_restart_age()
94 :
95 : use ice_fileunits, only: nu_restart_age
96 : use ice_state, only: trcrn
97 :
98 : ! local variables
99 :
100 : logical (kind=log_kind) :: &
101 : diag
102 : integer (kind=int_kind) :: nt_iage
103 : character(len=*),parameter :: subname='(read_restart_age)'
104 :
105 12 : call icepack_query_tracer_indices(nt_iage_out=nt_iage)
106 12 : call icepack_warnings_flush(nu_diag)
107 12 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
108 0 : file=__FILE__, line=__LINE__)
109 :
110 12 : diag = .true.
111 :
112 12 : if (my_task == master_task) write(nu_diag,*) subname,'min/max age (s)'
113 :
114 0 : call read_restart_field(nu_restart_age,0,trcrn(:,:,nt_iage,:,:),'ruf8', &
115 12 : 'iage',ncat,diag,field_loc_center,field_type_scalar)
116 :
117 12 : end subroutine read_restart_age
118 :
119 : !=======================================================================
120 :
121 : ! Dumps all values needed for restarting
122 : ! author Elizabeth C. Hunke, LANL
123 :
124 49 : subroutine write_restart_FY()
125 :
126 : use ice_fileunits, only: nu_dump_FY
127 : use ice_flux, only: frz_onset
128 : use ice_state, only: trcrn
129 :
130 : ! local variables
131 :
132 : logical (kind=log_kind) :: diag
133 : integer (kind=int_kind) :: nt_FY
134 : character(len=*),parameter :: subname='(write_restart_FY)'
135 :
136 49 : call icepack_query_tracer_indices(nt_FY_out=nt_FY)
137 49 : call icepack_warnings_flush(nu_diag)
138 49 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
139 0 : file=__FILE__, line=__LINE__)
140 :
141 49 : diag = .true.
142 :
143 : !-----------------------------------------------------------------
144 :
145 0 : call write_restart_field(nu_dump_FY,0,trcrn(:,:,nt_FY,:,:),'ruf8', &
146 49 : 'FY',ncat,diag)
147 : call write_restart_field(nu_dump_FY,0,frz_onset,'ruf8', &
148 49 : 'frz_onset',1,diag)
149 :
150 49 : end subroutine write_restart_FY
151 :
152 : !=======================================================================
153 :
154 : ! Reads all values needed for an ice FY restart
155 : ! author Elizabeth C. Hunke, LANL
156 :
157 12 : subroutine read_restart_FY()
158 :
159 : use ice_fileunits, only: nu_restart_FY
160 : use ice_flux, only: frz_onset
161 : use ice_state, only: trcrn
162 :
163 : ! local variables
164 :
165 : logical (kind=log_kind) :: &
166 : diag
167 : integer (kind=int_kind) :: nt_FY
168 : character(len=*),parameter :: subname='(read_restart_FY)'
169 :
170 12 : call icepack_query_tracer_indices(nt_FY_out=nt_FY)
171 12 : call icepack_warnings_flush(nu_diag)
172 12 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
173 0 : file=__FILE__, line=__LINE__)
174 :
175 12 : diag = .true.
176 :
177 12 : if (my_task == master_task) write(nu_diag,*) subname,'min/max first-year ice area'
178 :
179 0 : call read_restart_field(nu_restart_FY,0,trcrn(:,:,nt_FY,:,:),'ruf8', &
180 12 : 'FY',ncat,diag,field_loc_center,field_type_scalar)
181 :
182 12 : if (my_task == master_task) write(nu_diag,*) subname,'min/max frz_onset'
183 :
184 : call read_restart_field(nu_restart_FY,0,frz_onset,'ruf8', &
185 12 : 'frz_onset',1,diag,field_loc_center,field_type_scalar)
186 :
187 12 : end subroutine read_restart_FY
188 :
189 : !=======================================================================
190 :
191 : ! Dumps all values needed for restarting
192 : !
193 : ! author Elizabeth C. Hunke, LANL
194 :
195 49 : subroutine write_restart_lvl()
196 :
197 : use ice_fileunits, only: nu_dump_lvl
198 : use ice_state, only: trcrn
199 :
200 : ! local variables
201 :
202 : logical (kind=log_kind) :: diag
203 : integer (kind=int_kind) :: nt_alvl, nt_vlvl
204 : character(len=*),parameter :: subname='(write_restart_lvl)'
205 :
206 49 : call icepack_query_tracer_indices(nt_alvl_out=nt_alvl, nt_vlvl_out=nt_vlvl)
207 49 : call icepack_warnings_flush(nu_diag)
208 49 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
209 0 : file=__FILE__, line=__LINE__)
210 :
211 49 : diag = .true.
212 :
213 : !-----------------------------------------------------------------
214 :
215 0 : call write_restart_field(nu_dump_lvl,0,trcrn(:,:,nt_alvl,:,:),'ruf8', &
216 49 : 'alvl',ncat,diag)
217 0 : call write_restart_field(nu_dump_lvl,0,trcrn(:,:,nt_vlvl,:,:),'ruf8', &
218 49 : 'vlvl',ncat,diag)
219 :
220 49 : end subroutine write_restart_lvl
221 :
222 : !=======================================================================
223 :
224 : ! Reads all values needed for an ice lvl restart
225 : !
226 : ! author Elizabeth C. Hunke, LANL
227 :
228 12 : subroutine read_restart_lvl()
229 :
230 : use ice_fileunits, only: nu_restart_lvl
231 : use ice_state, only: trcrn
232 :
233 : ! local variables
234 :
235 : logical (kind=log_kind) :: &
236 : diag
237 : integer (kind=int_kind) :: nt_alvl, nt_vlvl
238 : character(len=*),parameter :: subname='(read_restart_lvl)'
239 :
240 12 : call icepack_query_tracer_indices(nt_alvl_out=nt_alvl, nt_vlvl_out=nt_vlvl)
241 12 : call icepack_warnings_flush(nu_diag)
242 12 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
243 0 : file=__FILE__, line=__LINE__)
244 :
245 12 : diag = .true.
246 :
247 12 : if (my_task == master_task) write(nu_diag,*) subname,'min/max level ice area, volume'
248 :
249 0 : call read_restart_field(nu_restart_lvl,0,trcrn(:,:,nt_alvl,:,:),'ruf8', &
250 12 : 'alvl',ncat,diag,field_loc_center,field_type_scalar)
251 0 : call read_restart_field(nu_restart_lvl,0,trcrn(:,:,nt_vlvl,:,:),'ruf8', &
252 12 : 'vlvl',ncat,diag,field_loc_center,field_type_scalar)
253 :
254 12 : end subroutine read_restart_lvl
255 :
256 : !=======================================================================
257 : !
258 : ! Dumps all values needed for restarting
259 : !
260 : ! authors Elizabeth C. Hunke, LANL
261 :
262 49 : subroutine write_restart_pond_lvl()
263 :
264 : use ice_arrays_column, only: dhsn, ffracn
265 : use ice_fileunits, only: nu_dump_pond
266 : use ice_flux, only: fsnow
267 : use ice_state, only: trcrn
268 :
269 : ! local variables
270 :
271 : logical (kind=log_kind) :: diag
272 : integer (kind=int_kind) :: nt_apnd, nt_hpnd, nt_ipnd
273 : character(len=*),parameter :: subname='(write_restart_pond_lvl)'
274 :
275 : call icepack_query_tracer_indices(nt_apnd_out=nt_apnd, nt_hpnd_out=nt_hpnd, &
276 49 : nt_ipnd_out=nt_ipnd)
277 49 : call icepack_warnings_flush(nu_diag)
278 49 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
279 0 : file=__FILE__, line=__LINE__)
280 :
281 49 : diag = .true.
282 :
283 0 : call write_restart_field(nu_dump_pond,0, trcrn(:,:,nt_apnd,:,:),'ruf8', &
284 49 : 'apnd',ncat,diag)
285 0 : call write_restart_field(nu_dump_pond,0, trcrn(:,:,nt_hpnd,:,:),'ruf8', &
286 49 : 'hpnd',ncat,diag)
287 0 : call write_restart_field(nu_dump_pond,0, trcrn(:,:,nt_ipnd,:,:),'ruf8', &
288 49 : 'ipnd',ncat,diag)
289 : call write_restart_field(nu_dump_pond,0, fsnow(:,:, :),'ruf8', &
290 49 : 'fsnow',1,diag)
291 : call write_restart_field(nu_dump_pond,0, dhsn(:,:, :,:),'ruf8', &
292 49 : 'dhs',ncat,diag)
293 : call write_restart_field(nu_dump_pond,0,ffracn(:,:, :,:),'ruf8', &
294 49 : 'ffrac',ncat,diag)
295 :
296 49 : end subroutine write_restart_pond_lvl
297 :
298 : !=======================================================================
299 :
300 : ! Reads all values needed for a meltpond volume restart
301 : !
302 : ! authors Elizabeth C. Hunke, LANL
303 :
304 12 : subroutine read_restart_pond_lvl()
305 :
306 : use ice_arrays_column, only: dhsn, ffracn
307 : use ice_fileunits, only: nu_restart_pond
308 : use ice_flux, only: fsnow
309 : use ice_state, only: trcrn
310 :
311 : ! local variables
312 :
313 : logical (kind=log_kind) :: &
314 : diag
315 : integer (kind=int_kind) :: nt_apnd, nt_hpnd, nt_ipnd
316 : character(len=*),parameter :: subname='(read_restart_pond_lvl)'
317 :
318 : call icepack_query_tracer_indices(nt_apnd_out=nt_apnd, nt_hpnd_out=nt_hpnd, &
319 12 : nt_ipnd_out=nt_ipnd)
320 12 : call icepack_warnings_flush(nu_diag)
321 12 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
322 0 : file=__FILE__, line=__LINE__)
323 :
324 12 : diag = .true.
325 :
326 12 : if (my_task == master_task) write(nu_diag,*) subname,'min/max level-ice ponds'
327 :
328 0 : call read_restart_field(nu_restart_pond,0, trcrn(:,:,nt_apnd,:,:),'ruf8', &
329 12 : 'apnd',ncat,diag,field_loc_center,field_type_scalar)
330 0 : call read_restart_field(nu_restart_pond,0, trcrn(:,:,nt_hpnd,:,:),'ruf8', &
331 12 : 'hpnd',ncat,diag,field_loc_center,field_type_scalar)
332 0 : call read_restart_field(nu_restart_pond,0, trcrn(:,:,nt_ipnd,:,:),'ruf8', &
333 12 : 'ipnd',ncat,diag,field_loc_center,field_type_scalar)
334 : call read_restart_field(nu_restart_pond,0, fsnow(:,:, :),'ruf8', &
335 12 : 'fsnow',1,diag,field_loc_center,field_type_scalar)
336 : call read_restart_field(nu_restart_pond,0, dhsn(:,:, :,:),'ruf8', &
337 12 : 'dhs',ncat,diag,field_loc_center,field_type_scalar)
338 : call read_restart_field(nu_restart_pond,0,ffracn(:,:, :,:),'ruf8', &
339 12 : 'ffrac',ncat,diag,field_loc_center,field_type_scalar)
340 :
341 12 : end subroutine read_restart_pond_lvl
342 :
343 : !=======================================================================
344 :
345 : ! Dumps all values needed for restarting
346 : !
347 : ! authors Elizabeth C. Hunke, LANL
348 : ! David A. Bailey, NCAR
349 :
350 0 : subroutine write_restart_pond_topo()
351 :
352 : use ice_fileunits, only: nu_dump_pond
353 : use ice_state, only: trcrn
354 :
355 : ! local variables
356 :
357 : logical (kind=log_kind) :: diag
358 : integer (kind=int_kind) :: nt_apnd, nt_hpnd, nt_ipnd
359 : character(len=*),parameter :: subname='(write_restart_pond_topo)'
360 :
361 : call icepack_query_tracer_indices(nt_apnd_out=nt_apnd, nt_hpnd_out=nt_hpnd, &
362 0 : nt_ipnd_out=nt_ipnd)
363 0 : call icepack_warnings_flush(nu_diag)
364 0 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
365 0 : file=__FILE__, line=__LINE__)
366 :
367 0 : diag = .true.
368 :
369 0 : call write_restart_field(nu_dump_pond,0,trcrn(:,:,nt_apnd,:,:),'ruf8', &
370 0 : 'apnd',ncat,diag)
371 0 : call write_restart_field(nu_dump_pond,0,trcrn(:,:,nt_hpnd,:,:),'ruf8', &
372 0 : 'hpnd',ncat,diag)
373 0 : call write_restart_field(nu_dump_pond,0,trcrn(:,:,nt_ipnd,:,:),'ruf8', &
374 0 : 'ipnd',ncat,diag)
375 :
376 0 : end subroutine write_restart_pond_topo
377 :
378 : !=======================================================================
379 :
380 : ! Reads all values needed for a meltpond volume restart
381 : !
382 : ! authors Elizabeth C. Hunke, LANL
383 : ! David A. Bailey, NCAR
384 :
385 0 : subroutine read_restart_pond_topo()
386 :
387 : use ice_fileunits, only: nu_restart_pond
388 : use ice_state, only: trcrn
389 :
390 : ! local variables
391 :
392 : logical (kind=log_kind) :: &
393 : diag
394 : integer (kind=int_kind) :: nt_apnd, nt_hpnd, nt_ipnd
395 : character(len=*),parameter :: subname='(read_restart_pond_topo)'
396 :
397 : call icepack_query_tracer_indices(nt_apnd_out=nt_apnd, nt_hpnd_out=nt_hpnd, &
398 0 : nt_ipnd_out=nt_ipnd)
399 0 : call icepack_warnings_flush(nu_diag)
400 0 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
401 0 : file=__FILE__, line=__LINE__)
402 :
403 0 : diag = .true.
404 :
405 0 : if (my_task == master_task) write(nu_diag,*) subname,'min/max topo ponds'
406 :
407 0 : call read_restart_field(nu_restart_pond,0,trcrn(:,:,nt_apnd,:,:),'ruf8', &
408 0 : 'apnd',ncat,diag,field_loc_center,field_type_scalar)
409 0 : call read_restart_field(nu_restart_pond,0,trcrn(:,:,nt_hpnd,:,:),'ruf8', &
410 0 : 'hpnd',ncat,diag,field_loc_center,field_type_scalar)
411 0 : call read_restart_field(nu_restart_pond,0,trcrn(:,:,nt_ipnd,:,:),'ruf8', &
412 0 : 'ipnd',ncat,diag,field_loc_center,field_type_scalar)
413 :
414 0 : end subroutine read_restart_pond_topo
415 :
416 : !=======================================================================
417 :
418 : ! Dumps all values needed for restarting snow redistribution/metamorphism
419 : ! author Elizabeth C. Hunke, LANL
420 :
421 0 : subroutine write_restart_snow()
422 :
423 : use ice_fileunits, only: nu_dump_snow
424 : use ice_state, only: trcrn
425 :
426 : ! local variables
427 :
428 : logical (kind=log_kind) :: diag
429 : integer (kind=int_kind) :: nt_smice, nt_smliq, nt_rhos, nt_rsnw, k
430 : character(len=3) :: ck
431 : character(len=*),parameter :: subname='(write_restart_snow)'
432 :
433 : call icepack_query_tracer_indices(nt_smice_out=nt_smice, &
434 0 : nt_smliq_out=nt_smliq, nt_rhos_out=nt_rhos, nt_rsnw_out=nt_rsnw)
435 0 : call icepack_warnings_flush(nu_diag)
436 0 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
437 0 : file=__FILE__, line=__LINE__)
438 :
439 0 : diag = .true.
440 :
441 : !-----------------------------------------------------------------
442 :
443 0 : do k = 1,nslyr
444 0 : write(ck,'(i3.3)') k
445 0 : call write_restart_field(nu_dump_snow,0, trcrn(:,:,nt_smice+k-1,:,:), &
446 0 : 'ruf8','smice'//trim(ck),ncat,diag)
447 0 : call write_restart_field(nu_dump_snow,0, trcrn(:,:,nt_smliq+k-1,:,:), &
448 0 : 'ruf8','smliq'//trim(ck),ncat,diag)
449 0 : call write_restart_field(nu_dump_snow,0, trcrn(:,:,nt_rhos+k-1,:,:), &
450 0 : 'ruf8','rhos'//trim(ck),ncat,diag)
451 0 : call write_restart_field(nu_dump_snow,0, trcrn(:,:,nt_rsnw+k-1,:,:), &
452 0 : 'ruf8','rsnw'//trim(ck),ncat,diag)
453 : enddo
454 :
455 0 : end subroutine write_restart_snow
456 :
457 : !=======================================================================
458 :
459 : ! Reads all values needed for a restart with snow redistribution/metamorphism
460 : ! author Elizabeth C. Hunke, LANL
461 :
462 0 : subroutine read_restart_snow()
463 :
464 : use ice_fileunits, only: nu_restart_snow
465 : use ice_state, only: trcrn
466 :
467 : ! local variables
468 :
469 : logical (kind=log_kind) :: &
470 : diag
471 : integer (kind=int_kind) :: nt_smice, nt_smliq, nt_rhos, nt_rsnw, k
472 : character(len=3) :: ck
473 : character(len=*),parameter :: subname='(read_restart_snow)'
474 :
475 : call icepack_query_tracer_indices(nt_smice_out=nt_smice, &
476 0 : nt_smliq_out=nt_smliq, nt_rhos_out=nt_rhos, nt_rsnw_out=nt_rsnw)
477 0 : call icepack_warnings_flush(nu_diag)
478 0 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
479 0 : file=__FILE__, line=__LINE__)
480 :
481 0 : diag = .true.
482 :
483 0 : if (my_task == master_task) write(nu_diag,*) subname,'min/max snow tracers'
484 :
485 0 : do k=1,nslyr
486 0 : write(ck,'(i3.3)') k
487 0 : call read_restart_field(nu_restart_snow,0,trcrn(:,:,nt_smice+k-1,:,:), &
488 : 'ruf8','smice'//trim(ck),ncat,diag, & ! LCOV_EXCL_LINE
489 0 : field_type=field_type_scalar,field_loc=field_loc_center)
490 0 : call read_restart_field(nu_restart_snow,0,trcrn(:,:,nt_smliq+k-1,:,:), &
491 : 'ruf8','smliq'//trim(ck),ncat,diag, & ! LCOV_EXCL_LINE
492 0 : field_type=field_type_scalar,field_loc=field_loc_center)
493 0 : call read_restart_field(nu_restart_snow,0,trcrn(:,:,nt_rhos+k-1,:,:), &
494 : 'ruf8','rhos'//trim(ck),ncat,diag, & ! LCOV_EXCL_LINE
495 0 : field_type=field_type_scalar,field_loc=field_loc_center)
496 0 : call read_restart_field(nu_restart_snow,0,trcrn(:,:,nt_rsnw+k-1,:,:), &
497 : 'ruf8','rsnw'//trim(ck),ncat,diag, & ! LCOV_EXCL_LINE
498 0 : field_type=field_type_scalar,field_loc=field_loc_center)
499 : enddo
500 :
501 0 : end subroutine read_restart_snow
502 :
503 : !=======================================================================
504 :
505 : ! Dumps all values needed for restarting
506 : ! author Elizabeth C. Hunke, LANL
507 :
508 0 : subroutine write_restart_fsd()
509 :
510 : use ice_fileunits, only: nu_dump_fsd
511 : use ice_state, only: trcrn
512 :
513 : ! local variables
514 :
515 : logical (kind=log_kind) :: diag
516 : integer (kind=int_kind) :: nt_fsd, k
517 : character(len=3) :: ck
518 : character(len=*),parameter :: subname='(write_restart_fsd)'
519 :
520 0 : call icepack_query_tracer_indices(nt_fsd_out=nt_fsd)
521 0 : call icepack_warnings_flush(nu_diag)
522 0 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
523 0 : file=__FILE__, line=__LINE__)
524 :
525 0 : diag = .true.
526 :
527 : !-----------------------------------------------------------------
528 :
529 0 : do k=1,nfsd
530 0 : write(ck,'(i3.3)') k
531 0 : call write_restart_field(nu_dump_fsd,0, trcrn(:,:,nt_fsd+k-1,:,:), &
532 0 : 'ruf8','fsd'//trim(ck),ncat,diag)
533 : enddo
534 :
535 0 : end subroutine write_restart_fsd
536 :
537 : !=======================================================================
538 :
539 : ! Reads all values needed for an ice fsd restart
540 : ! author Elizabeth C. Hunke, LANL
541 :
542 0 : subroutine read_restart_fsd()
543 :
544 : use ice_fileunits, only: nu_restart_fsd
545 : use ice_state, only: trcrn
546 :
547 : ! local variables
548 :
549 : logical (kind=log_kind) :: &
550 : diag
551 : integer (kind=int_kind) :: nt_fsd, k
552 : character(len=3) :: ck
553 : character(len=*),parameter :: subname='(read_restart_fsd)'
554 :
555 0 : call icepack_query_tracer_indices(nt_fsd_out=nt_fsd)
556 0 : call icepack_warnings_flush(nu_diag)
557 0 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
558 0 : file=__FILE__, line=__LINE__)
559 :
560 0 : diag = .true.
561 :
562 0 : if (my_task == master_task) write(nu_diag,*) subname,'min/max fsd (s)'
563 :
564 0 : do k=1,nfsd
565 0 : write(ck,'(i3.3)') k
566 0 : call read_restart_field(nu_restart_fsd,0,trcrn(:,:,nt_fsd+k-1,:,:), &
567 : 'ruf8','fsd'//trim(ck),ncat,diag, & ! LCOV_EXCL_LINE
568 0 : field_type=field_type_scalar,field_loc=field_loc_center)
569 : enddo
570 :
571 0 : end subroutine read_restart_fsd
572 :
573 : !=======================================================================
574 :
575 : ! Dumps all values needed for restarting
576 : ! author Elizabeth C. Hunke, LANL
577 :
578 0 : subroutine write_restart_iso()
579 :
580 : use ice_domain_size, only: n_iso
581 : use ice_fileunits, only: nu_dump_iso
582 : use ice_state, only: trcrn
583 :
584 : ! local variables
585 :
586 : logical (kind=log_kind) :: diag
587 : integer (kind=int_kind) :: nt_isosno, nt_isoice, k
588 : character(len=3) :: ck
589 : character(len=*),parameter :: subname='(write_restart_iso)'
590 :
591 0 : call icepack_query_tracer_indices(nt_isosno_out=nt_isosno, nt_isoice_out=nt_isoice)
592 0 : call icepack_warnings_flush(nu_diag)
593 0 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
594 0 : file=__FILE__, line=__LINE__)
595 :
596 0 : diag = .true.
597 :
598 : !-----------------------------------------------------------------
599 :
600 0 : do k = 1, n_iso
601 0 : write(ck,'(i3.3)') k
602 0 : call write_restart_field(nu_dump_iso,0, trcrn(:,:,nt_isosno+k-1,:,:), &
603 0 : 'ruf8','isosno'//trim(ck),ncat,diag)
604 : enddo
605 :
606 0 : do k = 1, n_iso
607 0 : write(ck,'(i3.3)') k
608 0 : call write_restart_field(nu_dump_iso,0, trcrn(:,:,nt_isoice+k-1,:,:), &
609 0 : 'ruf8','isoice'//trim(ck),ncat,diag)
610 : enddo
611 :
612 0 : end subroutine write_restart_iso
613 :
614 : !=======================================================================
615 :
616 : ! Reads all values needed to restart isotope tracers
617 : ! author Elizabeth C. Hunke, LANL
618 :
619 0 : subroutine read_restart_iso()
620 :
621 : use ice_domain_size, only: n_iso
622 : use ice_fileunits, only: nu_restart_iso
623 : use ice_state, only: trcrn
624 :
625 : ! local variables
626 :
627 : logical (kind=log_kind) :: &
628 : diag
629 : integer (kind=int_kind) :: nt_isosno, nt_isoice, k
630 : character(len=3) :: ck
631 : character(len=*),parameter :: subname='(read_restart_iso)'
632 :
633 0 : call icepack_query_tracer_indices(nt_isosno_out=nt_isosno, nt_isoice_out=nt_isoice)
634 0 : call icepack_warnings_flush(nu_diag)
635 0 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
636 0 : file=__FILE__, line=__LINE__)
637 :
638 0 : diag = .true.
639 :
640 0 : do k = 1, n_iso
641 0 : write(ck,'(i3.3)') k
642 0 : call read_restart_field(nu_restart_iso,0,trcrn(:,:,nt_isosno+k-1,:,:), &
643 : 'ruf8','isosno'//trim(ck),ncat,diag, & ! LCOV_EXCL_LINE
644 0 : field_type=field_type_scalar,field_loc=field_loc_center)
645 : enddo
646 :
647 0 : do k = 1, n_iso
648 0 : write(ck,'(i3.3)') k
649 0 : call read_restart_field(nu_restart_iso,0,trcrn(:,:,nt_isoice+k-1,:,:), &
650 : 'ruf8','isoice'//trim(ck),ncat,diag, & ! LCOV_EXCL_LINE
651 0 : field_type=field_type_scalar,field_loc=field_loc_center)
652 : enddo
653 :
654 0 : end subroutine read_restart_iso
655 :
656 : !=======================================================================
657 :
658 : ! Dumps all values needed for restarting
659 : !
660 : ! authors Elizabeth Hunke, LANL (original version)
661 : ! David Bailey, NCAR
662 : ! Marika Holland, NCAR
663 :
664 0 : subroutine write_restart_aero()
665 :
666 : use ice_domain_size, only: n_aero
667 : use ice_state, only: trcrn
668 : use ice_fileunits, only: nu_dump_aero
669 :
670 : ! local variables
671 :
672 : integer (kind=int_kind) :: &
673 : k ! loop indices
674 :
675 : logical (kind=log_kind) :: diag
676 :
677 : character (len=3) :: nchar
678 : integer (kind=int_kind) :: nt_aero
679 : character(len=*),parameter :: subname='(write_restart_aero)'
680 :
681 : !-----------------------------------------------------------------
682 :
683 0 : if (my_task == master_task) write(nu_diag,*) subname,'aerosols'
684 :
685 0 : call icepack_query_tracer_indices(nt_aero_out=nt_aero)
686 0 : call icepack_warnings_flush(nu_diag)
687 0 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
688 0 : file=__FILE__, line=__LINE__)
689 :
690 0 : diag = .true.
691 :
692 0 : do k = 1, n_aero
693 0 : write(nchar,'(i3.3)') k
694 : call write_restart_field(nu_dump_aero,0, &
695 : trcrn(:,:,nt_aero +(k-1)*4,:,:),'ruf8','aerosnossl'//nchar, & ! LCOV_EXCL_LINE
696 0 : ncat,diag)
697 : call write_restart_field(nu_dump_aero,0, &
698 : trcrn(:,:,nt_aero+1+(k-1)*4,:,:),'ruf8','aerosnoint'//nchar, & ! LCOV_EXCL_LINE
699 0 : ncat,diag)
700 : call write_restart_field(nu_dump_aero,0, &
701 : trcrn(:,:,nt_aero+2+(k-1)*4,:,:),'ruf8','aeroicessl'//nchar, & ! LCOV_EXCL_LINE
702 0 : ncat,diag)
703 : call write_restart_field(nu_dump_aero,0, &
704 : trcrn(:,:,nt_aero+3+(k-1)*4,:,:),'ruf8','aeroiceint'//nchar, & ! LCOV_EXCL_LINE
705 0 : ncat,diag)
706 : enddo
707 :
708 0 : end subroutine write_restart_aero
709 :
710 : !=======================================================================
711 :
712 : ! Reads all values needed for an ice aerosol restart
713 : !
714 : ! authors Elizabeth Hunke, LANL (original version)
715 : ! David Bailey, NCAR
716 : ! Marika Holland, NCAR
717 :
718 0 : subroutine read_restart_aero()
719 :
720 : use ice_domain_size, only: n_aero
721 : use ice_state, only: trcrn
722 : use ice_fileunits, only: nu_restart_aero
723 :
724 : ! local variables
725 :
726 : integer (kind=int_kind) :: &
727 : k ! loop indices
728 :
729 : logical (kind=log_kind) :: &
730 : diag
731 : integer (kind=int_kind) :: nt_aero
732 :
733 : character (len=3) :: nchar
734 : character(len=*),parameter :: subname='(read_restart_aero)'
735 :
736 : !-----------------------------------------------------------------
737 :
738 0 : if (my_task == master_task) write(nu_diag,*) subname,'aerosols'
739 :
740 0 : call icepack_query_tracer_indices(nt_aero_out=nt_aero)
741 0 : call icepack_warnings_flush(nu_diag)
742 0 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
743 0 : file=__FILE__, line=__LINE__)
744 :
745 0 : diag = .true.
746 :
747 0 : do k = 1, n_aero
748 0 : write(nchar,'(i3.3)') k
749 : call read_restart_field(nu_restart_aero,0, &
750 : trcrn(:,:,nt_aero +(k-1)*4,:,:),'ruf8','aerosnossl'//trim(nchar), & ! LCOV_EXCL_LINE
751 0 : ncat,diag,field_type=field_type_scalar,field_loc=field_loc_center)
752 : call read_restart_field(nu_restart_aero,0, &
753 : trcrn(:,:,nt_aero+1+(k-1)*4,:,:),'ruf8','aerosnoint'//trim(nchar), & ! LCOV_EXCL_LINE
754 0 : ncat,diag,field_type=field_type_scalar,field_loc=field_loc_center)
755 : call read_restart_field(nu_restart_aero,0, &
756 : trcrn(:,:,nt_aero+2+(k-1)*4,:,:),'ruf8','aeroicessl'//trim(nchar), & ! LCOV_EXCL_LINE
757 0 : ncat,diag,field_type=field_type_scalar,field_loc=field_loc_center)
758 : call read_restart_field(nu_restart_aero,0, &
759 : trcrn(:,:,nt_aero+3+(k-1)*4,:,:),'ruf8','aeroiceint'//trim(nchar), & ! LCOV_EXCL_LINE
760 0 : ncat,diag,field_type=field_type_scalar,field_loc=field_loc_center)
761 : enddo
762 :
763 0 : end subroutine read_restart_aero
764 :
765 : !=======================================================================
766 :
767 0 : subroutine read_restart_hbrine()
768 :
769 : ! Reads all values needed for hbrine
770 : ! author Elizabeth C. Hunke, LANL
771 :
772 : use ice_arrays_column, only: first_ice_real, first_ice
773 : use ice_blocks, only: block, get_block
774 : use ice_communicate, only: my_task, master_task
775 : use ice_domain, only: nblocks, blocks_ice
776 : use ice_fileunits, only: nu_restart_hbrine
777 : use ice_state, only: trcrn
778 : use ice_restart,only: read_restart_field
779 :
780 : ! local variables
781 :
782 : integer (kind=int_kind) :: &
783 : i, j, n, iblk , & ! horizontal indices ! LCOV_EXCL_LINE
784 : ilo,ihi,jlo,jhi ! beginning and end of physical domain
785 :
786 : type (block) :: &
787 : this_block ! block information for current block
788 :
789 : logical (kind=log_kind) :: &
790 : diag
791 :
792 : integer (kind=int_kind) :: nt_fbri
793 :
794 : character(len=*),parameter :: subname='(read_restart_hbrine)'
795 :
796 0 : call icepack_query_tracer_indices(nt_fbri_out=nt_fbri)
797 0 : call icepack_warnings_flush(nu_diag)
798 0 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
799 0 : file=__FILE__, line=__LINE__)
800 :
801 0 : diag = .true.
802 :
803 0 : if (my_task == master_task) write(nu_diag,*) subname,'brine restart'
804 :
805 0 : call read_restart_field(nu_restart_hbrine,0,trcrn(:,:,nt_fbri,:,:),'ruf8', &
806 0 : 'fbrn',ncat,diag,field_loc_center,field_type_scalar)
807 : call read_restart_field(nu_restart_hbrine,0,first_ice_real(:,:,:,:),'ruf8', &
808 0 : 'first_ice',ncat,diag,field_loc_center,field_type_scalar)
809 :
810 0 : !$OMP PARALLEL DO PRIVATE(iblk,i,j,n,ilo,ihi,jlo,jhi,this_block)
811 0 : do iblk = 1, nblocks
812 :
813 0 : this_block = get_block(blocks_ice(iblk),iblk)
814 0 : ilo = this_block%ilo
815 0 : ihi = this_block%ihi
816 0 : jlo = this_block%jlo
817 0 : jhi = this_block%jhi
818 :
819 0 : do j = jlo, jhi
820 0 : do i = ilo, ihi
821 0 : do n = 1, ncat
822 0 : if (first_ice_real(i,j,n,iblk) >= p5) then
823 0 : first_ice (i,j,n,iblk) = .true.
824 : else
825 0 : first_ice (i,j,n,iblk) = .false.
826 : endif
827 : enddo ! ncat
828 : enddo ! i
829 : enddo ! j
830 : enddo ! iblk
831 : !$OMP END PARALLEL DO
832 :
833 0 : end subroutine read_restart_hbrine
834 :
835 : !=======================================================================
836 :
837 0 : subroutine write_restart_hbrine()
838 :
839 : ! Dumps all values needed for a hbrine restart
840 : ! author Elizabeth C. Hunke, LANL
841 :
842 : use ice_arrays_column, only: first_ice, first_ice_real
843 : use ice_blocks, only: block, get_block
844 : use ice_domain, only: nblocks, blocks_ice
845 : use ice_fileunits, only: nu_dump_hbrine
846 : use ice_grid, only: tmask
847 : use ice_state, only: trcrn
848 : use ice_restart,only: write_restart_field
849 :
850 : ! local variables
851 :
852 : integer (kind=int_kind) :: &
853 : i, j, n, iblk , & ! horizontal indices ! LCOV_EXCL_LINE
854 : ilo,ihi,jlo,jhi ! beginning and end of physical domain
855 :
856 : logical (kind=log_kind) :: diag
857 :
858 : integer (kind=int_kind) :: nt_fbri
859 :
860 : type (block) :: &
861 : this_block ! block information for current block
862 : character(len=*),parameter :: subname='(write_restart_hbrine)'
863 :
864 0 : call icepack_query_tracer_indices(nt_fbri_out=nt_fbri)
865 0 : call icepack_warnings_flush(nu_diag)
866 0 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
867 0 : file=__FILE__, line=__LINE__)
868 :
869 0 : diag = .true.
870 :
871 0 : !$OMP PARALLEL DO PRIVATE(iblk,i,j,n,ilo,ihi,jlo,jhi,this_block)
872 0 : do iblk = 1, nblocks
873 :
874 0 : this_block = get_block(blocks_ice(iblk),iblk)
875 0 : ilo = this_block%ilo
876 0 : ihi = this_block%ihi
877 0 : jlo = this_block%jlo
878 0 : jhi = this_block%jhi
879 :
880 0 : do j = jlo, jhi
881 0 : do i = ilo, ihi
882 0 : do n = 1, ncat
883 : ! zero out first_ice over land
884 0 : if (tmask(i,j,iblk) .and. first_ice (i,j,n,iblk)) then
885 0 : first_ice_real(i,j,n,iblk) = c1
886 : else
887 0 : first_ice_real(i,j,n,iblk) = c0
888 : endif
889 : enddo ! n
890 : enddo ! i
891 : enddo ! j
892 : enddo ! iblk
893 : !$OMP END PARALLEL DO
894 :
895 0 : call write_restart_field(nu_dump_hbrine,0,trcrn(:,:,nt_fbri,:,:),'ruf8', &
896 0 : 'fbrn',ncat,diag)
897 : call write_restart_field(nu_dump_hbrine,0,first_ice_real(:,:,:,:),'ruf8', &
898 0 : 'first_ice',ncat,diag)
899 :
900 0 : end subroutine write_restart_hbrine
901 :
902 : !=======================================================================
903 : !
904 : ! Dumps all values needed for a bgc restart
905 : !
906 : ! author Elizabeth C. Hunke, LANL
907 :
908 0 : subroutine write_restart_bgc()
909 :
910 : use ice_blocks, only: block, get_block
911 : use ice_domain, only: nblocks, blocks_ice
912 : use ice_domain_size, only: ncat, n_algae, n_doc, n_dic, &
913 : n_don, n_zaero, n_fed, n_fep
914 : use ice_fileunits, only: nu_dump_bgc
915 : use ice_flux_bgc, only: nit, amm, sil, dmsp, dms, algalN, &
916 : doc, don, dic, fed, fep, zaeros, hum
917 : use ice_grid, only: tmask
918 : use ice_state, only: trcrn
919 : use ice_flux, only: sss
920 : use ice_restart, only: write_restart_field
921 :
922 : ! local variables
923 :
924 : integer (kind=int_kind) :: &
925 : i, j, k, iblk , & ! horizontal, vertical and block indices ! LCOV_EXCL_LINE
926 : mm , & ! n_algae ! LCOV_EXCL_LINE
927 : ilo,ihi,jlo,jhi ! beginning and end of physical domain
928 :
929 : logical (kind=log_kind) :: diag
930 :
931 : character (len=3) :: nchar, ncharb
932 :
933 : integer (kind=int_kind) :: nt_bgc_Am, &
934 : nt_bgc_DMS, nt_bgc_DMSPd, & ! LCOV_EXCL_LINE
935 : nt_bgc_DMSPp, nt_bgc_Nit, nt_bgc_Sil, & ! LCOV_EXCL_LINE
936 : nt_bgc_PON, nt_zbgc_frac, nt_bgc_hum, nbtrcr
937 :
938 : integer (kind=int_kind), dimension(icepack_max_algae) :: &
939 : nt_bgc_N , & ! diatoms, phaeocystis, pico/small ! LCOV_EXCL_LINE
940 : nt_bgc_C , & ! diatoms, phaeocystis, pico/small ! LCOV_EXCL_LINE
941 : nt_bgc_chl ! diatoms, phaeocystis, pico/small
942 :
943 : integer (kind=int_kind), dimension(icepack_max_doc) :: &
944 : nt_bgc_DOC ! dissolved organic carbon
945 :
946 : integer (kind=int_kind), dimension(icepack_max_don) :: &
947 : nt_bgc_DON ! dissolved organic nitrogen
948 :
949 : integer (kind=int_kind), dimension(icepack_max_dic) :: &
950 : nt_bgc_DIC ! dissolved inorganic carbon
951 :
952 : integer (kind=int_kind), dimension(icepack_max_fe) :: &
953 : nt_bgc_Fed, & ! dissolved iron ! LCOV_EXCL_LINE
954 : nt_bgc_Fep ! particulate iron
955 :
956 : integer (kind=int_kind), dimension(icepack_max_aero) :: &
957 : nt_zaero ! black carbon and other aerosols
958 :
959 : logical (kind=log_kind) :: tr_bgc_Nit, tr_bgc_Am, tr_bgc_Sil,&
960 : tr_bgc_DMS, tr_bgc_PON, tr_bgc_N, tr_bgc_C, & ! LCOV_EXCL_LINE
961 : tr_bgc_DON, tr_bgc_Fe, tr_zaero , tr_bgc_chl, & ! LCOV_EXCL_LINE
962 : tr_bgc_hum
963 :
964 : logical (kind=log_kind) :: skl_bgc
965 :
966 : type (block) :: &
967 : this_block ! block information for current block
968 :
969 : character(len=*),parameter :: subname='(write_restart_bgc)'
970 :
971 0 : call icepack_query_parameters(skl_bgc_out=skl_bgc)
972 0 : call icepack_query_tracer_sizes(nbtrcr_out=nbtrcr)
973 : call icepack_query_tracer_flags(tr_bgc_Nit_out=tr_bgc_Nit, &
974 : tr_bgc_Am_out=tr_bgc_Am, tr_bgc_Sil_out=tr_bgc_Sil, & ! LCOV_EXCL_LINE
975 : tr_bgc_DMS_out=tr_bgc_DMS, tr_bgc_PON_out=tr_bgc_PON, & ! LCOV_EXCL_LINE
976 : tr_bgc_N_out=tr_bgc_N, tr_bgc_C_out=tr_bgc_C, & ! LCOV_EXCL_LINE
977 : tr_bgc_DON_out=tr_bgc_DON, tr_bgc_Fe_out=tr_bgc_Fe, tr_zaero_out=tr_zaero, & ! LCOV_EXCL_LINE
978 0 : tr_bgc_chl_out=tr_bgc_chl, tr_bgc_hum_out=tr_bgc_hum)
979 : call icepack_query_tracer_indices(nt_bgc_Am_out=nt_bgc_Am, &
980 : nt_bgc_DMS_out=nt_bgc_DMS, nt_bgc_DMSPd_out=nt_bgc_DMSPd, & ! LCOV_EXCL_LINE
981 : nt_bgc_C_out=nt_bgc_C, nt_bgc_chl_out=nt_bgc_chl, & ! LCOV_EXCL_LINE
982 : nt_bgc_DMSPp_out=nt_bgc_DMSPp, nt_bgc_Nit_out=nt_bgc_Nit, & ! LCOV_EXCL_LINE
983 : nt_bgc_Sil_out=nt_bgc_Sil, nt_bgc_PON_out=nt_bgc_PON, & ! LCOV_EXCL_LINE
984 : nt_bgc_DON_out=nt_bgc_DON, nt_bgc_DOC_out=nt_bgc_DOC, nt_bgc_DIC_out=nt_bgc_DIC, & ! LCOV_EXCL_LINE
985 : nt_bgc_N_out=nt_bgc_N, nt_zaero_out=nt_zaero, nt_bgc_Fed_out=nt_bgc_Fed, & ! LCOV_EXCL_LINE
986 0 : nt_bgc_hum_out=nt_bgc_hum, nt_bgc_Fep_out=nt_bgc_Fep, nt_zbgc_frac_out=nt_zbgc_frac)
987 0 : call icepack_warnings_flush(nu_diag)
988 0 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
989 0 : file=__FILE__, line=__LINE__)
990 :
991 0 : diag = .true.
992 :
993 : !-----------------------------------------------------------------
994 : ! Zero out tracers over land
995 : !-----------------------------------------------------------------
996 :
997 0 : !$OMP PARALLEL DO PRIVATE(iblk,i,j,ilo,ihi,jlo,jhi,this_block)
998 0 : do iblk = 1, nblocks
999 0 : this_block = get_block(blocks_ice(iblk),iblk)
1000 0 : ilo = this_block%ilo
1001 0 : ihi = this_block%ihi
1002 0 : jlo = this_block%jlo
1003 0 : jhi = this_block%jhi
1004 0 : do j = jlo, jhi
1005 0 : do i = ilo, ihi
1006 0 : if (.not. tmask(i,j,iblk)) then
1007 0 : if (tr_bgc_N ) algalN(i,j,:,iblk) = c0
1008 0 : if (tr_bgc_C ) doc (i,j,:,iblk) = c0
1009 0 : if (tr_bgc_C ) dic (i,j,:,iblk) = c0
1010 0 : if (tr_bgc_Nit) nit (i,j ,iblk) = c0
1011 0 : if (tr_bgc_Am ) amm (i,j ,iblk) = c0
1012 0 : if (tr_bgc_Sil) sil (i,j ,iblk) = c0
1013 0 : if (tr_bgc_hum) hum (i,j ,iblk) = c0
1014 0 : if (tr_bgc_DMS) dms (i,j ,iblk) = c0
1015 0 : if (tr_bgc_DMS) dmsp (i,j ,iblk) = c0
1016 0 : if (tr_bgc_DON) don (i,j,:,iblk) = c0
1017 0 : if (tr_bgc_Fe ) fed (i,j,:,iblk) = c0
1018 0 : if (tr_bgc_Fe ) fep (i,j,:,iblk) = c0
1019 : endif
1020 : enddo
1021 : enddo
1022 : enddo
1023 : !$OMP END PARALLEL DO
1024 :
1025 : !-----------------------------------------------------------------
1026 : ! Skeletal layer BGC
1027 : !-----------------------------------------------------------------
1028 :
1029 0 : if (skl_bgc) then
1030 0 : do k = 1, n_algae
1031 0 : write(nchar,'(i3.3)') k
1032 0 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_N(k),:,:), &
1033 0 : 'ruf8','bgc_N'//trim(nchar),ncat,diag)
1034 0 : if (tr_bgc_chl) &
1035 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_chl(k),:,:), & ! LCOV_EXCL_LINE
1036 0 : 'ruf8','bgc_chl'//trim(nchar),ncat,diag)
1037 : enddo
1038 0 : if (tr_bgc_C) then
1039 : ! do k = 1, n_algae
1040 : ! write(nchar,'(i3.3)') k
1041 : ! call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_C(k),:,:), &
1042 : ! 'ruf8','bgc_C'//trim(nchar),ncat,diag)
1043 : ! enddo
1044 0 : do k = 1, n_doc
1045 0 : write(nchar,'(i3.3)') k
1046 0 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_DOC(k),:,:), &
1047 0 : 'ruf8','bgc_DOC'//trim(nchar),ncat,diag)
1048 : enddo
1049 0 : do k = 1, n_dic
1050 0 : write(nchar,'(i3.3)') k
1051 0 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_DIC(k),:,:), &
1052 0 : 'ruf8','bgc_DIC'//trim(nchar),ncat,diag)
1053 : enddo
1054 : endif
1055 0 : if (tr_bgc_Nit) &
1056 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_Nit,:,:), & ! LCOV_EXCL_LINE
1057 0 : 'ruf8','bgc_Nit',ncat,diag)
1058 0 : if (tr_bgc_Am) &
1059 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_Am,:,:), & ! LCOV_EXCL_LINE
1060 0 : 'ruf8','bgc_Am',ncat,diag)
1061 0 : if (tr_bgc_Sil) &
1062 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_Sil,:,:), & ! LCOV_EXCL_LINE
1063 0 : 'ruf8','bgc_Sil',ncat,diag)
1064 0 : if (tr_bgc_hum) &
1065 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_hum,:,:), & ! LCOV_EXCL_LINE
1066 0 : 'ruf8','bgc_hum',ncat,diag)
1067 0 : if (tr_bgc_DMS) then
1068 0 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_DMSPp,:,:), &
1069 0 : 'ruf8','bgc_DMSPp',ncat,diag)
1070 0 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_DMSPd,:,:), &
1071 0 : 'ruf8','bgc_DMSPd',ncat,diag)
1072 0 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_DMS,:,:), &
1073 0 : 'ruf8','bgc_DMS',ncat,diag)
1074 : endif
1075 0 : if (tr_bgc_PON) &
1076 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_PON,:,:), & ! LCOV_EXCL_LINE
1077 0 : 'ruf8','bgc_PON',ncat,diag)
1078 :
1079 0 : if (tr_bgc_DON) then
1080 0 : do k = 1, n_don
1081 0 : write(nchar,'(i3.3)') k
1082 0 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_DON(k),:,:), &
1083 0 : 'ruf8','bgc_DON'//trim(nchar),ncat,diag)
1084 : enddo
1085 : endif
1086 0 : if (tr_bgc_Fe ) then
1087 0 : do k = 1, n_fed
1088 0 : write(nchar,'(i3.3)') k
1089 0 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_Fed (k),:,:), &
1090 0 : 'ruf8','bgc_Fed'//trim(nchar),ncat,diag)
1091 : enddo
1092 0 : do k = 1, n_fep
1093 0 : write(nchar,'(i3.3)') k
1094 0 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_Fep (k),:,:), &
1095 0 : 'ruf8','bgc_Fep'//trim(nchar),ncat,diag)
1096 : enddo
1097 : endif
1098 :
1099 : else
1100 :
1101 : !-----------------------------------------------------------------
1102 : ! Z layer BGC
1103 : !-----------------------------------------------------------------
1104 :
1105 0 : if (tr_bgc_Nit) then
1106 0 : do k = 1,nblyr+3
1107 0 : write(nchar,'(i3.3)') k
1108 0 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_Nit+k-1,:,:),'ruf8', &
1109 0 : 'bgc_Nit'//trim(nchar),ncat,diag)
1110 : enddo
1111 : endif
1112 0 : if (tr_bgc_N) then
1113 0 : do mm = 1,n_algae
1114 0 : write(ncharb, '(i3.3)') mm
1115 0 : do k = 1,nblyr+3
1116 0 : write(nchar,'(i3.3)') k
1117 : call write_restart_field(nu_dump_bgc,0, &
1118 : trcrn(:,:,nt_bgc_N(mm)+k-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1119 0 : 'bgc_N'//trim(ncharb)//trim(nchar),ncat,diag)
1120 : enddo
1121 0 : if (tr_bgc_chl) then
1122 0 : do k = 1,nblyr+3
1123 0 : write(nchar,'(i3.3)') k
1124 : call write_restart_field(nu_dump_bgc,0, &
1125 : trcrn(:,:,nt_bgc_chl(mm)+k-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1126 0 : 'bgc_chl'//trim(ncharb)//trim(nchar),ncat,diag)
1127 : enddo
1128 : endif
1129 : enddo !n_algae
1130 : endif ! tr_bgc_N
1131 0 : if (tr_bgc_C) then
1132 : ! do mm = 1,n_algae
1133 : ! write(ncharb, '(i3.3)') mm
1134 : ! do k = 1,nblyr+3
1135 : ! write(nchar,'(i3.3)') k
1136 : ! call write_restart_field(nu_dump_bgc,0, &
1137 : ! trcrn(:,:,nt_bgc_C(mm)+k-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1138 : ! 'bgc_C'//trim(ncharb)//trim(nchar),ncat,diag)
1139 : ! enddo
1140 : ! enddo
1141 0 : do mm = 1,n_doc
1142 0 : write(ncharb, '(i3.3)') mm
1143 0 : do k = 1,nblyr+3
1144 0 : write(nchar,'(i3.3)') k
1145 : call write_restart_field(nu_dump_bgc,0, &
1146 : trcrn(:,:,nt_bgc_DOC(mm)+k-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1147 0 : 'bgc_DOC'//trim(ncharb)//trim(nchar),ncat,diag)
1148 : enddo
1149 : enddo
1150 0 : do mm = 1,n_dic
1151 0 : write(ncharb, '(i3.3)') mm
1152 0 : do k = 1,nblyr+3
1153 0 : write(nchar,'(i3.3)') k
1154 : call write_restart_field(nu_dump_bgc,0, &
1155 : trcrn(:,:,nt_bgc_DIC(mm)+k-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1156 0 : 'bgc_DIC'//trim(ncharb)//trim(nchar),ncat,diag)
1157 : enddo
1158 : enddo
1159 : endif !tr_bgc_C
1160 0 : if (tr_bgc_Am) then
1161 0 : do k = 1,nblyr+3
1162 0 : write(nchar,'(i3.3)') k
1163 0 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_Am+k-1,:,:),'ruf8', &
1164 0 : 'bgc_Am'//trim(nchar),ncat,diag)
1165 : enddo
1166 : endif
1167 0 : if (tr_bgc_Sil) then
1168 0 : do k = 1,nblyr+3
1169 0 : write(nchar,'(i3.3)') k
1170 0 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_Sil+k-1,:,:),'ruf8', &
1171 0 : 'bgc_Sil'//trim(nchar),ncat,diag)
1172 : enddo
1173 : endif
1174 0 : if (tr_bgc_hum) then
1175 0 : do k = 1,nblyr+3
1176 0 : write(nchar,'(i3.3)') k
1177 0 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_hum+k-1,:,:),'ruf8', &
1178 0 : 'bgc_hum'//trim(nchar),ncat,diag)
1179 : enddo
1180 : endif
1181 0 : if (tr_bgc_DMS) then
1182 0 : do k = 1,nblyr+3
1183 0 : write(nchar,'(i3.3)') k
1184 0 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_DMSPp+k-1,:,:),'ruf8', &
1185 0 : 'bgc_DMSPp'//trim(nchar),ncat,diag)
1186 0 : write(nchar,'(i3.3)') k
1187 0 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_DMSPd+k-1,:,:),'ruf8', &
1188 0 : 'bgc_DMSPd'//trim(nchar),ncat,diag)
1189 0 : write(nchar,'(i3.3)') k
1190 0 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_DMS+k-1,:,:),'ruf8', &
1191 0 : 'bgc_DMS'//trim(nchar),ncat,diag)
1192 : enddo
1193 : endif
1194 0 : if (tr_bgc_PON) then
1195 0 : do k = 1,nblyr+3
1196 0 : write(nchar,'(i3.3)') k
1197 0 : call write_restart_field(nu_dump_bgc,0,trcrn(:,:,nt_bgc_PON+k-1,:,:),'ruf8', &
1198 0 : 'bgc_PON'//trim(nchar),ncat,diag)
1199 : enddo
1200 : endif
1201 0 : if (tr_bgc_DON) then
1202 0 : do mm = 1,n_don
1203 0 : write(ncharb, '(i3.3)') mm
1204 0 : do k = 1,nblyr+3
1205 0 : write(nchar,'(i3.3)') k
1206 : call write_restart_field(nu_dump_bgc,0, &
1207 : trcrn(:,:,nt_bgc_DON(mm)+k-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1208 0 : 'bgc_DON'//trim(ncharb)//trim(nchar),ncat,diag)
1209 : enddo
1210 : enddo
1211 : endif
1212 0 : if (tr_bgc_Fe ) then
1213 0 : do mm = 1,n_fed
1214 0 : write(ncharb, '(i3.3)') mm
1215 0 : do k = 1,nblyr+3
1216 0 : write(nchar,'(i3.3)') k
1217 : call write_restart_field(nu_dump_bgc,0, &
1218 : trcrn(:,:,nt_bgc_Fed(mm)+k-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1219 0 : 'bgc_Fed'//trim(ncharb)//trim(nchar),ncat,diag)
1220 : enddo
1221 : enddo
1222 0 : do mm = 1,n_fep
1223 0 : write(ncharb, '(i3.3)') mm
1224 0 : do k = 1,nblyr+3
1225 0 : write(nchar,'(i3.3)') k
1226 : call write_restart_field(nu_dump_bgc,0, &
1227 : trcrn(:,:,nt_bgc_Fep(mm)+k-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1228 0 : 'bgc_Fep'//trim(ncharb)//trim(nchar),ncat,diag)
1229 : enddo
1230 : enddo
1231 : endif
1232 0 : if (tr_zaero) then
1233 0 : do mm = 1,n_zaero
1234 0 : write(ncharb, '(i3.3)') mm
1235 0 : do k = 1,nblyr+3
1236 0 : write(nchar,'(i3.3)') k
1237 : call write_restart_field(nu_dump_bgc,0, &
1238 : trcrn(:,:,nt_zaero(mm)+k-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1239 0 : 'zaero'//trim(ncharb)//trim(nchar),ncat,diag)
1240 : enddo
1241 : enddo
1242 : endif
1243 0 : do mm = 1,nbtrcr
1244 0 : write(nchar,'(i3.3)') mm
1245 : call write_restart_field(nu_dump_bgc,0, &
1246 : trcrn(:,:,nt_zbgc_frac+mm-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1247 0 : 'zbgc_frac'//trim(nchar),ncat,diag)
1248 : enddo
1249 : endif
1250 :
1251 : !-----------------------------------------------------------------
1252 : ! Ocean BGC
1253 : !-----------------------------------------------------------------
1254 :
1255 0 : if (tr_bgc_N) then
1256 0 : do k = 1,n_algae
1257 0 : write(nchar,'(i3.3)') k
1258 0 : call write_restart_field(nu_dump_bgc,0,algalN(:,:,k,:),'ruf8','algalN'//trim(nchar),1,diag)
1259 : enddo !k
1260 : endif
1261 0 : if (tr_bgc_C) then
1262 0 : do k = 1,n_doc
1263 0 : write(nchar,'(i3.3)') k
1264 0 : call write_restart_field(nu_dump_bgc,0,doc(:,:,k,:),'ruf8','doc'//trim(nchar),1,diag)
1265 : enddo !k
1266 0 : do k = 1,n_dic
1267 0 : write(nchar,'(i3.3)') k
1268 0 : call write_restart_field(nu_dump_bgc,0,dic(:,:,k,:),'ruf8','dic'//trim(nchar),1,diag)
1269 : enddo !k
1270 : endif
1271 0 : if (tr_bgc_Nit) &
1272 0 : call write_restart_field(nu_dump_bgc,0,nit, 'ruf8','nit', 1,diag)
1273 0 : if (tr_bgc_Am) &
1274 0 : call write_restart_field(nu_dump_bgc,0,amm, 'ruf8','amm', 1,diag)
1275 0 : if (tr_bgc_Sil) &
1276 0 : call write_restart_field(nu_dump_bgc,0,sil, 'ruf8','sil', 1,diag)
1277 0 : if (tr_bgc_hum) &
1278 0 : call write_restart_field(nu_dump_bgc,0,hum, 'ruf8','hum', 1,diag)
1279 0 : if (tr_bgc_DMS) then
1280 0 : call write_restart_field(nu_dump_bgc,0,dmsp, 'ruf8','dmsp', 1,diag)
1281 0 : call write_restart_field(nu_dump_bgc,0,dms, 'ruf8','dms', 1,diag)
1282 : endif
1283 0 : if (tr_bgc_DON) then
1284 0 : do k = 1,n_don
1285 0 : write(nchar,'(i3.3)') k
1286 0 : call write_restart_field(nu_dump_bgc,0,don(:,:,k,:),'ruf8','don'//trim(nchar),1,diag)
1287 : enddo !k
1288 : endif
1289 0 : if (tr_bgc_Fe ) then
1290 0 : do k = 1,n_fed
1291 0 : write(nchar,'(i3.3)') k
1292 0 : call write_restart_field(nu_dump_bgc,0,fed(:,:,k,:),'ruf8','fed'//trim(nchar),1,diag)
1293 : enddo !k
1294 0 : do k = 1,n_fep
1295 0 : write(nchar,'(i3.3)') k
1296 0 : call write_restart_field(nu_dump_bgc,0,fep(:,:,k,:),'ruf8','fep'//trim(nchar),1,diag)
1297 : enddo !k
1298 : endif
1299 0 : if (tr_zaero) then
1300 0 : do k = 1,n_zaero
1301 0 : write(nchar,'(i3.3)') k
1302 0 : call write_restart_field(nu_dump_bgc,0,zaeros(:,:,k,:),'ruf8','zaeros'//trim(nchar),1,diag)
1303 : enddo !k
1304 : endif
1305 :
1306 0 : end subroutine write_restart_bgc
1307 :
1308 : !=======================================================================
1309 : !
1310 : ! Reads all values needed for a bgc restart
1311 : !
1312 : ! author Elizabeth C. Hunke, LANL
1313 :
1314 0 : subroutine read_restart_bgc()
1315 :
1316 : use ice_blocks, only: block, get_block
1317 : use ice_communicate, only: my_task, master_task
1318 : use ice_domain, only: nblocks, blocks_ice
1319 : use ice_domain_size, only: ncat, n_algae, n_doc, n_dic,&
1320 : n_don, n_zaero, n_fed, n_fep
1321 : use ice_fileunits, only: nu_restart_bgc
1322 : use ice_flux, only: sss
1323 : use ice_flux_bgc, only: nit, amm, sil, dmsp, dms, algalN, &
1324 : doc, don, dic, fed, fep, zaeros, hum
1325 : use ice_state, only: trcrn
1326 : use ice_restart, only: read_restart_field
1327 :
1328 : ! local variables
1329 :
1330 : type (block) :: &
1331 : this_block ! block information for current block
1332 :
1333 : integer (kind=int_kind) :: &
1334 : i, j, k, iblk, & ! indices ! LCOV_EXCL_LINE
1335 : mm , & ! n_algae ! LCOV_EXCL_LINE
1336 : ilo,ihi,jlo,jhi ! beginning and end of physical domain
1337 :
1338 : logical (kind=log_kind) :: diag
1339 :
1340 : integer (kind=int_kind) :: nt_bgc_Am, &
1341 : nt_bgc_DMS, nt_bgc_DMSPd, & ! LCOV_EXCL_LINE
1342 : nt_bgc_DMSPp, nt_bgc_Nit, nt_bgc_Sil, & ! LCOV_EXCL_LINE
1343 : nt_bgc_PON, nt_zbgc_frac, nt_bgc_hum, nbtrcr
1344 :
1345 : integer (kind=int_kind), dimension(icepack_max_algae) :: &
1346 : nt_bgc_N , & ! diatoms, phaeocystis, pico/small ! LCOV_EXCL_LINE
1347 : nt_bgc_C , & ! diatoms, phaeocystis, pico/small ! LCOV_EXCL_LINE
1348 : nt_bgc_chl ! diatoms, phaeocystis, pico/small
1349 :
1350 : integer (kind=int_kind), dimension(icepack_max_doc) :: &
1351 : nt_bgc_DOC ! dissolved organic carbon
1352 :
1353 : integer (kind=int_kind), dimension(icepack_max_don) :: &
1354 : nt_bgc_DON ! dissolved organic nitrogen
1355 :
1356 : integer (kind=int_kind), dimension(icepack_max_dic) :: &
1357 : nt_bgc_DIC ! dissolved inorganic carbon
1358 :
1359 : integer (kind=int_kind), dimension(icepack_max_fe) :: &
1360 : nt_bgc_Fed, & ! dissolved iron ! LCOV_EXCL_LINE
1361 : nt_bgc_Fep ! particulate iron
1362 :
1363 : integer (kind=int_kind), dimension(icepack_max_aero) :: &
1364 : nt_zaero ! black carbon and other aerosols
1365 :
1366 : logical (kind=log_kind) :: tr_bgc_Nit, tr_bgc_Am, tr_bgc_Sil,&
1367 : tr_bgc_DMS, tr_bgc_PON, tr_bgc_N, tr_bgc_C, & ! LCOV_EXCL_LINE
1368 : tr_bgc_DON, tr_bgc_Fe, tr_zaero , tr_bgc_chl, & ! LCOV_EXCL_LINE
1369 : tr_bgc_hum
1370 :
1371 : logical (kind=log_kind) :: skl_bgc
1372 :
1373 : character (len=3) :: nchar, ncharb
1374 :
1375 : character(len=*),parameter :: subname='(read_restart_bgc)'
1376 :
1377 0 : call icepack_query_parameters(skl_bgc_out=skl_bgc)
1378 0 : call icepack_query_tracer_sizes(nbtrcr_out=nbtrcr)
1379 : call icepack_query_tracer_flags(tr_bgc_Nit_out=tr_bgc_Nit, &
1380 : tr_bgc_Am_out=tr_bgc_Am, tr_bgc_Sil_out=tr_bgc_Sil, & ! LCOV_EXCL_LINE
1381 : tr_bgc_DMS_out=tr_bgc_DMS, tr_bgc_PON_out=tr_bgc_PON, & ! LCOV_EXCL_LINE
1382 : tr_bgc_N_out=tr_bgc_N, tr_bgc_C_out=tr_bgc_C, & ! LCOV_EXCL_LINE
1383 : tr_bgc_DON_out=tr_bgc_DON, tr_bgc_Fe_out=tr_bgc_Fe, tr_zaero_out=tr_zaero, & ! LCOV_EXCL_LINE
1384 0 : tr_bgc_chl_out=tr_bgc_chl, tr_bgc_hum_out=tr_bgc_hum)
1385 : call icepack_query_tracer_indices(nt_bgc_Am_out=nt_bgc_Am, &
1386 : nt_bgc_DMS_out=nt_bgc_DMS, nt_bgc_DMSPd_out=nt_bgc_DMSPd, & ! LCOV_EXCL_LINE
1387 : nt_bgc_C_out=nt_bgc_C, nt_bgc_chl_out=nt_bgc_chl, & ! LCOV_EXCL_LINE
1388 : nt_bgc_DMSPp_out=nt_bgc_DMSPp, nt_bgc_Nit_out=nt_bgc_Nit, & ! LCOV_EXCL_LINE
1389 : nt_bgc_Sil_out=nt_bgc_Sil, nt_bgc_PON_out=nt_bgc_PON, & ! LCOV_EXCL_LINE
1390 : nt_bgc_DON_out=nt_bgc_DON, nt_bgc_DOC_out=nt_bgc_DOC, nt_bgc_DIC_out=nt_bgc_DIC, & ! LCOV_EXCL_LINE
1391 : nt_bgc_N_out=nt_bgc_N, nt_zaero_out=nt_zaero, nt_bgc_Fed_out=nt_bgc_Fed, & ! LCOV_EXCL_LINE
1392 0 : nt_bgc_Fep_out=nt_bgc_Fep, nt_zbgc_frac_out=nt_zbgc_frac, nt_bgc_hum_out=nt_bgc_hum)
1393 0 : call icepack_warnings_flush(nu_diag)
1394 0 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
1395 0 : file=__FILE__, line=__LINE__)
1396 :
1397 0 : diag = .true.
1398 :
1399 : !-----------------------------------------------------------------
1400 : ! Skeletal Layer BGC
1401 : !-----------------------------------------------------------------
1402 0 : if (restart_bgc) then
1403 :
1404 0 : if (skl_bgc) then
1405 0 : if (my_task == master_task) write(nu_diag,*) subname,'skl bgc restart'
1406 :
1407 0 : do k = 1, n_algae
1408 0 : write(nchar,'(i3.3)') k
1409 : call read_restart_field(nu_restart_bgc,0, &
1410 : trcrn(:,:,nt_bgc_N(k),:,:), & ! LCOV_EXCL_LINE
1411 0 : 'ruf8','bgc_N'//trim(nchar),ncat,diag)
1412 0 : if (tr_bgc_chl) &
1413 : call read_restart_field(nu_restart_bgc,0, & ! LCOV_EXCL_LINE
1414 : trcrn(:,:,nt_bgc_chl(k),:,:), & ! LCOV_EXCL_LINE
1415 0 : 'ruf8','bgc_chl'//trim(nchar),ncat,diag)
1416 : enddo !k
1417 0 : if (tr_bgc_C) then
1418 : ! do k = 1, n_algae
1419 : ! write(nchar,'(i3.3)') k
1420 : ! call read_restart_field(nu_restart_bgc,0, &
1421 : ! trcrn(:,:,nt_bgc_C(k),:,:), & ! LCOV_EXCL_LINE
1422 : ! 'ruf8','bgc_C'//trim(nchar),ncat,diag)
1423 : ! enddo
1424 0 : do k = 1, n_doc
1425 0 : write(nchar,'(i3.3)') k
1426 : call read_restart_field(nu_restart_bgc,0, &
1427 : trcrn(:,:,nt_bgc_DOC(k),:,:), & ! LCOV_EXCL_LINE
1428 0 : 'ruf8','bgc_DOC'//trim(nchar),ncat,diag)
1429 : enddo
1430 0 : do k = 1, n_dic
1431 0 : write(nchar,'(i3.3)') k
1432 : call read_restart_field(nu_restart_bgc,0, &
1433 : trcrn(:,:,nt_bgc_DIC(k),:,:), & ! LCOV_EXCL_LINE
1434 0 : 'ruf8','bgc_DIC'//trim(nchar),ncat,diag)
1435 : enddo
1436 : endif
1437 0 : if (tr_bgc_Nit) &
1438 : call read_restart_field(nu_restart_bgc,0,trcrn(:,:,nt_bgc_Nit,:,:), & ! LCOV_EXCL_LINE
1439 0 : 'ruf8','bgc_Nit',ncat,diag)
1440 0 : if (tr_bgc_Am) &
1441 : call read_restart_field(nu_restart_bgc,0,trcrn(:,:,nt_bgc_Am,:,:), & ! LCOV_EXCL_LINE
1442 0 : 'ruf8','bgc_Am',ncat,diag)
1443 0 : if (tr_bgc_Sil) &
1444 : call read_restart_field(nu_restart_bgc,0,trcrn(:,:,nt_bgc_Sil,:,:), & ! LCOV_EXCL_LINE
1445 0 : 'ruf8','bgc_Sil',ncat,diag)
1446 0 : if (tr_bgc_hum) &
1447 : call read_restart_field(nu_restart_bgc,0,trcrn(:,:,nt_bgc_hum,:,:), & ! LCOV_EXCL_LINE
1448 0 : 'ruf8','bgc_hum',ncat,diag)
1449 0 : if(tr_bgc_DMS) then
1450 0 : call read_restart_field(nu_restart_bgc,0,trcrn(:,:,nt_bgc_DMSPp,:,:), &
1451 0 : 'ruf8','bgc_DMSPp',ncat,diag)
1452 0 : call read_restart_field(nu_restart_bgc,0,trcrn(:,:,nt_bgc_DMSPd,:,:), &
1453 0 : 'ruf8','bgc_DMSPd',ncat,diag)
1454 0 : call read_restart_field(nu_restart_bgc,0,trcrn(:,:,nt_bgc_DMS,:,:), &
1455 0 : 'ruf8','bgc_DMS',ncat,diag)
1456 : endif
1457 0 : if (tr_bgc_PON) &
1458 : call read_restart_field(nu_restart_bgc,0,trcrn(:,:,nt_bgc_PON,:,:), & ! LCOV_EXCL_LINE
1459 0 : 'ruf8','bgc_PON',ncat,diag)
1460 0 : if (tr_bgc_DON) then
1461 0 : do k = 1, n_don
1462 0 : write(nchar,'(i3.3)') k
1463 : call read_restart_field(nu_restart_bgc,0, &
1464 : trcrn(:,:,nt_bgc_DON(k),:,:), & ! LCOV_EXCL_LINE
1465 0 : 'ruf8','bgc_DON'//trim(nchar),ncat,diag)
1466 : enddo
1467 : endif
1468 0 : if (tr_bgc_Fe) then
1469 0 : do k = 1, n_fed
1470 0 : write(nchar,'(i3.3)') k
1471 : call read_restart_field(nu_restart_bgc,0, &
1472 : trcrn(:,:,nt_bgc_Fed (k),:,:), & ! LCOV_EXCL_LINE
1473 0 : 'ruf8','bgc_Fed'//trim(nchar),ncat,diag)
1474 : enddo
1475 0 : do k = 1, n_fep
1476 0 : write(nchar,'(i3.3)') k
1477 : call read_restart_field(nu_restart_bgc,0, &
1478 : trcrn(:,:,nt_bgc_Fep (k),:,:), & ! LCOV_EXCL_LINE
1479 0 : 'ruf8','bgc_Fep'//trim(nchar),ncat,diag)
1480 : enddo
1481 : endif
1482 :
1483 : else
1484 :
1485 : !-----------------------------------------------------------------
1486 : ! Z Layer BGC
1487 : !-----------------------------------------------------------------
1488 :
1489 0 : if (tr_bgc_Nit) then
1490 0 : if (my_task == master_task) write(nu_diag,*) subname,'z bgc restart: min/max Nitrate'
1491 0 : do k=1,nblyr+3
1492 0 : write(nchar,'(i3.3)') k
1493 0 : call read_restart_field(nu_restart_bgc,0,trcrn(:,:,nt_bgc_Nit+k-1,:,:),'ruf8', &
1494 0 : 'bgc_Nit'//trim(nchar),ncat,diag,field_loc_center,field_type_scalar)
1495 : enddo
1496 : endif !Nit
1497 0 : if (tr_bgc_N) then
1498 0 : do mm = 1,n_algae
1499 0 : write(ncharb,'(i3.3)') mm
1500 0 : if (my_task == master_task) write(nu_diag,*) subname,' min/max Algal N'
1501 0 : do k=1,nblyr+3
1502 0 : write(nchar,'(i3.3)') k
1503 : call read_restart_field(nu_restart_bgc,0, &
1504 : trcrn(:,:,nt_bgc_N(mm)+k-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1505 0 : 'bgc_N'//trim(ncharb)//trim(nchar),ncat,diag,field_loc_center,field_type_scalar)
1506 : enddo
1507 0 : if (tr_bgc_chl) then
1508 0 : if (my_task == master_task) write(nu_diag,*) subname,' min/max Algal chla'
1509 0 : do k=1,nblyr+3
1510 0 : write(nchar,'(i3.3)') k
1511 : call read_restart_field(nu_restart_bgc,0, &
1512 : trcrn(:,:,nt_bgc_chl(mm)+k-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1513 0 : 'bgc_chl'//trim(ncharb)//trim(nchar),ncat,diag,field_loc_center,field_type_scalar)
1514 : enddo
1515 : endif ! tr_bgc_chl
1516 : enddo !n_algae
1517 : endif ! tr_bgc_N
1518 0 : if (tr_bgc_C) then
1519 : ! do mm = 1,n_algae
1520 : ! write(ncharb,'(i3.3)') mm
1521 : ! if (my_task == master_task) write(nu_diag,*) subname,' min/max Algal C'
1522 : ! do k=1,nblyr+3
1523 : ! write(nchar,'(i3.3)') k
1524 : ! call read_restart_field(nu_restart_bgc,0, &
1525 : ! trcrn(:,:,nt_bgc_C(mm)+k-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1526 : ! 'bgc_C'//trim(ncharb)//trim(nchar),ncat,diag,field_loc_center,field_type_scalar)
1527 : ! enddo
1528 : ! enddo !mm
1529 0 : do mm = 1,n_doc
1530 0 : write(ncharb,'(i3.3)') mm
1531 0 : if (my_task == master_task) write(nu_diag,*) subname,' min/max DOC'
1532 0 : do k=1,nblyr+3
1533 0 : write(nchar,'(i3.3)') k
1534 : call read_restart_field(nu_restart_bgc,0, &
1535 : trcrn(:,:,nt_bgc_DOC(mm)+k-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1536 0 : 'bgc_DOC'//trim(ncharb)//trim(nchar),ncat,diag,field_loc_center,field_type_scalar)
1537 : enddo
1538 : enddo !mm
1539 0 : do mm = 1,n_dic
1540 0 : write(ncharb,'(i3.3)') mm
1541 0 : if (my_task == master_task) write(nu_diag,*) subname,' min/max DIC'
1542 0 : do k=1,nblyr+3
1543 0 : write(nchar,'(i3.3)') k
1544 : call read_restart_field(nu_restart_bgc,0, &
1545 : trcrn(:,:,nt_bgc_DIC(mm)+k-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1546 0 : 'bgc_DIC'//trim(ncharb)//trim(nchar),ncat,diag,field_loc_center,field_type_scalar)
1547 : enddo
1548 : enddo !mm
1549 : endif ! tr_bgc_C
1550 0 : if (tr_bgc_Am) then
1551 0 : if (my_task == master_task) write(nu_diag,*) subname,' min/max ammonium'
1552 0 : do k=1,nblyr+3
1553 0 : write(nchar,'(i3.3)') k
1554 0 : call read_restart_field(nu_restart_bgc,0,trcrn(:,:,nt_bgc_Am+k-1,:,:),'ruf8', &
1555 0 : 'bgc_Am'//trim(nchar),ncat,diag,field_loc_center,field_type_scalar)
1556 : enddo
1557 : endif
1558 0 : if (tr_bgc_Sil) then
1559 0 : if (my_task == master_task) write(nu_diag,*) subname,' min/max silicate'
1560 0 : do k=1,nblyr+3
1561 0 : write(nchar,'(i3.3)') k
1562 0 : call read_restart_field(nu_restart_bgc,0,trcrn(:,:,nt_bgc_Sil+k-1,:,:),'ruf8', &
1563 0 : 'bgc_Sil'//trim(nchar),ncat,diag,field_loc_center,field_type_scalar)
1564 : enddo
1565 : endif
1566 0 : if (tr_bgc_hum) then
1567 0 : if (my_task == master_task) write(nu_diag,*) subname,' min/max humic material'
1568 0 : do k=1,nblyr+3
1569 0 : write(nchar,'(i3.3)') k
1570 0 : call read_restart_field(nu_restart_bgc,0,trcrn(:,:,nt_bgc_hum+k-1,:,:),'ruf8', &
1571 0 : 'bgc_hum'//trim(nchar),ncat,diag,field_loc_center,field_type_scalar)
1572 : enddo
1573 : endif
1574 0 : if (tr_bgc_DMS) then
1575 0 : do k=1,nblyr+3
1576 0 : if (my_task == master_task) write(nu_diag,*) subname,' min/max DMSPp'
1577 0 : write(nchar,'(i3.3)') k
1578 0 : call read_restart_field(nu_restart_bgc,0,trcrn(:,:,nt_bgc_DMSPp+k-1,:,:),'ruf8', &
1579 0 : 'bgc_DMSPp'//trim(nchar),ncat,diag,field_loc_center,field_type_scalar)
1580 0 : if (my_task == master_task) write(nu_diag,*) subname,' min/max DMSPd'
1581 0 : write(nchar,'(i3.3)') k
1582 0 : call read_restart_field(nu_restart_bgc,0,trcrn(:,:,nt_bgc_DMSPd+k-1,:,:),'ruf8', &
1583 0 : 'bgc_DMSPd'//trim(nchar),ncat,diag,field_loc_center,field_type_scalar)
1584 0 : if (my_task == master_task) write(nu_diag,*) subname,' min/max DMS'
1585 0 : write(nchar,'(i3.3)') k
1586 0 : call read_restart_field(nu_restart_bgc,0,trcrn(:,:,nt_bgc_DMS+k-1,:,:),'ruf8', &
1587 0 : 'bgc_DMS'//trim(nchar),ncat,diag,field_loc_center,field_type_scalar)
1588 : enddo
1589 : endif
1590 0 : if (tr_bgc_PON) then
1591 0 : if (my_task == master_task) write(nu_diag,*) subname,' min/max PON'
1592 0 : do k=1,nblyr+3
1593 0 : write(nchar,'(i3.3)') k
1594 0 : call read_restart_field(nu_restart_bgc,0,trcrn(:,:,nt_bgc_PON+k-1,:,:),'ruf8', &
1595 0 : 'bgc_PON'//trim(nchar),ncat,diag,field_loc_center,field_type_scalar)
1596 : enddo
1597 : endif
1598 0 : if (tr_bgc_DON) then
1599 0 : do mm = 1,n_don
1600 0 : write(ncharb,'(i3.3)') mm
1601 0 : if (my_task == master_task) write(nu_diag,*) subname,' min/max DON'
1602 0 : do k=1,nblyr+3
1603 0 : write(nchar,'(i3.3)') k
1604 : call read_restart_field(nu_restart_bgc,0, &
1605 : trcrn(:,:,nt_bgc_DON(mm)+k-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1606 0 : 'bgc_DON'//trim(ncharb)//trim(nchar),ncat,diag,field_loc_center,field_type_scalar)
1607 : enddo
1608 : enddo !mm
1609 : endif
1610 0 : if (tr_bgc_Fe) then
1611 0 : do mm = 1,n_fed
1612 0 : write(ncharb,'(i3.3)') mm
1613 0 : if (my_task == master_task) write(nu_diag,*) subname,' min/max dFe '
1614 0 : do k=1,nblyr+3
1615 0 : write(nchar,'(i3.3)') k
1616 : call read_restart_field(nu_restart_bgc,0, &
1617 : trcrn(:,:,nt_bgc_Fed (mm)+k-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1618 0 : 'bgc_Fed'//trim(ncharb)//trim(nchar),ncat,diag,field_loc_center,field_type_scalar)
1619 : enddo
1620 : enddo !mm
1621 0 : do mm = 1,n_fep
1622 0 : write(ncharb,'(i3.3)') mm
1623 0 : if (my_task == master_task) write(nu_diag,*) subname,' min/max pFe '
1624 0 : do k=1,nblyr+3
1625 0 : write(nchar,'(i3.3)') k
1626 : call read_restart_field(nu_restart_bgc,0, &
1627 : trcrn(:,:,nt_bgc_Fep (mm)+k-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1628 0 : 'bgc_Fep'//trim(ncharb)//trim(nchar),ncat,diag,field_loc_center,field_type_scalar)
1629 : enddo
1630 : enddo !mm
1631 : endif
1632 0 : if (tr_zaero) then
1633 0 : do mm = 1,n_zaero
1634 0 : write(ncharb,'(i3.3)') mm
1635 0 : if (my_task == master_task) write(nu_diag,*) subname,' min/max z aerosols'
1636 0 : do k=1,nblyr+3
1637 0 : write(nchar,'(i3.3)') k
1638 : call read_restart_field(nu_restart_bgc,0, &
1639 : trcrn(:,:,nt_zaero(mm)+k-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1640 0 : 'zaero'//trim(ncharb)//trim(nchar),ncat,diag,field_loc_center,field_type_scalar)
1641 : enddo
1642 : enddo !mm
1643 : endif
1644 0 : do mm = 1,nbtrcr
1645 0 : write(nchar,'(i3.3)') mm
1646 : call read_restart_field(nu_restart_bgc,0, &
1647 : trcrn(:,:,nt_zbgc_frac+mm-1,:,:),'ruf8', & ! LCOV_EXCL_LINE
1648 0 : 'zbgc_frac'//trim(nchar),ncat,diag)
1649 : enddo
1650 : endif
1651 :
1652 : !-----------------------------------------------------------------
1653 : ! Ocean BGC
1654 : !-----------------------------------------------------------------
1655 :
1656 0 : if (my_task == master_task) write(nu_diag,*) subname,'mixed layer ocean bgc restart'
1657 0 : if (tr_bgc_N) then
1658 0 : do k = 1,n_algae
1659 0 : write(nchar,'(i3.3)') k
1660 0 : call read_restart_field(nu_restart_bgc,0,algalN(:,:,k,:),'ruf8','algalN'//trim(nchar),1,diag)
1661 : enddo !k
1662 : endif
1663 0 : if (tr_bgc_C) then
1664 0 : do k = 1,n_doc
1665 0 : write(nchar,'(i3.3)') k
1666 0 : call read_restart_field(nu_restart_bgc,0,doc(:,:,k,:),'ruf8','doc'//trim(nchar),1,diag)
1667 : enddo !k
1668 0 : do k = 1,n_dic
1669 0 : write(nchar,'(i3.3)') k
1670 0 : call read_restart_field(nu_restart_bgc,0,dic(:,:,k,:),'ruf8','dic'//trim(nchar),1,diag)
1671 : enddo !k
1672 : endif !tr_bgc_C
1673 :
1674 0 : if (tr_bgc_Nit) &
1675 : call read_restart_field(nu_restart_bgc,0,nit ,'ruf8','nit' ,& ! LCOV_EXCL_LINE
1676 0 : 1,diag,field_loc_center,field_type_scalar)
1677 0 : if (tr_bgc_Am) &
1678 : call read_restart_field(nu_restart_bgc,0,amm ,'ruf8','amm' ,& ! LCOV_EXCL_LINE
1679 0 : 1,diag,field_loc_center,field_type_scalar)
1680 0 : if (tr_bgc_Sil) &
1681 : call read_restart_field(nu_restart_bgc,0,sil ,'ruf8','sil' ,& ! LCOV_EXCL_LINE
1682 0 : 1,diag,field_loc_center,field_type_scalar)
1683 0 : if (tr_bgc_hum) &
1684 : call read_restart_field(nu_restart_bgc,0,hum ,'ruf8','hum' ,& ! LCOV_EXCL_LINE
1685 0 : 1,diag,field_loc_center,field_type_scalar)
1686 0 : if (tr_bgc_DMS) then
1687 0 : call read_restart_field(nu_restart_bgc,0,dmsp ,'ruf8','dmsp' ,1,diag)
1688 0 : call read_restart_field(nu_restart_bgc,0,dms ,'ruf8','dms' ,1,diag)
1689 : endif
1690 0 : if (tr_bgc_DON) then
1691 0 : do k = 1,n_don
1692 0 : write(nchar,'(i3.3)') k
1693 0 : call read_restart_field(nu_restart_bgc,0,don(:,:,k,:),'ruf8','don'//trim(nchar),1,diag)
1694 : enddo !k
1695 : endif
1696 0 : if (tr_bgc_Fe ) then
1697 0 : do k = 1,n_fed
1698 0 : write(nchar,'(i3.3)') k
1699 0 : call read_restart_field(nu_restart_bgc,0,fed(:,:,k,:),'ruf8','fed'//trim(nchar),1,diag)
1700 : enddo !k
1701 0 : do k = 1,n_fep
1702 0 : write(nchar,'(i3.3)') k
1703 0 : call read_restart_field(nu_restart_bgc,0,fep(:,:,k,:),'ruf8','fep'//trim(nchar),1,diag)
1704 : enddo !k
1705 : endif
1706 0 : if (tr_zaero) then
1707 0 : do k = 1,n_zaero
1708 0 : write(nchar,'(i3.3)') k
1709 0 : call read_restart_field(nu_restart_bgc,0,zaeros(:,:,k,:),'ruf8','zaeros'//trim(nchar),1,diag)
1710 : enddo !k
1711 : endif
1712 : endif ! restart_bgc
1713 :
1714 0 : end subroutine read_restart_bgc
1715 :
1716 : !=======================================================================
1717 :
1718 : end module ice_restart_column
1719 :
1720 : !=======================================================================
|