Line data Source code
1 : !=========================================================================
2 : !
3 : ! Restart routines for the column package.
4 : !
5 : ! author: Elizabeth C. Hunke, LANL
6 : !
7 : module icedrv_restart_bgc
8 :
9 : use icedrv_kinds
10 : use icedrv_constants
11 : use icedrv_domain_size, only: ncat, nblyr, nx
12 : use icepack_intfc, only: icepack_max_algae, icepack_max_doc, icepack_max_don
13 : use icepack_intfc, only: icepack_max_dic, icepack_max_aero, icepack_max_fe
14 : use icepack_intfc, only: icepack_warnings_flush, icepack_warnings_aborted
15 : use icepack_intfc, only: icepack_query_parameters, icepack_query_tracer_sizes
16 : use icepack_intfc, only: icepack_query_tracer_flags, icepack_query_tracer_indices
17 : use icedrv_system, only: icedrv_system_abort
18 :
19 : implicit none
20 :
21 : private
22 : public :: write_restart_bgc, &
23 : read_restart_bgc
24 :
25 : !=======================================================================
26 :
27 : contains
28 :
29 : !=======================================================================
30 : !
31 : ! Dumps all values needed for a bgc restart
32 : !
33 : ! author Elizabeth C. Hunke, LANL
34 :
35 25 : subroutine write_restart_bgc()
36 :
37 : use icedrv_arrays_column, only: Rayleigh_criteria, Rayleigh_real
38 : use icedrv_domain_size, only: n_algae, n_doc, n_dic
39 : use icedrv_domain_size, only: n_don, n_zaero, n_fed, n_fep
40 : use icedrv_flux, only: sss, nit, amm, sil, dmsp, dms, algalN
41 : use icedrv_flux, only: doc, don, dic, fed, fep, zaeros, hum
42 : use icedrv_state, only: trcrn
43 : use icedrv_restart, only: write_restart_field
44 :
45 : ! local variables
46 :
47 : integer (kind=int_kind) :: &
48 : i, k , & ! horizontal, vertical indices
49 : mm ! n_algae
50 :
51 : logical (kind=log_kind) :: skl_bgc, solve_zsal, z_tracers
52 :
53 : integer (kind=int_kind) :: nbtrcr
54 : logical (kind=log_kind) :: tr_bgc_Nit, tr_bgc_Am, tr_bgc_Sil, tr_bgc_hum
55 : logical (kind=log_kind) :: tr_bgc_DMS, tr_bgc_PON, tr_bgc_N, tr_bgc_C
56 : logical (kind=log_kind) :: tr_bgc_DON, tr_bgc_Fe, tr_zaero , tr_bgc_chl
57 : integer (kind=int_kind) :: nt_bgc_S, nt_bgc_Nit, nt_bgc_AM, nt_bgc_Sil
58 : integer (kind=int_kind) :: nt_bgc_hum, nt_bgc_PON
59 : integer (kind=int_kind) :: nt_bgc_DMSPp, nt_bgc_DMSPd, nt_bgc_DMS
60 : integer (kind=int_kind) :: nt_zbgc_frac
61 : integer (kind=int_kind), dimension(icepack_max_algae) :: nt_bgc_N
62 : integer (kind=int_kind), dimension(icepack_max_algae) :: nt_bgc_chl
63 : integer (kind=int_kind), dimension(icepack_max_algae) :: nt_bgc_C
64 : integer (kind=int_kind), dimension(icepack_max_doc) :: nt_bgc_DOC
65 : integer (kind=int_kind), dimension(icepack_max_don) :: nt_bgc_DON
66 : integer (kind=int_kind), dimension(icepack_max_dic) :: nt_bgc_DIC
67 : integer (kind=int_kind), dimension(icepack_max_aero) :: nt_zaero
68 : integer (kind=int_kind), dimension(icepack_max_fe) :: nt_bgc_Fed
69 : integer (kind=int_kind), dimension(icepack_max_fe) :: nt_bgc_Fep
70 :
71 : character(len=*), parameter :: subname='(write_restart_bgc)'
72 :
73 : !-----------------------------------------------------------------
74 : ! query Icepack values
75 : !-----------------------------------------------------------------
76 :
77 25 : call icepack_query_parameters(skl_bgc_out=skl_bgc)
78 25 : call icepack_query_parameters(solve_zsal_out=solve_zsal)
79 25 : call icepack_query_parameters(z_tracers_out=z_tracers)
80 25 : call icepack_query_tracer_sizes(nbtrcr_out=nbtrcr)
81 : call icepack_query_tracer_flags(tr_bgc_Nit_out=tr_bgc_Nit, &
82 : tr_bgc_Am_out=tr_bgc_Am, &
83 : tr_bgc_Sil_out=tr_bgc_Sil, tr_bgc_hum_out=tr_bgc_hum, &
84 : tr_bgc_DMS_out=tr_bgc_DMS, tr_bgc_PON_out=tr_bgc_PON, &
85 : tr_bgc_N_out=tr_bgc_N, tr_bgc_C_out=tr_bgc_C, &
86 : tr_bgc_DON_out=tr_bgc_DON, tr_bgc_Fe_out=tr_bgc_Fe, &
87 25 : tr_zaero_out=tr_zaero, tr_bgc_chl_out=tr_bgc_chl)
88 : call icepack_query_tracer_indices( nt_bgc_S_out=nt_bgc_S, &
89 : nt_bgc_N_out=nt_bgc_N, nt_bgc_AM_out=nt_bgc_AM, &
90 : nt_bgc_chl_out=nt_bgc_chl, nt_bgc_C_out=nt_bgc_C, &
91 : nt_bgc_DOC_out=nt_bgc_DOC, &
92 : nt_bgc_DIC_out=nt_bgc_DIC, nt_bgc_Nit_out=nt_bgc_Nit, &
93 : nt_bgc_Sil_out=nt_bgc_Sil, nt_bgc_hum_out=nt_bgc_hum, &
94 : nt_bgc_DMSPp_out=nt_bgc_DMSPp, nt_bgc_DMSPd_out=nt_bgc_DMSPd, &
95 : nt_bgc_DMS_out=nt_bgc_DMS, nt_bgc_PON_out=nt_bgc_PON, &
96 : nt_bgc_DON_out=nt_bgc_DON, nt_bgc_Fed_out=nt_bgc_Fed, &
97 : nt_bgc_Fep_out=nt_bgc_Fep, nt_zbgc_frac_out=nt_zbgc_frac, &
98 25 : nt_zaero_out=nt_zaero)
99 25 : call icepack_warnings_flush(nu_diag)
100 25 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
101 0 : file=__FILE__,line= __LINE__)
102 :
103 : !-----------------------------------------------------------------
104 : ! Salinity and extras
105 : !-----------------------------------------------------------------
106 25 : if (solve_zsal) then
107 :
108 0 : do k = 1, nblyr
109 0 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_S+k-1,:),ncat)
110 : enddo
111 :
112 0 : call write_restart_field(nu_dump,sss,1)
113 :
114 0 : do i = 1, nx
115 0 : if (Rayleigh_criteria(i)) then
116 0 : Rayleigh_real (i) = c1
117 0 : elseif (.NOT. Rayleigh_criteria(i)) then
118 0 : Rayleigh_real (i) = c0
119 : endif
120 : enddo
121 :
122 0 : call write_restart_field(nu_dump,Rayleigh_real,1)
123 :
124 : endif ! solve_zsal
125 :
126 : !-----------------------------------------------------------------
127 : ! Skeletal layer BGC
128 : !-----------------------------------------------------------------
129 :
130 25 : if (skl_bgc) then
131 4 : do k = 1, n_algae
132 3 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_N(k),:),ncat)
133 3 : if (tr_bgc_chl) &
134 1 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_chl(k),:),ncat)
135 : enddo
136 1 : if (tr_bgc_C) then
137 : ! do k = 1, n_algae
138 : ! call write_restart_field(nu_dump,trcrn(:,nt_bgc_C(k),:),ncat)
139 : ! enddo
140 3 : do k = 1, n_doc
141 3 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_DOC(k),:),ncat)
142 : enddo
143 1 : do k = 1, n_dic
144 1 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_DIC(k),:),ncat)
145 : enddo
146 : endif
147 1 : if (tr_bgc_Nit) &
148 1 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_Nit,:),ncat)
149 1 : if (tr_bgc_Am) &
150 1 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_Am,:),ncat)
151 1 : if (tr_bgc_Sil) &
152 1 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_Sil,:),ncat)
153 1 : if (tr_bgc_hum) &
154 1 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_hum,:),ncat)
155 1 : if (tr_bgc_DMS) then
156 1 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_DMSPp,:),ncat)
157 1 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_DMSPd,:),ncat)
158 1 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_DMS,:),ncat)
159 : endif
160 1 : if (tr_bgc_PON) &
161 1 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_PON,:),ncat)
162 1 : if (tr_bgc_DON) then
163 2 : do k = 1, n_don
164 2 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_DON(k),:),ncat)
165 : enddo
166 : endif
167 1 : if (tr_bgc_Fe ) then
168 2 : do k = 1, n_fed
169 2 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_Fed(k),:),ncat)
170 : enddo
171 2 : do k = 1, n_fep
172 2 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_Fep(k),:),ncat)
173 : enddo
174 : endif
175 :
176 24 : elseif (z_tracers) then
177 :
178 : !-----------------------------------------------------------------
179 : ! Z layer BGC
180 : !-----------------------------------------------------------------
181 :
182 24 : if (tr_bgc_Nit) then
183 24 : write(nu_diag,*) 'z bgc restart: min/max Nitrate'
184 264 : do k = 1, nblyr+3
185 264 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_Nit+k-1,:),ncat)
186 : enddo
187 : endif
188 24 : if (tr_bgc_N) then
189 96 : do mm = 1, n_algae
190 792 : do k = 1, nblyr+3
191 792 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_N(mm)+k-1,:),ncat)
192 : enddo
193 96 : if (tr_bgc_chl) then
194 0 : do k = 1, nblyr+3
195 0 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_chl(mm)+k-1,:),ncat)
196 : enddo
197 : endif
198 : enddo ! n_algae
199 : endif ! tr_bgc_N
200 24 : if (tr_bgc_C) then
201 : ! algal C is not yet distinct from algal N
202 : ! do mm = 1, n_algae
203 : ! do k = 1, nblyr+3
204 : ! call write_restart_field(nu_dump,trcrn(:,nt_bgc_C(mm)+k-1,:),ncat)
205 : ! enddo
206 : ! enddo
207 72 : do mm = 1, n_doc
208 552 : do k = 1, nblyr+3
209 528 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_DOC(mm)+k-1,:),ncat)
210 : enddo
211 : enddo
212 24 : do mm = 1, n_dic
213 24 : do k = 1, nblyr+3
214 0 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_DIC(mm)+k-1,:),ncat)
215 : enddo
216 : enddo
217 : endif !tr_bgc_C
218 24 : if (tr_bgc_Am) then
219 264 : do k = 1, nblyr+3
220 264 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_Am+k-1,:),ncat)
221 : enddo
222 : endif
223 24 : if (tr_bgc_Sil) then
224 264 : do k = 1, nblyr+3
225 264 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_Sil+k-1,:),ncat)
226 : enddo
227 : endif
228 24 : if (tr_bgc_hum) then
229 264 : do k = 1, nblyr+3
230 264 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_hum+k-1,:),ncat)
231 : enddo
232 : endif
233 24 : if (tr_bgc_DMS) then
234 264 : do k = 1, nblyr+3
235 240 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_DMSPp+k-1,:),ncat)
236 240 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_DMSPd+k-1,:),ncat)
237 264 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_DMS+k-1,:),ncat)
238 : enddo
239 : endif
240 24 : if (tr_bgc_PON) then
241 264 : do k = 1, nblyr+3
242 264 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_PON+k-1,:),ncat)
243 : enddo
244 : endif
245 24 : if (tr_bgc_DON) then
246 48 : do mm = 1, n_don
247 288 : do k = 1, nblyr+3
248 264 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_DON(mm)+k-1,:),ncat)
249 : enddo
250 : enddo
251 : endif
252 24 : if (tr_bgc_Fe ) then
253 2 : do mm = 1, n_fed
254 12 : do k = 1, nblyr+3
255 11 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_Fed(mm)+k-1,:),ncat)
256 : enddo
257 : enddo
258 2 : do mm = 1, n_fep
259 12 : do k = 1, nblyr+3
260 11 : call write_restart_field(nu_dump,trcrn(:,nt_bgc_Fep(mm)+k-1,:),ncat)
261 : enddo
262 : enddo
263 : endif
264 24 : if (tr_zaero) then
265 161 : do mm = 1, n_zaero
266 1541 : do k = 1, nblyr+3
267 1518 : call write_restart_field(nu_dump,trcrn(:,nt_zaero(mm)+k-1,:),ncat)
268 : enddo
269 : enddo
270 : endif
271 24 : if (nbtrcr > 0) then
272 500 : do mm = 1, nbtrcr
273 500 : call write_restart_field(nu_dump,trcrn(:,nt_zbgc_frac+mm-1,:),ncat)
274 : enddo
275 : endif
276 :
277 : !-----------------------------------------------------------------
278 : ! Ocean BGC
279 : !-----------------------------------------------------------------
280 :
281 24 : if (tr_bgc_N) then
282 96 : do k = 1, n_algae
283 96 : call write_restart_field(nu_dump,algalN(:,k),1)
284 : enddo ! k
285 : endif
286 24 : if (tr_bgc_C) then
287 72 : do k = 1, n_doc
288 72 : call write_restart_field(nu_dump,doc(:,k),1)
289 : enddo ! k
290 24 : do k = 1, n_dic
291 24 : call write_restart_field(nu_dump,dic(:,k),1)
292 : enddo ! k
293 : endif
294 24 : if (tr_bgc_Nit) &
295 24 : call write_restart_field(nu_dump,nit,1)
296 24 : if (tr_bgc_Am) &
297 24 : call write_restart_field(nu_dump,amm,1)
298 24 : if (tr_bgc_Sil) &
299 24 : call write_restart_field(nu_dump,sil,1)
300 24 : if (tr_bgc_hum) &
301 24 : call write_restart_field(nu_dump,hum,1)
302 24 : if (tr_bgc_DMS) then
303 24 : call write_restart_field(nu_dump,dmsp,1)
304 24 : call write_restart_field(nu_dump,dms,1)
305 : endif
306 24 : if (tr_bgc_DON) then
307 48 : do k = 1, n_don
308 48 : call write_restart_field(nu_dump,don(:,k),1)
309 : enddo ! k
310 : endif
311 24 : if (tr_bgc_Fe ) then
312 2 : do k = 1, n_fed
313 2 : call write_restart_field(nu_dump,fed(:,k),1)
314 : enddo ! k
315 2 : do k = 1, n_fep
316 2 : call write_restart_field(nu_dump,fep(:,k),1)
317 : enddo ! k
318 : endif
319 24 : if (tr_zaero) then
320 161 : do k = 1, n_zaero
321 161 : call write_restart_field(nu_dump,zaeros(:,k),1)
322 : enddo ! k
323 : endif
324 :
325 : endif ! skl_bgc or z_tracers
326 :
327 25 : end subroutine write_restart_bgc
328 :
329 : !=======================================================================
330 : !
331 : ! Reads all values needed for a bgc restart
332 : !
333 : ! author Elizabeth C. Hunke, LANL
334 :
335 1 : subroutine read_restart_bgc()
336 :
337 : use icedrv_arrays_column, only: Rayleigh_real, Rayleigh_criteria
338 : use icedrv_domain_size, only: n_algae, n_doc, n_dic
339 : use icedrv_domain_size, only: n_don, n_zaero, n_fed, n_fep
340 : use icedrv_flux, only: sss, nit, amm, sil, dmsp, dms, algalN
341 : use icedrv_flux, only: doc, don, dic, fed, fep, zaeros, hum
342 : use icedrv_state, only: trcrn
343 : use icedrv_restart, only: read_restart_field
344 :
345 : ! local variables
346 :
347 : integer (kind=int_kind) :: &
348 : i, k, & ! indices
349 : mm ! n_algae
350 :
351 : logical (kind=log_kind) :: skl_bgc, solve_zsal, z_tracers
352 :
353 : integer (kind=int_kind) :: nbtrcr
354 : logical (kind=log_kind) :: tr_bgc_Nit, tr_bgc_Am, tr_bgc_Sil, tr_bgc_hum
355 : logical (kind=log_kind) :: tr_bgc_DMS, tr_bgc_PON, tr_bgc_N, tr_bgc_C
356 : logical (kind=log_kind) :: tr_bgc_DON, tr_bgc_Fe, tr_zaero , tr_bgc_chl
357 : integer (kind=int_kind) :: nt_bgc_S, nt_bgc_Nit, nt_bgc_AM, nt_bgc_Sil
358 : integer (kind=int_kind) :: nt_bgc_hum, nt_bgc_PON
359 : integer (kind=int_kind) :: nt_bgc_DMSPp, nt_bgc_DMSPd, nt_bgc_DMS
360 : integer (kind=int_kind) :: nt_zbgc_frac
361 : integer (kind=int_kind), dimension(icepack_max_algae) :: nt_bgc_N
362 : integer (kind=int_kind), dimension(icepack_max_algae) :: nt_bgc_chl
363 : integer (kind=int_kind), dimension(icepack_max_algae) :: nt_bgc_C
364 : integer (kind=int_kind), dimension(icepack_max_doc) :: nt_bgc_DOC
365 : integer (kind=int_kind), dimension(icepack_max_don) :: nt_bgc_DON
366 : integer (kind=int_kind), dimension(icepack_max_dic) :: nt_bgc_DIC
367 : integer (kind=int_kind), dimension(icepack_max_aero) :: nt_zaero
368 : integer (kind=int_kind), dimension(icepack_max_fe) :: nt_bgc_Fed
369 : integer (kind=int_kind), dimension(icepack_max_fe) :: nt_bgc_Fep
370 :
371 : character(len=*), parameter :: subname='(read_restart_bgc)'
372 :
373 : !-----------------------------------------------------------------
374 : ! query Icepack values
375 : !-----------------------------------------------------------------
376 :
377 1 : call icepack_query_parameters(skl_bgc_out=skl_bgc)
378 1 : call icepack_query_parameters(solve_zsal_out=solve_zsal)
379 1 : call icepack_query_parameters(z_tracers_out=z_tracers)
380 1 : call icepack_query_tracer_sizes(nbtrcr_out=nbtrcr)
381 :
382 :
383 : call icepack_query_tracer_flags(tr_bgc_Nit_out=tr_bgc_Nit, &
384 : tr_bgc_Am_out=tr_bgc_Am, &
385 : tr_bgc_Sil_out=tr_bgc_Sil, tr_bgc_hum_out=tr_bgc_hum, &
386 : tr_bgc_DMS_out=tr_bgc_DMS, tr_bgc_PON_out=tr_bgc_PON, &
387 : tr_bgc_N_out=tr_bgc_N, tr_bgc_C_out=tr_bgc_C, &
388 : tr_bgc_DON_out=tr_bgc_DON, tr_bgc_Fe_out=tr_bgc_Fe, &
389 1 : tr_zaero_out=tr_zaero, tr_bgc_chl_out=tr_bgc_chl)
390 : call icepack_query_tracer_indices( nt_bgc_S_out=nt_bgc_S, &
391 : nt_bgc_N_out=nt_bgc_N, nt_bgc_AM_out=nt_bgc_AM, &
392 : nt_bgc_chl_out=nt_bgc_chl, nt_bgc_C_out=nt_bgc_C, &
393 : nt_bgc_DOC_out=nt_bgc_DOC, &
394 : nt_bgc_DIC_out=nt_bgc_DIC, nt_bgc_Nit_out=nt_bgc_Nit, &
395 : nt_bgc_Sil_out=nt_bgc_Sil, nt_bgc_hum_out=nt_bgc_hum, &
396 : nt_bgc_DMSPp_out=nt_bgc_DMSPp, nt_bgc_DMSPd_out=nt_bgc_DMSPd, &
397 : nt_bgc_DMS_out=nt_bgc_DMS, nt_bgc_PON_out=nt_bgc_PON, &
398 : nt_bgc_DON_out=nt_bgc_DON, nt_bgc_Fed_out=nt_bgc_Fed, &
399 : nt_bgc_Fep_out=nt_bgc_Fep, nt_zbgc_frac_out=nt_zbgc_frac, &
400 1 : nt_zaero_out=nt_zaero)
401 1 : call icepack_warnings_flush(nu_diag)
402 1 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
403 0 : file=__FILE__,line= __LINE__)
404 :
405 : !-----------------------------------------------------------------
406 : ! Salinity and extras
407 : !-----------------------------------------------------------------
408 :
409 1 : if (solve_zsal) then
410 :
411 0 : write(nu_diag,*)'zSalinity restart'
412 0 : do k = 1, nblyr
413 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_S+k-1,:),ncat)
414 : enddo
415 :
416 0 : write(nu_diag,*) 'min/max sea surface salinity'
417 0 : call read_restart_field(nu_restart,sss,1)
418 0 : write(nu_diag,*) 'min/max Rayleigh'
419 0 : call read_restart_field(nu_restart,Rayleigh_real,1)
420 :
421 0 : do i = 1, nx
422 0 : if (Rayleigh_real (i) .GE. c1) then
423 0 : Rayleigh_criteria (i) = .true.
424 0 : elseif (Rayleigh_real (i) < c1) then
425 0 : Rayleigh_criteria (i) = .false.
426 : endif
427 : enddo
428 :
429 : endif ! solve_zsal
430 :
431 : !-----------------------------------------------------------------
432 : ! Skeletal Layer BGC
433 : !-----------------------------------------------------------------
434 :
435 1 : if (skl_bgc) then
436 0 : write(nu_diag,*) 'skl bgc restart'
437 :
438 0 : write (nu_diag,*) 'min/max algal N, chl'
439 0 : do k = 1, n_algae
440 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_N(k),:),ncat)
441 0 : if (tr_bgc_chl) &
442 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_chl(k),:),ncat)
443 : enddo ! k
444 0 : if (tr_bgc_C) then
445 : ! algal C is not yet distinct from algal N
446 : ! write (nu_diag,*) 'min/max algal C, DOC, DIC'
447 : ! do k = 1, n_algae
448 : ! call read_restart_field(nu_restart,trcrn(:,nt_bgc_C(k),:),ncat)
449 : ! enddo
450 0 : write (nu_diag,*) 'min/max DOC, DIC'
451 0 : do k = 1, n_doc
452 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_DOC(k),:),ncat)
453 : enddo
454 0 : do k = 1, n_dic
455 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_DIC(k),:),ncat)
456 : enddo
457 : endif
458 0 : write (nu_diag,*) 'min/max Nit, Am, Sil, hum'
459 0 : if (tr_bgc_Nit) &
460 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_Nit,:),ncat)
461 0 : if (tr_bgc_Am) &
462 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_Am,:),ncat)
463 0 : if (tr_bgc_Sil) &
464 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_Sil,:),ncat)
465 0 : if (tr_bgc_hum) &
466 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_hum,:),ncat)
467 0 : if (tr_bgc_DMS) then
468 0 : write (nu_diag,*) 'min/max pDMSP, dDMSP, DMS'
469 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_DMSPp,:),ncat)
470 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_DMSPd,:),ncat)
471 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_DMS,:),ncat)
472 : endif
473 0 : write (nu_diag,*) 'min/max PON'
474 0 : if (tr_bgc_PON) &
475 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_PON,:),ncat)
476 0 : if (tr_bgc_DON) then
477 0 : write (nu_diag,*) 'min/max DON'
478 0 : do k = 1, n_don
479 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_DON(k),:),ncat)
480 : enddo
481 : endif
482 0 : if (tr_bgc_Fe) then
483 0 : write (nu_diag,*) 'min/max dFe, pFe'
484 0 : do k = 1, n_fed
485 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_Fed (k),:),ncat)
486 : enddo
487 0 : do k = 1, n_fep
488 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_Fep (k),:),ncat)
489 : enddo
490 : endif
491 :
492 1 : elseif (z_tracers) then
493 :
494 : !-----------------------------------------------------------------
495 : ! Z Layer BGC
496 : !-----------------------------------------------------------------
497 :
498 1 : if (tr_bgc_Nit) then
499 1 : write(nu_diag,*) 'z bgc restart: min/max Nitrate'
500 11 : do k=1, nblyr+3
501 11 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_Nit+k-1,:),ncat)
502 : enddo
503 : endif ! Nit
504 1 : if (tr_bgc_N) then
505 4 : do mm = 1, n_algae
506 3 : write(nu_diag,*) ' min/max Algal N', n_algae
507 33 : do k=1, nblyr+3
508 33 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_N(mm)+k-1,:),ncat)
509 : enddo
510 4 : if (tr_bgc_chl) then
511 0 : write(nu_diag,*) ' min/max Algal chla'
512 0 : do k=1, nblyr+3
513 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_chl(mm)+k-1,:),ncat)
514 : enddo
515 : endif ! tr_bgc_chl
516 : enddo ! n_algae
517 : endif ! tr_bgc_N
518 1 : if (tr_bgc_C) then
519 : ! algal C is not yet distinct from algal N
520 : ! do mm = 1, n_algae
521 : ! write(nu_diag,*) ' min/max algalC'
522 : ! do k=1, nblyr+3
523 : ! call read_restart_field(nu_restart,trcrn(:,nt_bgc_C(mm)+k-1,:),ncat)
524 : ! enddo
525 : ! enddo !mm
526 3 : do mm = 1, n_doc
527 2 : write(nu_diag,*) ' min/max DOC'
528 23 : do k=1, nblyr+3
529 22 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_DOC(mm)+k-1,:),ncat)
530 : enddo
531 : enddo !mm
532 1 : do mm = 1, n_dic
533 0 : write(nu_diag,*) ' min/max DIC'
534 1 : do k=1, nblyr+3
535 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_DIC(mm)+k-1,:),ncat)
536 : enddo
537 : enddo !mm
538 : endif ! tr_bgc_C
539 1 : if (tr_bgc_Am) then
540 1 : write(nu_diag,*) ' min/max ammonium'
541 11 : do k=1, nblyr+3
542 11 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_Am+k-1,:),ncat)
543 : enddo
544 : endif
545 1 : if (tr_bgc_Sil) then
546 1 : write(nu_diag,*) ' min/max silicate'
547 11 : do k=1, nblyr+3
548 11 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_Sil+k-1,:),ncat)
549 : enddo
550 : endif
551 1 : if (tr_bgc_hum) then
552 1 : write(nu_diag,*) ' min/max humic material'
553 11 : do k=1, nblyr+3
554 11 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_hum+k-1,:),ncat)
555 : enddo
556 : endif
557 1 : if (tr_bgc_DMS) then
558 1 : write (nu_diag,*) 'min/max pDMSP, dDMSP, DMS'
559 11 : do k=1, nblyr+3
560 10 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_DMSPp+k-1,:),ncat)
561 10 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_DMSPd+k-1,:),ncat)
562 11 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_DMS+k-1,:),ncat)
563 : enddo
564 : endif
565 1 : if (tr_bgc_PON) then
566 1 : write(nu_diag,*) ' min/max PON'
567 11 : do k=1, nblyr+3
568 11 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_PON+k-1,:),ncat)
569 : enddo
570 : endif
571 1 : if (tr_bgc_DON) then
572 2 : do mm = 1, n_don
573 1 : write(nu_diag,*) ' min/max DON'
574 12 : do k=1, nblyr+3
575 11 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_DON(mm)+k-1,:),ncat)
576 : enddo
577 : enddo !mm
578 : endif
579 1 : if (tr_bgc_Fe) then
580 0 : do mm = 1, n_fed
581 0 : write(nu_diag,*) ' min/max dFe '
582 0 : do k=1, nblyr+3
583 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_Fed (mm)+k-1,:),ncat)
584 : enddo
585 : enddo !mm
586 0 : do mm = 1, n_fep
587 0 : write(nu_diag,*) ' min/max pFe '
588 0 : do k=1, nblyr+3
589 0 : call read_restart_field(nu_restart,trcrn(:,nt_bgc_Fep (mm)+k-1,:),ncat)
590 : enddo
591 : enddo !mm
592 : endif
593 1 : if (tr_zaero) then
594 7 : do mm = 1, n_zaero
595 6 : write(nu_diag,*) ' min/max z aerosols'
596 67 : do k=1, nblyr+3
597 66 : call read_restart_field(nu_restart,trcrn(:,nt_zaero(mm)+k-1,:),ncat)
598 : enddo
599 : enddo !mm
600 : endif
601 1 : if (nbtrcr > 0) then
602 1 : write (nu_diag,*) 'min/max zbgc_frac'
603 21 : do mm = 1, nbtrcr
604 21 : call read_restart_field(nu_restart,trcrn(:,nt_zbgc_frac+mm-1,:),ncat)
605 : enddo
606 : endif
607 :
608 : !-----------------------------------------------------------------
609 : ! Ocean BGC
610 : !-----------------------------------------------------------------
611 :
612 1 : write(nu_diag,*) 'mixed layer ocean bgc restart'
613 1 : if (tr_bgc_N) then
614 1 : write (nu_diag,*) 'min/max algalN'
615 4 : do k = 1, n_algae
616 4 : call read_restart_field(nu_restart,algalN(:,k),1)
617 : enddo ! k
618 : endif
619 1 : if (tr_bgc_C) then
620 1 : write (nu_diag,*) 'min/max DOC, DIC'
621 3 : do k = 1, n_doc
622 3 : call read_restart_field(nu_restart,doc(:,k),1)
623 : enddo ! k
624 1 : do k = 1, n_dic
625 1 : call read_restart_field(nu_restart,dic(:,k),1)
626 : enddo ! k
627 : endif !tr_bgc_C
628 :
629 1 : write (nu_diag,*) 'min/max Nit, Am, Sil, hum, DMSP, DMS'
630 1 : if (tr_bgc_Nit) &
631 1 : call read_restart_field(nu_restart,nit,1)
632 1 : if (tr_bgc_Am) &
633 1 : call read_restart_field(nu_restart,amm ,1)
634 1 : if (tr_bgc_Sil) &
635 1 : call read_restart_field(nu_restart,sil,1)
636 1 : if (tr_bgc_hum) &
637 1 : call read_restart_field(nu_restart,hum,1)
638 1 : if (tr_bgc_DMS) then
639 1 : call read_restart_field(nu_restart,dmsp,1)
640 1 : call read_restart_field(nu_restart,dms,1)
641 : endif
642 1 : if (tr_bgc_DON) then
643 1 : write (nu_diag,*) 'min/max DON'
644 2 : do k = 1, n_don
645 2 : call read_restart_field(nu_restart,don(:,k),1)
646 : enddo ! k
647 : endif
648 1 : if (tr_bgc_Fe ) then
649 0 : write (nu_diag,*) 'min/max dFe, pFe'
650 0 : do k = 1, n_fed
651 0 : call read_restart_field(nu_restart,fed(:,k),1)
652 : enddo ! k
653 0 : do k = 1, n_fep
654 0 : call read_restart_field(nu_restart,fep(:,k),1)
655 : enddo ! k
656 : endif
657 1 : if (tr_zaero) then
658 1 : write (nu_diag,*) 'min/max zaeros'
659 7 : do k = 1, n_zaero
660 7 : call read_restart_field(nu_restart,zaeros(:,k),1)
661 : enddo ! k
662 : endif
663 :
664 : endif ! skl_bgc or z_tracers
665 :
666 1 : end subroutine read_restart_bgc
667 :
668 : !=======================================================================
669 :
670 : end module icedrv_restart_bgc
671 :
672 : !=======================================================================
|