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