Line data Source code
1 : !=======================================================================
2 : !
3 : ! Biogeochemistry driver
4 : !
5 : ! authors: Nicole Jeffery, LANL
6 : ! Scott Elliot, LANL
7 : ! Elizabeth C. Hunke, LANL
8 : !
9 : module icepack_zbgc
10 :
11 : use icepack_kinds
12 : use icepack_parameters, only: c0, c1, c2, p001, p1, p5, puny
13 : use icepack_parameters, only: depressT, rhosi, min_salin, salt_loss
14 : use icepack_parameters, only: fr_resp, algal_vel, R_dFe2dust, dustFe_sol, T_max
15 : use icepack_parameters, only: op_dep_min, fr_graze_s, fr_graze_e, fr_mort2min, fr_dFe
16 : use icepack_parameters, only: k_nitrif, t_iron_conv, max_loss, max_dfe_doc1
17 : use icepack_parameters, only: fr_resp_s, y_sk_DMS, t_sk_conv, t_sk_ox
18 : use icepack_parameters, only: scale_bgc, ktherm, skl_bgc
19 : use icepack_parameters, only: z_tracers, fsal, conserv_check
20 :
21 : use icepack_tracers, only: nt_sice, bio_index
22 : use icepack_tracers, only: tr_brine, nt_fbri, nt_qice, nt_Tsfc
23 : use icepack_tracers, only: nt_zbgc_frac
24 : use icepack_tracers, only: bio_index_o, bio_index
25 :
26 : use icepack_zbgc_shared, only: zbgc_init_frac
27 : use icepack_zbgc_shared, only: zbgc_frac_init
28 : use icepack_zbgc_shared, only: bgc_tracer_type, remap_zbgc
29 : use icepack_zbgc_shared, only: regrid_stationary
30 : use icepack_zbgc_shared, only: R_S2N, R_Si2N, R_Fe2C, R_Fe2N, R_Fe2DON, R_Fe2DOC
31 : use icepack_zbgc_shared, only: chlabs, alpha2max_low, beta2max
32 : use icepack_zbgc_shared, only: mu_max, grow_Tdep, fr_graze
33 : use icepack_zbgc_shared, only: mort_pre, mort_Tdep, k_exude
34 : use icepack_zbgc_shared, only: K_Nit, K_Am, K_Sil, K_Fe
35 : use icepack_zbgc_shared, only: f_don, kn_bac, f_don_Am, f_doc
36 : use icepack_zbgc_shared, only: f_exude, k_bac
37 : use icepack_zbgc_shared, only: tau_ret, tau_rel
38 : use icepack_zbgc_shared, only: R_C2N, R_CHL2N, f_abs_chl, R_C2N_DON
39 :
40 : use icepack_warnings, only: warnstr, icepack_warnings_add
41 : use icepack_warnings, only: icepack_warnings_setabort, icepack_warnings_aborted
42 :
43 : use icepack_brine, only: preflushing_changes, compute_microS_mushy
44 : use icepack_brine, only: update_hbrine
45 : use icepack_algae, only: zbio, sklbio
46 : use icepack_therm_shared, only: calculate_Tin_from_qin
47 : use icepack_itd, only: column_sum, column_conservation_check
48 :
49 : implicit none
50 :
51 : private
52 : public :: add_new_ice_bgc, &
53 : lateral_melt_bgc, & ! LCOV_EXCL_LINE
54 : icepack_init_bgc, & ! LCOV_EXCL_LINE
55 : icepack_init_zbgc, & ! LCOV_EXCL_LINE
56 : icepack_biogeochemistry, & ! LCOV_EXCL_LINE
57 : icepack_load_ocean_bio_array, & ! LCOV_EXCL_LINE
58 : icepack_init_ocean_bio
59 :
60 : !=======================================================================
61 :
62 : contains
63 :
64 : !=======================================================================
65 :
66 : ! Adjust biogeochemical tracers when new frazil ice forms
67 :
68 0 : subroutine add_new_ice_bgc (dt, nblyr, &
69 : ncat, nilyr, nltrcr, & ! LCOV_EXCL_LINE
70 : bgrid, cgrid, igrid, & ! LCOV_EXCL_LINE
71 : aicen_init, vicen_init, vi0_init, & ! LCOV_EXCL_LINE
72 : aicen, vicen, vsnon1, & ! LCOV_EXCL_LINE
73 : vi0new, & ! LCOV_EXCL_LINE
74 : ntrcr, trcrn, nbtrcr, & ! LCOV_EXCL_LINE
75 : sss, ocean_bio, flux_bio, & ! LCOV_EXCL_LINE
76 : hsurp)
77 :
78 : integer (kind=int_kind), intent(in) :: &
79 : nblyr , & ! number of bio layers ! LCOV_EXCL_LINE
80 : ncat , & ! number of thickness categories ! LCOV_EXCL_LINE
81 : nilyr , & ! number of ice layers ! LCOV_EXCL_LINE
82 : nltrcr , & ! number of zbgc tracers ! LCOV_EXCL_LINE
83 : nbtrcr , & ! number of biology tracers ! LCOV_EXCL_LINE
84 : ntrcr ! number of tracers in use
85 :
86 : real (kind=dbl_kind), dimension (nblyr+2), intent(in) :: &
87 : bgrid ! biology nondimensional vertical grid points
88 :
89 : real (kind=dbl_kind), dimension (nblyr+1), intent(in) :: &
90 : igrid ! biology vertical interface points
91 :
92 : real (kind=dbl_kind), dimension (nilyr+1), intent(in) :: &
93 : cgrid ! CICE vertical coordinate
94 :
95 : real (kind=dbl_kind), intent(in) :: &
96 : dt ! time step (s)
97 :
98 : real (kind=dbl_kind), dimension (:), intent(in) :: &
99 : aicen_init , & ! initial concentration of ice ! LCOV_EXCL_LINE
100 : vicen_init , & ! intiial volume per unit area of ice (m) ! LCOV_EXCL_LINE
101 : aicen , & ! concentration of ice ! LCOV_EXCL_LINE
102 : vicen ! volume per unit area of ice (m)
103 :
104 : real (kind=dbl_kind), intent(in) :: &
105 : vsnon1 ! category 1 snow volume per unit area (m)
106 :
107 : real (kind=dbl_kind), dimension (:,:), intent(inout) :: &
108 : trcrn ! ice tracers
109 :
110 : real (kind=dbl_kind), intent(in) :: &
111 : sss !sea surface salinity (ppt)
112 :
113 : real (kind=dbl_kind), intent(in) :: &
114 : vi0_init , & ! volume of new ice added to cat 1 (intial) ! LCOV_EXCL_LINE
115 : vi0new ! volume of new ice added to cat 1
116 :
117 : real (kind=dbl_kind), intent(in) :: &
118 : hsurp ! thickness of new ice added to each cat
119 :
120 : real (kind=dbl_kind), dimension (:), intent(inout) :: &
121 : flux_bio ! tracer flux to ocean from biology (mmol/m^2/s)
122 :
123 : real (kind=dbl_kind), dimension (:), intent(in) :: &
124 : ocean_bio ! ocean concentration of biological tracer
125 :
126 : ! local
127 :
128 : integer (kind=int_kind) :: &
129 : location , & ! 1 (add frazil to bottom), 0 (add frazil throughout) ! LCOV_EXCL_LINE
130 : n , & ! ice category index ! LCOV_EXCL_LINE
131 : k ! ice layer index
132 :
133 : real (kind=dbl_kind) :: &
134 : vbri1 , & ! starting volume of existing brine ! LCOV_EXCL_LINE
135 : vbri_init , & ! brine volume summed over categories ! LCOV_EXCL_LINE
136 0 : vbri_final ! brine volume summed over categories
137 :
138 : real (kind=dbl_kind) :: &
139 : vsurp , & ! volume of new ice added to each cat ! LCOV_EXCL_LINE
140 0 : vtmp ! total volume of new and old ice
141 :
142 : real (kind=dbl_kind), dimension (ncat) :: &
143 0 : vbrin ! trcrn(nt_fbri,n)*vicen(n)
144 :
145 : real (kind=dbl_kind) :: &
146 0 : vice_new ! vicen_init + vsurp
147 :
148 : real (kind=dbl_kind) :: &
149 : Tmlts ! melting temperature (oC)
150 :
151 : character (len=char_len) :: &
152 : fieldid ! field identifier
153 :
154 : character(len=*),parameter :: subname='(add_new_ice_bgc)'
155 :
156 : !-----------------------------------------------------------------
157 : ! brine
158 : !-----------------------------------------------------------------
159 0 : vbrin(:) = c0
160 0 : do n = 1, ncat
161 0 : vbrin(n) = vicen_init(n)
162 0 : if (tr_brine) vbrin(n) = trcrn(nt_fbri,n)*vicen_init(n)
163 : enddo
164 :
165 0 : call column_sum (ncat, vbrin, vbri_init)
166 0 : if (icepack_warnings_aborted(subname)) return
167 :
168 0 : vbri_init = vbri_init + vi0_init
169 0 : do k = 1, nbtrcr
170 0 : flux_bio(k) = flux_bio(k) &
171 0 : - vi0_init/dt*ocean_bio(k)*zbgc_init_frac(k)
172 : enddo
173 : !-----------------------------------------------------------------
174 : ! Distribute bgc in new ice volume among all ice categories by
175 : ! increasing ice thickness, leaving ice area unchanged.
176 : !-----------------------------------------------------------------
177 :
178 : ! Diffuse_bio handles concentration changes from ice growth/melt
179 : ! ice area does not change
180 : ! add salt to the bottom , location = 1
181 :
182 0 : vsurp = c0
183 0 : vtmp = c0
184 :
185 0 : do n = 1,ncat
186 :
187 0 : if (hsurp > c0) then
188 :
189 0 : vtmp = vbrin(n)
190 0 : vsurp = hsurp * aicen_init(n)
191 0 : vbrin(n) = vbrin(n) + vsurp
192 0 : vice_new = vicen_init(n) + vsurp
193 0 : if (tr_brine .and. vicen(n) > c0) then
194 0 : trcrn(nt_fbri,n) = vbrin(n)/vicen(n)
195 0 : elseif (tr_brine .and. vicen(n) <= c0) then
196 0 : trcrn(nt_fbri,n) = c1
197 : endif
198 :
199 0 : if (nltrcr > 0) then
200 0 : location = 1
201 : call adjust_tracer_profile(nbtrcr, dt, ntrcr, &
202 : aicen_init(n), & ! LCOV_EXCL_LINE
203 : vbrin(n), & ! LCOV_EXCL_LINE
204 : vice_new, & ! LCOV_EXCL_LINE
205 : trcrn(:,n), & ! LCOV_EXCL_LINE
206 : vtmp, & ! LCOV_EXCL_LINE
207 : vsurp, sss, & ! LCOV_EXCL_LINE
208 : nilyr, nblyr, & ! LCOV_EXCL_LINE
209 : bgrid, & ! LCOV_EXCL_LINE
210 : cgrid, & ! LCOV_EXCL_LINE
211 : ocean_bio, igrid, & ! LCOV_EXCL_LINE
212 0 : location)
213 0 : if (icepack_warnings_aborted(subname)) return
214 : endif ! nltrcr
215 : endif ! hsurp > 0
216 : enddo ! n
217 :
218 : !-----------------------------------------------------------------
219 : ! Combine bgc in new ice grown in open water with category 1 ice.
220 : !-----------------------------------------------------------------
221 :
222 0 : if (vi0new > c0) then
223 :
224 0 : vbri1 = vbrin(1)
225 0 : vbrin(1) = vbrin(1) + vi0new
226 0 : if (tr_brine .and. vicen(1) > c0) then
227 0 : trcrn(nt_fbri,1) = vbrin(1)/vicen(1)
228 0 : elseif (tr_brine .and. vicen(1) <= c0) then
229 0 : trcrn(nt_fbri,1) = c1
230 : endif
231 :
232 : ! Diffuse_bio handles concentration changes from ice growth/melt
233 : ! ice area changes
234 : ! add salt throughout, location = 0
235 :
236 0 : if (nltrcr > 0) then
237 0 : location = 0
238 : call adjust_tracer_profile(nbtrcr, dt, ntrcr, &
239 : aicen(1), & ! LCOV_EXCL_LINE
240 : vbrin(1), & ! LCOV_EXCL_LINE
241 : vicen(1), & ! LCOV_EXCL_LINE
242 : trcrn(:,1), & ! LCOV_EXCL_LINE
243 : vbri1, & ! LCOV_EXCL_LINE
244 : vi0new, sss, & ! LCOV_EXCL_LINE
245 : nilyr, nblyr, & ! LCOV_EXCL_LINE
246 : bgrid, & ! LCOV_EXCL_LINE
247 : cgrid, & ! LCOV_EXCL_LINE
248 : ocean_bio, igrid, & ! LCOV_EXCL_LINE
249 0 : location)
250 0 : if (icepack_warnings_aborted(subname)) return
251 :
252 : endif ! nltrcr > 0
253 : endif ! vi0new > 0
254 :
255 0 : if (tr_brine .and. conserv_check) then
256 0 : call column_sum (ncat, vbrin, vbri_final)
257 0 : if (icepack_warnings_aborted(subname)) return
258 :
259 0 : fieldid = subname//':vbrin'
260 : call column_conservation_check (fieldid, &
261 : vbri_init, vbri_final, & ! LCOV_EXCL_LINE
262 0 : puny)
263 0 : if (icepack_warnings_aborted(subname)) return
264 :
265 : endif ! conserv_check
266 :
267 : end subroutine add_new_ice_bgc
268 :
269 : !=======================================================================
270 :
271 : ! When sea ice melts laterally, flux bgc to ocean
272 :
273 0 : subroutine lateral_melt_bgc (dt, &
274 : ncat, nblyr, & ! LCOV_EXCL_LINE
275 : rside, vicen, & ! LCOV_EXCL_LINE
276 : trcrn, & ! LCOV_EXCL_LINE
277 0 : flux_bio, nbltrcr)
278 :
279 : integer (kind=int_kind), intent(in) :: &
280 : ncat , & ! number of thickness categories ! LCOV_EXCL_LINE
281 : nblyr , & ! number of bio layers ! LCOV_EXCL_LINE
282 : nbltrcr ! number of biology tracers
283 :
284 : real (kind=dbl_kind), intent(in) :: &
285 : dt ! time step (s)
286 :
287 : real (kind=dbl_kind), dimension(:), intent(in) :: &
288 : vicen ! volume per unit area of ice (m)
289 :
290 : real (kind=dbl_kind), dimension (:,:), intent(in) :: &
291 : trcrn ! tracer array
292 :
293 : real (kind=dbl_kind), intent(in) :: &
294 : rside ! fraction of ice that melts laterally
295 :
296 : real (kind=dbl_kind), dimension(:), intent(inout) :: &
297 : flux_bio ! biology tracer flux from layer bgc (mmol/m^2/s)
298 :
299 : ! local variables
300 :
301 : integer (kind=int_kind) :: &
302 : k , & ! layer index ! LCOV_EXCL_LINE
303 : m , & ! ! LCOV_EXCL_LINE
304 : n ! category index
305 :
306 : real (kind=dbl_kind) :: &
307 0 : zspace ! bio grid spacing
308 :
309 : character(len=*),parameter :: subname='(lateral_melt_bgc)'
310 :
311 0 : zspace = c1/(real(nblyr,kind=dbl_kind))
312 :
313 0 : do m = 1, nbltrcr
314 0 : do n = 1, ncat
315 0 : do k = 1, nblyr+1
316 0 : flux_bio(m) = flux_bio(m) + trcrn(nt_fbri,n) &
317 : * vicen(n)*zspace*trcrn(bio_index(m)+k-1,n) & ! LCOV_EXCL_LINE
318 0 : * rside/dt
319 : enddo
320 : enddo
321 : enddo
322 :
323 0 : end subroutine lateral_melt_bgc
324 :
325 : !=======================================================================
326 : !
327 : ! Add new ice tracers to the ice bottom and adjust the vertical profile
328 : !
329 : ! author: Nicole Jeffery, LANL
330 :
331 0 : subroutine adjust_tracer_profile (nbtrcr, dt, ntrcr, &
332 : aicen, vbrin, & ! LCOV_EXCL_LINE
333 : vicen, trcrn, & ! LCOV_EXCL_LINE
334 : vtmp, & ! LCOV_EXCL_LINE
335 : vsurp, sss, & ! LCOV_EXCL_LINE
336 : nilyr, nblyr, & ! LCOV_EXCL_LINE
337 : bgrid, & ! LCOV_EXCL_LINE
338 : cgrid, ocean_bio, & ! LCOV_EXCL_LINE
339 0 : igrid, location)
340 :
341 : integer (kind=int_kind), intent(in) :: &
342 : location , & ! 1 (add frazil to bottom), 0 (add frazil throughout) ! LCOV_EXCL_LINE
343 : ntrcr , & ! number of tracers in use ! LCOV_EXCL_LINE
344 : nilyr , & ! number of ice layers ! LCOV_EXCL_LINE
345 : nbtrcr , & ! number of biology tracers ! LCOV_EXCL_LINE
346 : nblyr ! number of biology layers
347 :
348 : real (kind=dbl_kind), intent(in) :: &
349 : dt ! timestep (s)
350 :
351 : real (kind=dbl_kind), intent(in) :: &
352 : aicen , & ! concentration of ice ! LCOV_EXCL_LINE
353 : vicen , & ! volume of ice ! LCOV_EXCL_LINE
354 : sss , & ! ocean salinity (ppt) ! LCOV_EXCL_LINE
355 : ! hsurp , & ! flags new ice added to each cat ! LCOV_EXCL_LINE
356 : vsurp , & ! volume of new ice added to each cat ! LCOV_EXCL_LINE
357 : vtmp ! total volume of new and old ice
358 :
359 : real (kind=dbl_kind), dimension (nbtrcr), intent(in) :: &
360 : ocean_bio
361 :
362 : real (kind=dbl_kind), intent(in) :: &
363 : vbrin ! fbri*volume per unit area of ice (m)
364 :
365 : real (kind=dbl_kind), dimension (nblyr+1), intent(in) :: &
366 : igrid ! zbio grid
367 :
368 : real (kind=dbl_kind), dimension (nblyr+2), intent(in) :: &
369 : bgrid ! zsal grid
370 :
371 : real (kind=dbl_kind), dimension (nilyr+1), intent(in) :: &
372 : cgrid ! CICE grid
373 :
374 : real (kind=dbl_kind), dimension (ntrcr), intent(inout) :: &
375 : trcrn ! ice tracers
376 :
377 : ! local variables
378 :
379 : real (kind=dbl_kind), dimension (ntrcr+2) :: &
380 : trtmp0, & ! temporary, remapped tracers ! LCOV_EXCL_LINE
381 0 : trtmp ! temporary, remapped tracers
382 :
383 : real (kind=dbl_kind) :: &
384 : hin , & ! ice height ! LCOV_EXCL_LINE
385 : hinS_new, & ! brine height ! LCOV_EXCL_LINE
386 : temp_S
387 :
388 : integer (kind=int_kind) :: &
389 : k, m
390 :
391 : real (kind=dbl_kind), dimension (nblyr+1) :: &
392 0 : C_stationary ! stationary bulk concentration*h (mmol/m^2)
393 :
394 : real (kind=dbl_kind), dimension (nblyr) :: &
395 : S_stationary ! stationary bulk concentration*h (ppt*m)
396 :
397 : real(kind=dbl_kind) :: &
398 : top_conc , & ! salinity or bgc ocean concentration of frazil ! LCOV_EXCL_LINE
399 : fluxb , & ! needed for regrid (set to zero here) ! LCOV_EXCL_LINE
400 : hbri_old , & ! previous timestep brine height ! LCOV_EXCL_LINE
401 0 : hbri ! brine height
402 :
403 : character(len=*),parameter :: subname='(adjust_tracer_profile)'
404 :
405 0 : trtmp0(:) = c0
406 0 : trtmp(:) = c0
407 0 : fluxb = c0
408 :
409 0 : if (location == 1 .and. vbrin > c0) then ! add frazil to bottom
410 :
411 0 : hbri = vbrin
412 0 : hbri_old = vtmp
413 :
414 0 : do m = 1, nbtrcr
415 0 : top_conc = ocean_bio(m)*zbgc_init_frac(m)
416 0 : do k = 1, nblyr+1
417 0 : C_stationary(k) = trcrn(bio_index(m) + k-1)* hbri_old
418 : enddo !k
419 : call regrid_stationary (C_stationary, hbri_old, &
420 : hbri, dt, & ! LCOV_EXCL_LINE
421 : ntrcr, & ! LCOV_EXCL_LINE
422 : nblyr, top_conc, & ! LCOV_EXCL_LINE
423 0 : igrid, fluxb )
424 0 : if (icepack_warnings_aborted(subname)) return
425 0 : do k = 1, nblyr+1
426 0 : trcrn(bio_index(m) + k-1) = C_stationary(k)/hbri
427 : enddo !k
428 : enddo !m
429 :
430 0 : elseif (vbrin > c0) then ! add frazil throughout location == 0 .and.
431 :
432 0 : do k = 1, nblyr+1
433 0 : do m = 1, nbtrcr
434 0 : trcrn(bio_index(m) + k-1) = (trcrn(bio_index(m) + k-1) * vtmp &
435 0 : + ocean_bio(m)*zbgc_init_frac(m) * vsurp) / vbrin
436 : enddo
437 : enddo
438 :
439 : endif ! location
440 :
441 : end subroutine adjust_tracer_profile
442 :
443 : !=======================================================================
444 : !autodocument_start icepack_init_bgc
445 : !
446 0 : subroutine icepack_init_bgc(ncat, nblyr, nilyr, ntrcr_o, &
447 : cgrid, igrid, ntrcr, nbtrcr, & ! LCOV_EXCL_LINE
448 0 : sicen, trcrn, sss, ocean_bio_all)
449 :
450 : integer (kind=int_kind), intent(in) :: &
451 : ncat , & ! number of thickness categories ! LCOV_EXCL_LINE
452 : nilyr , & ! number of ice layers ! LCOV_EXCL_LINE
453 : nblyr , & ! number of bio layers ! LCOV_EXCL_LINE
454 : ntrcr_o,& ! number of tracers not including bgc ! LCOV_EXCL_LINE
455 : ntrcr , & ! number of tracers in use ! LCOV_EXCL_LINE
456 : nbtrcr ! number of bio tracers in use
457 :
458 : real (kind=dbl_kind), dimension (nblyr+1), intent(inout) :: &
459 : igrid ! biology vertical interface points
460 :
461 : real (kind=dbl_kind), dimension (nilyr+1), intent(inout) :: &
462 : cgrid ! CICE vertical coordinate
463 :
464 : real (kind=dbl_kind), dimension(nilyr, ncat), intent(in) :: &
465 : sicen ! salinity on the cice grid
466 :
467 : real (kind=dbl_kind), dimension (:,:), intent(inout) :: &
468 : trcrn ! subset of tracer array (only bgc)
469 :
470 : real (kind=dbl_kind), intent(in) :: &
471 : sss ! sea surface salinity (ppt)
472 :
473 : real (kind=dbl_kind), dimension (:), intent(inout) :: &
474 : ocean_bio_all ! fixed order, all values even for tracers false
475 :
476 : !autodocument_end
477 :
478 : ! local variables
479 :
480 : integer (kind=int_kind) :: &
481 : k , & ! vertical index ! LCOV_EXCL_LINE
482 : n , & ! category index ! LCOV_EXCL_LINE
483 : mm ! bio tracer index
484 :
485 : real (kind=dbl_kind), dimension (ntrcr+2) :: &
486 0 : trtmp ! temporary, remapped tracers
487 :
488 : character(len=*),parameter :: subname='(icepack_init_bgc)'
489 :
490 : !-----------------------------------------------------------------------------
491 : ! Skeletal Layer Model
492 : ! All bgc tracers are Bulk quantities in units of mmol or mg per m^3
493 : ! The skeletal layer model assumes a constant
494 : ! layer depth (sk_l) and porosity (phi_sk)
495 : !-----------------------------------------------------------------------------
496 0 : if (skl_bgc) then
497 :
498 0 : do n = 1,ncat
499 0 : do mm = 1,nbtrcr
500 : ! bulk concentration (mmol or mg per m^3, or 10^-3 mmol/m^3)
501 0 : trcrn(bio_index(mm)-ntrcr_o, n) = ocean_bio_all(bio_index_o(mm))
502 : enddo ! nbtrcr
503 : enddo ! n
504 :
505 : !-----------------------------------------------------------------------------
506 : ! zbgc Model
507 : ! All bgc tracers are Bulk quantities in units of mmol or mg per m^3
508 : ! The vertical layer model uses prognosed porosity and layer depth
509 : !-----------------------------------------------------------------------------
510 :
511 : else ! not skl_bgc
512 :
513 0 : if (scale_bgc .and. ktherm == 2) then
514 0 : trtmp(:) = c0
515 0 : do n = 1,ncat
516 : call remap_zbgc(nilyr, &
517 : 1, & ! LCOV_EXCL_LINE
518 : sicen(:,n), trtmp, & ! LCOV_EXCL_LINE
519 : 0, nblyr+1, & ! LCOV_EXCL_LINE
520 : c1, c1, & ! LCOV_EXCL_LINE
521 : cgrid(2:nilyr+1), & ! LCOV_EXCL_LINE
522 : igrid(1:nblyr+1), & ! LCOV_EXCL_LINE
523 0 : sicen(1,n) )
524 0 : if (icepack_warnings_aborted(subname)) return
525 :
526 0 : do mm = 1,nbtrcr
527 0 : do k = 1, nblyr + 1
528 0 : trcrn(bio_index(mm)+k-1-ntrcr_o,n) = &
529 0 : (trtmp(k)/sss*ocean_bio_all(bio_index_o(mm)))
530 0 : trcrn(bio_index(mm)+nblyr+1-ntrcr_o:bio_index(mm)+nblyr+2-ntrcr_o,n) = c0 ! snow
531 : enddo ! k
532 : enddo ! mm
533 : enddo ! n
534 :
535 0 : elseif (nbtrcr > 0 .and. nt_fbri > 0) then ! not scale_bgc
536 :
537 0 : do n = 1,ncat
538 0 : do mm = 1,nbtrcr
539 0 : do k = 1, nblyr+1
540 0 : trcrn(bio_index(mm)+k-1-ntrcr_o,n) = ocean_bio_all(bio_index_o(mm)) &
541 0 : * zbgc_init_frac(mm)
542 0 : trcrn(bio_index(mm)+nblyr+1-ntrcr_o:bio_index(mm)+nblyr+2-ntrcr_o,n) = c0 ! snow
543 : enddo ! k
544 0 : trcrn(nt_zbgc_frac-1+mm-ntrcr_o,n) = zbgc_frac_init(mm)
545 : enddo ! mm
546 : enddo ! n
547 :
548 : endif ! scale_bgc
549 : endif ! skl_bgc
550 :
551 : end subroutine icepack_init_bgc
552 :
553 : !=======================================================================
554 : !autodocument_start icepack_init_zbgc
555 : !
556 :
557 74 : subroutine icepack_init_zbgc ( &
558 : R_Si2N_in, R_S2N_in, R_Fe2C_in, R_Fe2N_in, R_C2N_in, R_C2N_DON_in, & ! LCOV_EXCL_LINE
559 : R_chl2N_in, F_abs_chl_in, R_Fe2DON_in, R_Fe2DOC_in, chlabs_in, & ! LCOV_EXCL_LINE
560 : alpha2max_low_in, beta2max_in, mu_max_in, fr_graze_in, mort_pre_in, & ! LCOV_EXCL_LINE
561 : mort_Tdep_in, k_exude_in, K_Nit_in, K_Am_in, K_sil_in, K_Fe_in, & ! LCOV_EXCL_LINE
562 : f_don_in, kn_bac_in, f_don_Am_in, f_doc_in, f_exude_in, k_bac_in, & ! LCOV_EXCL_LINE
563 : grow_Tdep_in, zbgc_frac_init_in, & ! LCOV_EXCL_LINE
564 : zbgc_init_frac_in, tau_ret_in, tau_rel_in, bgc_tracer_type_in, & ! LCOV_EXCL_LINE
565 : fr_resp_in, algal_vel_in, R_dFe2dust_in, dustFe_sol_in, T_max_in, & ! LCOV_EXCL_LINE
566 : op_dep_min_in, fr_graze_s_in, fr_graze_e_in, fr_mort2min_in, fr_dFe_in, & ! LCOV_EXCL_LINE
567 : k_nitrif_in, t_iron_conv_in, max_loss_in, max_dfe_doc1_in, & ! LCOV_EXCL_LINE
568 : fr_resp_s_in, y_sk_DMS_in, t_sk_conv_in, t_sk_ox_in, fsal_in)
569 :
570 : real (kind=dbl_kind), optional :: R_C2N_in(:) ! algal C to N (mole/mole)
571 : real (kind=dbl_kind), optional :: R_chl2N_in(:) ! 3 algal chlorophyll to N (mg/mmol)
572 : real (kind=dbl_kind), optional :: F_abs_chl_in(:) ! to scale absorption in Dedd
573 : real (kind=dbl_kind), optional :: R_C2N_DON_in(:) ! increase compare to algal R_Fe2C
574 : real (kind=dbl_kind), optional :: R_Si2N_in(:) ! algal Sil to N (mole/mole)
575 : real (kind=dbl_kind), optional :: R_S2N_in(:) ! algal S to N (mole/mole)
576 : real (kind=dbl_kind), optional :: R_Fe2C_in(:) ! algal Fe to carbon (umol/mmol)
577 : real (kind=dbl_kind), optional :: R_Fe2N_in(:) ! algal Fe to N (umol/mmol)
578 : real (kind=dbl_kind), optional :: R_Fe2DON_in(:) ! Fe to N of DON (nmol/umol)
579 : real (kind=dbl_kind), optional :: R_Fe2DOC_in(:) ! Fe to C of DOC (nmol/umol)
580 :
581 : real (kind=dbl_kind), optional :: fr_resp_in ! frac of algal growth lost due to respiration
582 : real (kind=dbl_kind), optional :: algal_vel_in ! 0.5 cm/d(m/s) Lavoie 2005 1.5 cm/day
583 : real (kind=dbl_kind), optional :: R_dFe2dust_in ! g/g (3.5% content) Tagliabue 2009
584 : real (kind=dbl_kind), optional :: dustFe_sol_in ! solubility fraction
585 : real (kind=dbl_kind), optional :: T_max_in ! maximum temperature (C)
586 : real (kind=dbl_kind), optional :: op_dep_min_in ! Light attenuates for optical depths exceeding min
587 : real (kind=dbl_kind), optional :: fr_graze_s_in ! fraction of grazing spilled or slopped
588 : real (kind=dbl_kind), optional :: fr_graze_e_in ! fraction of assimilation excreted
589 : real (kind=dbl_kind), optional :: fr_mort2min_in ! fractionation of mortality to Am
590 : real (kind=dbl_kind), optional :: fr_dFe_in ! fraction of remineralized nitrogen
591 : ! (in units of algal iron)
592 : real (kind=dbl_kind), optional :: k_nitrif_in ! nitrification rate (1/day)
593 : real (kind=dbl_kind), optional :: t_iron_conv_in ! desorption loss pFe to dFe (day)
594 : real (kind=dbl_kind), optional :: max_loss_in ! restrict uptake to % of remaining value
595 : real (kind=dbl_kind), optional :: max_dfe_doc1_in ! max ratio of dFe to saccharides in the ice (nM Fe/muM C)
596 : real (kind=dbl_kind), optional :: fr_resp_s_in ! DMSPd fraction of respiration loss as DMSPd
597 : real (kind=dbl_kind), optional :: y_sk_DMS_in ! fraction conversion given high yield
598 : real (kind=dbl_kind), optional :: t_sk_conv_in ! Stefels conversion time (d)
599 : real (kind=dbl_kind), optional :: t_sk_ox_in ! DMS oxidation time (d)
600 : real (kind=dbl_kind), optional :: fsal_in ! salinity limitation factor (1)
601 :
602 : real (kind=dbl_kind), optional :: chlabs_in(:) ! chla absorption 1/m/(mg/m^3)
603 : real (kind=dbl_kind), optional :: alpha2max_low_in(:) ! light limitation (1/(W/m^2))
604 : real (kind=dbl_kind), optional :: beta2max_in(:) ! light inhibition (1/(W/m^2))
605 : real (kind=dbl_kind), optional :: mu_max_in(:) ! maximum growth rate (1/d)
606 : real (kind=dbl_kind), optional :: grow_Tdep_in(:) ! T dependence of growth (1/C)
607 : real (kind=dbl_kind), optional :: fr_graze_in(:) ! fraction of algae grazed
608 : real (kind=dbl_kind), optional :: mort_pre_in(:) ! mortality (1/day)
609 : real (kind=dbl_kind), optional :: mort_Tdep_in(:) ! T dependence of mortality (1/C)
610 : real (kind=dbl_kind), optional :: k_exude_in(:) ! algal carbon exudation rate (1/d)
611 : real (kind=dbl_kind), optional :: K_Nit_in(:) ! nitrate half saturation (mmol/m^3)
612 : real (kind=dbl_kind), optional :: K_Am_in(:) ! ammonium half saturation (mmol/m^3)
613 : real (kind=dbl_kind), optional :: K_Sil_in(:) ! silicon half saturation (mmol/m^3)
614 : real (kind=dbl_kind), optional :: K_Fe_in(:) ! iron half saturation or micromol/m^3
615 : real (kind=dbl_kind), optional :: f_don_in(:) ! fraction of spilled grazing to DON
616 : real (kind=dbl_kind), optional :: kn_bac_in(:) ! Bacterial degredation of DON (1/d)
617 : real (kind=dbl_kind), optional :: f_don_Am_in(:) ! fraction of remineralized DON to Am
618 : real (kind=dbl_kind), optional :: f_doc_in(:) ! fraction of mort_N that goes to each doc pool
619 : real (kind=dbl_kind), optional :: f_exude_in(:) ! fraction of exuded carbon to each DOC pool
620 : real (kind=dbl_kind), optional :: k_bac_in(:) ! Bacterial degredation of DOC (1/d)
621 :
622 : real (kind=dbl_kind), optional :: zbgc_frac_init_in(:) ! initializes mobile fraction
623 : real (kind=dbl_kind), optional :: bgc_tracer_type_in(:) ! described tracer in mobile or stationary phases
624 : real (kind=dbl_kind), optional :: zbgc_init_frac_in(:) ! fraction of ocean tracer concentration in new ice
625 : real (kind=dbl_kind), optional :: tau_ret_in(:) ! retention timescale (s), mobile to stationary phase
626 : real (kind=dbl_kind), optional :: tau_rel_in(:) ! release timescale (s), stationary to mobile phase
627 :
628 : !autodocument_end
629 :
630 : character(len=*),parameter :: subname='(icepack_init_zbgc)'
631 :
632 : !--------
633 :
634 185 : if (present(R_C2N_in)) R_C2N(:) = R_C2N_in(:)
635 185 : if (present(R_chl2N_in)) R_chl2N(:) = R_chl2N_in(:)
636 185 : if (present(F_abs_chl_in)) F_abs_chl(:) = F_abs_chl_in(:)
637 111 : if (present(R_C2N_DON_in)) R_C2N_DON(:) = R_C2N_DON_in(:)
638 185 : if (present(R_Si2N_in)) R_Si2N(:) = R_Si2N_in(:)
639 185 : if (present(R_S2N_in)) R_S2N(:) = R_S2N_in(:)
640 185 : if (present(R_Fe2C_in)) R_Fe2C(:) = R_Fe2C_in(:)
641 185 : if (present(R_Fe2N_in)) R_Fe2N(:) = R_Fe2N_in(:)
642 111 : if (present(R_Fe2DON_in)) R_Fe2DON(:) = R_Fe2DON_in(:)
643 185 : if (present(R_Fe2DOC_in)) R_Fe2DOC(:) = R_Fe2DOC_in(:)
644 :
645 74 : if (present(fr_resp_in)) fr_resp = fr_resp_in
646 74 : if (present(algal_vel_in)) algal_vel = algal_vel_in
647 74 : if (present(R_dFe2dust_in)) R_dFe2dust = R_dFe2dust_in
648 74 : if (present(dustFe_sol_in)) dustFe_sol = dustFe_sol_in
649 74 : if (present(T_max_in)) T_max = T_max_in
650 74 : if (present(op_dep_min_in)) op_dep_min = op_dep_min_in
651 74 : if (present(fr_graze_s_in)) fr_graze_s = fr_graze_s_in
652 74 : if (present(fr_graze_e_in)) fr_graze_e = fr_graze_e_in
653 74 : if (present(fr_mort2min_in)) fr_mort2min = fr_mort2min_in
654 74 : if (present(fr_dFe_in)) fr_dFe = fr_dFe_in
655 74 : if (present(k_nitrif_in)) k_nitrif = k_nitrif_in
656 74 : if (present(t_iron_conv_in)) t_iron_conv = t_iron_conv_in
657 74 : if (present(max_loss_in)) max_loss = max_loss_in
658 74 : if (present(max_dfe_doc1_in)) max_dfe_doc1 = max_dfe_doc1_in
659 74 : if (present(fr_resp_s_in)) fr_resp_s = fr_resp_s_in
660 74 : if (present(y_sk_DMS_in)) y_sk_DMS = y_sk_DMS_in
661 74 : if (present(t_sk_conv_in)) t_sk_conv = t_sk_conv_in
662 74 : if (present(t_sk_ox_in)) t_sk_ox = t_sk_ox_in
663 74 : if (present(fsal_in)) fsal = fsal_in
664 :
665 185 : if (present(chlabs_in)) chlabs(:) = chlabs_in(:)
666 185 : if (present(alpha2max_low_in)) alpha2max_low(:) = alpha2max_low_in(:)
667 185 : if (present(beta2max_in)) beta2max(:) = beta2max_in(:)
668 185 : if (present(mu_max_in)) mu_max(:) = mu_max_in(:)
669 185 : if (present(grow_Tdep_in)) grow_Tdep(:) = grow_Tdep_in(:)
670 185 : if (present(fr_graze_in)) fr_graze(:) = fr_graze_in(:)
671 185 : if (present(mort_pre_in)) mort_pre(:) = mort_pre_in(:)
672 185 : if (present(mort_Tdep_in)) mort_Tdep(:) = mort_Tdep_in(:)
673 185 : if (present(k_exude_in)) k_exude(:) = k_exude_in(:)
674 185 : if (present(K_Nit_in)) K_Nit(:) = K_Nit_in(:)
675 185 : if (present(K_Am_in)) K_Am(:) = K_Am_in(:)
676 185 : if (present(K_Sil_in)) K_Sil(:) = K_Sil_in(:)
677 185 : if (present(K_Fe_in)) K_Fe(:) = K_Fe_in(:)
678 111 : if (present(f_don_in)) f_don(:) = f_don_in(:)
679 111 : if (present(kn_bac_in)) kn_bac(:) = kn_bac_in(:)
680 111 : if (present(f_don_Am_in)) f_don_Am(:) = f_don_Am_in(:)
681 185 : if (present(f_doc_in)) f_doc(:) = f_doc_in(:)
682 185 : if (present(f_exude_in)) f_exude(:) = f_exude_in(:)
683 185 : if (present(k_bac_in)) k_bac(:) = k_bac_in(:)
684 :
685 1147 : if (present(zbgc_frac_init_in)) zbgc_frac_init(:) = zbgc_frac_init_in(:)
686 1147 : if (present(bgc_tracer_type_in)) bgc_tracer_type(:) = bgc_tracer_type_in(:)
687 1147 : if (present(zbgc_init_frac_in)) zbgc_init_frac(:) = zbgc_init_frac_in(:)
688 1147 : if (present(tau_ret_in)) tau_ret(:) = tau_ret_in(:)
689 1147 : if (present(tau_rel_in)) tau_rel(:) = tau_rel_in(:)
690 :
691 :
692 74 : end subroutine icepack_init_zbgc
693 :
694 : !=======================================================================
695 : !autodocument_start icepack_biogeochemistry
696 : !
697 :
698 0 : subroutine icepack_biogeochemistry(dt, &
699 : ntrcr, nbtrcr, & ! LCOV_EXCL_LINE
700 : upNO, upNH, iDi, iki, zfswin, & ! LCOV_EXCL_LINE
701 : zsal_tot, darcy_V, grow_net, & ! LCOV_EXCL_LINE
702 : PP_net, hbri,dhbr_bot, dhbr_top, Zoo,& ! LCOV_EXCL_LINE
703 : fbio_snoice, fbio_atmice, ocean_bio, & ! LCOV_EXCL_LINE
704 : first_ice, fswpenln, bphi, bTiz, ice_bio_net, & ! LCOV_EXCL_LINE
705 : snow_bio_net, fswthrun, Rayleigh_criteria, & ! LCOV_EXCL_LINE
706 : sice_rho, fzsal, fzsal_g, & ! LCOV_EXCL_LINE
707 : bgrid, igrid, icgrid, cgrid, & ! LCOV_EXCL_LINE
708 : nblyr, nilyr, nslyr, n_algae, n_zaero, ncat, & ! LCOV_EXCL_LINE
709 : n_doc, n_dic, n_don, n_fed, n_fep, & ! LCOV_EXCL_LINE
710 : meltbn, melttn, congeln, snoicen, & ! LCOV_EXCL_LINE
711 : sst, sss, fsnow, meltsn, & ! LCOV_EXCL_LINE
712 : hin_old, flux_bio, flux_bio_atm, & ! LCOV_EXCL_LINE
713 : aicen_init, vicen_init, aicen, vicen, vsnon, & ! LCOV_EXCL_LINE
714 0 : aice0, trcrn, vsnon_init, skl_bgc)
715 :
716 : real (kind=dbl_kind), intent(in) :: &
717 : dt ! time step
718 :
719 : integer (kind=int_kind), intent(in) :: &
720 : ncat, & ! LCOV_EXCL_LINE
721 : nilyr, & ! LCOV_EXCL_LINE
722 : nslyr, & ! LCOV_EXCL_LINE
723 : nblyr, & ! LCOV_EXCL_LINE
724 : ntrcr, & ! LCOV_EXCL_LINE
725 : nbtrcr, & ! LCOV_EXCL_LINE
726 : n_algae, n_zaero, & ! LCOV_EXCL_LINE
727 : n_doc, n_dic, n_don, n_fed, n_fep
728 :
729 : real (kind=dbl_kind), dimension (:), intent(inout) :: &
730 : bgrid , & ! biology nondimensional vertical grid points ! LCOV_EXCL_LINE
731 : igrid , & ! biology vertical interface points ! LCOV_EXCL_LINE
732 : cgrid , & ! CICE vertical coordinate ! LCOV_EXCL_LINE
733 : icgrid , & ! interface grid for CICE (shortwave variable) ! LCOV_EXCL_LINE
734 : ocean_bio , & ! contains all the ocean bgc tracer concentrations ! LCOV_EXCL_LINE
735 : fbio_snoice , & ! fluxes from snow to ice ! LCOV_EXCL_LINE
736 : fbio_atmice , & ! fluxes from atm to ice ! LCOV_EXCL_LINE
737 : dhbr_top , & ! brine top change ! LCOV_EXCL_LINE
738 : dhbr_bot , & ! brine bottom change ! LCOV_EXCL_LINE
739 : darcy_V , & ! darcy velocity positive up (m/s) ! LCOV_EXCL_LINE
740 : hin_old , & ! old ice thickness ! LCOV_EXCL_LINE
741 : sice_rho , & ! avg sea ice density (kg/m^3) ! LCOV_EXCL_LINE
742 : ice_bio_net , & ! depth integrated tracer (mmol/m^2) ! LCOV_EXCL_LINE
743 : snow_bio_net , & ! depth integrated snow tracer (mmol/m^2) ! LCOV_EXCL_LINE
744 : flux_bio ! all bio fluxes to ocean
745 :
746 : logical (kind=log_kind), dimension (:), intent(inout) :: &
747 : first_ice ! distinguishes ice that disappears (e.g. melts)
748 : ! and reappears (e.g. transport) in a grid cell
749 : ! during a single time step from ice that was
750 : ! there the entire time step (true until ice forms)
751 :
752 : real (kind=dbl_kind), dimension (:,:), intent(inout) :: &
753 : Zoo , & ! N losses accumulated in timestep (ie. zooplankton/bacteria) ! LCOV_EXCL_LINE
754 : ! mmol/m^3
755 : bphi , & ! porosity of layers
756 : bTiz , & ! layer temperatures interpolated on bio grid (C) ! LCOV_EXCL_LINE
757 : zfswin , & ! Shortwave flux into layers interpolated on bio grid (W/m^2) ! LCOV_EXCL_LINE
758 : iDi , & ! igrid Diffusivity (m^2/s) ! LCOV_EXCL_LINE
759 : iki , & ! Ice permeability (m^2) ! LCOV_EXCL_LINE
760 : trcrn ! tracers
761 :
762 : real (kind=dbl_kind), intent(inout) :: &
763 : grow_net , & ! Specific growth rate (/s) per grid cell ! LCOV_EXCL_LINE
764 : PP_net , & ! Total production (mg C/m^2/s) per grid cell ! LCOV_EXCL_LINE
765 : hbri , & ! brine height, area-averaged for comparison with hi (m) ! LCOV_EXCL_LINE
766 : upNO , & ! nitrate uptake rate (mmol/m^2/d) times aice ! LCOV_EXCL_LINE
767 : upNH ! ammonium uptake rate (mmol/m^2/d) times aice
768 :
769 : real (kind=dbl_kind), intent(inout), optional :: &
770 : zsal_tot ! Total ice salinity in per grid cell (g/m^2) (deprecated)
771 :
772 : real (kind=dbl_kind), intent(inout), optional :: &
773 : fzsal , & ! Total flux of salt to ocean at time step for conservation (deprecated) ! LCOV_EXCL_LINE
774 : fzsal_g ! Total gravity drainage flux (deprecated)
775 :
776 : logical (kind=log_kind), intent(inout), optional :: &
777 : Rayleigh_criteria ! .true. means Ra_c was reached (deprecated)
778 :
779 : real (kind=dbl_kind), dimension (:,:), intent(in) :: &
780 : fswpenln ! visible SW entering ice layers (W m-2)
781 :
782 : real (kind=dbl_kind), dimension (:), intent(in) :: &
783 : fswthrun , & ! SW through ice to ocean (W/m^2) ! LCOV_EXCL_LINE
784 : meltsn , & ! snow melt in category n (m) ! LCOV_EXCL_LINE
785 : melttn , & ! top melt in category n (m) ! LCOV_EXCL_LINE
786 : meltbn , & ! bottom melt in category n (m) ! LCOV_EXCL_LINE
787 : congeln , & ! congelation ice formation in category n (m) ! LCOV_EXCL_LINE
788 : snoicen , & ! snow-ice formation in category n (m) ! LCOV_EXCL_LINE
789 : flux_bio_atm, & ! all bio fluxes to ice from atmosphere ! LCOV_EXCL_LINE
790 : aicen_init , & ! initial ice concentration, for linear ITD ! LCOV_EXCL_LINE
791 : vicen_init , & ! initial ice volume (m), for linear ITD ! LCOV_EXCL_LINE
792 : vsnon_init , & ! initial snow volume (m), for aerosol ! LCOV_EXCL_LINE
793 : aicen , & ! concentration of ice ! LCOV_EXCL_LINE
794 : vicen , & ! volume per unit area of ice (m) ! LCOV_EXCL_LINE
795 : vsnon ! volume per unit area of snow (m)
796 :
797 : real (kind=dbl_kind), intent(in) :: &
798 : aice0 , & ! open water area fraction ! LCOV_EXCL_LINE
799 : sss , & ! sea surface salinity (ppt) ! LCOV_EXCL_LINE
800 : sst , & ! sea surface temperature (C) ! LCOV_EXCL_LINE
801 : fsnow ! snowfall rate (kg/m^2 s)
802 :
803 : logical (kind=log_kind), intent(in) :: &
804 : skl_bgc ! if true, solve skeletal biochemistry
805 :
806 : !autodocument_end
807 :
808 : ! local variables
809 :
810 : integer (kind=int_kind) :: &
811 : n, mm ! thickness category index
812 :
813 : real (kind=dbl_kind) :: &
814 : hin , & ! new ice thickness ! LCOV_EXCL_LINE
815 : hsn , & ! snow thickness (m) ! LCOV_EXCL_LINE
816 : hbr_old , & ! old brine thickness before growh/melt ! LCOV_EXCL_LINE
817 : dhice , & ! change due to sublimation/condensation (m) ! LCOV_EXCL_LINE
818 : kavg , & ! average ice permeability (m^2) ! LCOV_EXCL_LINE
819 : bphi_o , & ! surface ice porosity ! LCOV_EXCL_LINE
820 : hbrin , & ! brine height ! LCOV_EXCL_LINE
821 0 : dh_direct ! surface flooding or runoff
822 :
823 : real (kind=dbl_kind), dimension (nblyr+2) :: &
824 : ! Defined on Bio Grid points
825 0 : bSin , & ! salinity on the bio grid (ppt)
826 : brine_sal , & ! brine salinity (ppt) ! LCOV_EXCL_LINE
827 0 : brine_rho ! brine_density (kg/m^3)
828 :
829 : real (kind=dbl_kind), dimension (nblyr+1) :: &
830 : ! Defined on Bio Grid interfaces
831 0 : iphin , & ! porosity
832 : ibrine_sal , & ! brine salinity (ppt) ! LCOV_EXCL_LINE
833 : ibrine_rho , & ! brine_density (kg/m^3) ! LCOV_EXCL_LINE
834 0 : iTin ! Temperature on the interface grid (oC)
835 :
836 : real (kind=dbl_kind) :: &
837 0 : sloss ! brine flux contribution from surface runoff (g/m^2)
838 :
839 : ! for bgc sk
840 : real (kind=dbl_kind) :: &
841 : dh_bot_chl , & ! Chlorophyll may or may not flush ! LCOV_EXCL_LINE
842 : dh_top_chl , & ! Chlorophyll may or may not flush ! LCOV_EXCL_LINE
843 0 : darcy_V_chl
844 :
845 : character(len=*),parameter :: subname='(icepack_biogeochemistry)'
846 :
847 :
848 0 : do n = 1, ncat
849 :
850 : !-----------------------------------------------------------------
851 : ! initialize
852 : !-----------------------------------------------------------------
853 0 : hin_old(n) = c0
854 0 : if (aicen_init(n) > puny) then
855 0 : hin_old(n) = vicen_init(n) &
856 0 : / aicen_init(n)
857 : else
858 0 : first_ice(n) = .true.
859 0 : if (tr_brine) trcrn(nt_fbri,n) = c1
860 0 : if (z_tracers) then
861 0 : do mm = 1,nbtrcr
862 0 : trcrn(nt_zbgc_frac-1+mm,n) = zbgc_frac_init(mm)
863 : enddo
864 : endif
865 : endif
866 :
867 0 : if (aicen(n) > puny) then
868 :
869 0 : dh_top_chl = c0
870 0 : dh_bot_chl = c0
871 0 : darcy_V_chl= c0
872 0 : bSin(:) = c0
873 0 : hsn = c0
874 0 : hin = c0
875 0 : hbrin = c0
876 0 : kavg = c0
877 0 : bphi_o = c0
878 0 : sloss = c0
879 :
880 : !-----------------------------------------------------------------
881 : ! brine dynamics
882 : !-----------------------------------------------------------------
883 :
884 0 : dhbr_top(n) = c0
885 0 : dhbr_bot(n) = c0
886 :
887 0 : if (tr_brine) then
888 0 : if (trcrn(nt_fbri,n) .le. c0) trcrn(nt_fbri,n) = c1
889 :
890 0 : dhice = c0
891 0 : call preflushing_changes (aicen (n), &
892 : vicen (n), vsnon (n), & ! LCOV_EXCL_LINE
893 : meltbn (n), melttn (n), & ! LCOV_EXCL_LINE
894 : congeln (n), snoicen(n), & ! LCOV_EXCL_LINE
895 : hin_old (n), dhice, & ! LCOV_EXCL_LINE
896 : trcrn(nt_fbri,n), & ! LCOV_EXCL_LINE
897 : dhbr_top(n), dhbr_bot(n), & ! LCOV_EXCL_LINE
898 : hbr_old, hin, & ! LCOV_EXCL_LINE
899 0 : hsn, first_ice(n) )
900 0 : if (icepack_warnings_aborted(subname)) return
901 :
902 : ! Requires the average ice permeability = kavg(:)
903 : ! and the surface ice porosity = zphi_o(:)
904 : ! computed in "compute_microS" or from "thermosaline_vertical"
905 :
906 0 : iDi(:,n) = c0
907 :
908 : call compute_microS_mushy (nilyr, nblyr, &
909 : bgrid, cgrid, igrid, & ! LCOV_EXCL_LINE
910 : trcrn(:,n), hin_old(n), hbr_old, & ! LCOV_EXCL_LINE
911 : sss, sst, bTiz(:,n), & ! LCOV_EXCL_LINE
912 : iTin(:), bphi(:,n), kavg, & ! LCOV_EXCL_LINE
913 : bphi_o, bSin(:), & ! LCOV_EXCL_LINE
914 : brine_sal(:), brine_rho(:), iphin(:), & ! LCOV_EXCL_LINE
915 : ibrine_rho(:), ibrine_sal(:), sice_rho(n), & ! LCOV_EXCL_LINE
916 0 : iDi(:,n) )
917 0 : if (icepack_warnings_aborted(subname)) return
918 :
919 0 : call update_hbrine (melttn(n), &
920 : meltsn (n), dt, & ! LCOV_EXCL_LINE
921 : hin, hsn, & ! LCOV_EXCL_LINE
922 : hin_old (n), hbrin, & ! LCOV_EXCL_LINE
923 :
924 : hbr_old, &
925 : trcrn(nt_fbri,n), & ! LCOV_EXCL_LINE
926 : dhbr_top(n), dhbr_bot(n), & ! LCOV_EXCL_LINE
927 : dh_top_chl, dh_bot_chl, & ! LCOV_EXCL_LINE
928 : kavg, bphi_o, & ! LCOV_EXCL_LINE
929 : darcy_V (n), darcy_V_chl, & ! LCOV_EXCL_LINE
930 : bphi(2,n), aice0, & ! LCOV_EXCL_LINE
931 0 : dh_direct)
932 0 : if (icepack_warnings_aborted(subname)) return
933 :
934 0 : hbri = hbri + hbrin * aicen(n)
935 :
936 : endif ! tr_brine
937 :
938 : !-----------------------------------------------------------------
939 : ! biogeochemistry
940 : !-----------------------------------------------------------------
941 :
942 0 : if (z_tracers) then
943 :
944 : call zbio (dt, nblyr, &
945 : nslyr, nilyr, & ! LCOV_EXCL_LINE
946 : melttn(n), & ! LCOV_EXCL_LINE
947 : meltsn(n), meltbn (n), & ! LCOV_EXCL_LINE
948 : congeln(n), snoicen(n), & ! LCOV_EXCL_LINE
949 : nbtrcr, fsnow, & ! LCOV_EXCL_LINE
950 : ntrcr, trcrn(1:ntrcr,n), & ! LCOV_EXCL_LINE
951 : bio_index(1:nbtrcr), aicen_init(n), & ! LCOV_EXCL_LINE
952 : vicen_init(n), vsnon_init(n), & ! LCOV_EXCL_LINE
953 : vicen(n), vsnon(n), & ! LCOV_EXCL_LINE
954 : aicen(n), flux_bio_atm(1:nbtrcr), & ! LCOV_EXCL_LINE
955 : n, n_algae, & ! LCOV_EXCL_LINE
956 : n_doc, n_dic, & ! LCOV_EXCL_LINE
957 : n_don, & ! LCOV_EXCL_LINE
958 : n_fed, n_fep, & ! LCOV_EXCL_LINE
959 : n_zaero, first_ice(n), & ! LCOV_EXCL_LINE
960 : hin_old(n), ocean_bio(1:nbtrcr), & ! LCOV_EXCL_LINE
961 : bphi(:,n), iphin, & ! LCOV_EXCL_LINE
962 : iDi(:,n), & ! LCOV_EXCL_LINE
963 : fswpenln(:,n), & ! LCOV_EXCL_LINE
964 : dhbr_top(n), dhbr_bot(n), & ! LCOV_EXCL_LINE
965 : zfswin(:,n), & ! LCOV_EXCL_LINE
966 : hbrin, hbr_old, & ! LCOV_EXCL_LINE
967 : ! darcy_V(n), darcy_V_chl, & ! LCOV_EXCL_LINE
968 : darcy_V(n), & ! LCOV_EXCL_LINE
969 : bgrid, & ! LCOV_EXCL_LINE
970 : igrid, icgrid, & ! LCOV_EXCL_LINE
971 : bphi_o, & ! LCOV_EXCL_LINE
972 : iTin, & ! LCOV_EXCL_LINE
973 : Zoo(:,n), & ! LCOV_EXCL_LINE
974 : flux_bio(1:nbtrcr), dh_direct, & ! LCOV_EXCL_LINE
975 : upNO, upNH, & ! LCOV_EXCL_LINE
976 : fbio_snoice, fbio_atmice, & ! LCOV_EXCL_LINE
977 : PP_net, ice_bio_net (1:nbtrcr), & ! LCOV_EXCL_LINE
978 0 : snow_bio_net(1:nbtrcr),grow_net )
979 0 : if (icepack_warnings_aborted(subname)) return
980 :
981 0 : elseif (skl_bgc) then
982 :
983 : call sklbio (dt, ntrcr, &
984 : nbtrcr, n_algae, & ! LCOV_EXCL_LINE
985 : n_doc, & ! LCOV_EXCL_LINE
986 : n_dic, n_don, & ! LCOV_EXCL_LINE
987 : n_fed, n_fep, & ! LCOV_EXCL_LINE
988 : flux_bio (1:nbtrcr), ocean_bio(1:nbtrcr), & ! LCOV_EXCL_LINE
989 : aicen (n), & ! LCOV_EXCL_LINE
990 : meltbn (n), congeln (n), & ! LCOV_EXCL_LINE
991 : fswthrun (n), first_ice(n), & ! LCOV_EXCL_LINE
992 : trcrn (1:ntrcr,n), & ! LCOV_EXCL_LINE
993 : PP_net, upNO, & ! LCOV_EXCL_LINE
994 0 : upNH, grow_net )
995 0 : if (icepack_warnings_aborted(subname)) return
996 :
997 : endif ! skl_bgc
998 :
999 0 : first_ice(n) = .false.
1000 :
1001 : endif ! aicen > puny
1002 : enddo ! ncat
1003 :
1004 : end subroutine icepack_biogeochemistry
1005 :
1006 : !=======================================================================
1007 : !autodocument_start icepack_load_ocean_bio_array
1008 : ! basic initialization for ocean_bio_all
1009 :
1010 0 : subroutine icepack_load_ocean_bio_array(max_nbtrcr, &
1011 : max_algae, max_don, max_doc, max_dic, max_aero, max_fe, & ! LCOV_EXCL_LINE
1012 : nit, amm, sil, dmsp, dms, algalN, & ! LCOV_EXCL_LINE
1013 0 : doc, don, dic, fed, fep, zaeros, ocean_bio_all, hum)
1014 :
1015 : integer (kind=int_kind), intent(in) :: &
1016 : max_algae , & ! maximum number of algal types ! LCOV_EXCL_LINE
1017 : max_dic , & ! maximum number of dissolved inorganic carbon types ! LCOV_EXCL_LINE
1018 : max_doc , & ! maximum number of dissolved organic carbon types ! LCOV_EXCL_LINE
1019 : max_don , & ! maximum number of dissolved organic nitrogen types ! LCOV_EXCL_LINE
1020 : max_fe , & ! maximum number of iron types ! LCOV_EXCL_LINE
1021 : max_aero , & ! maximum number of aerosols ! LCOV_EXCL_LINE
1022 : max_nbtrcr ! maximum number of bio tracers
1023 :
1024 : real (kind=dbl_kind), intent(in) :: &
1025 : nit , & ! ocean nitrate (mmol/m^3) ! LCOV_EXCL_LINE
1026 : amm , & ! ammonia/um (mmol/m^3) ! LCOV_EXCL_LINE
1027 : sil , & ! silicate (mmol/m^3) ! LCOV_EXCL_LINE
1028 : dmsp , & ! dmsp (mmol/m^3) ! LCOV_EXCL_LINE
1029 : dms , & ! dms (mmol/m^3) ! LCOV_EXCL_LINE
1030 : hum ! humic material (mmol/m^3)
1031 :
1032 : real (kind=dbl_kind), dimension (max_algae), intent(in) :: &
1033 : algalN ! ocean algal nitrogen (mmol/m^3) (diatoms, phaeo, pico)
1034 :
1035 : real (kind=dbl_kind), dimension (max_doc), intent(in) :: &
1036 : doc ! ocean doc (mmol/m^3) (proteins, EPS, lipid)
1037 :
1038 : real (kind=dbl_kind), dimension (max_don), intent(in) :: &
1039 : don ! ocean don (mmol/m^3)
1040 :
1041 : real (kind=dbl_kind), dimension (max_dic), intent(in) :: &
1042 : dic ! ocean dic (mmol/m^3)
1043 :
1044 : real (kind=dbl_kind), dimension (max_fe), intent(in) :: &
1045 : fed, fep ! ocean disolved and particulate fe (nM)
1046 :
1047 : real (kind=dbl_kind), dimension (max_aero), intent(in) :: &
1048 : zaeros ! ocean aerosols (mmol/m^3)
1049 :
1050 : real (kind=dbl_kind), dimension (max_nbtrcr), intent(inout) :: &
1051 : ocean_bio_all ! fixed order, all values even for tracers false
1052 :
1053 : !autodocument_end
1054 :
1055 : ! local variables
1056 :
1057 : integer (kind=int_kind) :: &
1058 : k, ks ! tracer indices
1059 :
1060 : character(len=*),parameter :: subname='(icepack_load_ocean_bio_array)'
1061 :
1062 0 : ocean_bio_all(:) = c0
1063 :
1064 0 : do k = 1, max_algae
1065 0 : ocean_bio_all(k) = algalN(k) ! N
1066 0 : ks = max_algae + max_doc + max_dic + 1
1067 0 : ocean_bio_all(ks + k) = R_chl2N(k)*algalN(k)!chl
1068 : enddo
1069 :
1070 0 : ks = max_algae + 1
1071 0 : do k = 1, max_doc
1072 0 : ocean_bio_all(ks + k) = doc(k) ! doc
1073 : enddo
1074 0 : ks = ks + max_doc
1075 0 : do k = 1, max_dic
1076 0 : ocean_bio_all(ks + k) = dic(k) ! dic
1077 : enddo
1078 :
1079 0 : ks = 2*max_algae + max_doc + max_dic + 7
1080 0 : do k = 1, max_don
1081 0 : ocean_bio_all(ks + k) = don(k) ! don
1082 : enddo
1083 :
1084 0 : ks = max_algae + 1
1085 0 : ocean_bio_all(ks) = nit ! nit
1086 :
1087 0 : ks = 2*max_algae + max_doc + 2 + max_dic
1088 0 : ocean_bio_all(ks) = amm ! Am
1089 0 : ks = ks + 1
1090 0 : ocean_bio_all(ks) = sil ! Sil
1091 0 : ks = ks + 1
1092 0 : ocean_bio_all(ks) = R_S2N(1)*algalN(1) & ! DMSPp
1093 : + R_S2N(2)*algalN(2) & ! LCOV_EXCL_LINE
1094 0 : + R_S2N(3)*algalN(3)
1095 0 : ks = ks + 1
1096 0 : ocean_bio_all(ks) = dmsp ! DMSPd
1097 0 : ks = ks + 1
1098 0 : ocean_bio_all(ks) = dms ! DMS
1099 0 : ks = ks + 1
1100 0 : ocean_bio_all(ks) = nit ! PON
1101 0 : ks = 2*max_algae + max_doc + 7 + max_dic + max_don
1102 0 : do k = 1, max_fe
1103 0 : ocean_bio_all(ks + k) = fed(k) ! fed
1104 : enddo
1105 0 : ks = ks + max_fe
1106 0 : do k = 1, max_fe
1107 0 : ocean_bio_all(ks + k) = fep(k) ! fep
1108 : enddo
1109 0 : ks = ks + max_fe
1110 0 : do k = 1, max_aero
1111 0 : ocean_bio_all(ks+k) = zaeros(k) ! zaero
1112 : enddo
1113 0 : ks = ks + max_aero + 1
1114 0 : ocean_bio_all(ks) = hum ! humics
1115 :
1116 0 : end subroutine icepack_load_ocean_bio_array
1117 :
1118 : !=======================================================================
1119 : !autodocument_start icepack_init_ocean_bio
1120 : ! Initialize ocean concentration
1121 :
1122 0 : subroutine icepack_init_ocean_bio (amm, dmsp, dms, algalN, doc, dic, don, &
1123 : fed, fep, hum, nit, sil, zaeros, max_dic, max_don, max_fe, max_aero,& ! LCOV_EXCL_LINE
1124 0 : CToN, CToN_DON)
1125 :
1126 : integer (kind=int_kind), intent(in) :: &
1127 : max_dic, & ! LCOV_EXCL_LINE
1128 : max_don, & ! LCOV_EXCL_LINE
1129 : max_fe, & ! LCOV_EXCL_LINE
1130 : max_aero
1131 :
1132 : real (kind=dbl_kind), intent(out):: &
1133 : amm , & ! ammonium ! LCOV_EXCL_LINE
1134 : dmsp , & ! DMSPp ! LCOV_EXCL_LINE
1135 : dms , & ! DMS ! LCOV_EXCL_LINE
1136 : hum , & ! humic material ! LCOV_EXCL_LINE
1137 : nit , & ! nitrate ! LCOV_EXCL_LINE
1138 : sil ! silicate
1139 :
1140 : real (kind=dbl_kind), dimension(:), intent(out):: &
1141 : algalN , & ! algae ! LCOV_EXCL_LINE
1142 : doc , & ! DOC ! LCOV_EXCL_LINE
1143 : dic , & ! DIC ! LCOV_EXCL_LINE
1144 : don , & ! DON ! LCOV_EXCL_LINE
1145 : fed , & ! Dissolved Iron ! LCOV_EXCL_LINE
1146 : fep , & ! Particulate Iron ! LCOV_EXCL_LINE
1147 : zaeros ! BC and dust
1148 :
1149 : real (kind=dbl_kind), dimension(:), intent(inout), optional :: &
1150 : CToN , & ! carbon to nitrogen ratio for algae ! LCOV_EXCL_LINE
1151 : CToN_DON ! nitrogen to carbon ratio for proteins
1152 :
1153 : !autodocument_end
1154 :
1155 : ! local variables
1156 :
1157 : integer (kind=int_kind) :: &
1158 : k
1159 :
1160 : character(len=*),parameter :: subname='(icepack_init_ocean_bio)'
1161 :
1162 0 : if (present(CToN)) then
1163 0 : CToN(1) = R_C2N(1)
1164 0 : CToN(2) = R_C2N(2)
1165 0 : CToN(3) = R_C2N(3)
1166 : endif
1167 :
1168 0 : if (present(CToN_DON)) then
1169 0 : CToN_DON(1) = R_C2N_DON(1)
1170 : endif
1171 :
1172 0 : amm = c1 ! ISPOL < 1 mmol/m^3
1173 0 : dmsp = p1
1174 0 : dms = p1
1175 0 : algalN(1) = c1 !0.0026_dbl_kind ! ISPOL, Lannuzel 2013(pennate)
1176 0 : algalN(2) = 0.0057_dbl_kind ! ISPOL, Lannuzel 2013(small plankton)
1177 0 : algalN(3) = 0.0027_dbl_kind ! ISPOL, Lannuzel 2013(Phaeocystis)
1178 : ! 0.024_dbl_kind ! 5% of 1 mgchl/m^3
1179 0 : doc(1) = 16.2_dbl_kind ! 18% saccharides
1180 0 : doc(2) = 9.0_dbl_kind ! lipids
1181 0 : doc(3) = c1 !
1182 0 : do k = 1, max_dic
1183 0 : dic(k) = c1
1184 : enddo
1185 0 : do k = 1, max_don
1186 0 : don(k) = 12.9_dbl_kind
1187 : ! 64.3_dbl_kind ! 72% Total DOC~90 mmolC/m^3 ISPOL with N:C of 0.2
1188 : enddo
1189 : !ki = 1
1190 : !if (trim(fe_data_type) == 'clim') ki = 2
1191 0 : do k = 1, max_fe ! ki, max_fe
1192 0 : fed(k) = 0.4_dbl_kind ! c1 (nM) Lannuzel2007 DFe,
1193 : ! range 0.14-2.6 (nM) van der Merwe 2011
1194 : ! Tagliabue 2012 (0.4 nM)
1195 0 : fep(k) = c2 ! (nM) van der Merwe 2011
1196 : ! (0.6 to 2.9 nM ocean)
1197 : enddo
1198 0 : hum = c1 ! mmol C/m^3
1199 0 : nit = 12.0_dbl_kind
1200 0 : sil = 25.0_dbl_kind
1201 0 : do k = 1, max_aero
1202 0 : zaeros(k) = c0
1203 : enddo
1204 :
1205 :
1206 0 : end subroutine icepack_init_ocean_bio
1207 :
1208 : !=======================================================================
1209 :
1210 : end module icepack_zbgc
1211 :
1212 : !=======================================================================
|