Line data Source code
1 : !=======================================================================
2 : !
3 : ! Biogeochemistry variables
4 : !
5 : ! authors: Nicole Jeffery, LANL
6 : ! Scott Elliot, LANL
7 : ! Elizabeth C. Hunke, LANL
8 : !
9 : module icepack_zbgc_shared
10 :
11 : use icepack_kinds
12 : use icepack_parameters, only: p5, c0, c1, secday, puny
13 : use icepack_parameters, only: hs_ssl, sk_l
14 : use icepack_parameters, only: rhoi, cp_ocn, cp_ice, Lfresh
15 : use icepack_parameters, only: solve_zbgc
16 : use icepack_parameters, only: fr_resp
17 : use icepack_tracers, only: max_nbtrcr, max_algae, max_doc
18 : use icepack_tracers, only: max_don
19 : use icepack_tracers, only: nt_bgc_N, nt_fbri
20 : use icepack_warnings, only: warnstr, icepack_warnings_add
21 : use icepack_warnings, only: icepack_warnings_setabort, icepack_warnings_aborted
22 :
23 : implicit none
24 :
25 : private
26 : public :: calculate_qin_from_Sin, &
27 : remap_zbgc, & ! LCOV_EXCL_LINE
28 : zap_small_bgc, & ! LCOV_EXCL_LINE
29 : regrid_stationary, & ! LCOV_EXCL_LINE
30 : merge_bgc_fluxes, & ! LCOV_EXCL_LINE
31 : merge_bgc_fluxes_skl
32 :
33 : !-----------------------------------------------------------------
34 : ! Transport type
35 : !-----------------------------------------------------------------
36 : ! In delta Eddington, algal particles are assumed to cause no
37 : ! significant scattering (Brieglib and Light), only absorption
38 : ! in the visible spectral band (200-700 nm)
39 : ! Algal types: Diatoms, flagellates, Phaeocycstis
40 : ! DOC : Proteins, EPS, Lipids
41 : !-----------------------------------------------------------------
42 : !------------------------------------------------------------
43 : ! Aerosol order and type should be consistent with order/type
44 : ! specified in delta Eddington: 1) hydrophobic black carbon;
45 : ! 2) hydrophilic black carbon; 3) dust (0.05-0.5 micron);
46 : ! 4) dust (0.5-1.25 micron); 5) dust (1.25-2.5 micron);
47 : ! 6) dust (2.5-5 micron)
48 : !-------------------------------------------------------------
49 :
50 : ! bio parameters for algal_dyn
51 :
52 : real (kind=dbl_kind), dimension(max_algae), public :: &
53 : R_C2N , & ! algal C to N (mole/mole) ! LCOV_EXCL_LINE
54 : R_chl2N , & ! 3 algal chlorophyll to N (mg/mmol) ! LCOV_EXCL_LINE
55 : F_abs_chl ! to scale absorption in Dedd
56 :
57 : real (kind=dbl_kind), dimension(max_don), public :: & ! increase compare to algal R_Fe2C
58 : R_C2N_DON
59 :
60 : real (kind=dbl_kind), dimension(max_algae), public :: &
61 : R_Si2N , & ! algal Sil to N (mole/mole) ! LCOV_EXCL_LINE
62 : R_S2N , & ! algal S to N (mole/mole) ! LCOV_EXCL_LINE
63 : ! Marchetti et al 2006, 3 umol Fe/mol C for iron limited Pseudo-nitzschia
64 : R_Fe2C , & ! algal Fe to carbon (umol/mmol)
65 : R_Fe2N ! algal Fe to N (umol/mmol)
66 :
67 : real (kind=dbl_kind), dimension(max_don), public :: &
68 : R_Fe2DON ! Fe to N of DON (nmol/umol)
69 :
70 : real (kind=dbl_kind), dimension(max_doc), public :: &
71 : R_Fe2DOC ! Fe to C of DOC (nmol/umol)
72 :
73 : real (kind=dbl_kind), parameter, public :: &
74 : R_gC2molC = 12.01_dbl_kind ! mg/mmol C
75 :
76 : ! scavenging coefficient for tracers in snow
77 : ! bottom to last 6 are from Flanner et al., 2007
78 : ! very last one is for humic material
79 : real (kind=dbl_kind), parameter, dimension(max_nbtrcr), public :: &
80 : kscavz = (/ 0.03_dbl_kind, 0.03_dbl_kind, 0.03_dbl_kind, & ! LCOV_EXCL_LINE
81 : 0.03_dbl_kind, 0.03_dbl_kind, 0.03_dbl_kind, & ! LCOV_EXCL_LINE
82 : 0.03_dbl_kind, 0.03_dbl_kind, 0.03_dbl_kind, & ! LCOV_EXCL_LINE
83 : 0.03_dbl_kind, 0.03_dbl_kind, 0.03_dbl_kind, & ! LCOV_EXCL_LINE
84 : 0.03_dbl_kind, 0.03_dbl_kind, 0.03_dbl_kind, & ! LCOV_EXCL_LINE
85 : 0.03_dbl_kind, 0.03_dbl_kind, 0.03_dbl_kind, & ! LCOV_EXCL_LINE
86 : 0.03_dbl_kind, 0.03_dbl_kind, 0.03_dbl_kind, & ! LCOV_EXCL_LINE
87 : 0.03_dbl_kind, & ! LCOV_EXCL_LINE
88 : 0.03_dbl_kind, 0.20_dbl_kind, 0.02_dbl_kind, & ! LCOV_EXCL_LINE
89 : 0.02_dbl_kind, 0.01_dbl_kind, 0.01_dbl_kind, & ! LCOV_EXCL_LINE
90 : 0.03_dbl_kind /)
91 :
92 : !-----------------------------------------------------------------
93 : ! skeletal layer biogeochemistry
94 : !-----------------------------------------------------------------
95 :
96 : real (kind=dbl_kind), parameter, public :: &
97 : phi_sk = 0.30_dbl_kind ! skeletal layer porosity
98 :
99 : !-----------------------------------------------------------------
100 : ! general biogeochemistry
101 : !-----------------------------------------------------------------
102 :
103 : real (kind=dbl_kind), dimension(max_nbtrcr), public :: &
104 : zbgc_frac_init,&! initializes mobile fraction ! LCOV_EXCL_LINE
105 : bgc_tracer_type ! described tracer in mobile or stationary phases
106 : ! < 0 is purely mobile (eg. nitrate)
107 : ! > 0 has timescales for transitions between
108 : ! phases based on whether the ice is melting or growing
109 :
110 : real (kind=dbl_kind), dimension(max_nbtrcr), public :: &
111 : zbgc_init_frac, & ! fraction of ocean tracer concentration in new ice ! LCOV_EXCL_LINE
112 : tau_ret, & ! retention timescale (s), mobile to stationary phase ! LCOV_EXCL_LINE
113 : tau_rel ! release timescale (s), stationary to mobile phase
114 :
115 : !-----------------------------------------------------------------
116 : ! From algal_dyn in icepack_algae.F90 but not in namelist
117 : !-----------------------------------------------------------------
118 :
119 : real (kind=dbl_kind), dimension(max_algae), public :: &
120 : chlabs , & ! chla absorption 1/m/(mg/m^3) ! LCOV_EXCL_LINE
121 : alpha2max_low , & ! light limitation (1/(W/m^2)) ! LCOV_EXCL_LINE
122 : beta2max , & ! light inhibition (1/(W/m^2)) ! LCOV_EXCL_LINE
123 : mu_max , & ! maximum growth rate (1/d) ! LCOV_EXCL_LINE
124 : grow_Tdep , & ! T dependence of growth (1/C) ! LCOV_EXCL_LINE
125 : fr_graze , & ! fraction of algae grazed ! LCOV_EXCL_LINE
126 : mort_pre , & ! mortality (1/day) ! LCOV_EXCL_LINE
127 : mort_Tdep , & ! T dependence of mortality (1/C) ! LCOV_EXCL_LINE
128 : k_exude , & ! algal carbon exudation rate (1/d) ! LCOV_EXCL_LINE
129 : K_Nit , & ! nitrate half saturation (mmol/m^3) ! LCOV_EXCL_LINE
130 : K_Am , & ! ammonium half saturation (mmol/m^3) ! LCOV_EXCL_LINE
131 : K_Sil , & ! silicon half saturation (mmol/m^3) ! LCOV_EXCL_LINE
132 : K_Fe ! iron half saturation or micromol/m^3
133 :
134 : real (kind=dbl_kind), dimension(max_DON), public :: &
135 : f_don , & ! fraction of spilled grazing to DON ! LCOV_EXCL_LINE
136 : kn_bac , & ! Bacterial degredation of DON (1/d) ! LCOV_EXCL_LINE
137 : f_don_Am ! fraction of remineralized DON to Am
138 :
139 : real (kind=dbl_kind), dimension(max_DOC), public :: &
140 : f_doc , & ! fraction of mort_N that goes to each doc pool ! LCOV_EXCL_LINE
141 : f_exude , & ! fraction of exuded carbon to each DOC pool ! LCOV_EXCL_LINE
142 : k_bac ! Bacterial degredation of DOC (1/d)
143 :
144 : !-----------------------------------------------------------------
145 : ! brine
146 : !-----------------------------------------------------------------
147 :
148 : integer (kind=int_kind), parameter, public :: &
149 : exp_h = 3 ! power law for hierarchical model
150 :
151 : real (kind=dbl_kind), parameter, public :: &
152 : k_o = 3.e-8_dbl_kind, & ! permeability scaling factor (m^2) ! LCOV_EXCL_LINE
153 : thinS = 0.05_dbl_kind ! minimum ice thickness for brine
154 :
155 : real (kind=dbl_kind), public :: &
156 : flood_frac ! fraction of ocean/meltwater that floods !*****
157 :
158 : real (kind=dbl_kind), parameter, public :: &
159 : bphimin = 0.03_dbl_kind ! minimum porosity for zbgc only
160 :
161 : real (kind=dbl_kind), parameter, public :: &
162 : viscos_dynamic = 2.2_dbl_kind , & ! 1.8e-3_dbl_kind (pure water at 0^oC) (kg/m/s) ! LCOV_EXCL_LINE
163 : Dm = 1.0e-9_dbl_kind, & ! molecular diffusion (m^2/s) ! LCOV_EXCL_LINE
164 : Ra_c = 0.05_dbl_kind ! critical Rayleigh number for bottom convection
165 :
166 : !=======================================================================
167 :
168 : contains
169 :
170 : !=======================================================================
171 : !
172 : ! Compute the internal ice enthalpy using new salinity and Tin
173 : !
174 :
175 0 : function calculate_qin_from_Sin (Tin, Tmltk) &
176 0 : result(qin)
177 :
178 : real (kind=dbl_kind), intent(in) :: &
179 : Tin ,& ! internal temperature ! LCOV_EXCL_LINE
180 : Tmltk ! melting temperature at one level
181 :
182 : ! local variables
183 :
184 : real (kind=dbl_kind) :: &
185 : qin ! melting temperature at one level
186 :
187 : character(len=*),parameter :: subname='(calculate_qin_from_Sin)'
188 :
189 0 : qin =-rhoi*(cp_ice*(Tmltk-Tin) + Lfresh*(c1-Tmltk/Tin) - cp_ocn*Tmltk)
190 :
191 0 : end function calculate_qin_from_Sin
192 :
193 : !=======================================================================
194 : !
195 : ! Remaps tracer fields in a given category from one set of layers to another.
196 : ! Grids can be very different and so can vertical spaces.
197 :
198 0 : subroutine remap_zbgc(nlyrn, &
199 : it, & ! LCOV_EXCL_LINE
200 : trcrn, trtmp, & ! LCOV_EXCL_LINE
201 : nr0, nbyrn, & ! LCOV_EXCL_LINE
202 : hice, hinS, & ! LCOV_EXCL_LINE
203 : ice_grid, bio_grid, & ! LCOV_EXCL_LINE
204 : S_min )
205 :
206 : integer (kind=int_kind), intent(in) :: &
207 : it , & ! tracer index in top layer ! LCOV_EXCL_LINE
208 : nr0 , & ! receiver category ! LCOV_EXCL_LINE
209 : nlyrn , & ! number of ice layers ! LCOV_EXCL_LINE
210 : nbyrn ! number of biology layers
211 :
212 : real (kind=dbl_kind), dimension (:), intent(in) :: &
213 : trcrn ! ice tracers
214 :
215 : real (kind=dbl_kind), dimension (:), intent(inout) :: &
216 : trtmp ! temporary, remapped ice tracers
217 :
218 : real (kind=dbl_kind), dimension (:), intent(in) :: &
219 : ice_grid ! CICE grid cgrid(2:nilyr+1)
220 :
221 : real (kind=dbl_kind), dimension (:), intent(in) :: &
222 : bio_grid ! CICE grid grid(2:nbyrn+1)
223 :
224 : real(kind=dbl_kind), intent(in) :: &
225 : hice , & ! CICE ice thickness ! LCOV_EXCL_LINE
226 : hinS , & ! brine height ! LCOV_EXCL_LINE
227 : S_min ! for salinity on CICE grid
228 :
229 : ! local variables
230 :
231 : integer (kind=int_kind) :: &
232 : kd, kr, kdr , & ! more indices ! LCOV_EXCL_LINE
233 : kdi , & ! more indices ! LCOV_EXCL_LINE
234 : n_nd , & ! number of layers in donor ! LCOV_EXCL_LINE
235 : n_nr, n_plus ! number of layers in receiver
236 :
237 : real (kind=dbl_kind), dimension (nbyrn+3+nlyrn) :: &
238 : trdr , & ! combined tracer ! LCOV_EXCL_LINE
239 0 : trgrid ! combined grid
240 :
241 : real (kind=dbl_kind), dimension (nbyrn+nlyrn+3) :: &
242 : tracer , & ! temporary, ice tracers values ! LCOV_EXCL_LINE
243 : dgrid , & ! temporary, donor grid dimensional ! LCOV_EXCL_LINE
244 0 : rgrid ! temporary, receiver grid dimensional
245 :
246 : character(len=*),parameter :: subname='(remap_zbgc)'
247 :
248 0 : if ((hinS < c0) .OR. (hice < c0)) then
249 0 : call icepack_warnings_setabort(.true.,__FILE__,__LINE__)
250 0 : call icepack_warnings_add(subname//' ice: remap_layers_bgc error')
251 0 : return
252 : endif
253 :
254 0 : if (nr0 == 0) then ! cice to bio
255 :
256 0 : n_nd = nlyrn
257 0 : n_nr = nbyrn
258 0 : n_plus = 2
259 0 : dgrid (1) = min(-hice+hinS, -hinS+hice, c0)
260 0 : dgrid (nlyrn+2) = min(hinS, hice)
261 0 : tracer(1) = trcrn(it)
262 0 : tracer(nlyrn+2) = trcrn(it+nlyrn-1)
263 0 : rgrid (nbyrn+2) = min(hinS, hice)
264 0 : if (hice > hinS) then
265 0 : rgrid(1) = c0
266 0 : do kr = 1,n_nr
267 0 : rgrid(kr+1) = bio_grid(kr)*hinS
268 : enddo
269 0 : do kd = 1,n_nd
270 0 : dgrid(kd+1) = (ice_grid(kd)-c1)*hice+hinS
271 0 : tracer(kd+1) = trcrn(it+kd-1)
272 : enddo
273 : else
274 0 : rgrid(1) = -hinS + hice
275 0 : do kr = 1,n_nr
276 0 : rgrid(kr+1) = (bio_grid(kr)-c1)*hinS + hice
277 : enddo
278 0 : do kd = 1,n_nd
279 0 : dgrid(kd+1) = ice_grid(kd)*hice
280 0 : tracer(kd+1) = trcrn(it+kd-1)
281 : enddo
282 : endif
283 :
284 : else ! bio to cice
285 :
286 0 : n_nd = nbyrn
287 0 : n_nr = nlyrn
288 0 : if (hice > hinS) then ! add S_min to top layer
289 0 : n_plus = 3
290 0 : tracer(1) = S_min
291 0 : tracer(2) = S_min
292 0 : rgrid (1) = -hice + hinS
293 0 : rgrid (nlyrn+n_plus-1) = hinS
294 0 : do kr = 1,n_nr
295 0 : rgrid(kr+1) = (ice_grid(kr)-c1)*hice+ hinS
296 : enddo
297 0 : dgrid (1) = -hice+hinS
298 0 : dgrid (2) = (hinS-hice)*p5
299 0 : dgrid (nbyrn+n_plus) = hinS
300 0 : tracer(nbyrn+n_plus) = trcrn(it+nbyrn-1)
301 0 : do kd = 1,n_nd
302 0 : dgrid(kd+2) = bio_grid(kd)*hinS
303 0 : tracer(kd+2) = trcrn(it+kd-1)
304 : enddo
305 0 : tracer(n_plus) = (S_min*(hice-hinS) + &
306 : tracer(n_plus)*p5*(dgrid(n_plus+1)-dgrid(n_plus)))/ & ! LCOV_EXCL_LINE
307 0 : (hice-hinS+ p5*(dgrid(n_plus+1)-dgrid(n_plus)))
308 0 : tracer(1) = tracer(n_plus)
309 0 : tracer(2) = tracer(n_plus)
310 : else
311 0 : n_plus = 2
312 0 : tracer(1) = trcrn(it)
313 0 : tracer(nbyrn+2) = trcrn(it+nbyrn-1)
314 0 : dgrid (1) = hice-hinS
315 0 : dgrid (nbyrn+2) = hice
316 0 : rgrid (nlyrn+2) = hice
317 0 : rgrid (1) = c0
318 0 : do kd = 1,n_nd
319 0 : dgrid(kd+1) = (bio_grid(kd)-c1)*hinS + hice
320 0 : tracer(kd+1) = trcrn(it+kd-1)
321 : enddo
322 0 : do kr = 1,n_nr
323 0 : rgrid(kr+1) = ice_grid(kr)*hice
324 : enddo
325 : endif
326 :
327 : endif
328 :
329 0 : kdr = 0 ! combined indices
330 0 : kdi = 1
331 :
332 0 : do kr = 1, n_nr
333 0 : do kd = kdi, n_nd+n_plus
334 0 : if (dgrid(kd) < rgrid(kr+1)) then
335 0 : kdr = kdr+1
336 0 : trgrid(kdr) = dgrid(kd)
337 0 : trdr (kdr) = tracer(kd)
338 0 : elseif (dgrid(kd) > rgrid(kr+1)) then
339 0 : kdr = kdr + 1
340 0 : kdi = kd
341 0 : trgrid(kdr) = rgrid(kr+1)
342 0 : trtmp (it+kr-1) = trdr(kdr-1) &
343 : + (rgrid(kr+1) - trgrid(kdr-1)) & ! LCOV_EXCL_LINE
344 : * (tracer(kd) - trdr(kdr-1)) & ! LCOV_EXCL_LINE
345 0 : / (dgrid(kd) - trgrid(kdr-1))
346 0 : trdr(kdr) = trtmp(it+kr-1)
347 0 : EXIT
348 : else
349 0 : kdr = kdr+1
350 0 : kdi = kd+1
351 0 : trgrid(kdr) = rgrid(kr+1)
352 0 : trtmp (it+kr-1) = tracer(kd)
353 0 : trdr (kdr) = tracer(kd)
354 0 : EXIT
355 : endif
356 : enddo
357 : enddo
358 :
359 : end subroutine remap_zbgc
360 :
361 : !=======================================================================
362 :
363 : ! remove tracer for very small fractional areas
364 :
365 0 : subroutine zap_small_bgc (zlevels, dflux_bio, &
366 0 : dt, zvol, btrcr)
367 :
368 : integer (kind=int_kind), intent(in) :: &
369 : zlevels ! number of vertical levels in ice
370 :
371 : real (kind=dbl_kind), intent(in) :: &
372 : dt ! time step (s)
373 :
374 : real (kind=dbl_kind), intent(inout) :: &
375 : dflux_bio ! zapped bio tracer flux from biology (mmol/m^2/s)
376 :
377 : real (kind=dbl_kind), dimension (zlevels), intent(in) :: &
378 : btrcr , & ! zapped bio tracer flux from biology (mmol/m^2/s) ! LCOV_EXCL_LINE
379 : zvol ! ice volume (m)
380 :
381 : ! local variables
382 :
383 : integer (kind=int_kind) :: &
384 : k ! layer index
385 :
386 : character(len=*),parameter :: subname='(zap_small_bgc)'
387 :
388 0 : do k = 1, zlevels
389 0 : dflux_bio = dflux_bio + btrcr(k)*zvol(k)/dt
390 : enddo
391 :
392 0 : end subroutine zap_small_bgc
393 :
394 : !=======================================================================
395 : !
396 : ! authors Nicole Jeffery, LANL
397 :
398 0 : subroutine regrid_stationary (C_stationary, hbri_old, &
399 : hbri, dt, & ! LCOV_EXCL_LINE
400 : ntrcr, nblyr, & ! LCOV_EXCL_LINE
401 : top_conc, igrid, & ! LCOV_EXCL_LINE
402 : flux_bio, & ! LCOV_EXCL_LINE
403 : melt_b, con_gel)
404 :
405 : integer (kind=int_kind), intent(in) :: &
406 : ntrcr, & ! number of tracers ! LCOV_EXCL_LINE
407 : nblyr ! number of bio layers
408 :
409 : real (kind=dbl_kind), intent(inout) :: &
410 : flux_bio ! ocean tracer flux (mmol/m^2/s) positive into ocean
411 :
412 : real (kind=dbl_kind), dimension (nblyr+1), intent(inout) :: &
413 : C_stationary ! stationary bulk concentration*h (mmol/m^2)
414 :
415 : real (kind=dbl_kind), dimension (nblyr+1), intent(in) :: &
416 : igrid ! CICE bio grid
417 :
418 : real(kind=dbl_kind), intent(in) :: &
419 : dt , & ! time step ! LCOV_EXCL_LINE
420 : top_conc , & ! c0 or frazil concentration ! LCOV_EXCL_LINE
421 : hbri_old , & ! previous timestep brine height ! LCOV_EXCL_LINE
422 : hbri ! brine height
423 :
424 : real(kind=dbl_kind), intent(in), optional :: &
425 : melt_b, & ! bottom melt (m) ! LCOV_EXCL_LINE
426 : con_gel ! bottom growth (m)
427 :
428 : ! local variables
429 :
430 : integer (kind=int_kind) :: k, nt, nr
431 :
432 : real (kind=dbl_kind), dimension (ntrcr+2) :: &
433 : trtmp0, & ! temporary, remapped tracers ! LCOV_EXCL_LINE
434 0 : trtmp
435 :
436 : real (kind=dbl_kind):: &
437 : meltb, & ! ice bottom melt (m) ! LCOV_EXCL_LINE
438 : congel, & ! ice bottom growth (m) ! LCOV_EXCL_LINE
439 : htemp, & ! ice thickness after melt (m) ! LCOV_EXCL_LINE
440 : dflux, & ! regrid flux correction (mmol/m^2) ! LCOV_EXCL_LINE
441 : sum_i, & ! total tracer before melt loss ! LCOV_EXCL_LINE
442 : sum_f, & ! total tracer after melt ! LCOV_EXCL_LINE
443 : hice, & ! LCOV_EXCL_LINE
444 0 : hbio
445 :
446 : real (kind=dbl_kind), dimension(nblyr+1):: &
447 0 : zspace
448 :
449 : character(len=*),parameter :: subname='(regrid_stationary)'
450 :
451 : ! initialize
452 :
453 0 : zspace(:) = c1/(real(nblyr,kind=dbl_kind))
454 0 : zspace(1) = p5*zspace(1)
455 0 : zspace(nblyr+1) = zspace(1)
456 0 : trtmp0(:) = c0
457 0 : trtmp(:) = c0
458 0 : meltb = c0
459 0 : nt = 1
460 0 : nr = 0
461 0 : sum_i = c0
462 0 : sum_f = c0
463 0 : meltb = c0
464 0 : congel = c0
465 0 : dflux = c0
466 :
467 : !---------------------
468 : ! compute initial sum
469 : !----------------------
470 :
471 0 : do k = 1, nblyr+1
472 0 : sum_i = sum_i + C_stationary(k)*zspace(k)
473 :
474 : enddo
475 :
476 0 : if (present(melt_b)) then
477 0 : meltb = melt_b
478 : endif
479 0 : if (present(con_gel)) then
480 0 : congel = con_gel
481 : endif
482 :
483 0 : if (hbri_old > c0) then
484 0 : do k = 1, nblyr+1
485 0 : trtmp0(nblyr+2-k) = C_stationary(k)/hbri_old ! reverse order
486 : enddo ! k
487 : endif
488 :
489 0 : htemp = c0
490 :
491 0 : if (meltb > c0) then
492 0 : htemp = hbri_old-meltb
493 0 : nr = 0
494 0 : hice = hbri_old
495 0 : hbio = htemp
496 0 : elseif (congel > c0) then
497 0 : htemp = hbri_old+congel
498 0 : nr = 1
499 0 : hice = htemp
500 0 : hbio = hbri_old
501 0 : elseif (hbri .gt. hbri_old) then
502 0 : htemp = hbri
503 0 : nr = 1
504 0 : hice = htemp
505 0 : hbio = hbri_old
506 : endif
507 :
508 : !-----------------------------------------------------------------
509 : ! Regrid C_stationary to add or remove bottom layer(s)
510 : !-----------------------------------------------------------------
511 0 : if (htemp > c0) then
512 : call remap_zbgc (nblyr+1, &
513 : nt, & ! LCOV_EXCL_LINE
514 : trtmp0(1:ntrcr), & ! LCOV_EXCL_LINE
515 : trtmp, & ! LCOV_EXCL_LINE
516 : nr, nblyr+1, & ! LCOV_EXCL_LINE
517 : hice, hbio, & ! LCOV_EXCL_LINE
518 : igrid(1:nblyr+1), & ! LCOV_EXCL_LINE
519 0 : igrid(1:nblyr+1), top_conc )
520 0 : if (icepack_warnings_aborted(subname)) return
521 :
522 0 : trtmp0(:) = c0
523 0 : do k = 1,nblyr+1
524 0 : trtmp0(nblyr+2-k) = trtmp(nt + k-1)
525 : enddo !k
526 :
527 0 : do k = 1, nblyr+1
528 0 : C_stationary(k) = trtmp0(k)*htemp
529 0 : sum_f = sum_f + C_stationary(k)*zspace(k)
530 : enddo ! k
531 :
532 0 : if (congel > c0 .and. top_conc .le. c0 .and. abs(sum_i-sum_f) > puny) then
533 0 : dflux = sum_i - sum_f
534 0 : sum_f = c0
535 0 : do k = 1,nblyr+1
536 0 : C_stationary(k) = max(c0,C_stationary(k) + dflux)
537 0 : sum_f = sum_f + C_stationary(k)*zspace(k)
538 : enddo
539 : endif
540 :
541 0 : flux_bio = flux_bio + (sum_i -sum_f)/dt
542 : endif
543 :
544 : end subroutine regrid_stationary
545 :
546 : !=======================================================================
547 : !
548 : ! Aggregate flux information from all ice thickness categories
549 : ! for z layer biogeochemistry
550 : !
551 0 : subroutine merge_bgc_fluxes (dt, nblyr, &
552 : nslyr, & ! LCOV_EXCL_LINE
553 : bio_index, n_algae, & ! LCOV_EXCL_LINE
554 : nbtrcr, aicen, & ! LCOV_EXCL_LINE
555 : vicen, vsnon, & ! LCOV_EXCL_LINE
556 : iphin, & ! LCOV_EXCL_LINE
557 : trcrn, & ! LCOV_EXCL_LINE
558 : flux_bion, flux_bio, & ! LCOV_EXCL_LINE
559 : upNOn, upNHn, & ! LCOV_EXCL_LINE
560 : upNO, upNH, & ! LCOV_EXCL_LINE
561 : zbgc_snown, zbgc_atmn, & ! LCOV_EXCL_LINE
562 : zbgc_snow, zbgc_atm, & ! LCOV_EXCL_LINE
563 : PP_net, ice_bio_net,& ! LCOV_EXCL_LINE
564 : snow_bio_net, grow_alg, & ! LCOV_EXCL_LINE
565 : grow_net)
566 :
567 : real (kind=dbl_kind), intent(in) :: &
568 : dt ! timestep (s)
569 :
570 : integer (kind=int_kind), intent(in) :: &
571 : nblyr , & ! number of bio layers ! LCOV_EXCL_LINE
572 : nslyr , & ! number of snow layers ! LCOV_EXCL_LINE
573 : n_algae , & ! number of algal tracers ! LCOV_EXCL_LINE
574 : nbtrcr ! number of biology tracer tracers
575 :
576 : integer (kind=int_kind), dimension(:), intent(in) :: &
577 : bio_index ! relates bio indices, ie. nlt_bgc_N to nt_bgc_N
578 :
579 : real (kind=dbl_kind), dimension (:), intent(in) :: &
580 : trcrn , & ! input tracer fields ! LCOV_EXCL_LINE
581 : iphin ! porosity
582 :
583 : real (kind=dbl_kind), intent(in):: &
584 : aicen , & ! concentration of ice ! LCOV_EXCL_LINE
585 : vicen , & ! volume of ice (m) ! LCOV_EXCL_LINE
586 : vsnon ! volume of snow(m)
587 :
588 : ! single category rates
589 : real (kind=dbl_kind), dimension(:), intent(in):: &
590 : zbgc_snown , & ! bio flux from snow to ice per cat (mmol/m^3*m) ! LCOV_EXCL_LINE
591 : zbgc_atmn , & ! bio flux from atm to ice per cat (mmol/m^3*m) ! LCOV_EXCL_LINE
592 : flux_bion
593 :
594 : ! single category rates
595 : real (kind=dbl_kind), dimension(:,:), intent(in):: &
596 : upNOn , & ! nitrate uptake rate per cat (mmol/m^3/s) ! LCOV_EXCL_LINE
597 : upNHn , & ! ammonium uptake rate per cat (mmol/m^3/s) ! LCOV_EXCL_LINE
598 : grow_alg ! algal growth rate per cat (mmolN/m^3/s)
599 :
600 : ! cumulative fluxes
601 : real (kind=dbl_kind), dimension(:), intent(inout):: &
602 : flux_bio , & ! ! LCOV_EXCL_LINE
603 : zbgc_snow , & ! bio flux from snow to ice per cat (mmol/m^2/s) ! LCOV_EXCL_LINE
604 : zbgc_atm , & ! bio flux from atm to ice per cat (mmol/m^2/s) ! LCOV_EXCL_LINE
605 : ice_bio_net, & ! integrated ice tracers mmol or mg/m^2) ! LCOV_EXCL_LINE
606 : snow_bio_net ! integrated snow tracers mmol or mg/m^2)
607 :
608 : ! cumulative variables and rates
609 : real (kind=dbl_kind), intent(inout):: &
610 : PP_net , & ! net PP (mg C/m^2/d) times aice ! LCOV_EXCL_LINE
611 : grow_net , & ! net specific growth (m/d) times vice ! LCOV_EXCL_LINE
612 : upNO , & ! tot nitrate uptake rate (mmol/m^2/d) times aice ! LCOV_EXCL_LINE
613 : upNH ! tot ammonium uptake rate (mmol/m^2/d) times aice
614 :
615 : ! local variables
616 :
617 : real (kind=dbl_kind) :: &
618 : tmp , & ! temporary ! LCOV_EXCL_LINE
619 : dvssl , & ! volume of snow surface layer (m) ! LCOV_EXCL_LINE
620 0 : dvint ! volume of snow interior (m)
621 :
622 : integer (kind=int_kind) :: &
623 : k, mm ! tracer indice
624 :
625 : real (kind=dbl_kind), dimension (nblyr+1) :: &
626 0 : zspace
627 :
628 : character(len=*),parameter :: subname='(merge_bgc_fluxes)'
629 :
630 : !-----------------------------------------------------------------
631 : ! Column summation
632 : !-----------------------------------------------------------------
633 0 : zspace(:) = c1/real(nblyr,kind=dbl_kind)
634 0 : zspace(1) = p5/real(nblyr,kind=dbl_kind)
635 0 : zspace(nblyr+1) = p5/real(nblyr,kind=dbl_kind)
636 :
637 0 : do mm = 1, nbtrcr
638 0 : do k = 1, nblyr+1
639 0 : ice_bio_net(mm) = ice_bio_net(mm) &
640 : + trcrn(bio_index(mm)+k-1) & ! LCOV_EXCL_LINE
641 : * trcrn(nt_fbri) & ! LCOV_EXCL_LINE
642 0 : * vicen*zspace(k)
643 : enddo ! k
644 :
645 : !-----------------------------------------------------------------
646 : ! Merge fluxes
647 : !-----------------------------------------------------------------
648 0 : dvssl = min(p5*vsnon/real(nslyr,kind=dbl_kind), hs_ssl*aicen) ! snow surface layer
649 0 : dvint = vsnon - dvssl ! snow interior
650 0 : snow_bio_net(mm) = snow_bio_net(mm) &
651 : + trcrn(bio_index(mm)+nblyr+1)*dvssl & ! LCOV_EXCL_LINE
652 0 : + trcrn(bio_index(mm)+nblyr+2)*dvint
653 0 : flux_bio (mm) = flux_bio (mm) + flux_bion (mm)*aicen
654 0 : zbgc_snow (mm) = zbgc_snow(mm) + zbgc_snown(mm)*aicen/dt
655 0 : zbgc_atm (mm) = zbgc_atm (mm) + zbgc_atmn (mm)*aicen/dt
656 : enddo ! mm
657 :
658 0 : if (solve_zbgc) then
659 0 : do mm = 1, n_algae
660 0 : do k = 1, nblyr+1
661 0 : tmp = iphin(k)*trcrn(nt_fbri)*vicen*zspace(k)*secday
662 0 : PP_net = PP_net + grow_alg(k,mm)*tmp &
663 0 : * (c1-fr_resp)* R_C2N(mm)*R_gC2molC
664 0 : grow_net = grow_net + grow_alg(k,mm)*tmp &
665 0 : / (trcrn(nt_bgc_N(mm)+k-1)+puny)
666 0 : upNO = upNO + upNOn (k,mm)*tmp
667 0 : upNH = upNH + upNHn (k,mm)*tmp
668 : enddo ! k
669 : enddo ! mm
670 : endif
671 :
672 0 : end subroutine merge_bgc_fluxes
673 :
674 : !=======================================================================
675 :
676 : ! Aggregate flux information from all ice thickness categories
677 : ! for skeletal layer biogeochemistry
678 : !
679 : ! author: Elizabeth C. Hunke and William H. Lipscomb, LANL
680 :
681 0 : subroutine merge_bgc_fluxes_skl (nbtrcr, n_algae, &
682 : aicen, trcrn, & ! LCOV_EXCL_LINE
683 : flux_bion, flux_bio, & ! LCOV_EXCL_LINE
684 : PP_net, upNOn, & ! LCOV_EXCL_LINE
685 : upNHn, upNO, & ! LCOV_EXCL_LINE
686 : upNH, grow_net, & ! LCOV_EXCL_LINE
687 0 : grow_alg)
688 :
689 : integer (kind=int_kind), intent(in) :: &
690 : nbtrcr , & ! number of bgc tracers ! LCOV_EXCL_LINE
691 : n_algae ! number of autotrophs
692 :
693 : ! single category fluxes
694 : real (kind=dbl_kind), intent(in):: &
695 : aicen ! category ice area fraction
696 :
697 : real (kind=dbl_kind), dimension (:), intent(in) :: &
698 : trcrn ! Bulk tracer concentration (mmol N or mg/m^3)
699 :
700 : real (kind=dbl_kind), dimension(:), intent(in):: &
701 : flux_bion ! all bio fluxes to ocean, on categories
702 :
703 : real (kind=dbl_kind), dimension(:), intent(inout):: &
704 : flux_bio ! all bio fluxes to ocean, aggregated
705 :
706 : real (kind=dbl_kind), dimension(:), intent(in):: &
707 : grow_alg, & ! algal growth rate (mmol/m^3/s) ! LCOV_EXCL_LINE
708 : upNOn , & ! nitrate uptake rate per cat (mmol/m^3/s) ! LCOV_EXCL_LINE
709 : upNHn ! ammonium uptake rate per cat (mmol/m^3/s)
710 :
711 : ! history output
712 : real (kind=dbl_kind), intent(inout):: &
713 : PP_net , & ! Bulk net PP (mg C/m^2/d) ! LCOV_EXCL_LINE
714 : grow_net, & ! net specific growth (/d) ! LCOV_EXCL_LINE
715 : upNO , & ! tot nitrate uptake rate (mmol/m^2/d) ! LCOV_EXCL_LINE
716 : upNH ! tot ammonium uptake rate (mmol/m^2/d)
717 :
718 : ! local variables
719 :
720 : integer (kind=int_kind) :: &
721 : k, mm ! tracer indices
722 :
723 : real (kind=dbl_kind) :: &
724 0 : tmp ! temporary
725 :
726 : character(len=*),parameter :: subname='(merge_bgc_fluxes_skl)'
727 :
728 : !-----------------------------------------------------------------
729 : ! Merge fluxes
730 : !-----------------------------------------------------------------
731 :
732 0 : do k = 1,nbtrcr
733 0 : flux_bio (k) = flux_bio(k) + flux_bion(k)*aicen
734 : enddo
735 :
736 0 : do mm = 1, n_algae
737 0 : tmp = phi_sk * sk_l * aicen * secday
738 : PP_net = PP_net &
739 : + grow_alg(mm) * tmp & ! LCOV_EXCL_LINE
740 0 : * R_C2N(mm) * R_gC2molC * (c1-fr_resp)
741 : grow_net = grow_net &
742 : + grow_alg(mm) * tmp & ! LCOV_EXCL_LINE
743 0 : / (trcrn(nt_bgc_N(mm))+puny)
744 0 : upNO = upNO + upNOn(mm) * tmp
745 0 : upNH = upNH + upNHn(mm) * tmp
746 : enddo
747 :
748 0 : end subroutine merge_bgc_fluxes_skl
749 :
750 : !=======================================================================
751 :
752 : end module icepack_zbgc_shared
753 :
754 : !=======================================================================
|