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