Line data Source code
1 : !=========================================================================
2 : !
3 : ! Initialization routines for the column package.
4 : !
5 : ! author: Elizabeth C. Hunke, LANL
6 : !
7 : module icedrv_init_column
8 :
9 : use icedrv_kinds
10 : use icedrv_domain_size, only: ncat, nx, ncat
11 : use icedrv_domain_size, only: max_ntrcr, nblyr, nilyr, nslyr
12 : use icedrv_domain_size, only: n_algae, n_zaero, n_doc, n_dic, n_don
13 : use icedrv_domain_size, only: n_fed, n_fep, max_nsw, n_bgc, n_aero
14 : use icedrv_constants, only: c1, c2, p5, c0, p1
15 : use icedrv_constants, only: nu_diag, nu_nml
16 : use icepack_intfc, only: icepack_max_don, icepack_max_doc, icepack_max_dic
17 : use icepack_intfc, only: icepack_max_algae, icepack_max_aero, icepack_max_fe
18 : use icepack_intfc, only: icepack_max_nbtrcr
19 : use icepack_intfc, only: icepack_warnings_flush, icepack_warnings_aborted
20 : use icepack_intfc, only: icepack_init_tracer_sizes, icepack_init_tracer_flags
21 : use icepack_intfc, only: icepack_init_tracer_indices
22 : use icepack_intfc, only: icepack_init_parameters
23 : use icepack_intfc, only: icepack_query_tracer_sizes, icepack_query_tracer_flags
24 : use icepack_intfc, only: icepack_query_tracer_indices
25 : use icepack_intfc, only: icepack_query_parameters
26 : use icepack_intfc, only: icepack_init_zbgc
27 : use icepack_intfc, only: icepack_init_thermo
28 : use icepack_intfc, only: icepack_step_radiation, icepack_init_orbit
29 : use icepack_intfc, only: icepack_init_bgc, icepack_init_zsalinity
30 : use icepack_intfc, only: icepack_init_ocean_bio, icepack_load_ocean_bio_array
31 : use icepack_intfc, only: icepack_init_hbrine
32 : use icedrv_system, only: icedrv_system_abort
33 :
34 : implicit none
35 :
36 : private
37 : public :: init_thermo_vertical, init_shortwave, &
38 : init_bgc, init_hbrine, init_zbgc
39 :
40 : !=======================================================================
41 :
42 : contains
43 :
44 : !=======================================================================
45 : !
46 : ! Initialize the vertical profile of ice salinity and melting temperature.
47 : !
48 : ! authors: C. M. Bitz, UW
49 : ! William H. Lipscomb, LANL
50 :
51 45 : subroutine init_thermo_vertical
52 :
53 : use icedrv_flux, only: salinz, Tmltz
54 :
55 : integer (kind=int_kind) :: &
56 : i, & ! horizontal indices
57 : k ! ice layer index
58 :
59 : real (kind=dbl_kind), dimension(nilyr+1) :: &
60 132 : sprofile ! vertical salinity profile
61 :
62 : real (kind=dbl_kind) :: &
63 16 : depressT
64 :
65 : character(len=*), parameter :: subname='(init_thermo_vertical)'
66 :
67 : !-----------------------------------------------------------------
68 : ! query Icepack values
69 : !-----------------------------------------------------------------
70 :
71 45 : call icepack_query_parameters(depressT_out=depressT)
72 45 : call icepack_init_thermo(nilyr=nilyr, sprofile=sprofile)
73 45 : call icepack_warnings_flush(nu_diag)
74 45 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
75 0 : file=__FILE__, line=__LINE__)
76 :
77 : !-----------------------------------------------------------------
78 : ! Prescibe vertical profile of salinity and melting temperature.
79 : ! Note this profile is only used for BL99 thermodynamics.
80 : !-----------------------------------------------------------------
81 :
82 225 : do i = 1, nx
83 1521 : do k = 1, nilyr+1
84 1296 : salinz(i,k) = sprofile(k)
85 1476 : Tmltz (i,k) = -salinz(i,k)*depressT
86 : enddo ! k
87 : enddo ! i
88 :
89 45 : end subroutine init_thermo_vertical
90 :
91 : !=======================================================================
92 : !
93 : ! Initialize shortwave
94 :
95 45 : subroutine init_shortwave
96 :
97 : use icedrv_arrays_column, only: fswpenln, Iswabsn, Sswabsn, albicen
98 : use icedrv_arrays_column, only: albsnon, alvdrn, alidrn, alvdfn, alidfn
99 : use icedrv_arrays_column, only: fswsfcn, fswthrun, ffracn, snowfracn
100 : use icedrv_arrays_column, only: fswintn, albpndn, apeffn, trcrn_sw, dhsn
101 : use icedrv_arrays_column, only: kaer_tab, waer_tab, gaer_tab
102 : use icedrv_arrays_column, only: kaer_bc_tab, waer_bc_tab, gaer_bc_tab
103 : use icedrv_arrays_column, only: swgrid, igrid, bcenh
104 : use icedrv_calendar, only: istep1, dt, calendar_type
105 : use icedrv_calendar, only: days_per_year, nextsw_cday, yday, sec
106 : use icedrv_system, only: icedrv_system_abort
107 : use icedrv_flux, only: alvdf, alidf, alvdr, alidr
108 : use icedrv_flux, only: alvdr_ai, alidr_ai, alvdf_ai, alidf_ai
109 : use icedrv_flux, only: swvdr, swvdf, swidr, swidf, scale_factor, snowfrac
110 : use icedrv_flux, only: albice, albsno, albpnd, apeff_ai, coszen, fsnow
111 : use icedrv_init, only: tlat, tlon, tmask
112 : use icedrv_restart_shared, only: restart
113 : use icedrv_state, only: aicen, vicen, vsnon, trcrn
114 :
115 : integer (kind=int_kind) :: &
116 : i, k , & ! horizontal indices
117 : n ! thickness category index
118 :
119 : real (kind=dbl_kind) :: &
120 16 : netsw ! flag for shortwave radiation presence
121 :
122 : logical (kind=log_kind) :: &
123 : l_print_point, & ! flag to print designated grid point diagnostics
124 : dEdd_algae, & ! from icepack
125 : modal_aero ! from icepack
126 :
127 : character (len=char_len) :: &
128 : shortwave ! from icepack
129 :
130 : real (kind=dbl_kind), dimension(ncat) :: &
131 92 : fbri ! brine height to ice thickness
132 :
133 : real (kind=dbl_kind), allocatable, dimension(:,:) :: &
134 45 : ztrcr_sw
135 :
136 : logical (kind=log_kind) :: tr_brine, tr_zaero, tr_bgc_N
137 : integer (kind=int_kind) :: nt_alvl, nt_apnd, nt_hpnd, nt_ipnd, nt_aero, &
138 : nt_fbri, nt_tsfc, ntrcr, nbtrcr_sw, nlt_chl_sw
139 : integer (kind=int_kind), dimension(icepack_max_aero) :: nlt_zaero_sw
140 : integer (kind=int_kind), dimension(icepack_max_aero) :: nt_zaero
141 : integer (kind=int_kind), dimension(icepack_max_algae) :: nt_bgc_N
142 16 : real (kind=dbl_kind) :: puny
143 :
144 : character(len=*), parameter :: subname='(init_shortwave)'
145 :
146 : !-----------------------------------------------------------------
147 : ! query Icepack values
148 : !-----------------------------------------------------------------
149 :
150 45 : call icepack_query_parameters(puny_out=puny)
151 45 : call icepack_query_parameters(shortwave_out=shortwave)
152 45 : call icepack_query_parameters(dEdd_algae_out=dEdd_algae)
153 45 : call icepack_query_parameters(modal_aero_out=modal_aero)
154 : call icepack_query_tracer_sizes(ntrcr_out=ntrcr, &
155 45 : nbtrcr_sw_out=nbtrcr_sw)
156 : call icepack_query_tracer_flags(tr_brine_out=tr_brine, &
157 45 : tr_zaero_out=tr_zaero, tr_bgc_N_out=tr_bgc_N)
158 : call icepack_query_tracer_indices(nt_alvl_out=nt_alvl, &
159 : nt_apnd_out=nt_apnd, nt_hpnd_out=nt_hpnd, &
160 : nt_ipnd_out=nt_ipnd, nt_aero_out=nt_aero, &
161 : nt_fbri_out=nt_fbri, nt_tsfc_out=nt_tsfc, &
162 : nt_bgc_N_out=nt_bgc_N, nt_zaero_out=nt_zaero, &
163 45 : nlt_chl_sw_out=nlt_chl_sw, nlt_zaero_sw_out=nlt_zaero_sw)
164 45 : call icepack_warnings_flush(nu_diag)
165 45 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
166 0 : file=__FILE__,line= __LINE__)
167 :
168 : !-----------------------------------------------------------------
169 : ! Initialize
170 : !-----------------------------------------------------------------
171 :
172 45 : allocate(ztrcr_sw(nbtrcr_sw, ncat))
173 :
174 45 : fswpenln(:,:,:) = c0
175 45 : Iswabsn(:,:,:) = c0
176 45 : Sswabsn(:,:,:) = c0
177 :
178 225 : do i = 1, nx
179 :
180 180 : l_print_point = .false.
181 :
182 180 : alvdf(i) = c0
183 180 : alidf(i) = c0
184 180 : alvdr(i) = c0
185 180 : alidr(i) = c0
186 180 : alvdr_ai(i) = c0
187 180 : alidr_ai(i) = c0
188 180 : alvdf_ai(i) = c0
189 180 : alidf_ai(i) = c0
190 180 : albice(i) = c0
191 180 : albsno(i) = c0
192 180 : albpnd(i) = c0
193 180 : snowfrac(i) = c0
194 180 : apeff_ai(i) = c0
195 :
196 1077 : do n = 1, ncat
197 852 : alvdrn(i,n) = c0
198 852 : alidrn(i,n) = c0
199 852 : alvdfn(i,n) = c0
200 852 : alidfn(i,n) = c0
201 852 : fswsfcn(i,n) = c0
202 852 : fswintn(i,n) = c0
203 1032 : fswthrun(i,n) = c0
204 : enddo ! ncat
205 :
206 : enddo
207 :
208 225 : do i = 1, nx
209 :
210 180 : if (trim(shortwave) == 'dEdd') then ! delta Eddington
211 :
212 : ! initialize orbital parameters
213 : ! These come from the driver in the coupled model.
214 156 : call icepack_warnings_flush(nu_diag)
215 156 : call icepack_init_orbit()
216 156 : call icepack_warnings_flush(nu_diag)
217 156 : if (icepack_warnings_aborted()) &
218 0 : call icedrv_system_abort(i, istep1, subname, __FILE__, __LINE__)
219 : endif
220 :
221 180 : fbri(:) = c0
222 1032 : ztrcr_sw(:,:) = c0
223 1032 : do n = 1, ncat
224 1032 : if (tr_brine) fbri(n) = trcrn(i,nt_fbri,n)
225 : enddo
226 :
227 180 : if (tmask(i)) then
228 : call icepack_step_radiation ( &
229 : dt=dt, ncat=ncat, &
230 : nblyr=nblyr, &
231 : nilyr=nilyr, nslyr=nslyr, &
232 : dEdd_algae=dEdd_algae, &
233 : swgrid=swgrid(:), &
234 : igrid=igrid(:), &
235 : fbri=fbri(:), &
236 : aicen=aicen(i,:), &
237 : vicen=vicen(i,:), &
238 : vsnon=vsnon(i,:), &
239 : Tsfcn=trcrn(i,nt_Tsfc,:), &
240 : alvln=trcrn(i,nt_alvl,:), &
241 : apndn=trcrn(i,nt_apnd,:), &
242 : hpndn=trcrn(i,nt_hpnd,:), &
243 : ipndn=trcrn(i,nt_ipnd,:), &
244 0 : aeron=trcrn(i,nt_aero:nt_aero+4*n_aero-1,:), &
245 0 : bgcNn=trcrn(i,nt_bgc_N(1):nt_bgc_N(1)+n_algae*(nblyr+3)-1,:), &
246 0 : zaeron=trcrn(i,nt_zaero(1):nt_zaero(1)+n_zaero*(nblyr+3)-1,:), &
247 : trcrn_bgcsw=ztrcr_sw, &
248 48 : TLAT=TLAT(i), TLON=TLON(i), &
249 : calendar_type=calendar_type, &
250 : days_per_year=days_per_year, &
251 : nextsw_cday=nextsw_cday, yday=yday, sec=sec, &
252 : kaer_tab=kaer_tab, kaer_bc_tab=kaer_bc_tab(:,:), &
253 : waer_tab=waer_tab, waer_bc_tab=waer_bc_tab(:,:), &
254 : gaer_tab=gaer_tab, gaer_bc_tab=gaer_bc_tab(:,:), &
255 : bcenh=bcenh(:,:,:), &
256 : modal_aero=modal_aero, &
257 48 : swvdr=swvdr(i), swvdf=swvdf(i), &
258 48 : swidr=swidr(i), swidf=swidf(i), &
259 48 : coszen=coszen(i), fsnow=fsnow(i), &
260 : alvdrn=alvdrn(i,:), alvdfn=alvdfn(i,:), &
261 : alidrn=alidrn(i,:), alidfn=alidfn(i,:), &
262 : fswsfcn=fswsfcn(i,:), fswintn=fswintn(i,:), &
263 : fswthrun=fswthrun(i,:), fswpenln=fswpenln(i,:,:), &
264 : Sswabsn=Sswabsn(i,:,:), Iswabsn=Iswabsn(i,:,:), &
265 : albicen=albicen(i,:), albsnon=albsnon(i,:), &
266 : albpndn=albpndn(i,:), apeffn=apeffn(i,:), &
267 : snowfracn=snowfracn(i,:), &
268 : dhsn=dhsn(i,:), ffracn=ffracn(i,:), &
269 : l_print_point=l_print_point, &
270 231 : initonly = .true.)
271 : endif
272 180 : call icepack_warnings_flush(nu_diag)
273 180 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
274 0 : file=__FILE__, line=__LINE__)
275 :
276 : !-----------------------------------------------------------------
277 : ! Define aerosol tracer on shortwave grid
278 : !-----------------------------------------------------------------
279 :
280 225 : if (dEdd_algae .and. (tr_zaero .or. tr_bgc_N)) then
281 0 : do n = 1, ncat
282 180 : do k = 1, nbtrcr_sw
283 0 : trcrn_sw(i,k,n) = ztrcr_sw(k,n)
284 : enddo
285 : enddo
286 : endif
287 :
288 : enddo
289 :
290 : !-----------------------------------------------------------------
291 : ! Aggregate albedos
292 : !-----------------------------------------------------------------
293 :
294 258 : do n = 1, ncat
295 1110 : do i = 1, nx
296 :
297 1065 : if (aicen(i,n) > puny) then
298 :
299 366 : alvdf(i) = alvdf(i) + alvdfn(i,n)*aicen(i,n)
300 366 : alidf(i) = alidf(i) + alidfn(i,n)*aicen(i,n)
301 366 : alvdr(i) = alvdr(i) + alvdrn(i,n)*aicen(i,n)
302 366 : alidr(i) = alidr(i) + alidrn(i,n)*aicen(i,n)
303 :
304 366 : netsw = swvdr(i) + swidr(i) + swvdf(i) + swidf(i)
305 366 : if (netsw > puny) then ! sun above horizon
306 227 : albice(i) = albice(i) + albicen(i,n)*aicen(i,n)
307 227 : albsno(i) = albsno(i) + albsnon(i,n)*aicen(i,n)
308 227 : albpnd(i) = albpnd(i) + albpndn(i,n)*aicen(i,n)
309 : endif
310 :
311 366 : apeff_ai(i) = apeff_ai(i) + apeffn(i,n)*aicen(i,n)
312 366 : snowfrac(i) = snowfrac(i) + snowfracn(i,n)*aicen(i,n)
313 :
314 : endif ! aicen > puny
315 : enddo ! i
316 : enddo ! ncat
317 :
318 225 : do i = 1, nx
319 :
320 : !----------------------------------------------------------------
321 : ! Store grid box mean albedos and fluxes before scaling by aice
322 : !----------------------------------------------------------------
323 :
324 180 : alvdf_ai (i) = alvdf (i)
325 180 : alidf_ai (i) = alidf (i)
326 180 : alvdr_ai (i) = alvdr (i)
327 180 : alidr_ai (i) = alidr (i)
328 :
329 : !----------------------------------------------------------------
330 : ! Save net shortwave for scaling factor in scale_factor
331 : !----------------------------------------------------------------
332 225 : if (.not. restart) then
333 60 : scale_factor(i) = &
334 60 : swvdr(i)*(c1 - alvdr_ai(i)) &
335 60 : + swvdf(i)*(c1 - alvdf_ai(i)) &
336 60 : + swidr(i)*(c1 - alidr_ai(i)) &
337 120 : + swidf(i)*(c1 - alidf_ai(i))
338 : endif
339 :
340 : enddo ! i
341 :
342 45 : deallocate(ztrcr_sw)
343 :
344 45 : end subroutine init_shortwave
345 :
346 : !=======================================================================
347 :
348 : ! Initialize vertical profile for biogeochemistry
349 :
350 5 : subroutine init_bgc()
351 :
352 : use icedrv_arrays_column, only: zfswin, trcrn_sw
353 : use icedrv_arrays_column, only: ocean_bio_all, ice_bio_net, snow_bio_net
354 : use icedrv_arrays_column, only: cgrid, igrid, bphi, iDi, bTiz, iki
355 : use icedrv_arrays_column, only: Rayleigh_criteria, Rayleigh_real
356 : use icedrv_calendar, only: istep1
357 : use icedrv_system, only: icedrv_system_abort
358 : use icedrv_flux, only: sss, nit, amm, sil, dmsp, dms, algalN, &
359 : doc, don, dic, fed, fep, zaeros, hum
360 : use icedrv_forcing_bgc, only: get_forcing_bgc
361 : use icedrv_state, only: trcrn
362 :
363 : ! local variables
364 :
365 : integer (kind=int_kind) :: &
366 : i , & ! horizontal indices
367 : k , & ! vertical index
368 : n ! category index
369 :
370 : integer (kind=int_kind) :: &
371 : max_nbtrcr, max_algae, max_don, max_doc, max_dic, max_aero, max_fe
372 :
373 : logical (kind=log_kind) :: &
374 : RayleighC, &
375 : solve_zsal
376 :
377 : real(kind=dbl_kind) :: &
378 3 : RayleighR
379 :
380 : real(kind=dbl_kind), dimension(max_ntrcr,ncat) :: &
381 2923 : trcrn_bgc
382 :
383 : real(kind=dbl_kind), dimension(nilyr,ncat) :: &
384 123 : sicen
385 :
386 : integer (kind=int_kind) :: &
387 : nbtrcr, ntrcr, ntrcr_o, &
388 : nt_sice, nt_bgc_S
389 :
390 : character(len=*), parameter :: subname='(init_bgc)'
391 :
392 : ! Initialize
393 :
394 5 : bphi(:,:,:) = c0 ! initial porosity for no ice
395 5 : iDi (:,:,:) = c0 ! interface diffusivity
396 5 : bTiz(:,:,:) = c0 ! initial bio grid ice temperature
397 5 : iki (:,:,:) = c0 ! permeability
398 :
399 5 : ocean_bio_all(:,:) = c0
400 5 : ice_bio_net (:,:) = c0 ! integrated ice tracer conc (mmol/m^2 or mg/m^2)
401 5 : snow_bio_net (:,:) = c0 ! integrated snow tracer conc (mmol/m^2 or mg/m^2)
402 5 : zfswin (:,:,:) = c0 ! shortwave flux on bio grid
403 5 : trcrn_sw (:,:,:) = c0 ! tracers active in the shortwave calculation
404 5 : trcrn_bgc (:,:) = c0
405 :
406 : !-----------------------------------------------------------------
407 : ! query Icepack values
408 : !-----------------------------------------------------------------
409 :
410 5 : call icepack_query_parameters(solve_zsal_out=solve_zsal)
411 : call icepack_query_tracer_sizes(max_nbtrcr_out=max_nbtrcr, &
412 : max_algae_out=max_algae, max_don_out=max_don, max_doc_out=max_doc, &
413 5 : max_dic_out=max_dic, max_aero_out=max_aero, max_fe_out=max_fe)
414 : call icepack_query_tracer_sizes(nbtrcr_out=nbtrcr, ntrcr_out=ntrcr, &
415 5 : ntrcr_o_out=ntrcr_o)
416 5 : call icepack_query_tracer_indices(nt_sice_out=nt_sice, nt_bgc_S_out=nt_bgc_S)
417 5 : call icepack_warnings_flush(nu_diag)
418 5 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
419 0 : file=__FILE__,line= __LINE__)
420 :
421 : !-----------------------------------------------------------------
422 : ! zsalinity initialization
423 : !-----------------------------------------------------------------
424 :
425 5 : if (solve_zsal) then
426 0 : do i = 1, nx
427 : call icepack_init_zsalinity(ncat=ncat, nblyr=nblyr, &
428 : ntrcr_o = ntrcr_o, &
429 : Rayleigh_criteria = RayleighC, &
430 : Rayleigh_real = RayleighR, &
431 : trcrn_bgc = trcrn_bgc, &
432 : nt_bgc_S = nt_bgc_S, &
433 0 : sss = sss(i))
434 0 : Rayleigh_real (i) = RayleighR
435 0 : Rayleigh_criteria(i) = RayleighC
436 0 : do n = 1,ncat
437 0 : do k = 1, nblyr
438 0 : trcrn(i,nt_bgc_S+k-1,n) = trcrn_bgc(nt_bgc_S-1+k-ntrcr_o,n)
439 : enddo
440 : enddo
441 : enddo ! i
442 0 : call icepack_warnings_flush(nu_diag)
443 0 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
444 0 : file=__FILE__, line=__LINE__)
445 : endif ! solve_zsal
446 :
447 : !-----------------------------------------------------------------
448 : ! biogeochemistry initialization
449 : !-----------------------------------------------------------------
450 :
451 : !-----------------------------------------------------------------
452 : ! Initial Ocean Values if not coupled to the ocean bgc
453 : !-----------------------------------------------------------------
454 25 : do i = 1, nx
455 : call icepack_init_ocean_bio ( &
456 12 : amm=amm(i), dmsp=dmsp(i), dms=dms(i), &
457 : doc=doc(i,:), dic =dic(i,:), don=don(i,:), &
458 12 : fed=fed(i,:), fep =fep(i,:), hum=hum(i), &
459 12 : nit=nit(i), sil =sil(i), &
460 : zaeros=zaeros(i,:), &
461 : algalN=algalN(i,:), &
462 : max_dic=max_dic, max_don =max_don, &
463 37 : max_fe =max_fe, max_aero=max_aero)
464 : enddo ! i
465 5 : call icepack_warnings_flush(nu_diag)
466 5 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
467 0 : file=__FILE__, line=__LINE__)
468 :
469 5 : call get_forcing_bgc ! defines nit and sil
470 :
471 25 : do i = 1, nx
472 :
473 120 : do n = 1, ncat
474 800 : do k = 1, nilyr
475 800 : sicen(k,n) = trcrn(i,nt_sice+k-1,n)
476 : enddo
477 17260 : do k = ntrcr_o+1, ntrcr
478 17240 : trcrn_bgc(k-ntrcr_o,n) = trcrn(i,k,n)
479 : enddo
480 : enddo
481 :
482 : call icepack_load_ocean_bio_array(max_nbtrcr=max_nbtrcr, &
483 : max_algae=max_algae, max_don=max_don, max_doc=max_doc, &
484 : max_aero =max_aero, max_dic=max_dic, max_fe =max_fe, &
485 12 : nit =nit(i), amm=amm(i), sil =sil(i), &
486 12 : dmsp=dmsp(i), dms=dms(i), algalN=algalN(i,:), &
487 : doc =doc(i,:), don=don(i,:), dic =dic(i,:), &
488 : fed =fed(i,:), fep=fep(i,:), zaeros=zaeros(i,:), &
489 20 : ocean_bio_all=ocean_bio_all(i,:), hum=hum(i))
490 20 : call icepack_warnings_flush(nu_diag)
491 20 : if (icepack_warnings_aborted()) call icedrv_system_abort(i, istep1, subname, &
492 5 : __FILE__, __LINE__)
493 :
494 : enddo ! i
495 :
496 25 : do i = 1, nx
497 : call icepack_init_bgc(ncat=ncat, nblyr=nblyr, nilyr=nilyr, ntrcr_o=ntrcr_o, &
498 : cgrid=cgrid, igrid=igrid, ntrcr=ntrcr, nbtrcr=nbtrcr, &
499 : sicen=sicen(:,:), trcrn=trcrn_bgc(:,:), &
500 20 : sss=sss(i), ocean_bio_all=ocean_bio_all(i,:))
501 20 : call icepack_warnings_flush(nu_diag)
502 20 : if (icepack_warnings_aborted()) call icedrv_system_abort(i, istep1, subname, &
503 5 : __FILE__, __LINE__)
504 : enddo ! i
505 :
506 5 : end subroutine init_bgc
507 :
508 : !=======================================================================
509 :
510 : ! Initialize brine height tracer
511 :
512 5 : subroutine init_hbrine()
513 :
514 : use icedrv_arrays_column, only: first_ice, bgrid, igrid, cgrid
515 : use icedrv_arrays_column, only: icgrid, swgrid
516 : use icedrv_state, only: trcrn
517 :
518 3 : real (kind=dbl_kind) :: phi_snow
519 : integer (kind=int_kind) :: nt_fbri
520 : logical (kind=log_kind) :: tr_brine
521 : character(len=*), parameter :: subname='(init_hbrine)'
522 :
523 : !-----------------------------------------------------------------
524 : ! query Icepack values
525 : !-----------------------------------------------------------------
526 :
527 5 : call icepack_query_parameters(phi_snow_out=phi_snow)
528 5 : call icepack_query_tracer_flags(tr_brine_out=tr_brine)
529 5 : call icepack_query_tracer_indices(nt_fbri_out=nt_fbri)
530 5 : call icepack_warnings_flush(nu_diag)
531 5 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
532 0 : file=__FILE__,line= __LINE__)
533 :
534 : !-----------------------------------------------------------------
535 :
536 : call icepack_init_hbrine(bgrid=bgrid, igrid=igrid, cgrid=cgrid, icgrid=icgrid, &
537 5 : swgrid=swgrid, nblyr=nblyr, nilyr=nilyr, phi_snow=phi_snow)
538 5 : call icepack_warnings_flush(nu_diag)
539 5 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
540 0 : file=__FILE__, line=__LINE__)
541 :
542 5 : call icepack_init_parameters(phi_snow_in=phi_snow)
543 5 : call icepack_warnings_flush(nu_diag)
544 5 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
545 0 : file=__FILE__,line= __LINE__)
546 :
547 130 : first_ice(:,:) = .true.
548 130 : if (tr_brine) trcrn(:,nt_fbri,:) = c1
549 :
550 5 : end subroutine init_hbrine
551 :
552 : !=======================================================================
553 :
554 : ! Namelist variables, set to default values; may be altered at run time
555 : !
556 : ! author Elizabeth C. Hunke, LANL
557 : ! Nicole Jeffery, LANL
558 :
559 45 : subroutine init_zbgc
560 :
561 : use icedrv_state, only: trcr_base, trcr_depend, n_trcr_strata
562 : use icedrv_state, only: nt_strata
563 : use icedrv_forcing, only: bgc_data_type
564 :
565 : character (len=char_len) :: &
566 : shortwave ! from icepack
567 :
568 : integer (kind=int_kind) :: &
569 : ntrcr, nbtrcr, nbtrcr_sw, &
570 : ntrcr_o, nt_fbri, &
571 : nt_bgc_Nit, nt_bgc_Am, nt_bgc_Sil, &
572 : nt_bgc_DMS, nt_bgc_PON, nt_bgc_S, &
573 : nt_bgc_DMSPp, nt_bgc_DMSPd, &
574 : nt_zbgc_frac, nlt_chl_sw, &
575 : nlt_bgc_Nit, nlt_bgc_Am, nlt_bgc_Sil, &
576 : nlt_bgc_DMS, nlt_bgc_DMSPp,nlt_bgc_DMSPd,&
577 : nlt_bgc_PON, nt_bgc_hum, nlt_bgc_hum
578 :
579 : integer (kind=int_kind), dimension(icepack_max_aero) :: &
580 : nlt_zaero_sw ! points to aerosol in trcrn_sw
581 :
582 : integer (kind=int_kind), dimension(icepack_max_algae) :: &
583 : nlt_bgc_N , & ! algae
584 : nlt_bgc_C , & !
585 : nlt_bgc_chl
586 :
587 : integer (kind=int_kind), dimension(icepack_max_doc) :: &
588 : nlt_bgc_DOC ! disolved organic carbon
589 :
590 : integer (kind=int_kind), dimension(icepack_max_don) :: &
591 : nlt_bgc_DON !
592 :
593 : integer (kind=int_kind), dimension(icepack_max_dic) :: &
594 : nlt_bgc_DIC ! disolved inorganic carbon
595 :
596 : integer (kind=int_kind), dimension(icepack_max_fe) :: &
597 : nlt_bgc_Fed , & !
598 : nlt_bgc_Fep !
599 :
600 : integer (kind=int_kind), dimension(icepack_max_aero) :: &
601 : nlt_zaero ! non-reacting layer aerosols
602 :
603 : integer (kind=int_kind), dimension(icepack_max_algae) :: &
604 : nt_bgc_N , & ! diatoms, phaeocystis, pico/small
605 : nt_bgc_C , & ! diatoms, phaeocystis, pico/small
606 : nt_bgc_chl ! diatoms, phaeocystis, pico/small
607 :
608 : integer (kind=int_kind), dimension(icepack_max_doc) :: &
609 : nt_bgc_DOC ! dissolved organic carbon
610 :
611 : integer (kind=int_kind), dimension(icepack_max_don) :: &
612 : nt_bgc_DON ! dissolved organic nitrogen
613 :
614 : integer (kind=int_kind), dimension(icepack_max_dic) :: &
615 : nt_bgc_DIC ! dissolved inorganic carbon
616 :
617 : integer (kind=int_kind), dimension(icepack_max_fe) :: &
618 : nt_bgc_Fed , & ! dissolved iron
619 : nt_bgc_Fep ! particulate iron
620 :
621 : integer (kind=int_kind), dimension(icepack_max_aero) :: &
622 : nt_zaero , & ! black carbon and other aerosols
623 : nt_zaero_sw ! for shortwave calculation
624 :
625 : integer (kind=int_kind), dimension(icepack_max_nbtrcr) :: &
626 : bio_index_o ! relates nlt_bgc_NO to ocean concentration index
627 :
628 : integer (kind=int_kind), dimension(icepack_max_nbtrcr) :: &
629 : bio_index ! relates bio indices nlt_bgc_N to nt_bgc_N
630 :
631 : logical (kind=log_kind) :: &
632 : tr_brine, &
633 : tr_bgc_Nit, tr_bgc_Am, tr_bgc_Sil, &
634 : tr_bgc_DMS, tr_bgc_PON, &
635 : tr_bgc_N, tr_bgc_C, tr_bgc_chl, &
636 : tr_bgc_DON, tr_bgc_Fe, tr_zaero, &
637 : tr_bgc_hum, tr_aero
638 :
639 : integer (kind=int_kind) :: &
640 : ktherm
641 :
642 : logical (kind=log_kind) :: &
643 : solve_zsal, skl_bgc, z_tracers, scale_bgc, solve_zbgc, dEdd_algae, &
644 : modal_aero, restore_bgc
645 :
646 : character (char_len) :: &
647 : bgc_flux_type
648 :
649 : real (kind=dbl_kind) :: &
650 48 : grid_o, l_sk, initbio_frac, &
651 48 : frazil_scav, grid_oS, l_skS, &
652 16 : phi_snow, &
653 16 : ratio_Si2N_diatoms , ratio_Si2N_sp , ratio_Si2N_phaeo , &
654 16 : ratio_S2N_diatoms , ratio_S2N_sp , ratio_S2N_phaeo , &
655 16 : ratio_Fe2C_diatoms , ratio_Fe2C_sp , ratio_Fe2C_phaeo , &
656 16 : ratio_Fe2N_diatoms , ratio_Fe2N_sp , ratio_Fe2N_phaeo , &
657 16 : ratio_Fe2DON , ratio_Fe2DOC_s , ratio_Fe2DOC_l , &
658 32 : fr_resp , tau_min , tau_max , &
659 48 : algal_vel , R_dFe2dust , dustFe_sol , &
660 16 : chlabs_diatoms , chlabs_sp , chlabs_phaeo , &
661 16 : alpha2max_low_diatoms,alpha2max_low_sp , alpha2max_low_phaeo, &
662 16 : beta2max_diatoms , beta2max_sp , beta2max_phaeo , &
663 16 : mu_max_diatoms , mu_max_sp , mu_max_phaeo , &
664 16 : grow_Tdep_diatoms , grow_Tdep_sp , grow_Tdep_phaeo , &
665 48 : fr_graze_diatoms , fr_graze_sp , fr_graze_phaeo , &
666 16 : mort_pre_diatoms , mort_pre_sp , mort_pre_phaeo , &
667 16 : mort_Tdep_diatoms , mort_Tdep_sp , mort_Tdep_phaeo , &
668 16 : k_exude_diatoms , k_exude_sp , k_exude_phaeo , &
669 16 : K_Nit_diatoms , K_Nit_sp , K_Nit_phaeo , &
670 16 : K_Am_diatoms , K_Am_sp , K_Am_phaeo , &
671 16 : K_Sil_diatoms , K_Sil_sp , K_Sil_phaeo , &
672 16 : K_Fe_diatoms , K_Fe_sp , K_Fe_phaeo , &
673 32 : f_don_protein , kn_bac_protein , f_don_Am_protein , &
674 32 : f_doc_s , f_doc_l , f_exude_s , &
675 32 : f_exude_l , k_bac_s , k_bac_l , &
676 48 : T_max , fsal , op_dep_min , &
677 48 : fr_graze_s , fr_graze_e , fr_mort2min , &
678 48 : fr_dFe , k_nitrif , t_iron_conv , &
679 32 : max_loss , max_dfe_doc1 , fr_resp_s , &
680 32 : y_sk_DMS , t_sk_conv , t_sk_ox , &
681 16 : algaltype_diatoms , algaltype_sp , algaltype_phaeo , &
682 48 : nitratetype , ammoniumtype , silicatetype , &
683 32 : dmspptype , dmspdtype , humtype , &
684 32 : doctype_s , doctype_l , dontype_protein , &
685 48 : fedtype_1 , feptype_1 , zaerotype_bc1 , &
686 16 : zaerotype_bc2 , zaerotype_dust1 , zaerotype_dust2 , &
687 32 : zaerotype_dust3 , zaerotype_dust4 , ratio_C2N_diatoms , &
688 32 : ratio_C2N_sp , ratio_C2N_phaeo , ratio_chl2N_diatoms, &
689 32 : ratio_chl2N_sp , ratio_chl2N_phaeo , F_abs_chl_diatoms , &
690 32 : F_abs_chl_sp , F_abs_chl_phaeo , ratio_C2N_proteins
691 :
692 : real (kind=dbl_kind), dimension(icepack_max_dic) :: &
693 32 : dictype
694 :
695 : real (kind=dbl_kind), dimension(icepack_max_algae) :: &
696 64 : algaltype ! tau_min for both retention and release
697 :
698 : real (kind=dbl_kind), dimension(icepack_max_doc) :: &
699 64 : doctype
700 :
701 : real (kind=dbl_kind), dimension(icepack_max_don) :: &
702 32 : dontype
703 :
704 : real (kind=dbl_kind), dimension(icepack_max_fe) :: &
705 48 : fedtype
706 :
707 : real (kind=dbl_kind), dimension(icepack_max_fe) :: &
708 48 : feptype
709 :
710 : real (kind=dbl_kind), dimension(icepack_max_aero) :: &
711 112 : zaerotype
712 :
713 : real (kind=dbl_kind), dimension(icepack_max_algae) :: &
714 64 : R_C2N , & ! algal C to N (mole/mole)
715 64 : R_chl2N , & ! 3 algal chlorophyll to N (mg/mmol)
716 64 : F_abs_chl ! to scale absorption in Dedd
717 :
718 : real (kind=dbl_kind), dimension(icepack_max_don) :: & ! increase compare to algal R_Fe2C
719 32 : R_C2N_DON
720 :
721 : real (kind=dbl_kind), dimension(icepack_max_algae) :: &
722 64 : R_Si2N , & ! algal Sil to N (mole/mole)
723 64 : R_S2N , & ! algal S to N (mole/mole)
724 : ! Marchetti et al 2006, 3 umol Fe/mol C for iron limited Pseudo-nitzschia
725 64 : R_Fe2C , & ! algal Fe to carbon (umol/mmol)
726 64 : R_Fe2N ! algal Fe to N (umol/mmol)
727 :
728 : real (kind=dbl_kind), dimension(icepack_max_don) :: &
729 32 : R_Fe2DON ! Fe to N of DON (nmol/umol)
730 :
731 : real (kind=dbl_kind), dimension(icepack_max_doc) :: &
732 64 : R_Fe2DOC ! Fe to C of DOC (nmol/umol)
733 :
734 : real (kind=dbl_kind), dimension(icepack_max_algae) :: &
735 64 : chlabs , & ! chla absorption 1/m/(mg/m^3)
736 64 : alpha2max_low , & ! light limitation (1/(W/m^2))
737 64 : beta2max , & ! light inhibition (1/(W/m^2))
738 64 : mu_max , & ! maximum growth rate (1/d)
739 64 : grow_Tdep , & ! T dependence of growth (1/C)
740 64 : fr_graze , & ! fraction of algae grazed
741 64 : mort_pre , & ! mortality (1/day)
742 64 : mort_Tdep , & ! T dependence of mortality (1/C)
743 64 : k_exude , & ! algal carbon exudation rate (1/d)
744 64 : K_Nit , & ! nitrate half saturation (mmol/m^3)
745 64 : K_Am , & ! ammonium half saturation (mmol/m^3)
746 64 : K_Sil , & ! silicon half saturation (mmol/m^3)
747 64 : K_Fe ! iron half saturation or micromol/m^3
748 :
749 : real (kind=dbl_kind), dimension(icepack_max_DON) :: &
750 32 : f_don , & ! fraction of spilled grazing to DON
751 32 : kn_bac , & ! Bacterial degredation of DON (1/d)
752 32 : f_don_Am ! fraction of remineralized DON to Am
753 :
754 : real (kind=dbl_kind), dimension(icepack_max_DOC) :: &
755 64 : f_exude , & ! fraction of exuded carbon to each DOC pool
756 64 : k_bac ! Bacterial degredation of DOC (1/d)
757 :
758 : real (kind=dbl_kind), dimension(icepack_max_nbtrcr) :: &
759 480 : zbgc_frac_init,&! initializes mobile fraction
760 480 : bgc_tracer_type ! described tracer in mobile or stationary phases
761 : ! < 0 is purely mobile (eg. nitrate)
762 : ! > 0 has timescales for transitions between
763 : ! phases based on whether the ice is melting or growing
764 :
765 : real (kind=dbl_kind), dimension(icepack_max_nbtrcr) :: &
766 480 : zbgc_init_frac, & ! fraction of ocean tracer concentration in new ice
767 480 : tau_ret, & ! retention timescale (s), mobile to stationary phase
768 480 : tau_rel ! release timescale (s), stationary to mobile phase
769 :
770 : character (32) :: &
771 : nml_filename = 'icepack_in' ! namelist input file name
772 :
773 : integer (kind=int_kind) :: &
774 : nml_error, & ! namelist i/o error flag
775 : k, mm , & ! loop index
776 : ntd , & ! for tracer dependency calculation
777 : nk , & !
778 : nt_depend
779 :
780 : character(len=*), parameter :: subname='(init_zbgc)'
781 :
782 : !------------------------------------------------------------
783 : ! Tracers have mobile and stationary phases.
784 : ! ice growth allows for retention, ice melt facilitates mobility
785 : ! bgc_tracer_type defines the exchange timescales between these phases
786 : ! -1 : entirely in the mobile phase, no exchange (this is the default)
787 : ! 0 : retention time scale is tau_min, release time scale is tau_max
788 : ! 1 : retention time scale is tau_max, release time scale is tau_min
789 : ! 0.5: retention time scale is tau_min, release time scale is tau_min
790 : ! 2 : retention time scale is tau_max, release time scale is tau_max
791 : !------------------------------------------------------------
792 :
793 : !-----------------------------------------------------------------
794 : ! namelist variables
795 : !-----------------------------------------------------------------
796 :
797 : namelist /zbgc_nml/ &
798 : tr_brine, tr_zaero, modal_aero, skl_bgc, &
799 : z_tracers, dEdd_algae, solve_zbgc, bgc_flux_type, &
800 : restore_bgc, scale_bgc, solve_zsal, bgc_data_type, &
801 : tr_bgc_Nit, tr_bgc_C, tr_bgc_chl, tr_bgc_Am, tr_bgc_Sil, &
802 : tr_bgc_DMS, tr_bgc_PON, tr_bgc_hum, tr_bgc_DON, tr_bgc_Fe, &
803 : grid_o, l_sk, grid_oS, &
804 : l_skS, phi_snow, initbio_frac, frazil_scav, &
805 : ratio_Si2N_diatoms , ratio_Si2N_sp , ratio_Si2N_phaeo , &
806 : ratio_S2N_diatoms , ratio_S2N_sp , ratio_S2N_phaeo , &
807 : ratio_Fe2C_diatoms , ratio_Fe2C_sp , ratio_Fe2C_phaeo , &
808 : ratio_Fe2N_diatoms , ratio_Fe2N_sp , ratio_Fe2N_phaeo , &
809 : ratio_Fe2DON , ratio_Fe2DOC_s , ratio_Fe2DOC_l , &
810 : fr_resp , tau_min , tau_max , &
811 : algal_vel , R_dFe2dust , dustFe_sol , &
812 : chlabs_diatoms , chlabs_sp , chlabs_phaeo , &
813 : alpha2max_low_diatoms,alpha2max_low_sp , alpha2max_low_phaeo, &
814 : beta2max_diatoms , beta2max_sp , beta2max_phaeo , &
815 : mu_max_diatoms , mu_max_sp , mu_max_phaeo , &
816 : grow_Tdep_diatoms , grow_Tdep_sp , grow_Tdep_phaeo , &
817 : fr_graze_diatoms , fr_graze_sp , fr_graze_phaeo , &
818 : mort_pre_diatoms , mort_pre_sp , mort_pre_phaeo , &
819 : mort_Tdep_diatoms , mort_Tdep_sp , mort_Tdep_phaeo , &
820 : k_exude_diatoms , k_exude_sp , k_exude_phaeo , &
821 : K_Nit_diatoms , K_Nit_sp , K_Nit_phaeo , &
822 : K_Am_diatoms , K_Am_sp , K_Am_phaeo , &
823 : K_Sil_diatoms , K_Sil_sp , K_Sil_phaeo , &
824 : K_Fe_diatoms , K_Fe_sp , K_Fe_phaeo , &
825 : f_don_protein , kn_bac_protein , f_don_Am_protein , &
826 : f_doc_s , f_doc_l , f_exude_s , &
827 : f_exude_l , k_bac_s , k_bac_l , &
828 : T_max , fsal , op_dep_min , &
829 : fr_graze_s , fr_graze_e , fr_mort2min , &
830 : fr_dFe , k_nitrif , t_iron_conv , &
831 : max_loss , max_dfe_doc1 , fr_resp_s , &
832 : y_sk_DMS , t_sk_conv , t_sk_ox , &
833 : algaltype_diatoms , algaltype_sp , algaltype_phaeo , &
834 : nitratetype , ammoniumtype , silicatetype , &
835 : dmspptype , dmspdtype , humtype , &
836 : doctype_s , doctype_l , dontype_protein , &
837 : fedtype_1 , feptype_1 , zaerotype_bc1 , &
838 : zaerotype_bc2 , zaerotype_dust1 , zaerotype_dust2 , &
839 : zaerotype_dust3 , zaerotype_dust4 , ratio_C2N_diatoms , &
840 : ratio_C2N_sp , ratio_C2N_phaeo , ratio_chl2N_diatoms, &
841 : ratio_chl2N_sp , ratio_chl2N_phaeo , F_abs_chl_diatoms , &
842 : F_abs_chl_sp , F_abs_chl_phaeo , ratio_C2N_proteins
843 :
844 : !-----------------------------------------------------------------
845 : ! query Icepack values
846 : !-----------------------------------------------------------------
847 :
848 45 : call icepack_query_tracer_sizes(ntrcr_out=ntrcr)
849 45 : call icepack_query_tracer_flags(tr_aero_out=tr_aero)
850 : call icepack_query_parameters(ktherm_out=ktherm, shortwave_out=shortwave, &
851 : scale_bgc_out=scale_bgc, skl_bgc_out=skl_bgc, z_tracers_out=z_tracers, &
852 : dEdd_algae_out=dEdd_algae, solve_zbgc_out=solve_zbgc, phi_snow_out=phi_snow, &
853 : bgc_flux_type_out=bgc_flux_type, grid_o_out=grid_o, l_sk_out=l_sk, &
854 : initbio_frac_out=initbio_frac, frazil_scav_out=frazil_scav, &
855 : algal_vel_out=algal_vel, R_dFe2dust_out=R_dFe2dust, &
856 : dustFe_sol_out=dustFe_sol, T_max_out=T_max, fsal_out=fsal, &
857 : op_dep_min_out=op_dep_min, fr_graze_s_out=fr_graze_s, &
858 : fr_graze_e_out=fr_graze_e, fr_mort2min_out=fr_mort2min, &
859 : fr_dFe_out=fr_dFe, k_nitrif_out=k_nitrif, t_iron_conv_out=t_iron_conv, &
860 : max_loss_out=max_loss, max_dfe_doc1_out=max_dfe_doc1, &
861 : fr_resp_s_out=fr_resp_s, y_sk_DMS_out=y_sk_DMS, t_sk_conv_out=t_sk_conv, &
862 45 : t_sk_ox_out=t_sk_ox)
863 45 : call icepack_warnings_flush(nu_diag)
864 45 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
865 0 : file=__FILE__, line=__LINE__)
866 :
867 : !-----------------------------------------------------------------
868 : ! default values
869 : !-----------------------------------------------------------------
870 45 : tr_brine = .false. ! brine height differs from ice height
871 45 : tr_zaero = .false. ! z aerosol tracers
872 45 : modal_aero = .false. ! use modal aerosol treatment of aerosols
873 45 : restore_bgc = .false. ! restore bgc if true
874 45 : solve_zsal = .false. ! update salinity tracer profile from solve_S_dt
875 45 : bgc_data_type = 'default'! source of bgc data
876 45 : tr_bgc_PON = .false. !---------------------------------------------
877 45 : tr_bgc_Nit = .false. ! biogeochemistry (skl or zbgc)
878 45 : tr_bgc_C = .false. ! if skl_bgc = .true. then skl
879 45 : tr_bgc_chl = .false. ! if z_tracers = .true. then vertically resolved
880 45 : tr_bgc_Sil = .false. ! if z_tracers + solve_zbgc = .true. then
881 45 : tr_bgc_Am = .false. ! vertically resolved with reactions
882 45 : tr_bgc_DMS = .false. !------------------------------------------------
883 45 : tr_bgc_DON = .false. !
884 45 : tr_bgc_hum = .false. !
885 45 : tr_bgc_Fe = .false. !
886 45 : tr_bgc_N = .true. !
887 :
888 : ! z biology parameters
889 45 : ratio_Si2N_diatoms = 1.8_dbl_kind ! algal Si to N (mol/mol)
890 45 : ratio_Si2N_sp = c0 ! diatoms, small plankton, phaeocystis
891 45 : ratio_Si2N_phaeo = c0
892 45 : ratio_S2N_diatoms = 0.03_dbl_kind ! algal S to N (mol/mol)
893 45 : ratio_S2N_sp = 0.03_dbl_kind
894 45 : ratio_S2N_phaeo = 0.03_dbl_kind
895 45 : ratio_Fe2C_diatoms = 0.0033_dbl_kind ! algal Fe to C (umol/mol)
896 45 : ratio_Fe2C_sp = 0.0033_dbl_kind
897 45 : ratio_Fe2C_phaeo = p1
898 45 : ratio_Fe2N_diatoms = 0.023_dbl_kind ! algal Fe to N (umol/mol)
899 45 : ratio_Fe2N_sp = 0.023_dbl_kind
900 45 : ratio_Fe2N_phaeo = 0.7_dbl_kind
901 45 : ratio_Fe2DON = 0.023_dbl_kind ! Fe to N of DON (nmol/umol)
902 45 : ratio_Fe2DOC_s = p1 ! Fe to C of DOC (nmol/umol) saccharids
903 45 : ratio_Fe2DOC_l = 0.033_dbl_kind ! Fe to C of DOC (nmol/umol) lipids
904 45 : tau_min = 5200.0_dbl_kind ! rapid mobile to stationary exchanges (s)
905 45 : tau_max = 1.73e5_dbl_kind ! long time mobile to stationary exchanges (s)
906 45 : chlabs_diatoms = 0.03_dbl_kind ! chl absorption (1/m/(mg/m^3))
907 45 : chlabs_sp = 0.01_dbl_kind
908 45 : chlabs_phaeo = 0.05_dbl_kind
909 45 : alpha2max_low_diatoms = 0.8_dbl_kind ! light limitation (1/(W/m^2))
910 45 : alpha2max_low_sp = 0.67_dbl_kind
911 45 : alpha2max_low_phaeo = 0.67_dbl_kind
912 45 : beta2max_diatoms = 0.018_dbl_kind ! light inhibition (1/(W/m^2))
913 45 : beta2max_sp = 0.0025_dbl_kind
914 45 : beta2max_phaeo = 0.01_dbl_kind
915 45 : mu_max_diatoms = 1.2_dbl_kind ! maximum growth rate (1/day)
916 45 : mu_max_sp = 0.851_dbl_kind
917 45 : mu_max_phaeo = 0.851_dbl_kind
918 45 : grow_Tdep_diatoms = 0.06_dbl_kind ! Temperature dependence of growth (1/C)
919 45 : grow_Tdep_sp = 0.06_dbl_kind
920 45 : grow_Tdep_phaeo = 0.06_dbl_kind
921 45 : fr_graze_diatoms = 0.01_dbl_kind ! Fraction grazed
922 45 : fr_graze_sp = p1
923 45 : fr_graze_phaeo = p1
924 45 : mort_pre_diatoms = 0.007_dbl_kind! Mortality (1/day)
925 45 : mort_pre_sp = 0.007_dbl_kind
926 45 : mort_pre_phaeo = 0.007_dbl_kind
927 45 : mort_Tdep_diatoms = 0.03_dbl_kind ! T dependence of mortality (1/C)
928 45 : mort_Tdep_sp = 0.03_dbl_kind
929 45 : mort_Tdep_phaeo = 0.03_dbl_kind
930 45 : k_exude_diatoms = c0 ! algal exudation (1/d)
931 45 : k_exude_sp = c0
932 45 : k_exude_phaeo = c0
933 45 : K_Nit_diatoms = c1 ! nitrate half saturation (mmol/m^3)
934 45 : K_Nit_sp = c1
935 45 : K_Nit_phaeo = c1
936 45 : K_Am_diatoms = 0.3_dbl_kind ! ammonium half saturation (mmol/m^3)
937 45 : K_Am_sp = 0.3_dbl_kind
938 45 : K_Am_phaeo = 0.3_dbl_kind
939 45 : K_Sil_diatoms = 4.0_dbl_kind ! silicate half saturation (mmol/m^3)
940 45 : K_Sil_sp = c0
941 45 : K_Sil_phaeo = c0
942 45 : K_Fe_diatoms = c1 ! iron half saturation (nM)
943 45 : K_Fe_sp = 0.2_dbl_kind
944 45 : K_Fe_phaeo = p1
945 45 : f_don_protein = 0.6_dbl_kind ! fraction of spilled grazing to proteins
946 45 : kn_bac_protein = 0.03_dbl_kind ! Bacterial degredation of DON (1/d)
947 45 : f_don_Am_protein = 0.25_dbl_kind ! fraction of remineralized DON to ammonium
948 45 : f_doc_s = 0.4_dbl_kind ! fraction of mortality to DOC
949 45 : f_doc_l = 0.4_dbl_kind
950 45 : f_exude_s = c1 ! fraction of exudation to DOC
951 45 : f_exude_l = c1
952 45 : k_bac_s = 0.03_dbl_kind ! Bacterial degredation of DOC (1/d)
953 45 : k_bac_l = 0.03_dbl_kind
954 45 : algaltype_diatoms = c0 ! ------------------
955 45 : algaltype_sp = p5 !
956 45 : algaltype_phaeo = p5 !
957 45 : nitratetype = -c1 ! mobility type between
958 45 : ammoniumtype = c1 ! stationary <--> mobile
959 45 : silicatetype = -c1 !
960 45 : dmspptype = p5 !
961 45 : dmspdtype = -c1 !
962 45 : humtype = c1 !
963 45 : doctype_s = p5 !
964 45 : doctype_l = p5 !
965 45 : dontype_protein = p5 !
966 45 : fedtype_1 = p5 !
967 45 : feptype_1 = p5 !
968 45 : zaerotype_bc1 = c1 !
969 45 : zaerotype_bc2 = c1 !
970 45 : zaerotype_dust1 = c1 !
971 45 : zaerotype_dust2 = c1 !
972 45 : zaerotype_dust3 = c1 !
973 45 : zaerotype_dust4 = c1 !--------------------
974 45 : ratio_C2N_diatoms = 7.0_dbl_kind ! algal C to N ratio (mol/mol)
975 45 : ratio_C2N_sp = 7.0_dbl_kind
976 45 : ratio_C2N_phaeo = 7.0_dbl_kind
977 45 : ratio_chl2N_diatoms= 2.1_dbl_kind ! algal chlorophyll to N ratio (mg/mmol)
978 45 : ratio_chl2N_sp = 1.1_dbl_kind
979 45 : ratio_chl2N_phaeo = 0.84_dbl_kind
980 45 : F_abs_chl_diatoms = 2.0_dbl_kind ! scales absorbed radiation for dEdd
981 45 : F_abs_chl_sp = 4.0_dbl_kind
982 45 : F_abs_chl_phaeo = 5.0_dbl_kind
983 45 : ratio_C2N_proteins = 7.0_dbl_kind ! ratio of C to N in proteins (mol/mol)
984 :
985 : !-----------------------------------------------------------------
986 : ! read from input file
987 : !-----------------------------------------------------------------
988 :
989 45 : open (nu_nml, file=trim(nml_filename), status='old',iostat=nml_error)
990 45 : if (nml_error /= 0) then
991 0 : nml_error = -1
992 : else
993 45 : nml_error = 1
994 : endif
995 :
996 45 : print*,'Reading zbgc_nml'
997 90 : do while (nml_error > 0)
998 90 : read(nu_nml, nml=zbgc_nml,iostat=nml_error)
999 : end do
1000 45 : if (nml_error == 0) close(nu_nml)
1001 45 : if (nml_error /= 0) then
1002 0 : print*,'error reading zbgc namelist'
1003 0 : call icedrv_system_abort(file=__FILE__,line=__LINE__)
1004 : endif
1005 :
1006 : !-----------------------------------------------------------------
1007 : ! resolve conflicts
1008 : !-----------------------------------------------------------------
1009 : ! zsalinity and brine
1010 : !-----------------------------------------------------------------
1011 45 : if (solve_zsal .and. TRZS == 0) then
1012 0 : write(nu_diag,*) 'WARNING: solve_zsal=T but 0 zsalinity tracers'
1013 0 : write(nu_diag,*) 'WARNING: setting solve_zsal = F'
1014 0 : solve_zsal = .false.
1015 : elseif (solve_zsal .and. nblyr < 1) then
1016 : write(nu_diag,*) 'WARNING: solve_zsal=T but 0 zsalinity tracers'
1017 : write(nu_diag,*) 'WARNING: setting solve_zsal = F'
1018 : solve_zsal = .false.
1019 : endif
1020 :
1021 45 : if (solve_zsal .and. ((.not. tr_brine) .or. (ktherm /= 1))) then
1022 0 : write(nu_diag,*) 'WARNING: solve_zsal needs tr_brine=T and ktherm=1'
1023 0 : write(nu_diag,*) 'WARNING: setting tr_brine=T and ktherm=1'
1024 0 : tr_brine = .true.
1025 0 : ktherm = 1
1026 : endif
1027 :
1028 40 : if (tr_brine .and. TRBRI == 0 ) then
1029 : write(nu_diag,*) &
1030 0 : 'WARNING: tr_brine=T but no brine height compiled'
1031 : write(nu_diag,*) &
1032 0 : 'WARNING: setting solve_zsal and tr_brine = F'
1033 0 : solve_zsal = .false.
1034 0 : tr_brine = .false.
1035 : elseif (tr_brine .and. nblyr < 1 ) then
1036 : write(nu_diag,*) &
1037 : 'WARNING: tr_brine=T but no biology layers compiled'
1038 : write(nu_diag,*) &
1039 : 'WARNING: setting solve_zsal and tr_brine = F'
1040 : solve_zsal = .false.
1041 : tr_brine = .false.
1042 : endif
1043 :
1044 45 : write(nu_diag,1010) ' tr_brine = ', tr_brine
1045 45 : if (tr_brine) then
1046 5 : write(nu_diag,1005) ' phi_snow = ', phi_snow
1047 : endif
1048 45 : if (solve_zsal) then
1049 0 : write(nu_diag,1010) ' solve_zsal = ', solve_zsal
1050 0 : write(nu_diag,1000) ' grid_oS = ', grid_oS
1051 0 : write(nu_diag,1005) ' l_skS = ', l_skS
1052 : endif
1053 :
1054 : !-----------------------------------------------------------------
1055 : ! biogeochemistry
1056 : !-----------------------------------------------------------------
1057 :
1058 45 : if (.not. tr_brine) then
1059 40 : if (solve_zbgc) then
1060 0 : write(nu_diag,*) 'WARNING: tr_brine = F and solve_zbgc = T'
1061 0 : write(nu_diag,*) 'WARNING: setting solve_zbgc = F'
1062 0 : solve_zbgc = .false.
1063 : endif
1064 40 : if (skl_bgc) then
1065 0 : write(nu_diag,*) 'WARNING: tr_brine = F and skl_bgc = T'
1066 0 : write(nu_diag,*) 'WARNING: setting skl_bgc = F'
1067 0 : skl_bgc = .false.
1068 : endif
1069 40 : if (tr_zaero) then
1070 0 : write(nu_diag,*) 'WARNING: tr_brine = F and tr_zaero = T'
1071 0 : write(nu_diag,*) 'WARNING: setting tr_zaero = F'
1072 0 : tr_zaero = .false.
1073 : endif
1074 : endif
1075 :
1076 45 : if ((skl_bgc .AND. solve_zbgc) .or. (skl_bgc .AND. z_tracers)) then
1077 0 : print*, 'ERROR: skl_bgc and (solve_zbgc or z_tracers) are both true'
1078 0 : call icedrv_system_abort(file=__FILE__,line=__LINE__)
1079 : endif
1080 :
1081 45 : if (skl_bgc .AND. tr_zaero) then
1082 0 : write(nu_diag,*) 'WARNING: skl bgc does not use vertical tracers'
1083 0 : write(nu_diag,*) 'WARNING: setting tr_zaero = F'
1084 0 : tr_zaero = .false.
1085 : endif
1086 :
1087 45 : if (dEdd_algae .AND. trim(shortwave) /= 'dEdd') then
1088 0 : write(nu_diag,*) 'WARNING: dEdd_algae = T but shortwave /= dEdd'
1089 0 : write(nu_diag,*) 'WARNING: setting dEdd_algae = F'
1090 0 : dEdd_algae = .false.
1091 : endif
1092 :
1093 45 : if (dEdd_algae .AND. (.NOT. tr_bgc_N) .AND. (.NOT. tr_zaero)) then
1094 0 : write(nu_diag,*) 'WARNING: need tr_bgc_N or tr_zaero for dEdd_algae'
1095 0 : write(nu_diag,*) 'WARNING: setting dEdd_algae = F'
1096 0 : dEdd_algae = .false.
1097 : endif
1098 :
1099 45 : if (modal_aero .AND. (.NOT. tr_zaero) .AND. (.NOT. tr_aero)) then
1100 0 : modal_aero = .false.
1101 : endif
1102 :
1103 45 : if (modal_aero .AND. trim(shortwave) /= 'dEdd') then
1104 0 : write(nu_diag,*) 'WARNING: modal_aero = T but shortwave /= dEdd'
1105 0 : write(nu_diag,*) 'WARNING: setting modal_aero = F'
1106 0 : modal_aero = .false.
1107 : endif
1108 : if (n_algae > icepack_max_algae) then
1109 : print*, 'error:number of algal types exceeds icepack_max_algae'
1110 : call icedrv_system_abort(file=__FILE__,line=__LINE__)
1111 : endif
1112 : if (n_doc > icepack_max_doc) then
1113 : print*, 'error:number of algal types exceeds icepack_max_doc'
1114 : call icedrv_system_abort(file=__FILE__,line=__LINE__)
1115 : endif
1116 : if (n_dic > icepack_max_dic) then
1117 : print*, 'error:number of dic types exceeds icepack_max_dic'
1118 : call icedrv_system_abort(file=__FILE__,line=__LINE__)
1119 : endif
1120 : if (n_don > icepack_max_don) then
1121 : print*, 'error:number of don types exceeds icepack_max_don'
1122 : call icedrv_system_abort(file=__FILE__,line=__LINE__)
1123 : endif
1124 : if (n_fed > icepack_max_fe) then
1125 : print*, 'error:number of dissolved fe types exceeds icepack_max_fe'
1126 : call icedrv_system_abort(file=__FILE__,line=__LINE__)
1127 : endif
1128 : if (n_fep > icepack_max_fe) then
1129 : print*, 'error:number of particulate fe types exceeds icepack_max_fe'
1130 : call icedrv_system_abort(file=__FILE__,line=__LINE__)
1131 : endif
1132 :
1133 44 : if ((TRBGCS == 0 .and. skl_bgc) .or. (TRALG == 0 .and. skl_bgc)) then
1134 : write(nu_diag,*) &
1135 0 : 'WARNING: skl_bgc=T but 0 bgc or algal tracers compiled'
1136 : write(nu_diag,*) &
1137 0 : 'WARNING: setting skl_bgc = F'
1138 0 : skl_bgc = .false.
1139 : endif
1140 :
1141 41 : if ((TRBGCZ == 0 .and. solve_zbgc) .or. (TRALG == 0 .and. solve_zbgc)) then
1142 : write(nu_diag,*) &
1143 0 : 'WARNING: solve_zbgc=T but 0 zbgc or algal tracers compiled'
1144 : write(nu_diag,*) &
1145 0 : 'WARNING: setting solve_zbgc = F'
1146 0 : solve_zbgc = .false.
1147 : endif
1148 :
1149 45 : if (solve_zbgc .and. .not. z_tracers) z_tracers = .true.
1150 45 : if (skl_bgc .or. solve_zbgc) then
1151 5 : tr_bgc_N = .true. ! minimum NP biogeochemistry
1152 5 : tr_bgc_Nit = .true.
1153 : else
1154 40 : tr_bgc_N = .false.
1155 40 : tr_bgc_C = .false.
1156 40 : tr_bgc_chl = .false.
1157 40 : tr_bgc_Nit = .false.
1158 40 : tr_bgc_Am = .false.
1159 40 : tr_bgc_Sil = .false.
1160 40 : tr_bgc_hum = .false.
1161 40 : tr_bgc_DMS = .false.
1162 40 : tr_bgc_PON = .false.
1163 40 : tr_bgc_DON = .false.
1164 40 : tr_bgc_Fe = .false.
1165 : endif
1166 :
1167 : !-----------------------------------------------------------------
1168 : ! z layer aerosols
1169 : !-----------------------------------------------------------------
1170 45 : if (tr_zaero .and. .not. z_tracers) z_tracers = .true.
1171 :
1172 : if (n_zaero > icepack_max_aero) then
1173 : print*, 'error:number of z aerosols exceeds icepack_max_aero'
1174 : call icedrv_system_abort(file=__FILE__,line=__LINE__)
1175 : endif
1176 :
1177 : if (skl_bgc .and. n_bgc < 2) then
1178 : write (nu_diag,*) ' '
1179 : write (nu_diag,*) 'comp_ice must have number of bgc tracers >= 2'
1180 : write (nu_diag,*) 'number of bgc tracers compiled:',n_bgc
1181 : call icedrv_system_abort(file=__FILE__,line=__LINE__)
1182 : endif
1183 :
1184 : if (solve_zbgc .and. n_bgc < 2) then
1185 : write (nu_diag,*) ' '
1186 : write (nu_diag,*) 'comp_ice must have number of zbgc tracers >= 2'
1187 : write (nu_diag,*) 'number of bgc tracers compiled:',n_bgc
1188 : call icedrv_system_abort(file=__FILE__,line=__LINE__)
1189 : endif
1190 :
1191 42 : if (tr_zaero .and. TRZAERO < 1) then
1192 0 : write (nu_diag,*) ' '
1193 0 : write (nu_diag,*) 'comp_ice must have number of TRZAERO > 0'
1194 0 : write (nu_diag,*) 'in order to solve z aerosols:',TRZAERO
1195 0 : call icedrv_system_abort(file=__FILE__,line=__LINE__)
1196 : endif
1197 :
1198 : !-----------------------------------------------------------------
1199 : ! end conflict resolution
1200 : !-----------------------------------------------------------------
1201 : !-----------------------------------------------------------------
1202 : ! set Icepack values
1203 : !-----------------------------------------------------------------
1204 :
1205 : call icepack_init_parameters(ktherm_in=ktherm, shortwave_in=shortwave, &
1206 : scale_bgc_in=scale_bgc, skl_bgc_in=skl_bgc, z_tracers_in=z_tracers, &
1207 : dEdd_algae_in=dEdd_algae, solve_zbgc_in=solve_zbgc, &
1208 : bgc_flux_type_in=bgc_flux_type, grid_o_in=grid_o, l_sk_in=l_sk, &
1209 : initbio_frac_in=initbio_frac, frazil_scav_in=frazil_scav, &
1210 : grid_oS_in=grid_oS, l_skS_in=l_skS, phi_snow_in=phi_snow, &
1211 : algal_vel_in=algal_vel, R_dFe2dust_in=R_dFe2dust, &
1212 : dustFe_sol_in=dustFe_sol, T_max_in=T_max, fsal_in=fsal, &
1213 : op_dep_min_in=op_dep_min, fr_graze_s_in=fr_graze_s, &
1214 : fr_graze_e_in=fr_graze_e, fr_mort2min_in=fr_mort2min, &
1215 : fr_dFe_in=fr_dFe, k_nitrif_in=k_nitrif, t_iron_conv_in=t_iron_conv, &
1216 : max_loss_in=max_loss, max_dfe_doc1_in=max_dfe_doc1, fr_resp_in=fr_resp, &
1217 : fr_resp_s_in=fr_resp_s, y_sk_DMS_in=y_sk_DMS, t_sk_conv_in=t_sk_conv, &
1218 45 : t_sk_ox_in=t_sk_ox, modal_aero_in=modal_aero)
1219 45 : call icepack_warnings_flush(nu_diag)
1220 45 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
1221 0 : file=__FILE__, line=__LINE__)
1222 :
1223 : !-----------------------------------------------------------------
1224 : ! initialize zbgc tracer indices
1225 : !-----------------------------------------------------------------
1226 :
1227 45 : ntrcr_o = ntrcr
1228 45 : nt_fbri = 0
1229 45 : if (tr_brine) then
1230 5 : nt_fbri = ntrcr + 1 ! ice volume fraction with salt
1231 5 : ntrcr = ntrcr + 1
1232 5 : trcr_depend(nt_fbri) = 1 ! volume-weighted
1233 5 : trcr_base (nt_fbri,1) = c0 ! volume-weighted
1234 5 : trcr_base (nt_fbri,2) = c1 ! volume-weighted
1235 5 : trcr_base (nt_fbri,3) = c0 ! volume-weighted
1236 5 : n_trcr_strata(nt_fbri) = 0
1237 5 : nt_strata (nt_fbri,1) = 0
1238 5 : nt_strata (nt_fbri,2) = 0
1239 : endif
1240 :
1241 45 : ntd = 0 ! if nt_fbri /= 0 then use fbri dependency
1242 45 : if (nt_fbri == 0) ntd = -1 ! otherwise make tracers depend on ice volume
1243 :
1244 45 : nt_bgc_S = 0
1245 45 : if (solve_zsal) then ! .true. only if tr_brine = .true.
1246 0 : nt_bgc_S = ntrcr + 1
1247 0 : ntrcr = ntrcr + nblyr
1248 0 : do k = 1,nblyr
1249 0 : trcr_depend(nt_bgc_S + k - 1) = 2 + nt_fbri + ntd
1250 0 : trcr_base (nt_bgc_S,1) = c0 ! default: ice area
1251 0 : trcr_base (nt_bgc_S,2) = c1
1252 0 : trcr_base (nt_bgc_S,3) = c0
1253 0 : n_trcr_strata(nt_bgc_S) = 1
1254 0 : nt_strata(nt_bgc_S,1) = nt_fbri
1255 0 : nt_strata(nt_bgc_S,2) = 0
1256 : enddo
1257 : endif
1258 :
1259 : !-----------------------------------------------------------------
1260 : ! biogeochemistry
1261 : !-----------------------------------------------------------------
1262 :
1263 45 : nbtrcr = 0
1264 45 : nbtrcr_sw = 0
1265 :
1266 : ! vectors of size icepack_max_algae
1267 45 : nlt_bgc_N(:) = 0
1268 45 : nlt_bgc_C(:) = 0
1269 45 : nlt_bgc_chl(:) = 0
1270 45 : nt_bgc_N(:) = 0
1271 45 : nt_bgc_C(:) = 0
1272 45 : nt_bgc_chl(:) = 0
1273 :
1274 : ! vectors of size icepack_max_dic
1275 45 : nlt_bgc_DIC(:) = 0
1276 45 : nt_bgc_DIC(:) = 0
1277 :
1278 : ! vectors of size icepack_max_doc
1279 45 : nlt_bgc_DOC(:) = 0
1280 45 : nt_bgc_DOC(:) = 0
1281 :
1282 : ! vectors of size icepack_max_don
1283 45 : nlt_bgc_DON(:) = 0
1284 45 : nt_bgc_DON(:) = 0
1285 :
1286 : ! vectors of size icepack_max_fe
1287 45 : nlt_bgc_Fed(:) = 0
1288 45 : nlt_bgc_Fep(:) = 0
1289 45 : nt_bgc_Fed(:) = 0
1290 45 : nt_bgc_Fep(:) = 0
1291 :
1292 : ! vectors of size icepack_max_aero
1293 45 : nlt_zaero(:) = 0
1294 45 : nlt_zaero_sw(:) = 0
1295 45 : nt_zaero(:) = 0
1296 45 : nt_zaero_sw(:) = 0
1297 :
1298 45 : nlt_bgc_Nit = 0
1299 45 : nlt_bgc_Am = 0
1300 45 : nlt_bgc_Sil = 0
1301 45 : nlt_bgc_DMSPp = 0
1302 45 : nlt_bgc_DMSPd = 0
1303 45 : nlt_bgc_DMS = 0
1304 45 : nlt_bgc_PON = 0
1305 45 : nlt_bgc_hum = 0
1306 45 : nlt_chl_sw = 0
1307 45 : bio_index(:) = 0
1308 45 : bio_index_o(:) = 0
1309 :
1310 45 : nt_bgc_Nit = 0
1311 45 : nt_bgc_Am = 0
1312 45 : nt_bgc_Sil = 0
1313 45 : nt_bgc_DMSPp = 0
1314 45 : nt_bgc_DMSPd = 0
1315 45 : nt_bgc_DMS = 0
1316 45 : nt_bgc_PON = 0
1317 45 : nt_bgc_hum = 0
1318 45 : nt_zbgc_frac = 0
1319 :
1320 : !-----------------------------------------------------------------
1321 : ! Define array parameters
1322 : !-----------------------------------------------------------------
1323 :
1324 45 : R_Si2N(1) = ratio_Si2N_diatoms
1325 45 : R_Si2N(2) = ratio_Si2N_sp
1326 45 : R_Si2N(3) = ratio_Si2N_phaeo
1327 :
1328 45 : R_S2N(1) = ratio_S2N_diatoms
1329 45 : R_S2N(2) = ratio_S2N_sp
1330 45 : R_S2N(3) = ratio_S2N_phaeo
1331 :
1332 45 : R_Fe2C(1) = ratio_Fe2C_diatoms
1333 45 : R_Fe2C(2) = ratio_Fe2C_sp
1334 45 : R_Fe2C(3) = ratio_Fe2C_phaeo
1335 :
1336 45 : R_Fe2N(1) = ratio_Fe2N_diatoms
1337 45 : R_Fe2N(2) = ratio_Fe2N_sp
1338 45 : R_Fe2N(3) = ratio_Fe2N_phaeo
1339 :
1340 45 : R_C2N(1) = ratio_C2N_diatoms
1341 45 : R_C2N(2) = ratio_C2N_sp
1342 45 : R_C2N(3) = ratio_C2N_phaeo
1343 :
1344 45 : R_chl2N(1) = ratio_chl2N_diatoms
1345 45 : R_chl2N(2) = ratio_chl2N_sp
1346 45 : R_chl2N(3) = ratio_chl2N_phaeo
1347 :
1348 45 : F_abs_chl(1) = F_abs_chl_diatoms
1349 45 : F_abs_chl(2) = F_abs_chl_sp
1350 45 : F_abs_chl(3) = F_abs_chl_phaeo
1351 :
1352 45 : R_Fe2DON(1) = ratio_Fe2DON
1353 45 : R_C2N_DON(1) = ratio_C2N_proteins
1354 :
1355 45 : R_Fe2DOC(1) = ratio_Fe2DOC_s
1356 45 : R_Fe2DOC(2) = ratio_Fe2DOC_l
1357 45 : R_Fe2DOC(3) = c0
1358 :
1359 45 : chlabs(1) = chlabs_diatoms
1360 45 : chlabs(2) = chlabs_sp
1361 45 : chlabs(3) = chlabs_phaeo
1362 :
1363 45 : alpha2max_low(1) = alpha2max_low_diatoms
1364 45 : alpha2max_low(2) = alpha2max_low_sp
1365 45 : alpha2max_low(3) = alpha2max_low_phaeo
1366 :
1367 45 : beta2max(1) = beta2max_diatoms
1368 45 : beta2max(2) = beta2max_sp
1369 45 : beta2max(3) = beta2max_phaeo
1370 :
1371 45 : mu_max(1) = mu_max_diatoms
1372 45 : mu_max(2) = mu_max_sp
1373 45 : mu_max(3) = mu_max_phaeo
1374 :
1375 45 : grow_Tdep(1) = grow_Tdep_diatoms
1376 45 : grow_Tdep(2) = grow_Tdep_sp
1377 45 : grow_Tdep(3) = grow_Tdep_phaeo
1378 :
1379 45 : fr_graze(1) = fr_graze_diatoms
1380 45 : fr_graze(2) = fr_graze_sp
1381 45 : fr_graze(3) = fr_graze_phaeo
1382 :
1383 45 : mort_pre(1) = mort_pre_diatoms
1384 45 : mort_pre(2) = mort_pre_sp
1385 45 : mort_pre(3) = mort_pre_phaeo
1386 :
1387 45 : mort_Tdep(1) = mort_Tdep_diatoms
1388 45 : mort_Tdep(2) = mort_Tdep_sp
1389 45 : mort_Tdep(3) = mort_Tdep_phaeo
1390 :
1391 45 : k_exude(1) = k_exude_diatoms
1392 45 : k_exude(2) = k_exude_sp
1393 45 : k_exude(3) = k_exude_phaeo
1394 :
1395 45 : K_Nit(1) = K_Nit_diatoms
1396 45 : K_Nit(2) = K_Nit_sp
1397 45 : K_Nit(3) = K_Nit_phaeo
1398 :
1399 45 : K_Am(1) = K_Am_diatoms
1400 45 : K_Am(2) = K_Am_sp
1401 45 : K_Am(3) = K_Am_phaeo
1402 :
1403 45 : K_Sil(1) = K_Sil_diatoms
1404 45 : K_Sil(2) = K_Sil_sp
1405 45 : K_Sil(3) = K_Sil_phaeo
1406 :
1407 45 : K_Fe(1) = K_Fe_diatoms
1408 45 : K_Fe(2) = K_Fe_sp
1409 45 : K_Fe(3) = K_Fe_phaeo
1410 :
1411 45 : f_don(1) = f_don_protein
1412 45 : kn_bac(1) = kn_bac_protein
1413 45 : f_don_Am(1) = f_don_Am_protein
1414 :
1415 45 : f_exude(1) = f_exude_s
1416 45 : f_exude(2) = f_exude_l
1417 45 : k_bac(1) = k_bac_s
1418 45 : k_bac(2) = k_bac_l
1419 :
1420 90 : dictype(:) = -c1
1421 :
1422 45 : algaltype(1) = algaltype_diatoms
1423 45 : algaltype(2) = algaltype_sp
1424 45 : algaltype(3) = algaltype_phaeo
1425 :
1426 45 : doctype(1) = doctype_s
1427 45 : doctype(2) = doctype_l
1428 :
1429 45 : dontype(1) = dontype_protein
1430 :
1431 45 : fedtype(1) = fedtype_1
1432 45 : feptype(1) = feptype_1
1433 :
1434 45 : zaerotype(1) = zaerotype_bc1
1435 45 : zaerotype(2) = zaerotype_bc2
1436 45 : zaerotype(3) = zaerotype_dust1
1437 45 : zaerotype(4) = zaerotype_dust2
1438 45 : zaerotype(5) = zaerotype_dust3
1439 45 : zaerotype(6) = zaerotype_dust4
1440 :
1441 : call icepack_init_zbgc ( &
1442 : R_Si2N_in=R_Si2N, &
1443 : R_S2N_in=R_S2N, R_Fe2C_in=R_Fe2C, R_Fe2N_in=R_Fe2N, R_C2N_in=R_C2N, &
1444 : R_chl2N_in=R_chl2N, F_abs_chl_in=F_abs_chl, R_Fe2DON_in=R_Fe2DON, &
1445 : R_C2N_DON_in=R_C2N_DON, &
1446 : R_Fe2DOC_in=R_Fe2DOC, &
1447 : chlabs_in=chlabs, alpha2max_low_in=alpha2max_low, beta2max_in=beta2max, &
1448 : mu_max_in=mu_max, grow_Tdep_in=grow_Tdep, fr_graze_in=fr_graze, &
1449 : mort_pre_in=mort_pre, &
1450 : mort_Tdep_in=mort_Tdep, k_exude_in=k_exude, &
1451 : K_Nit_in=K_Nit, K_Am_in=K_Am, K_sil_in=K_Sil, K_Fe_in=K_Fe, &
1452 : f_don_in=f_don, kn_bac_in=kn_bac, f_don_Am_in=f_don_Am, f_exude_in=f_exude, &
1453 : k_bac_in=k_bac, &
1454 : fr_resp_in=fr_resp, algal_vel_in=algal_vel, R_dFe2dust_in=R_dFe2dust, &
1455 : dustFe_sol_in=dustFe_sol, T_max_in=T_max, fr_mort2min_in=fr_mort2min, &
1456 : fr_dFe_in=fr_dFe, op_dep_min_in=op_dep_min, &
1457 : fr_graze_s_in=fr_graze_s, fr_graze_e_in=fr_graze_e, &
1458 : k_nitrif_in=k_nitrif, t_iron_conv_in=t_iron_conv, &
1459 : max_loss_in=max_loss, max_dfe_doc1_in=max_dfe_doc1, &
1460 : fr_resp_s_in=fr_resp_s, y_sk_DMS_in=y_sk_DMS, &
1461 45 : t_sk_conv_in=t_sk_conv, t_sk_ox_in=t_sk_ox, fsal_in=fsal)
1462 45 : call icepack_warnings_flush(nu_diag)
1463 45 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
1464 0 : file=__FILE__, line=__LINE__)
1465 :
1466 45 : if (skl_bgc) then
1467 :
1468 1 : nk = 1
1469 1 : nt_depend = 0
1470 :
1471 1 : if (dEdd_algae) then
1472 0 : nlt_chl_sw = 1
1473 0 : nbtrcr_sw = nilyr+nslyr+2 ! only the bottom layer will be nonzero
1474 : endif
1475 :
1476 44 : elseif (z_tracers) then ! defined on nblyr+1 in ice
1477 : ! and 2 snow layers (snow surface + interior)
1478 :
1479 4 : nk = nblyr + 1
1480 4 : nt_depend = 2 + nt_fbri + ntd
1481 :
1482 4 : if (tr_bgc_N) then
1483 4 : if (dEdd_algae) then
1484 0 : nlt_chl_sw = 1
1485 0 : nbtrcr_sw = nilyr+nslyr+2
1486 : endif
1487 : endif ! tr_bgc_N
1488 :
1489 : endif ! skl_bgc or z_tracers
1490 :
1491 45 : if (skl_bgc .or. z_tracers) then
1492 :
1493 : !-----------------------------------------------------------------
1494 : ! assign tracer indices and dependencies
1495 : ! bgc_tracer_type: < 0 purely mobile , >= 0 stationary
1496 : !------------------------------------------------------------------
1497 :
1498 5 : if (tr_bgc_N) then
1499 20 : do mm = 1, n_algae
1500 : call init_bgc_trcr(nk, nt_fbri, &
1501 9 : nt_bgc_N(mm), nlt_bgc_N(mm), &
1502 9 : algaltype(mm), nt_depend, &
1503 : ntrcr, nbtrcr, &
1504 : bgc_tracer_type, trcr_depend, &
1505 : trcr_base, n_trcr_strata, &
1506 15 : nt_strata, bio_index)
1507 20 : bio_index_o(nlt_bgc_N(mm)) = mm
1508 : enddo ! mm
1509 : endif ! tr_bgc_N
1510 :
1511 5 : if (tr_bgc_Nit) then
1512 : call init_bgc_trcr(nk, nt_fbri, &
1513 : nt_bgc_Nit, nlt_bgc_Nit, &
1514 : nitratetype, nt_depend, &
1515 : ntrcr, nbtrcr, &
1516 : bgc_tracer_type, trcr_depend, &
1517 : trcr_base, n_trcr_strata, &
1518 5 : nt_strata, bio_index)
1519 5 : bio_index_o(nlt_bgc_Nit) = icepack_max_algae + 1
1520 : endif ! tr_bgc_Nit
1521 :
1522 5 : if (tr_bgc_C) then
1523 : !
1524 : ! Algal C is not yet distinct from algal N
1525 : ! * Reqires exudation and/or changing C:N ratios
1526 : ! for implementation
1527 : !
1528 : ! do mm = 1,n_algae
1529 : ! call init_bgc_trcr(nk, nt_fbri, &
1530 : ! nt_bgc_C(mm), nlt_bgc_C(mm), &
1531 : ! algaltype(mm), nt_depend, &
1532 : ! ntrcr, nbtrcr, &
1533 : ! bgc_tracer_type, trcr_depend, &
1534 : ! trcr_base, n_trcr_strata, &
1535 : ! nt_strata, bio_index)
1536 : ! bio_index_o(nlt_bgc_C(mm)) = icepack_max_algae + 1 + mm
1537 : ! enddo ! mm
1538 :
1539 15 : do mm = 1, n_doc
1540 : call init_bgc_trcr(nk, nt_fbri, &
1541 6 : nt_bgc_DOC(mm), nlt_bgc_DOC(mm), &
1542 6 : doctype(mm), nt_depend, &
1543 : ntrcr, nbtrcr, &
1544 : bgc_tracer_type, trcr_depend, &
1545 : trcr_base, n_trcr_strata, &
1546 10 : nt_strata, bio_index)
1547 15 : bio_index_o(nlt_bgc_DOC(mm)) = icepack_max_algae + 1 + mm
1548 : enddo ! mm
1549 5 : do mm = 1, n_dic
1550 : call init_bgc_trcr(nk, nt_fbri, &
1551 0 : nt_bgc_DIC(mm), nlt_bgc_DIC(mm), &
1552 0 : dictype(mm), nt_depend, &
1553 : ntrcr, nbtrcr, &
1554 : bgc_tracer_type, trcr_depend, &
1555 : trcr_base, n_trcr_strata, &
1556 0 : nt_strata, bio_index)
1557 5 : bio_index_o(nlt_bgc_DIC(mm)) = icepack_max_algae + icepack_max_doc + 1 + mm
1558 : enddo ! mm
1559 : endif ! tr_bgc_C
1560 :
1561 5 : if (tr_bgc_chl) then
1562 0 : do mm = 1, n_algae
1563 : call init_bgc_trcr(nk, nt_fbri, &
1564 0 : nt_bgc_chl(mm), nlt_bgc_chl(mm), &
1565 0 : algaltype(mm), nt_depend, &
1566 : ntrcr, nbtrcr, &
1567 : bgc_tracer_type, trcr_depend, &
1568 : trcr_base, n_trcr_strata, &
1569 0 : nt_strata, bio_index)
1570 0 : bio_index_o(nlt_bgc_chl(mm)) = icepack_max_algae + 1 + icepack_max_doc + icepack_max_dic + mm
1571 : enddo ! mm
1572 : endif ! tr_bgc_chl
1573 :
1574 5 : if (tr_bgc_Am) then
1575 : call init_bgc_trcr(nk, nt_fbri, &
1576 : nt_bgc_Am, nlt_bgc_Am, &
1577 : ammoniumtype, nt_depend, &
1578 : ntrcr, nbtrcr, &
1579 : bgc_tracer_type, trcr_depend, &
1580 : trcr_base, n_trcr_strata, &
1581 5 : nt_strata, bio_index)
1582 5 : bio_index_o(nlt_bgc_Am) = 2*icepack_max_algae + icepack_max_doc + icepack_max_dic + 2
1583 : endif
1584 5 : if (tr_bgc_Sil) then
1585 : call init_bgc_trcr(nk, nt_fbri, &
1586 : nt_bgc_Sil, nlt_bgc_Sil, &
1587 : silicatetype, nt_depend, &
1588 : ntrcr, nbtrcr, &
1589 : bgc_tracer_type, trcr_depend, &
1590 : trcr_base, n_trcr_strata, &
1591 5 : nt_strata, bio_index)
1592 5 : bio_index_o(nlt_bgc_Sil) = 2*icepack_max_algae + icepack_max_doc + icepack_max_dic + 3
1593 : endif
1594 5 : if (tr_bgc_DMS) then ! all together
1595 : call init_bgc_trcr(nk, nt_fbri, &
1596 : nt_bgc_DMSPp, nlt_bgc_DMSPp, &
1597 : dmspptype, nt_depend, &
1598 : ntrcr, nbtrcr, &
1599 : bgc_tracer_type, trcr_depend, &
1600 : trcr_base, n_trcr_strata, &
1601 5 : nt_strata, bio_index)
1602 5 : bio_index_o(nlt_bgc_DMSPp) = 2*icepack_max_algae + icepack_max_doc + icepack_max_dic + 4
1603 :
1604 : call init_bgc_trcr(nk, nt_fbri, &
1605 : nt_bgc_DMSPd, nlt_bgc_DMSPd, &
1606 : dmspdtype, nt_depend, &
1607 : ntrcr, nbtrcr, &
1608 : bgc_tracer_type, trcr_depend, &
1609 : trcr_base, n_trcr_strata, &
1610 5 : nt_strata, bio_index)
1611 5 : bio_index_o(nlt_bgc_DMSPd) = 2*icepack_max_algae + icepack_max_doc + icepack_max_dic + 5
1612 :
1613 : call init_bgc_trcr(nk, nt_fbri, &
1614 : nt_bgc_DMS, nlt_bgc_DMS, &
1615 : dmspdtype, nt_depend, &
1616 : ntrcr, nbtrcr, &
1617 : bgc_tracer_type, trcr_depend, &
1618 : trcr_base, n_trcr_strata, &
1619 5 : nt_strata, bio_index)
1620 5 : bio_index_o(nlt_bgc_DMS) = 2*icepack_max_algae + icepack_max_doc + icepack_max_dic + 6
1621 : endif
1622 5 : if (tr_bgc_PON) then
1623 : call init_bgc_trcr(nk, nt_fbri, &
1624 : nt_bgc_PON, nlt_bgc_PON, &
1625 : nitratetype, nt_depend, &
1626 : ntrcr, nbtrcr, &
1627 : bgc_tracer_type, trcr_depend, &
1628 : trcr_base, n_trcr_strata, &
1629 5 : nt_strata, bio_index)
1630 5 : bio_index_o(nlt_bgc_PON) = 2*icepack_max_algae + icepack_max_doc + icepack_max_dic + 7
1631 : endif
1632 5 : if (tr_bgc_DON) then
1633 10 : do mm = 1, n_don
1634 : call init_bgc_trcr(nk, nt_fbri, &
1635 3 : nt_bgc_DON(mm), nlt_bgc_DON(mm), &
1636 3 : dontype(mm), nt_depend, &
1637 : ntrcr, nbtrcr, &
1638 : bgc_tracer_type, trcr_depend, &
1639 : trcr_base, n_trcr_strata, &
1640 5 : nt_strata, bio_index)
1641 10 : bio_index_o(nlt_bgc_DON(mm)) = 2*icepack_max_algae + icepack_max_doc + icepack_max_dic + 7 + mm
1642 : enddo ! mm
1643 : endif ! tr_bgc_DON
1644 5 : if (tr_bgc_Fe) then
1645 4 : do mm = 1, n_fed
1646 : call init_bgc_trcr(nk, nt_fbri, &
1647 2 : nt_bgc_Fed(mm), nlt_bgc_Fed(mm), &
1648 2 : fedtype(mm), nt_depend, &
1649 : ntrcr, nbtrcr, &
1650 : bgc_tracer_type, trcr_depend, &
1651 : trcr_base, n_trcr_strata, &
1652 2 : nt_strata, bio_index)
1653 4 : bio_index_o(nlt_bgc_Fed(mm)) = 2*icepack_max_algae + icepack_max_doc + icepack_max_dic &
1654 6 : + icepack_max_don + 7 + mm
1655 : enddo ! mm
1656 4 : do mm = 1, n_fep
1657 : call init_bgc_trcr(nk, nt_fbri, &
1658 2 : nt_bgc_Fep(mm), nlt_bgc_Fep(mm), &
1659 2 : feptype(mm), nt_depend, &
1660 : ntrcr, nbtrcr, &
1661 : bgc_tracer_type, trcr_depend, &
1662 : trcr_base, n_trcr_strata, &
1663 2 : nt_strata, bio_index)
1664 4 : bio_index_o(nlt_bgc_Fep(mm)) = 2*icepack_max_algae + icepack_max_doc + icepack_max_dic &
1665 6 : + icepack_max_don + icepack_max_fe + 7 + mm
1666 : enddo ! mm
1667 : endif ! tr_bgc_Fe
1668 :
1669 5 : if (tr_bgc_hum) then
1670 : call init_bgc_trcr(nk, nt_fbri, &
1671 : nt_bgc_hum, nlt_bgc_hum, &
1672 : humtype, nt_depend, &
1673 : ntrcr, nbtrcr, &
1674 : bgc_tracer_type, trcr_depend, &
1675 : trcr_base, n_trcr_strata, &
1676 5 : nt_strata, bio_index)
1677 3 : bio_index_o(nlt_bgc_hum) = 2*icepack_max_algae + icepack_max_doc + 8 + icepack_max_dic &
1678 5 : + icepack_max_don + 2*icepack_max_fe + icepack_max_aero
1679 : endif
1680 : endif ! skl_bgc or z_tracers
1681 :
1682 45 : if (z_tracers) then ! defined on nblyr+1 in ice
1683 : ! and 2 snow layers (snow surface + interior)
1684 :
1685 4 : nk = nblyr + 1
1686 4 : nt_depend = 2 + nt_fbri + ntd
1687 :
1688 : ! z layer aerosols
1689 4 : if (tr_zaero) then
1690 21 : do mm = 1, n_zaero
1691 18 : if (dEdd_algae) then
1692 0 : nlt_zaero_sw(mm) = nbtrcr_sw + 1
1693 0 : nbtrcr_sw = nbtrcr_sw + nilyr + nslyr+2
1694 : endif
1695 : call init_bgc_trcr(nk, nt_fbri, &
1696 6 : nt_zaero(mm), nlt_zaero(mm), &
1697 6 : zaerotype(mm), nt_depend, &
1698 : ntrcr, nbtrcr, &
1699 : bgc_tracer_type, trcr_depend, &
1700 : trcr_base, n_trcr_strata, &
1701 18 : nt_strata, bio_index)
1702 12 : bio_index_o(nlt_zaero(mm)) = 2*icepack_max_algae + icepack_max_doc + icepack_max_dic &
1703 27 : + icepack_max_don + 2*icepack_max_fe + 7 + mm
1704 : enddo ! mm
1705 : endif ! tr_zaero
1706 :
1707 : !echmod keep trcr indices etc here but move zbgc_frac_init, zbgc_init_frac, tau_ret, tau_rel to icepack
1708 4 : if (nbtrcr > 0) then
1709 4 : nt_zbgc_frac = ntrcr + 1
1710 4 : ntrcr = ntrcr + nbtrcr
1711 80 : do k = 1,nbtrcr
1712 76 : zbgc_frac_init(k) = c1
1713 76 : trcr_depend(nt_zbgc_frac+k-1) = 2+nt_fbri
1714 76 : trcr_base(nt_zbgc_frac+ k - 1,1) = c0
1715 76 : trcr_base(nt_zbgc_frac+ k - 1,2) = c1
1716 76 : trcr_base(nt_zbgc_frac+ k - 1,3) = c0
1717 76 : n_trcr_strata(nt_zbgc_frac+ k - 1)= 1
1718 76 : nt_strata(nt_zbgc_frac+ k - 1,1) = nt_fbri
1719 76 : nt_strata(nt_zbgc_frac+ k - 1,2) = 0
1720 76 : tau_ret(k) = c1
1721 76 : tau_rel(k) = c1
1722 80 : if (bgc_tracer_type(k) >= c0 .and. bgc_tracer_type(k) < p5) then
1723 4 : tau_ret(k) = tau_min
1724 4 : tau_rel(k) = tau_max
1725 4 : zbgc_frac_init(k) = c1
1726 72 : elseif (bgc_tracer_type(k) >= p5 .and. bgc_tracer_type(k) < c1) then
1727 26 : tau_ret(k) = tau_min
1728 26 : tau_rel(k) = tau_min
1729 26 : zbgc_frac_init(k) = c1
1730 46 : elseif (bgc_tracer_type(k) >= c1 .and. bgc_tracer_type(k) < c2) then
1731 26 : tau_ret(k) = tau_max
1732 26 : tau_rel(k) = tau_min
1733 26 : zbgc_frac_init(k) = c1
1734 20 : elseif (bgc_tracer_type(k) >= c2 ) then
1735 0 : tau_ret(k) = tau_max
1736 0 : tau_rel(k) = tau_max
1737 0 : zbgc_frac_init(k) = c1
1738 : endif
1739 : enddo
1740 : endif
1741 :
1742 : endif ! z_tracers
1743 :
1744 137 : do k = 1, nbtrcr
1745 92 : zbgc_init_frac(k) = frazil_scav
1746 137 : if (bgc_tracer_type(k) < c0) zbgc_init_frac(k) = initbio_frac
1747 : enddo
1748 :
1749 : !-----------------------------------------------------------------
1750 : ! set values in icepack
1751 : !-----------------------------------------------------------------
1752 :
1753 : call icepack_init_zbgc( &
1754 : zbgc_init_frac_in=zbgc_init_frac, tau_ret_in=tau_ret, tau_rel_in=tau_rel, &
1755 45 : zbgc_frac_init_in=zbgc_frac_init, bgc_tracer_type_in=bgc_tracer_type)
1756 45 : call icepack_warnings_flush(nu_diag)
1757 45 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
1758 0 : file=__FILE__, line=__LINE__)
1759 :
1760 : call icepack_init_tracer_sizes( &
1761 : n_algae_in=n_algae, &
1762 : n_DOC_in=n_DOC, n_DON_in=n_DON, n_DIC_in=n_DIC, &
1763 : n_fed_in=n_fed, n_fep_in=n_fep, n_zaero_in=n_zaero, &
1764 45 : ntrcr_in=ntrcr, ntrcr_o_in=ntrcr_o, nbtrcr_in=nbtrcr, nbtrcr_sw_in=nbtrcr_sw)
1765 45 : call icepack_warnings_flush(nu_diag)
1766 45 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
1767 0 : file=__FILE__, line=__LINE__)
1768 :
1769 : call icepack_init_tracer_flags( &
1770 : tr_brine_in =tr_brine, &
1771 : tr_bgc_Nit_in=tr_bgc_Nit, tr_bgc_Am_in =tr_bgc_Am, tr_bgc_Sil_in=tr_bgc_Sil, &
1772 : tr_bgc_DMS_in=tr_bgc_DMS, tr_bgc_PON_in=tr_bgc_PON, &
1773 : tr_bgc_N_in =tr_bgc_N, tr_bgc_C_in =tr_bgc_C, tr_bgc_chl_in=tr_bgc_chl, &
1774 : tr_bgc_DON_in=tr_bgc_DON, tr_bgc_Fe_in =tr_bgc_Fe, tr_zaero_in =tr_zaero, &
1775 45 : tr_bgc_hum_in=tr_bgc_hum)
1776 45 : call icepack_warnings_flush(nu_diag)
1777 45 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
1778 0 : file=__FILE__, line=__LINE__)
1779 :
1780 : call icepack_init_tracer_indices( &
1781 : nt_fbri_in=nt_fbri, &
1782 : nt_bgc_Nit_in=nt_bgc_Nit, nt_bgc_Am_in=nt_bgc_Am, nt_bgc_Sil_in=nt_bgc_Sil, &
1783 : nt_bgc_DMS_in=nt_bgc_DMS, nt_bgc_PON_in=nt_bgc_PON, nt_bgc_S_in=nt_bgc_S, &
1784 : nt_bgc_N_in=nt_bgc_N, nt_bgc_chl_in=nt_bgc_chl, nt_bgc_hum_in=nt_bgc_hum, &
1785 : nt_bgc_DOC_in=nt_bgc_DOC, nt_bgc_DON_in=nt_bgc_DON, nt_bgc_DIC_in=nt_bgc_DIC, &
1786 : nt_zaero_in=nt_zaero, nt_bgc_DMSPp_in=nt_bgc_DMSPp, nt_bgc_DMSPd_in=nt_bgc_DMSPd, &
1787 : nt_bgc_Fed_in=nt_bgc_Fed, nt_bgc_Fep_in=nt_bgc_Fep, nt_zbgc_frac_in=nt_zbgc_frac, &
1788 : nlt_chl_sw_in=nlt_chl_sw, nlt_bgc_Sil_in=nlt_bgc_Sil, nlt_zaero_sw_in=nlt_zaero_sw, &
1789 : nlt_bgc_N_in=nlt_bgc_N, nlt_bgc_Nit_in=nlt_bgc_Nit, nlt_bgc_Am_in=nlt_bgc_Am, &
1790 : nlt_bgc_DMS_in=nlt_bgc_DMS, nlt_bgc_DMSPp_in=nlt_bgc_DMSPp, &
1791 : nlt_bgc_DMSPd_in=nlt_bgc_DMSPd, &
1792 : nlt_bgc_DIC_in=nlt_bgc_DIC, nlt_bgc_DOC_in=nlt_bgc_DOC, nlt_bgc_PON_in=nlt_bgc_PON, &
1793 : nlt_bgc_DON_in=nlt_bgc_DON, nlt_bgc_Fed_in=nlt_bgc_Fed, nlt_bgc_Fep_in=nlt_bgc_Fep, &
1794 : nlt_bgc_chl_in=nlt_bgc_chl, nlt_bgc_hum_in=nlt_bgc_hum, nlt_zaero_in=nlt_zaero, &
1795 45 : bio_index_o_in=bio_index_o, bio_index_in=bio_index)
1796 45 : call icepack_warnings_flush(nu_diag)
1797 45 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
1798 0 : file=__FILE__, line=__LINE__)
1799 :
1800 : !-----------------------------------------------------------------
1801 : ! final consistency checks
1802 : !-----------------------------------------------------------------
1803 45 : if (nbtrcr > icepack_max_nbtrcr) then
1804 0 : write (nu_diag,*) ' '
1805 0 : write (nu_diag,*) 'nbtrcr > icepack_max_nbtrcr'
1806 0 : write (nu_diag,*) 'nbtrcr, icepack_max_nbtrcr:',nbtrcr, icepack_max_nbtrcr
1807 0 : call icedrv_system_abort(file=__FILE__,line=__LINE__)
1808 : endif
1809 45 : if (.NOT. dEdd_algae) nbtrcr_sw = 1
1810 :
1811 45 : if (nbtrcr_sw > max_nsw) then
1812 0 : write (nu_diag,*) ' '
1813 0 : write (nu_diag,*) 'nbtrcr_sw > max_nsw'
1814 0 : write (nu_diag,*) 'nbtrcr_sw, max_nsw:',nbtrcr_sw, max_nsw
1815 0 : call icedrv_system_abort(file=__FILE__,line=__LINE__)
1816 : endif
1817 :
1818 45 : if (ntrcr > max_ntrcr) then
1819 0 : write(nu_diag,*) 'max_ntrcr < number of namelist tracers'
1820 0 : write(nu_diag,*) 'max_ntrcr = ',max_ntrcr,' ntrcr = ',ntrcr
1821 0 : call icedrv_system_abort(file=__FILE__,line=__LINE__)
1822 : endif
1823 :
1824 : !-----------------------------------------------------------------
1825 : ! spew
1826 : !-----------------------------------------------------------------
1827 45 : if (skl_bgc) then
1828 :
1829 1 : write(nu_diag,1010) ' skl_bgc = ', skl_bgc
1830 1 : write(nu_diag,1030) ' bgc_flux_type = ', bgc_flux_type
1831 1 : write(nu_diag,1010) ' restore_bgc = ', restore_bgc
1832 1 : write(nu_diag,*) ' bgc_data_type = ', &
1833 2 : trim(bgc_data_type)
1834 1 : write(nu_diag,1020) ' number of bio tracers = ', nbtrcr
1835 1 : write(nu_diag,1020) ' number of Isw tracers = ', nbtrcr_sw
1836 1 : write(nu_diag,1020) ' number of autotrophs = ', n_algae
1837 1 : write(nu_diag,1020) ' number of doc = ', n_doc
1838 1 : write(nu_diag,1020) ' number of dic = ', n_dic
1839 1 : write(nu_diag,1020) ' number of don = ', n_don
1840 1 : write(nu_diag,1020) ' number of fed = ', n_fed
1841 1 : write(nu_diag,1020) ' number of fep = ', n_fep
1842 1 : write(nu_diag,1010) ' tr_bgc_N = ', tr_bgc_N
1843 1 : write(nu_diag,1010) ' tr_bgc_C = ', tr_bgc_C
1844 1 : write(nu_diag,1010) ' tr_bgc_chl = ', tr_bgc_chl
1845 1 : write(nu_diag,1010) ' tr_bgc_Nit = ', tr_bgc_Nit
1846 1 : write(nu_diag,1010) ' tr_bgc_Am = ', tr_bgc_Am
1847 1 : write(nu_diag,1010) ' tr_bgc_Sil = ', tr_bgc_Sil
1848 1 : write(nu_diag,1010) ' tr_bgc_hum = ', tr_bgc_hum
1849 1 : write(nu_diag,1010) ' tr_bgc_DMS = ', tr_bgc_DMS
1850 1 : write(nu_diag,1010) ' tr_bgc_PON = ', tr_bgc_PON
1851 1 : write(nu_diag,1010) ' tr_bgc_DON = ', tr_bgc_DON
1852 1 : write(nu_diag,1010) ' tr_bgc_Fe = ', tr_bgc_Fe
1853 :
1854 44 : elseif (z_tracers) then
1855 :
1856 4 : write(nu_diag,*) ' bgc_data_type = ', &
1857 8 : trim(bgc_data_type)
1858 4 : write(nu_diag,1010) ' dEdd_algae = ', dEdd_algae
1859 4 : write(nu_diag,1010) ' modal_aero = ', modal_aero
1860 4 : write(nu_diag,1010) ' scale_bgc = ', scale_bgc
1861 4 : write(nu_diag,1010) ' solve_zbgc = ', solve_zbgc
1862 4 : write(nu_diag,1020) ' number of ztracers = ', nbtrcr
1863 4 : write(nu_diag,1020) ' number of Isw tracers = ', nbtrcr_sw
1864 4 : write(nu_diag,1020) ' number of autotrophs = ', n_algae
1865 4 : write(nu_diag,1020) ' number of doc = ', n_doc
1866 4 : write(nu_diag,1020) ' number of dic = ', n_dic
1867 4 : write(nu_diag,1020) ' number of fed = ', n_fed
1868 4 : write(nu_diag,1020) ' number of fep = ', n_fep
1869 4 : write(nu_diag,1020) ' number of aerosols = ', n_zaero
1870 4 : write(nu_diag,1010) ' tr_zaero = ', tr_zaero
1871 4 : write(nu_diag,1010) ' tr_bgc_Nit = ', tr_bgc_Nit
1872 4 : write(nu_diag,1010) ' tr_bgc_N = ', tr_bgc_N
1873 4 : write(nu_diag,1010) ' tr_bgc_Am = ', tr_bgc_Am
1874 4 : write(nu_diag,1010) ' tr_bgc_C = ', tr_bgc_C
1875 4 : write(nu_diag,1010) ' tr_bgc_Sil = ', tr_bgc_Sil
1876 4 : write(nu_diag,1010) ' tr_bgc_hum = ', tr_bgc_hum
1877 4 : write(nu_diag,1010) ' tr_bgc_chl = ', tr_bgc_chl
1878 4 : write(nu_diag,1010) ' tr_bgc_DMS = ', tr_bgc_DMS
1879 4 : write(nu_diag,1010) ' tr_bgc_PON = ', tr_bgc_PON
1880 4 : write(nu_diag,1010) ' tr_bgc_DON = ', tr_bgc_DON
1881 4 : write(nu_diag,1010) ' tr_bgc_Fe = ', tr_bgc_Fe
1882 4 : write(nu_diag,1000) ' grid_o = ', grid_o
1883 4 : write(nu_diag,1005) ' l_sk = ', l_sk
1884 4 : write(nu_diag,1000) ' initbio_frac = ', initbio_frac
1885 4 : write(nu_diag,1000) ' frazil_scav = ', frazil_scav
1886 :
1887 : endif ! skl_bgc or solve_bgc
1888 :
1889 : 1000 format (a30,2x,f9.2) ! a30 to align formatted, unformatted statements
1890 : 1005 format (a30,2x,f9.6) ! float
1891 : 1010 format (a30,2x,l6) ! logical
1892 : 1020 format (a30,2x,i6) ! integer
1893 : 1030 format (a30, a8) ! character
1894 :
1895 45 : end subroutine init_zbgc
1896 :
1897 : !=======================================================================
1898 :
1899 92 : subroutine init_bgc_trcr(nk, nt_fbri, &
1900 : nt_bgc, nlt_bgc, &
1901 : bgctype, nt_depend, &
1902 : ntrcr, nbtrcr, &
1903 92 : bgc_tracer_type, trcr_depend, &
1904 92 : trcr_base, n_trcr_strata, &
1905 92 : nt_strata, bio_index)
1906 :
1907 : integer (kind=int_kind), intent(in) :: &
1908 : nk , & ! counter
1909 : nt_depend , & ! tracer dependency index
1910 : nt_fbri
1911 :
1912 : integer (kind=int_kind), intent(inout) :: &
1913 : ntrcr , & ! number of tracers
1914 : nbtrcr , & ! number of bio tracers
1915 : nt_bgc , & ! tracer index
1916 : nlt_bgc ! bio tracer index
1917 :
1918 : integer (kind=int_kind), dimension(:), intent(inout) :: &
1919 : trcr_depend , & ! tracer dependencies
1920 : n_trcr_strata, & ! number of underlying tracer layers
1921 : bio_index !
1922 :
1923 : integer (kind=int_kind), dimension(:,:), intent(inout) :: &
1924 : nt_strata ! indices of underlying tracer layers
1925 :
1926 : real (kind=dbl_kind), dimension(:,:), intent(inout) :: &
1927 : trcr_base ! = 0 or 1 depending on tracer dependency
1928 : ! argument 2: (1) aice, (2) vice, (3) vsno
1929 :
1930 : real (kind=dbl_kind), intent(in) :: &
1931 : bgctype ! bio tracer transport type (mobile vs stationary)
1932 :
1933 : real (kind=dbl_kind), dimension(:), intent(inout) :: &
1934 : bgc_tracer_type ! bio tracer transport type array
1935 :
1936 : ! local variables
1937 :
1938 : integer (kind=int_kind) :: &
1939 : k , & ! loop index
1940 : n_strata , & ! temporary values
1941 : nt_strata1, & !
1942 : nt_strata2
1943 :
1944 : real (kind=dbl_kind) :: &
1945 52 : trcr_base1, & ! temporary values
1946 52 : trcr_base2, &
1947 52 : trcr_base3
1948 :
1949 : character(len=*), parameter :: subname='(init_bgc_trcr)'
1950 :
1951 : !--------
1952 :
1953 92 : nt_bgc = ntrcr + 1
1954 92 : nbtrcr = nbtrcr + 1
1955 92 : nlt_bgc = nbtrcr
1956 92 : bgc_tracer_type(nbtrcr) = bgctype
1957 :
1958 92 : if (nk > 1) then
1959 : ! include vertical bgc in snow
1960 228 : do k = nk, nk+1
1961 152 : ntrcr = ntrcr + 1
1962 152 : trcr_depend (nt_bgc + k ) = 2 ! snow volume
1963 152 : trcr_base (nt_bgc + k,1) = c0
1964 152 : trcr_base (nt_bgc + k,2) = c0
1965 152 : trcr_base (nt_bgc + k,3) = c1
1966 152 : n_trcr_strata(nt_bgc + k ) = 0
1967 152 : nt_strata (nt_bgc + k,1) = 0
1968 228 : nt_strata (nt_bgc + k,2) = 0
1969 : enddo
1970 :
1971 76 : trcr_base1 = c0
1972 76 : trcr_base2 = c1
1973 76 : trcr_base3 = c0
1974 76 : n_strata = 1
1975 76 : nt_strata1 = nt_fbri
1976 76 : nt_strata2 = 0
1977 : else ! nk = 1
1978 16 : trcr_base1 = c1
1979 16 : trcr_base2 = c0
1980 16 : trcr_base3 = c0
1981 16 : n_strata = 0
1982 16 : nt_strata1 = 0
1983 16 : nt_strata2 = 0
1984 : endif ! nk
1985 :
1986 716 : do k = 1, nk !in ice
1987 624 : ntrcr = ntrcr + 1
1988 624 : trcr_depend (nt_bgc + k - 1 ) = nt_depend
1989 624 : trcr_base (nt_bgc + k - 1,1) = trcr_base1
1990 624 : trcr_base (nt_bgc + k - 1,2) = trcr_base2
1991 624 : trcr_base (nt_bgc + k - 1,3) = trcr_base3
1992 624 : n_trcr_strata(nt_bgc + k - 1 ) = n_strata
1993 624 : nt_strata (nt_bgc + k - 1,1) = nt_strata1
1994 716 : nt_strata (nt_bgc + k - 1,2) = nt_strata2
1995 : enddo
1996 :
1997 92 : bio_index (nlt_bgc) = nt_bgc
1998 :
1999 92 : end subroutine init_bgc_trcr
2000 :
2001 : !=======================================================================
2002 :
2003 : end module icedrv_init_column
2004 :
2005 : !=======================================================================
|