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
13 :
14 : !use icepack_parameters, only: c0, c1, c2, p001, p1, p5, puny
15 : !use icepack_parameters, only: depressT, rhosi, min_salin, salt_loss
16 : !use icepack_parameters, only: fr_resp, algal_vel, R_dFe2dust, dustFe_sol, T_max
17 : !use icepack_parameters, only: op_dep_min, fr_graze_s, fr_graze_e, fr_mort2min, fr_dFe
18 : !use icepack_parameters, only: k_nitrif, t_iron_conv, max_loss, max_dfe_doc1
19 : !use icepack_parameters, only: fr_resp_s, y_sk_DMS, t_sk_conv, t_sk_ox
20 : !use icepack_parameters, only: scale_bgc, ktherm, skl_bgc
21 : !use icepack_parameters, only: z_tracers, fsal, conserv_check
22 :
23 : use icepack_tracers, only: max_algae, max_dic, max_doc, max_don, max_fe
24 : use icepack_tracers, only: max_aero, max_nbtrcr
25 : use icepack_tracers, only: ncat, nilyr, nslyr, nblyr, nbtrcr, ntrcr, ntrcr_o
26 : use icepack_tracers, only: tr_brine, nt_fbri, nt_qice, nt_Tsfc
27 : use icepack_tracers, only: tr_zaero, tr_bgc_Nit, tr_bgc_N
28 : use icepack_tracers, only: tr_bgc_DON, tr_bgc_C, tr_bgc_chl
29 : use icepack_tracers, only: tr_bgc_Am, tr_bgc_Sil, tr_bgc_DMS
30 : use icepack_tracers, only: tr_bgc_Fe, tr_bgc_hum, tr_bgc_PON
31 : use icepack_tracers, only: nt_bgc_Nit, nlt_bgc_Nit
32 : use icepack_tracers, only: nt_bgc_N, nlt_bgc_N, nt_bgc_Am, nlt_bgc_Am
33 : use icepack_tracers, only: nt_bgc_DMSPp, nlt_bgc_DMSPp, nt_bgc_Sil, nlt_bgc_Sil
34 : use icepack_tracers, only: nt_bgc_DMSPd, nlt_bgc_DMSPd, nt_bgc_DMS, nlt_bgc_DMS
35 : use icepack_tracers, only: nt_bgc_hum, nlt_bgc_hum, nt_bgc_PON, nlt_bgc_PON
36 : use icepack_tracers, only: nt_bgc_C, nlt_bgc_C, nt_bgc_chl, nlt_bgc_chl
37 : use icepack_tracers, only: nt_bgc_DOC, nlt_bgc_DOC, nt_bgc_DON, nlt_bgc_DON
38 : use icepack_tracers, only: nt_bgc_DIC, nlt_bgc_DIC, nt_bgc_Fed, nlt_bgc_Fed
39 : use icepack_tracers, only: nt_sice, bio_index, bio_index_o
40 : use icepack_tracers, only: nt_zaero, nlt_zaero, nt_bgc_Fep, nlt_bgc_Fep
41 : use icepack_tracers, only: nlt_zaero_sw, nlt_chl_sw, nt_zbgc_frac
42 : use icepack_tracers, only: n_algae, n_doc, n_dic, n_don, n_fed, n_fep, n_zaero
43 :
44 : use icepack_zbgc_shared, only: zbgc_init_frac
45 : use icepack_zbgc_shared, only: zbgc_frac_init
46 : use icepack_zbgc_shared, only: bgrid, cgrid, igrid, icgrid
47 : use icepack_zbgc_shared, only: bgc_tracer_type, remap_zbgc
48 : use icepack_zbgc_shared, only: R_S2N, R_Si2N, R_Fe2C, R_Fe2N, R_Fe2DON, R_Fe2DOC
49 : use icepack_zbgc_shared, only: chlabs, alpha2max_low, beta2max
50 : use icepack_zbgc_shared, only: mu_max, grow_Tdep, fr_graze
51 : use icepack_zbgc_shared, only: mort_pre, mort_Tdep, k_exude
52 : use icepack_zbgc_shared, only: K_Nit, K_Am, K_Sil, K_Fe
53 : use icepack_zbgc_shared, only: f_don, kn_bac, f_don_Am, f_doc
54 : use icepack_zbgc_shared, only: f_exude, k_bac
55 : use icepack_zbgc_shared, only: tau_ret, tau_rel
56 : use icepack_zbgc_shared, only: R_C2N, R_chl2N, f_abs_chl, R_C2N_DON
57 : use icepack_zbgc_shared, only: doc_pool_fractions
58 : use icepack_zbgc_shared, only: algaltype, doctype, dictype
59 : use icepack_zbgc_shared, only: dontype, fedtype, feptype, zaerotype
60 :
61 : use icepack_warnings, only: warnstr, icepack_warnings_add
62 : use icepack_warnings, only: icepack_warnings_setabort, icepack_warnings_aborted
63 :
64 : use icepack_brine, only: preflushing_changes, compute_microS_mushy
65 : use icepack_brine, only: update_hbrine
66 : use icepack_algae, only: zbio, sklbio
67 : use icepack_therm_shared, only: calculate_Tin_from_qin
68 : use icepack_itd, only: column_sum, column_conservation_check
69 :
70 : implicit none
71 :
72 : private
73 : public :: add_new_ice_bgc, &
74 : lateral_melt_bgc, &
75 : icepack_init_bgc, &
76 : icepack_init_zbgc, &
77 : icepack_biogeochemistry, &
78 : icepack_load_ocean_bio_array, &
79 : icepack_init_ocean_bio
80 :
81 : !=======================================================================
82 :
83 : contains
84 :
85 : !=======================================================================
86 :
87 : ! Adjust biogeochemical tracers when new frazil ice forms
88 :
89 182052 : subroutine add_new_ice_bgc (dt, ncats, &
90 364104 : aicen_init, vicen_init, vi0_init, &
91 182052 : aicen, vicen, vin0new, &
92 182052 : trcrn, &
93 182052 : ocean_bio, flux_bio, hsurp, &
94 182052 : d_an_tot)
95 :
96 : integer (kind=int_kind), intent(in) :: &
97 : ncats ! 1 without floe size distribution or ncat
98 :
99 : real (kind=dbl_kind), intent(in) :: &
100 : dt ! time step (s)
101 :
102 : real (kind=dbl_kind), dimension (:), intent(in) :: &
103 : aicen_init , & ! initial concentration of ice
104 : vicen_init , & ! initial volume per unit area of ice (m)
105 : aicen , & ! concentration of ice
106 : vicen , & ! volume per unit area of ice (m)
107 : d_an_tot
108 :
109 : real (kind=dbl_kind), dimension (:,:), intent(inout) :: &
110 : trcrn ! ice tracers
111 :
112 : real (kind=dbl_kind), dimension (:), intent(in) :: &
113 : vin0new ! ice tracers
114 :
115 : real (kind=dbl_kind), intent(in) :: &
116 : vi0_init ! volume of new ice added to cat 1 (intial)
117 :
118 : real (kind=dbl_kind), intent(in) :: &
119 : hsurp ! thickness of new ice added to each cat
120 :
121 : real (kind=dbl_kind), dimension (:), intent(inout) :: &
122 : flux_bio ! tracer flux to ocean from biology (mmol/m^2/s)
123 :
124 : real (kind=dbl_kind), dimension (:), intent(in) :: &
125 : ocean_bio ! ocean concentration of biological tracer
126 :
127 : ! local
128 :
129 : integer (kind=int_kind) :: &
130 : location , & ! 1 (add frazil to bottom), 0 (add frazil throughout)
131 : m , & ! bio index
132 : n , & ! ice category index
133 : k , & ! ice layer index
134 : nbiolayer
135 :
136 : real (kind=dbl_kind) :: &
137 : vbri1 , & ! starting volume of existing brine
138 : vbri_init , & ! brine volume summed over categories
139 : vbri_final ! brine volume summed over categories
140 :
141 : real (kind=dbl_kind) :: &
142 : vsurp , & ! volume of new ice added to each cat
143 : vtmp ! total volume of new and old ice
144 :
145 : real (kind=dbl_kind), dimension (ncat) :: &
146 182052 : vbrin ! trcrn(nt_fbri,n)*vicen(n)
147 :
148 : real (kind=dbl_kind) :: &
149 : vice_new , & ! vicen_init + vsurp
150 : bio0new ! ocean_bio * zbgc_init_fac
151 :
152 : real (kind=dbl_kind) :: &
153 : Tmlts ! melting temperature (oC)
154 :
155 : character (len=char_len) :: &
156 : fieldid ! field identifier
157 :
158 : character(len=*),parameter :: subname='(add_new_ice_bgc)'
159 :
160 : real (kind=dbl_kind), dimension (nblyr+1) :: &
161 182052 : zspace ! vertical grid spacing
162 : !-----------------------------------------------------------------
163 : ! brine
164 : !-----------------------------------------------------------------
165 :
166 1638468 : zspace(:) = c1/real(nblyr,kind=dbl_kind)
167 182052 : zspace(1) = p5*zspace(2)
168 182052 : zspace(nblyr+1) = zspace(1)
169 1092312 : vbrin(:) = c0
170 1092312 : do n = 1, ncat
171 910260 : vbrin(n) = vicen_init(n)
172 1092312 : if (tr_brine) vbrin(n) = trcrn(nt_fbri,n)*vicen_init(n)
173 : enddo
174 :
175 182052 : call column_sum (ncat, vbrin, vbri_init)
176 182052 : if (icepack_warnings_aborted(subname)) return
177 :
178 182052 : vbri_init = vbri_init + vi0_init
179 2773692 : do k = 1, nbtrcr
180 : flux_bio(k) = flux_bio(k) &
181 2773692 : - vi0_init/dt*ocean_bio(k)*zbgc_init_frac(k)
182 : enddo
183 : !-----------------------------------------------------------------
184 : ! Distribute bgc in new ice volume among all ice categories by
185 : ! increasing ice thickness, leaving ice area unchanged.
186 : !-----------------------------------------------------------------
187 :
188 : ! Diffuse_bio handles concentration changes from ice growth/melt
189 : ! ice area does not change
190 : ! add salt to the bottom , location = 1
191 :
192 182052 : vsurp = c0
193 182052 : vtmp = c0
194 :
195 1092312 : do n = 1,ncat
196 :
197 1092312 : if (hsurp > c0) then
198 :
199 11305 : vtmp = vbrin(n)
200 11305 : vsurp = hsurp * aicen_init(n)
201 11305 : vbrin(n) = vbrin(n) + vsurp
202 11305 : vice_new = vicen_init(n) + vsurp
203 :
204 11305 : if (tr_brine .and. vice_new > puny) then !c0) then
205 11193 : trcrn(nt_fbri,n) = vbrin(n)/vice_new
206 112 : elseif (tr_brine .and. vicen(n) <= c0) then
207 112 : trcrn(nt_fbri,n) = c1
208 : endif
209 :
210 11305 : if (nbtrcr > 0) then
211 48220 : do m = 1, nbtrcr
212 36915 : bio0new = ocean_bio(m)*zbgc_init_frac(m)
213 36915 : nbiolayer = nblyr+1
214 : call update_vertical_bio_tracers(nbiolayer, trcrn(bio_index(m):bio_index(m) + nblyr,n), &
215 48220 : vtmp, vbrin(n), bio0new,zspace(:))
216 : enddo !nbtrcr
217 11305 : if (icepack_warnings_aborted(subname)) return
218 : endif ! nbtrcr
219 : endif ! hsurp > 0
220 : enddo ! n
221 :
222 : !-----------------------------------------------------------------
223 : ! Combine bgc in new ice grown in open water with category 1 ice.
224 : !-----------------------------------------------------------------
225 364104 : do n = 1, ncats
226 364104 : if (vin0new(n) > c0 .and. d_an_tot(n) > c0) then
227 :
228 29074 : vbri1 = vbrin(n)
229 29074 : vbrin(n) = vbrin(n) + vin0new(n)
230 29074 : if (tr_brine .and. vicen(n) > puny) then
231 29074 : trcrn(nt_fbri,n) = vbrin(n)/vicen(n)
232 0 : elseif (tr_brine .and. vicen(n) <= puny) then
233 0 : trcrn(nt_fbri,n) = c1
234 : endif
235 :
236 : ! Diffuse_bio handles concentration changes from ice growth/melt
237 : ! ice area changes
238 : ! add salt throughout, location = 0
239 :
240 29074 : if (nbtrcr > 0 .and. vbrin(n) > puny) then
241 129178 : do m = 1, nbtrcr
242 100104 : bio0new = ocean_bio(m)*zbgc_init_frac(m)
243 930010 : do k = 1, nblyr+1
244 : trcrn(bio_index(m) + k-1,n) = &
245 900936 : (trcrn(bio_index(m) + k-1,n)*vbri1 + bio0new * vin0new(n))/vbrin(n)
246 : enddo
247 : enddo
248 :
249 29074 : if (icepack_warnings_aborted(subname)) return
250 :
251 : endif ! nbtrcr > 0
252 : endif ! vin0new(n) > 0
253 : enddo ! n = 1,ncats
254 :
255 182052 : 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, &
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 141858 : subroutine lateral_melt_bgc (dt, &
274 141858 : rsiden, vicen_init,&
275 141858 : trcrn, flux_bio)
276 :
277 : real (kind=dbl_kind), intent(in) :: &
278 : dt ! time step (s)
279 :
280 : real (kind=dbl_kind), dimension(:), intent(in) :: &
281 : vicen_init! volume per unit area of ice (m)
282 :
283 : real (kind=dbl_kind), dimension (:,:), intent(in) :: &
284 : trcrn ! tracer array
285 :
286 : real (kind=dbl_kind), dimension(:), intent(in) :: &
287 : rsiden ! fraction of ice that melts laterally
288 :
289 : real (kind=dbl_kind), dimension(:), intent(inout) :: &
290 : flux_bio ! biology tracer flux from layer bgc (mmol/m^2/s)
291 :
292 : ! local variables
293 :
294 : real (kind=dbl_kind) :: &
295 : total_bio_initial, & ! initial column tracer concentration (mmol/m2)
296 : total_bio_final ! final column tracer concentration (mmol/m20
297 :
298 : integer (kind=int_kind) :: &
299 : k , & ! layer index
300 : m , & !
301 : n ! category index
302 :
303 : real (kind=dbl_kind), dimension (nblyr+1) :: &
304 283716 : zspace ! vertical grid spacing
305 :
306 : character(len=*),parameter :: subname='(lateral_melt_bgc)'
307 :
308 1276722 : zspace(:) = c1/real(nblyr,kind=dbl_kind)
309 141858 : zspace(1) = p5*zspace(2)
310 141858 : zspace(nblyr+1) = zspace(1)
311 :
312 2411094 : do m = 1, nbtrcr
313 13757274 : do n = 1, ncat
314 104384856 : do k = 1, nblyr+1
315 : flux_bio(m) = flux_bio(m) + trcrn(nt_fbri,n) &
316 : * vicen_init(n)*zspace(k)*trcrn(bio_index(m)+k-1,n) &
317 102115620 : * rsiden(n)/dt
318 : enddo
319 : enddo
320 : enddo
321 :
322 141858 : end subroutine lateral_melt_bgc
323 :
324 : !=======================================================================
325 : !autodocument_start icepack_init_bgc
326 : !
327 36 : subroutine icepack_init_bgc( &
328 36 : sicen, trcrn, sss, ocean_bio_all, DOCPoolFractions)
329 :
330 : real (kind=dbl_kind), dimension(nilyr, ncat), intent(in) :: &
331 : sicen ! salinity on the cice grid
332 :
333 : real (kind=dbl_kind), dimension (:,:), intent(inout) :: &
334 : trcrn ! subset of tracer array (only bgc)
335 :
336 : real (kind=dbl_kind), intent(in) :: &
337 : sss ! sea surface salinity (ppt)
338 :
339 : real (kind=dbl_kind), dimension (:), intent(in) :: &
340 : ocean_bio_all ! fixed order, all values even for tracers false
341 :
342 : real (kind=dbl_kind), dimension (:), optional, intent(out) :: &
343 : DOCPoolFractions ! Fraction of DOC in polysacharids, lipids, and proteins
344 :
345 : !autodocument_end
346 :
347 : ! local variables
348 :
349 : integer (kind=int_kind) :: &
350 : k , & ! vertical index
351 : n , & ! category index
352 : mm ! bio tracer index
353 :
354 : real (kind=dbl_kind), dimension (ntrcr+2) :: &
355 36 : trtmp ! temporary, remapped tracers
356 :
357 : character(len=*),parameter :: subname='(icepack_init_bgc)'
358 :
359 : !-----------------------------------------------------------------------------
360 : ! Skeletal Layer Model
361 : ! All bgc tracers are Bulk quantities in units of mmol or mg per m^3
362 : ! The skeletal layer model assumes a constant
363 : ! layer depth (sk_l) and porosity (phi_sk)
364 : !-----------------------------------------------------------------------------
365 36 : if (.not. restartbgc) then
366 :
367 36 : if (skl_bgc) then
368 :
369 0 : do n = 1,ncat
370 0 : do mm = 1,nbtrcr
371 : ! bulk concentration (mmol or mg per m^3, or 10^-3 mmol/m^3)
372 0 : trcrn(bio_index(mm)-ntrcr_o, n) = ocean_bio_all(bio_index_o(mm))
373 : enddo ! nbtrcr
374 : enddo ! n
375 :
376 : !-----------------------------------------------------------------------------
377 : ! zbgc Model
378 : ! All bgc tracers are Bulk quantities in units of mmol or mg per m^3
379 : ! The vertical layer model uses prognosed porosity and layer depth
380 : !-----------------------------------------------------------------------------
381 :
382 : else ! not skl_bgc
383 :
384 36 : if (scale_bgc .and. ktherm == 2) then
385 0 : trtmp(:) = c0
386 0 : do n = 1,ncat
387 : call remap_zbgc(nilyr, &
388 : 1, &
389 : sicen(:,n), trtmp, &
390 : 0, nblyr+1, &
391 : c1, c1, &
392 : cgrid(2:nilyr+1), &
393 : igrid(1:nblyr+1), &
394 0 : sicen(1,n) )
395 0 : if (icepack_warnings_aborted(subname)) return
396 :
397 0 : do mm = 1,nbtrcr
398 0 : do k = 1, nblyr + 1
399 : trcrn(bio_index(mm)+k-1-ntrcr_o,n) = &
400 0 : (trtmp(k)/sss*ocean_bio_all(bio_index_o(mm)))
401 0 : trcrn(bio_index(mm)+nblyr+1-ntrcr_o:bio_index(mm)+nblyr+2-ntrcr_o,n) = c0 ! snow
402 : enddo ! k
403 : enddo ! mm
404 : enddo ! n
405 :
406 36 : elseif (nbtrcr > 0 .and. nt_fbri > 0) then ! not scale_bgc
407 :
408 216 : do n = 1,ncat
409 2976 : do mm = 1,nbtrcr
410 24840 : do k = 1, nblyr+1
411 : trcrn(bio_index(mm)+k-1-ntrcr_o,n) = ocean_bio_all(bio_index_o(mm)) &
412 22080 : * zbgc_init_frac(mm)
413 69000 : trcrn(bio_index(mm)+nblyr+1-ntrcr_o:bio_index(mm)+nblyr+2-ntrcr_o,n) = c0 ! snow
414 : enddo ! k
415 2940 : trcrn(nt_zbgc_frac-1+mm-ntrcr_o,n) = zbgc_frac_init(mm)
416 : enddo ! mm
417 : enddo ! n
418 :
419 : endif ! scale_bgc
420 : endif ! skl_bgc
421 : endif ! restart
422 :
423 36 : if (present(DOCPoolFractions)) then
424 0 : DOCPoolFractions(:) = c1
425 0 : if (.not. use_macromolecules) then
426 0 : do mm = 1,max_doc
427 0 : DOCPoolFractions(mm) = doc_pool_fractions(mm)
428 : end do
429 : end if
430 : endif
431 :
432 : end subroutine icepack_init_bgc
433 :
434 : !=======================================================================
435 : !autodocument_start icepack_init_zbgc
436 : !
437 :
438 83 : subroutine icepack_init_zbgc (&
439 83 : zbgc_frac_init_in, zbgc_init_frac_in, tau_ret_in, tau_rel_in, bgc_tracer_type_in)
440 :
441 : real (kind=dbl_kind), dimension (:), intent(in), optional :: zbgc_frac_init_in(:) ! initializes mobile fraction
442 : real (kind=dbl_kind), dimension (:), intent(in), optional :: bgc_tracer_type_in(:) ! described tracer in mobile or stationary phases
443 : real (kind=dbl_kind), dimension (:), intent(in), optional :: zbgc_init_frac_in(:) ! fraction of ocean tracer concentration in new ice
444 : real (kind=dbl_kind), dimension (:), intent(in), optional :: tau_ret_in(:) ! retention timescale (s), mobile to stationary phase
445 : real (kind=dbl_kind), dimension (:), intent(in), optional :: tau_rel_in(:) ! release timescale (s), stationary to mobile phase
446 :
447 : !autodocument_end
448 :
449 : character(len=*),parameter :: subname='(icepack_init_zbgc)'
450 :
451 : !--------
452 :
453 2490 : if (present(zbgc_frac_init_in)) zbgc_frac_init(:) = zbgc_frac_init_in(:)
454 2490 : if (present(bgc_tracer_type_in)) bgc_tracer_type(:) = bgc_tracer_type_in(:)
455 2490 : if (present(zbgc_init_frac_in)) zbgc_init_frac(:) = zbgc_init_frac_in(:)
456 2490 : if (present(tau_ret_in)) tau_ret(:) = tau_ret_in(:)
457 2490 : if (present(tau_rel_in)) tau_rel(:) = tau_rel_in(:)
458 :
459 83 : R_Si2N(1) = ratio_Si2N_diatoms
460 83 : R_Si2N(2) = ratio_Si2N_sp
461 83 : R_Si2N(3) = ratio_Si2N_phaeo
462 :
463 83 : R_S2N(1) = ratio_S2N_diatoms
464 83 : R_S2N(2) = ratio_S2N_sp
465 83 : R_S2N(3) = ratio_S2N_phaeo
466 :
467 83 : R_Fe2C(1) = ratio_Fe2C_diatoms
468 83 : R_Fe2C(2) = ratio_Fe2C_sp
469 83 : R_Fe2C(3) = ratio_Fe2C_phaeo
470 :
471 83 : R_Fe2N(1) = ratio_Fe2N_diatoms
472 83 : R_Fe2N(2) = ratio_Fe2N_sp
473 83 : R_Fe2N(3) = ratio_Fe2N_phaeo
474 :
475 83 : R_C2N(1) = ratio_C2N_diatoms
476 83 : R_C2N(2) = ratio_C2N_sp
477 83 : R_C2N(3) = ratio_C2N_phaeo
478 :
479 83 : R_chl2N(1) = ratio_chl2N_diatoms
480 83 : R_chl2N(2) = ratio_chl2N_sp
481 83 : R_chl2N(3) = ratio_chl2N_phaeo
482 :
483 83 : F_abs_chl(1) = F_abs_chl_diatoms
484 83 : F_abs_chl(2) = F_abs_chl_sp
485 83 : F_abs_chl(3) = F_abs_chl_phaeo
486 :
487 83 : R_Fe2DON(1) = ratio_Fe2DON
488 83 : R_C2N_DON(1) = ratio_C2N_proteins
489 :
490 83 : R_Fe2DOC(1) = ratio_Fe2DOC_s
491 83 : R_Fe2DOC(2) = ratio_Fe2DOC_l
492 83 : R_Fe2DOC(3) = c0
493 :
494 83 : chlabs(1) = chlabs_diatoms
495 83 : chlabs(2) = chlabs_sp
496 83 : chlabs(3) = chlabs_phaeo
497 :
498 83 : alpha2max_low(1) = alpha2max_low_diatoms
499 83 : alpha2max_low(2) = alpha2max_low_sp
500 83 : alpha2max_low(3) = alpha2max_low_phaeo
501 :
502 83 : beta2max(1) = beta2max_diatoms
503 83 : beta2max(2) = beta2max_sp
504 83 : beta2max(3) = beta2max_phaeo
505 :
506 83 : mu_max(1) = mu_max_diatoms
507 83 : mu_max(2) = mu_max_sp
508 83 : mu_max(3) = mu_max_phaeo
509 :
510 83 : grow_Tdep(1) = grow_Tdep_diatoms
511 83 : grow_Tdep(2) = grow_Tdep_sp
512 83 : grow_Tdep(3) = grow_Tdep_phaeo
513 :
514 83 : fr_graze(1) = fr_graze_diatoms
515 83 : fr_graze(2) = fr_graze_sp
516 83 : fr_graze(3) = fr_graze_phaeo
517 :
518 83 : mort_pre(1) = mort_pre_diatoms
519 83 : mort_pre(2) = mort_pre_sp
520 83 : mort_pre(3) = mort_pre_phaeo
521 :
522 83 : mort_Tdep(1) = mort_Tdep_diatoms
523 83 : mort_Tdep(2) = mort_Tdep_sp
524 83 : mort_Tdep(3) = mort_Tdep_phaeo
525 :
526 83 : k_exude(1) = k_exude_diatoms
527 83 : k_exude(2) = k_exude_sp
528 83 : k_exude(3) = k_exude_phaeo
529 :
530 83 : K_Nit(1) = K_Nit_diatoms
531 83 : K_Nit(2) = K_Nit_sp
532 83 : K_Nit(3) = K_Nit_phaeo
533 :
534 83 : K_Am(1) = K_Am_diatoms
535 83 : K_Am(2) = K_Am_sp
536 83 : K_Am(3) = K_Am_phaeo
537 :
538 83 : K_Sil(1) = K_Sil_diatoms
539 83 : K_Sil(2) = K_Sil_sp
540 83 : K_Sil(3) = K_Sil_phaeo
541 :
542 83 : K_Fe(1) = K_Fe_diatoms
543 83 : K_Fe(2) = K_Fe_sp
544 83 : K_Fe(3) = K_Fe_phaeo
545 :
546 83 : f_doc(:) = c0
547 83 : f_doc(1) = f_doc_s
548 83 : f_doc(2) = f_doc_l
549 :
550 83 : f_don(1) = f_don_protein
551 83 : kn_bac(1) = kn_bac_protein
552 83 : f_don_Am(1) = f_don_Am_protein
553 :
554 83 : f_exude(:) = c0
555 83 : f_exude(1) = f_exude_s
556 83 : f_exude(2) = f_exude_l
557 :
558 83 : k_bac(:) = c0
559 83 : k_bac(1) = k_bac_s
560 83 : k_bac(2) = k_bac_l
561 :
562 83 : algaltype(1) = algaltype_diatoms
563 83 : algaltype(2) = algaltype_sp
564 83 : algaltype(3) = algaltype_phaeo
565 :
566 83 : doctype(:) = c0
567 83 : doctype(1) = doctype_s
568 83 : doctype(2) = doctype_l
569 :
570 166 : dictype(:) = dictype_1
571 :
572 166 : dontype(:) = dontype_protein
573 :
574 249 : fedtype(:) = fedtype_1
575 249 : feptype(:) = feptype_1
576 :
577 83 : zaerotype(1) = zaerotype_bc1
578 83 : zaerotype(2) = zaerotype_bc2
579 83 : zaerotype(3) = zaerotype_dust1
580 83 : zaerotype(4) = zaerotype_dust2
581 83 : zaerotype(5) = zaerotype_dust3
582 83 : zaerotype(6) = zaerotype_dust4
583 :
584 83 : end subroutine icepack_init_zbgc
585 :
586 : !=======================================================================
587 : !autodocument_start icepack_biogeochemistry
588 : !
589 242736 : subroutine icepack_biogeochemistry(dt, &
590 0 : upNO, upNH, iDi, iki, zfswin, &
591 0 : darcy_V, grow_net, &
592 485472 : PP_net, hbri,dhbr_bot, dhbr_top, Zoo,&
593 485472 : fbio_snoice, fbio_atmice, ocean_bio_dh, ocean_bio, &
594 970944 : first_ice, fswpenln, bphi, bTiz, ice_bio_net, &
595 485472 : snow_bio_net, totalChla, fswthrun, &
596 242736 : meltbn, melttn, congeln, snoicen, &
597 242736 : sst, sss, Tf, fsnow, meltsn, & !hmix, &
598 728208 : hin_old, flux_bio, flux_bio_atm, &
599 728208 : aicen_init, vicen_init, aicen, vicen, vsnon, &
600 485472 : aice0, trcrn, vsnon_init, &
601 485472 : flux_bion, bioPorosityIceCell, &
602 242736 : bioSalinityIceCell, bioTemperatureIceCell)
603 :
604 : real (kind=dbl_kind), intent(in) :: &
605 : dt ! time step
606 :
607 : real (kind=dbl_kind), dimension (:), intent(inout) :: &
608 : fbio_snoice , & ! fluxes from snow to ice
609 : fbio_atmice , & ! fluxes from atm to ice
610 : dhbr_top , & ! brine top change
611 : dhbr_bot , & ! brine bottom change
612 : darcy_V , & ! darcy velocity positive up (m/s)
613 : hin_old , & ! old ice thickness
614 : ice_bio_net , & ! depth integrated tracer (mmol/m^2)
615 : snow_bio_net , & ! depth integrated snow tracer (mmol/m^2)
616 : flux_bio ! all bio fluxes to ocean
617 :
618 : real (kind=dbl_kind), dimension (:), intent(in) :: &
619 : ocean_bio ! contains the ocean bgc tracer concentrations in use (mmol/m^3)
620 :
621 : real (kind=dbl_kind), optional, dimension (:), intent(out) :: &
622 : ocean_bio_dh ! The ocean bgc tracer concentrations in use * brine thickness * phi (mmol/m^2)
623 :
624 : logical (kind=log_kind), dimension (:), intent(inout) :: &
625 : first_ice ! distinguishes ice that disappears (e.g. melts)
626 : ! and reappears (e.g. transport) in a grid cell
627 : ! during a single time step from ice that was
628 : ! there the entire time step (true until ice forms)
629 :
630 : real (kind=dbl_kind), optional, dimension (:,:), intent(out) :: &
631 : flux_bion ! per categeory ice to ocean biogeochemistry flux (mmol/m2/s)
632 :
633 : real (kind=dbl_kind), optional, dimension (:), intent(inout) :: &
634 : bioPorosityIceCell, & ! category average porosity on the interface bio grid
635 : bioSalinityIceCell, & ! (ppt) category average porosity on the interface bio grid
636 : bioTemperatureIceCell ! (oC) category average porosity on the interface bio grid
637 :
638 : real (kind=dbl_kind), dimension (:,:), intent(inout) :: &
639 : Zoo , & ! N losses accumulated in timestep (ie. zooplankton/bacteria)
640 : ! mmol/m^3
641 : bphi , & ! porosity of layers
642 : bTiz , & ! layer temperatures interpolated on bio grid (C)
643 : zfswin , & ! Shortwave flux into layers interpolated on bio grid (W/m^2)
644 : iDi , & ! igrid Diffusivity (m^2/s)
645 : iki , & ! Ice permeability (m^2)
646 : trcrn ! tracers
647 :
648 : real (kind=dbl_kind), intent(inout) :: &
649 : grow_net , & ! Specific growth rate (/s) per grid cell
650 : PP_net , & ! Total production (mg C/m^2/s) per grid cell
651 : hbri , & ! brine height, area-averaged for comparison with hi (m)
652 : upNO , & ! nitrate uptake rate (mmol/m^2/d) times aice
653 : upNH ! ammonium uptake rate (mmol/m^2/d) times aice
654 :
655 : real (kind=dbl_kind), optional, intent(inout) :: &
656 : totalChla ! ice integrated chla and summed over all algal groups (mg/m^2)
657 :
658 : real (kind=dbl_kind), dimension (:,:), intent(in) :: &
659 : fswpenln ! visible SW entering ice layers (W m-2)
660 :
661 : real (kind=dbl_kind), dimension (:), intent(in) :: &
662 : fswthrun , & ! SW through ice to ocean (W/m^2)
663 : meltsn , & ! snow melt in category n (m)
664 : melttn , & ! top melt in category n (m)
665 : meltbn , & ! bottom melt in category n (m)
666 : congeln , & ! congelation ice formation in category n (m)
667 : snoicen , & ! snow-ice formation in category n (m)
668 : flux_bio_atm, & ! all bio fluxes to ice from atmosphere
669 : aicen_init , & ! initial ice concentration, for linear ITD
670 : vicen_init , & ! initial ice volume (m), for linear ITD
671 : vsnon_init , & ! initial snow volume (m), for aerosol
672 : aicen , & ! concentration of ice
673 : vicen , & ! volume per unit area of ice (m)
674 : vsnon ! volume per unit area of snow (m)
675 :
676 : real (kind=dbl_kind), intent(in) :: &
677 : aice0 , & ! open water area fraction
678 : sss , & ! sea surface salinity (ppt)
679 : sst , & ! sea surface temperature (C)
680 : !hmix , & ! mixed layer depth (m)
681 : Tf , & ! basal freezing temperature (C)
682 : fsnow ! snowfall rate (kg/m^2 s)
683 :
684 : !autodocument_end
685 :
686 : ! local variables
687 :
688 : integer (kind=int_kind) :: &
689 : k , & ! vertical index
690 : n, mm ! thickness category index
691 :
692 : real (kind=dbl_kind) :: &
693 : hin , & ! new ice thickness
694 : hsn , & ! snow thickness (m)
695 : hbr_old , & ! old brine thickness before growh/melt
696 : dhice , & ! change due to sublimation/condensation (m)
697 : kavg , & ! average ice permeability (m^2)
698 : bphi_o , & ! surface ice porosity
699 : hbrin , & ! brine height
700 : dh_direct ! surface flooding or runoff
701 :
702 : real (kind=dbl_kind), dimension (nblyr+2) :: &
703 : ! Defined on Bio Grid points
704 485472 : bSin , & ! salinity on the bio grid (ppt)
705 485472 : brine_sal , & ! brine salinity (ppt)
706 242736 : brine_rho ! brine_density (kg/m^3)
707 :
708 : real (kind=dbl_kind), dimension (nblyr+1) :: &
709 : ! Defined on Bio Grid interfaces
710 485472 : iphin , & ! porosity
711 485472 : ibrine_sal , & ! brine salinity (ppt)
712 485472 : ibrine_rho , & ! brine_density (kg/m^3)
713 485472 : iSin , & ! Salinity on the interface grid (ppt)
714 485472 : iTin ! Temperature on the interface grid (oC)
715 :
716 : real (kind=dbl_kind) :: &
717 : sloss ! brine flux contribution from surface runoff (g/m^2)
718 :
719 : real (kind=dbl_kind), dimension (nbtrcr) :: &
720 485472 : flux_bion_n ! temporary
721 :
722 : ! for bgc sk
723 : real (kind=dbl_kind) :: &
724 : dh_bot_chl , & ! Chlorophyll may or may not flush
725 : dh_top_chl , & ! Chlorophyll may or may not flush
726 : darcy_V_chl
727 :
728 : real (kind=dbl_kind), dimension (nblyr+1) :: &
729 242736 : zspace ! vertical grid spacing
730 :
731 : character(len=*),parameter :: subname='(icepack_biogeochemistry)'
732 :
733 2184624 : zspace(:) = c1/real(nblyr,kind=dbl_kind)
734 242736 : zspace(1) = p5*zspace(2)
735 242736 : zspace(nblyr+1) = zspace(1)
736 :
737 242736 : if (present(bioPorosityIceCell)) bioPorosityIceCell(:) = c0
738 242736 : if (present(bioSalinityIceCell)) bioSalinityIceCell(:) = c0
739 242736 : if (present(bioTemperatureIceCell)) bioTemperatureIceCell(:) = c0
740 :
741 1456416 : do n = 1, ncat
742 :
743 : !-----------------------------------------------------------------
744 : ! initialize
745 : !-----------------------------------------------------------------
746 18491280 : flux_bion_n(:) = c0
747 1213680 : hin_old(n) = c0
748 1213680 : if (aicen_init(n) > puny) then
749 : hin_old(n) = vicen_init(n) &
750 660884 : / aicen_init(n)
751 : else
752 552796 : first_ice(n) = .true.
753 552796 : if (tr_brine) trcrn(nt_fbri,n) = c1
754 552796 : if (z_tracers) then
755 10002791 : do mm = 1,nbtrcr
756 10002791 : trcrn(nt_zbgc_frac-1+mm,n) = zbgc_frac_init(mm)
757 : enddo
758 : endif
759 : endif
760 :
761 1456416 : if (aicen(n) > puny) then
762 :
763 660884 : dh_top_chl = c0
764 660884 : dh_bot_chl = c0
765 660884 : darcy_V_chl= c0
766 6608840 : bSin(:) = c0
767 660884 : hsn = c0
768 660884 : hin = c0
769 660884 : hbrin = c0
770 660884 : kavg = c0
771 660884 : bphi_o = c0
772 660884 : sloss = c0
773 :
774 : !-----------------------------------------------------------------
775 : ! brine dynamics
776 : !-----------------------------------------------------------------
777 :
778 660884 : dhbr_top(n) = c0
779 660884 : dhbr_bot(n) = c0
780 :
781 660884 : if (tr_brine) then
782 660884 : if (trcrn(nt_fbri,n) .le. c0) trcrn(nt_fbri,n) = c1
783 :
784 660884 : dhice = c0
785 : call preflushing_changes (aicen (n), &
786 : vicen (n), vsnon (n), &
787 : meltbn (n), melttn (n), &
788 : congeln (n), snoicen(n), &
789 : hin_old (n), dhice, &
790 : trcrn(nt_fbri,n), &
791 : dhbr_top(n), dhbr_bot(n), &
792 : hbr_old, hin, &
793 660884 : hsn)
794 660884 : if (icepack_warnings_aborted(subname)) return
795 :
796 : ! Requires the average ice permeability = kavg(:)
797 : ! and the surface ice porosity = zphi_o(:)
798 : ! computed in "compute_microS" or from "thermosaline_vertical"
799 :
800 5947956 : iDi(:,n) = c0
801 :
802 : call compute_microS_mushy ( &
803 : trcrn(:,n), hin_old(n), hbr_old, &
804 : sss, sst, bTiz(:,n), &
805 : iTin(:), bphi(:,n), kavg, &
806 : bphi_o, bSin(:), &
807 : brine_sal(:), brine_rho(:), iphin(:), &
808 : ibrine_rho(:), ibrine_sal(:), &
809 660884 : iDi(:,n) , iSin(:))
810 660884 : if (icepack_warnings_aborted(subname)) return
811 :
812 : call update_hbrine (melttn(n), &
813 : meltsn (n), dt, &
814 : hin, hsn, &
815 : hin_old (n), hbrin, &
816 : hbr_old, &
817 : trcrn(nt_fbri,n), &
818 : dhbr_top(n), dhbr_bot(n), &
819 : dh_top_chl, dh_bot_chl, &
820 : kavg, bphi_o, &
821 : darcy_V (n), darcy_V_chl, &
822 : bphi(2,n), aice0, &
823 660884 : dh_direct)
824 660884 : if (icepack_warnings_aborted(subname)) return
825 :
826 660884 : hbri = hbri + hbrin * aicen(n)
827 :
828 : endif ! tr_brine
829 :
830 : !-----------------------------------------------------------------
831 : ! biogeochemistry
832 : !-----------------------------------------------------------------
833 :
834 660884 : if (z_tracers) then
835 :
836 : call zbio (dt, &
837 : melttn(n), &
838 : meltsn(n), meltbn (n), &
839 : congeln(n), snoicen(n), &
840 : fsnow, trcrn(1:ntrcr,n), &
841 : bio_index(1:nbtrcr), bio_index_o(:), &
842 : aicen_init(n), &
843 : vicen_init(n), vsnon_init(n), &
844 : vicen(n), vsnon(n), &
845 : aicen(n), flux_bio_atm(1:nbtrcr), &
846 : n, first_ice(n), &
847 : hin_old(n), ocean_bio(1:nbtrcr), &
848 : ocean_bio_dh, &
849 : bphi(:,n), iphin, &
850 : iDi(:,n), &
851 : fswpenln(:,n), &
852 : dhbr_top(n), dhbr_bot(n), &
853 : zfswin(:,n), &
854 : hbrin, hbr_old, &
855 : darcy_V(n), &
856 : bphi_o, &
857 : iTin, &
858 : Zoo(:,n), &
859 : flux_bio(1:nbtrcr), dh_direct, &
860 : upNO, upNH, &
861 : fbio_snoice, fbio_atmice, &
862 : PP_net, ice_bio_net (1:nbtrcr), &
863 : snow_bio_net(1:nbtrcr),grow_net, &
864 : totalChla, &
865 : flux_bion_n(1:nbtrcr), iSin, &
866 : bioPorosityIceCell, bioSalinityIceCell, &
867 660884 : bioTemperatureIceCell )
868 660884 : if (icepack_warnings_aborted(subname)) return
869 :
870 660884 : if (present(flux_bion)) then
871 0 : do mm = 1, nbtrcr
872 0 : flux_bion(mm,n) = flux_bion_n(mm)
873 : enddo
874 : endif
875 :
876 0 : elseif (skl_bgc) then
877 :
878 : call sklbio (dt, Tf, &
879 : flux_bio (1:nbtrcr), ocean_bio(1:nbtrcr), &
880 : aicen (n), &
881 : meltbn (n), congeln (n), &
882 : fswthrun (n), first_ice(n), &
883 : trcrn (1:ntrcr,n), &
884 : PP_net, upNO, &
885 0 : upNH, grow_net )
886 0 : if (icepack_warnings_aborted(subname)) return
887 :
888 : endif ! skl_bgc
889 :
890 660884 : first_ice(n) = .false.
891 :
892 : else
893 552796 : if (z_tracers) then
894 10002791 : do mm = 1, nbtrcr
895 85602751 : do k = 1, nblyr+1
896 75599960 : if (present(flux_bion)) then
897 : flux_bion(mm,n) = flux_bion(mm,n) + trcrn(bio_index(mm) + k-1,n) * &
898 0 : hin_old(n) * zspace(k)/dt * trcrn(nt_fbri,n)
899 : endif
900 : flux_bio(mm) = flux_bio(mm) + trcrn(bio_index(mm) + k-1,n) * &
901 75599960 : vicen_init(n) * zspace(k)/dt * trcrn(nt_fbri,n)
902 85049955 : trcrn(bio_index(mm) + k-1,n) = c0
903 : enddo
904 : enddo
905 : endif
906 : endif ! aicen > puny
907 : enddo ! ncat
908 :
909 : end subroutine icepack_biogeochemistry
910 :
911 : !=======================================================================
912 : !autodocument_start icepack_load_ocean_bio_array
913 : ! basic initialization for ocean_bio_all
914 :
915 242772 : subroutine icepack_load_ocean_bio_array(&
916 0 : nit, amm, sil, dmsp, dms, algalN, &
917 485544 : doc, don, dic, fed, fep, zaeros, ocean_bio_all, hum)
918 :
919 : real (kind=dbl_kind), intent(in) :: &
920 : nit , & ! ocean nitrate (mmol/m^3)
921 : amm , & ! ammonia/um (mmol/m^3)
922 : sil , & ! silicate (mmol/m^3)
923 : dmsp , & ! dmsp (mmol/m^3)
924 : dms , & ! dms (mmol/m^3)
925 : hum ! humic material (mmol/m^3)
926 :
927 : real (kind=dbl_kind), dimension (:), intent(in) :: &
928 : algalN ! ocean algal nitrogen (mmol/m^3) (diatoms, phaeo, pico)
929 :
930 : real (kind=dbl_kind), dimension (:), intent(in) :: &
931 : doc ! ocean doc (mmol/m^3) (proteins, EPS, lipid)
932 :
933 : real (kind=dbl_kind), dimension (:), intent(in) :: &
934 : don ! ocean don (mmol/m^3)
935 :
936 : real (kind=dbl_kind), dimension (:), intent(in) :: &
937 : dic ! ocean dic (mmol/m^3)
938 :
939 : real (kind=dbl_kind), dimension (:), intent(in) :: &
940 : fed, fep ! ocean disolved and particulate fe (nM)
941 :
942 : real (kind=dbl_kind), dimension (:), intent(in) :: &
943 : zaeros ! ocean aerosols (mmol/m^3)
944 :
945 : real (kind=dbl_kind), dimension (:), intent(out) :: &
946 : ocean_bio_all ! fixed order, all values even for tracers false
947 :
948 : !autodocument_end
949 :
950 : ! local variables
951 :
952 : integer (kind=int_kind) :: &
953 : k, ks ! tracer indices
954 :
955 : character(len=*),parameter :: subname='(icepack_load_ocean_bio_array)'
956 :
957 7283160 : ocean_bio_all(:) = c0
958 :
959 971088 : do k = 1, max_algae
960 728316 : ocean_bio_all(k) = algalN(k) ! N
961 728316 : ks = max_algae + max_doc + max_dic + 1
962 971088 : ocean_bio_all(ks + k) = R_chl2N(k)*algalN(k)!chl
963 : enddo
964 :
965 242772 : ks = max_algae + 1
966 971088 : do k = 1, max_doc
967 971088 : ocean_bio_all(ks + k) = doc(k) ! doc
968 : enddo
969 242772 : ks = ks + max_doc
970 485544 : do k = 1, max_dic
971 485544 : ocean_bio_all(ks + k) = dic(k) ! dic
972 : enddo
973 :
974 242772 : ks = 2*max_algae + max_doc + max_dic + 7
975 485544 : do k = 1, max_don
976 485544 : ocean_bio_all(ks + k) = don(k) ! don
977 : enddo
978 :
979 242772 : ks = max_algae + 1
980 242772 : ocean_bio_all(ks) = nit ! nit
981 :
982 242772 : ks = 2*max_algae + max_doc + 2 + max_dic
983 242772 : ocean_bio_all(ks) = amm ! Am
984 242772 : ks = ks + 1
985 242772 : ocean_bio_all(ks) = sil ! Sil
986 242772 : ks = ks + 1
987 : ocean_bio_all(ks) = R_S2N(1)*algalN(1) & ! DMSPp
988 : + R_S2N(2)*algalN(2) &
989 242772 : + R_S2N(3)*algalN(3)
990 242772 : ks = ks + 1
991 242772 : ocean_bio_all(ks) = dmsp ! DMSPd
992 242772 : ks = ks + 1
993 242772 : ocean_bio_all(ks) = dms ! DMS
994 242772 : ks = ks + 1
995 242772 : ocean_bio_all(ks) = nit ! PON
996 242772 : ks = 2*max_algae + max_doc + 7 + max_dic + max_don
997 728316 : do k = 1, max_fe
998 728316 : ocean_bio_all(ks + k) = fed(k) ! fed
999 : enddo
1000 242772 : ks = ks + max_fe
1001 728316 : do k = 1, max_fe
1002 728316 : ocean_bio_all(ks + k) = fep(k) ! fep
1003 : enddo
1004 242772 : ks = ks + max_fe
1005 1699404 : do k = 1, max_aero
1006 1699404 : ocean_bio_all(ks+k) = zaeros(k) ! zaero
1007 : enddo
1008 242772 : ks = ks + max_aero + 1
1009 242772 : ocean_bio_all(ks) = hum ! humics
1010 :
1011 242772 : end subroutine icepack_load_ocean_bio_array
1012 :
1013 : !=======================================================================
1014 : !autodocument_start icepack_init_ocean_bio
1015 : ! Initialize ocean concentration
1016 :
1017 36 : subroutine icepack_init_ocean_bio (amm, dmsp, dms, algalN, doc, dic, don, &
1018 72 : fed, fep, hum, nit, sil, zaeros,CToN, CToN_DON)
1019 :
1020 : real (kind=dbl_kind), intent(out), optional:: &
1021 : amm , & ! ammonium
1022 : dmsp , & ! DMSPp
1023 : dms , & ! DMS
1024 : hum , & ! humic material
1025 : nit , & ! nitrate
1026 : sil ! silicate
1027 :
1028 : real (kind=dbl_kind), dimension(:), intent(out), optional:: &
1029 : algalN , & ! algae
1030 : doc , & ! DOC
1031 : dic , & ! DIC
1032 : don , & ! DON
1033 : fed , & ! Dissolved Iron
1034 : fep , & ! Particulate Iron
1035 : zaeros ! BC and dust
1036 :
1037 : real (kind=dbl_kind), dimension(:), intent(out), optional :: &
1038 : CToN , & ! carbon to nitrogen ratio for algae
1039 : CToN_DON ! nitrogen to carbon ratio for proteins
1040 :
1041 : !autodocument_end
1042 :
1043 : ! local variables
1044 :
1045 : integer (kind=int_kind) :: &
1046 : k
1047 :
1048 : character(len=*),parameter :: subname='(icepack_init_ocean_bio)'
1049 :
1050 36 : if (present(CToN)) then
1051 0 : CToN(:) = c0
1052 0 : CToN(1) = R_C2N(1)
1053 0 : CToN(2) = R_C2N(2)
1054 0 : CToN(3) = R_C2N(3)
1055 : endif
1056 :
1057 36 : if (present(CToN_DON)) then
1058 0 : CToN_DON(:) = c0
1059 0 : CToN_DON(1) = R_C2N_DON(1)
1060 : endif
1061 :
1062 36 : if (present(amm)) &
1063 36 : amm = c1 ! ISPOL < 1 mmol/m^3
1064 36 : if (present(dmsp)) &
1065 36 : dmsp = p1
1066 36 : if (present(dms)) &
1067 36 : dms = p1
1068 36 : if (present(algalN)) then
1069 144 : algalN(:) = c0
1070 36 : algalN(1) = c1 !0.0026_dbl_kind ! ISPOL, Lannuzel 2013(pennate)
1071 36 : algalN(2) = 0.0057_dbl_kind ! ISPOL, Lannuzel 2013(small plankton)
1072 36 : algalN(3) = 0.0027_dbl_kind ! ISPOL, Lannuzel 2013(Phaeocystis)
1073 : endif
1074 36 : if (present(doc))then
1075 144 : doc(:) = c0
1076 36 : doc(1) = 16.2_dbl_kind ! 18% saccharides
1077 36 : doc(2) = 9.0_dbl_kind ! lipids
1078 36 : doc(3) = c1 !
1079 : endif
1080 36 : if (present(dic)) then
1081 72 : dic(:) = c0
1082 36 : dic(1) = 1950.0_dbl_kind ! 1950-2260 mmol C/m3 (Tynan et al. 2015)
1083 : endif
1084 36 : if (present(don)) then
1085 72 : don(:) = c0
1086 36 : don(1) = 12.9_dbl_kind ! 64.3_dbl_kind ! 72% Total DOC~90 mmolC/m^3 ISPOL with N:C of 0.2
1087 : endif
1088 36 : if (present(fed)) then
1089 108 : fed(:) = c0
1090 36 : fed(1) = 0.4_dbl_kind ! c1 (nM) Lannuzel2007 DFe,! range 0.14-2.6 (nM) van der Merwe 2011
1091 : ! Tagliabue 2012 (0.4 nM)
1092 : endif
1093 36 : if (present(fep)) then
1094 108 : fep(:) = c0
1095 36 : fep(1) = c2 ! (nM) van der Merwe 2011
1096 : ! (0.6 to 2.9 nM ocean)
1097 : endif
1098 36 : if (present(hum)) &
1099 36 : hum = c1 ! mmol C/m^3
1100 36 : if (present(nit)) &
1101 36 : nit = 12.0_dbl_kind
1102 36 : if (present(sil)) &
1103 36 : sil = 25.0_dbl_kind
1104 36 : if (present(zaeros)) &
1105 252 : zaeros(:) = c0
1106 :
1107 36 : end subroutine icepack_init_ocean_bio
1108 : !
1109 : !=======================================================================
1110 : !
1111 : ! Given some added new ice to the base of the existing ice, recalculate
1112 : ! vertical bio tracer so that new grid cells are all the same size.
1113 : !
1114 : ! author: N. Jeffery, LANL
1115 : ! date : Mar 3, 2024
1116 : !
1117 : ! Based on update_vertical_tracers modified for vertical biogeochemistry
1118 : !
1119 36915 : subroutine update_vertical_bio_tracers(nbiolyr, trc, h1, h2, trc0, zspace)
1120 :
1121 : integer (kind=int_kind), intent(in) :: &
1122 : nbiolyr ! number of bio layers nblyr+1
1123 :
1124 : real (kind=dbl_kind), dimension(:), intent(inout) :: &
1125 : trc ! vertical tracer
1126 :
1127 : real (kind=dbl_kind), intent(in) :: &
1128 : h1, & ! old thickness
1129 : h2, & ! new thickness
1130 : trc0 ! tracer value of added ice on ice bottom
1131 :
1132 : real (kind=dbl_kind), dimension(nbiolyr), intent(in) :: &
1133 : zspace
1134 :
1135 : ! local variables
1136 :
1137 73830 : real(kind=dbl_kind), dimension(nbiolyr) :: trc2 ! updated tracer temporary
1138 :
1139 : ! vertical indices for old and new grid
1140 : integer :: k1, k2
1141 :
1142 : real (kind=dbl_kind) :: &
1143 : z1a, z1b, & ! upper, lower boundary of old cell/added new ice at bottom
1144 : z2a, z2b, & ! upper, lower boundary of new cell
1145 : overlap , & ! overlap between old and new cell
1146 : rnilyr
1147 :
1148 : !rnilyr = real(nilyr,dbl_kind)
1149 36915 : z2a = c0
1150 36915 : z2b = c0
1151 36915 : if (h2 > puny) then
1152 : ! loop over new grid cells
1153 309051 : do k2 = 1, nbiolyr
1154 :
1155 : ! initialize new tracer
1156 274712 : trc2(k2) = c0
1157 :
1158 : ! calculate upper and lower boundary of new cell
1159 274712 : z2a = z2b !((k2 - 1) * h2) * zspace(k2)+z2b ! / rnilyr
1160 274712 : z2b = z2b + h2 * zspace(k2) !(k2 * h2) * zspace(k2)+z2a !/ rnilyr
1161 :
1162 274712 : z1a = c0
1163 274712 : z1b = c0
1164 : ! loop over old grid cells
1165 2472408 : do k1 = 1, nbiolyr
1166 :
1167 : ! calculate upper and lower boundary of old cell
1168 2197696 : z1a = z1b !((k1 - 1) * h1) * zspace(k1)+z1b !/ rnilyr
1169 2197696 : z1b = z1b + h1 * zspace(k1) !(k1 * h1) * zspace(k1)+z1a !/ rnilyr
1170 :
1171 : ! calculate overlap between old and new cell
1172 2197696 : overlap = max(min(z1b, z2b) - max(z1a, z2a), c0)
1173 :
1174 : ! aggregate old grid cell contribution to new cell
1175 2472408 : trc2(k2) = trc2(k2) + overlap * trc(k1)
1176 :
1177 : enddo ! k1
1178 :
1179 : ! calculate upper and lower boundary of added new ice at bottom
1180 274712 : z1a = h1
1181 274712 : z1b = h2
1182 :
1183 : ! calculate overlap between added ice and new cell
1184 274712 : overlap = max(min(z1b, z2b) - max(z1a, z2a), c0)
1185 : ! aggregate added ice contribution to new cell
1186 274712 : trc2(k2) = trc2(k2) + overlap * trc0
1187 : ! renormalize new grid cell
1188 309051 : trc2(k2) = trc2(k2)/zspace(k2)/h2 !(rnilyr * trc2(k2)) / h2
1189 :
1190 : enddo ! k2
1191 : else
1192 23184 : trc2 = trc
1193 : endif
1194 : ! update vertical tracer array with the adjusted tracer
1195 332235 : trc = trc2
1196 :
1197 36915 : end subroutine update_vertical_bio_tracers
1198 :
1199 : !=======================================================================
1200 :
1201 : end module icepack_zbgc
1202 :
1203 : !=======================================================================
|