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