Line data Source code
1 : !=======================================================================
2 :
3 : ! Diagnostic information output during run for biogeochemistry
4 : !
5 : ! authors: Elizabeth C. Hunke, LANL
6 : ! Nicole Jeffery, LANL
7 :
8 : module icedrv_diagnostics_bgc
9 :
10 : use icedrv_kinds
11 : use icedrv_constants, only: nu_diag, nu_diag_out
12 : use icedrv_constants, only: c0, mps_to_cmpdy, c100, p5, c1
13 : use icedrv_domain_size, only: nx
14 : use icepack_intfc, only: icepack_warnings_flush, icepack_warnings_aborted
15 : use icepack_intfc, only: icepack_query_parameters
16 : use icepack_intfc, only: icepack_query_tracer_flags, icepack_query_tracer_indices
17 : use icepack_intfc, only: icepack_max_algae, icepack_max_aero, icepack_max_fe
18 : use icepack_intfc, only: icepack_max_doc, icepack_max_don
19 : use icedrv_system, only: icedrv_system_abort
20 :
21 : implicit none
22 : private
23 : public :: hbrine_diags, bgc_diags, zsal_diags
24 :
25 : !=======================================================================
26 :
27 : contains
28 :
29 : !=======================================================================
30 :
31 : !
32 : ! Writes diagnostic info (max, min, global sums, etc) to standard out
33 : !
34 : ! authors: Nicole Jeffery, LANL
35 :
36 925 : subroutine hbrine_diags ()
37 :
38 : use icedrv_arrays_column, only: darcy_V
39 : use icedrv_diagnostics, only: nx_names
40 : use icedrv_domain_size, only: nilyr
41 : use icedrv_state, only: aice, aicen, vicen, vice, trcr, trcrn
42 :
43 : ! local variables
44 :
45 : integer (kind=int_kind) :: &
46 : k, n
47 :
48 : integer (kind=int_kind) :: &
49 : ktherm
50 :
51 : integer (kind=int_kind) :: &
52 : nt_sice, &
53 : nt_fbri
54 :
55 : ! fields at diagnostic points
56 : real (kind=dbl_kind) :: &
57 285 : phinS, phinS1, pdarcy_V, pfbri
58 :
59 : real (kind=dbl_kind), dimension(nilyr) :: &
60 4275 : pSin, pSin1
61 :
62 : character(len=*), parameter :: subname='(hbrine_diags)'
63 :
64 : !-----------------------------------------------------------------
65 : ! query Icepack values
66 : !-----------------------------------------------------------------
67 :
68 925 : call icepack_query_parameters(ktherm_out=ktherm)
69 925 : call icepack_query_tracer_indices(nt_sice_out=nt_sice, nt_fbri_out=nt_fbri)
70 925 : call icepack_warnings_flush(nu_diag)
71 925 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
72 0 : file=__FILE__,line= __LINE__)
73 :
74 : !-----------------------------------------------------------------
75 : ! Dynamic brine height
76 : !-----------------------------------------------------------------
77 : ! NOTE these are computed for the last timestep only (not avg)
78 :
79 4625 : do n = 1, nx
80 3700 : phinS1 = c0
81 3700 : phinS = c0
82 3700 : pfbri = trcrn(n,nt_fbri,1)
83 3700 : pdarcy_V = darcy_V(n,1)
84 3700 : if (aice (n) > c0) phinS = trcr (n,nt_fbri )*vice (n )/aice (n )
85 3700 : if (aicen(n,1)> c0) phinS1 = trcrn(n,nt_fbri,1)*vicen(n,1)/aicen(n,1)
86 29600 : do k = 1 , nilyr
87 25900 : pSin (k) = trcr (n,nt_sice+k-1 )
88 29600 : pSin1(k) = trcrn(n,nt_sice+k-1,1)
89 : enddo
90 :
91 : !-----------------------------------------------------------------
92 : ! start spewing
93 : !-----------------------------------------------------------------
94 :
95 3700 : write(nu_diag_out+n-1,899) nx_names(n)
96 :
97 3700 : write(nu_diag_out+n-1,*) '------ hbrine ------'
98 3700 : write(nu_diag_out+n-1,900) 'hbrine, (m) = ',phinS
99 3700 : write(nu_diag_out+n-1,900) 'fbri, cat1 (m) = ',pfbri
100 3700 : write(nu_diag_out+n-1,900) 'hbrine cat1, (m) = ',phinS1
101 3700 : write(nu_diag_out+n-1,900) 'darcy_V cat1, (m/s)= ',pdarcy_V
102 4625 : if (ktherm == 2) then
103 3700 : write(nu_diag_out+n-1,*) ' '
104 3700 : write(nu_diag_out+n-1,*) '------ Thermosaline Salinity ------'
105 3700 : write(nu_diag_out+n-1,803) 'Sice1 cat1 S (ppt)'
106 3700 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
107 29600 : write(nu_diag_out+n-1,802) (pSin1(k), k = 1,nilyr)
108 3700 : write(nu_diag_out+n-1,*) ' '
109 3700 : write(nu_diag_out+n-1,803) 'Sice bulk S (ppt) '
110 3700 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
111 29600 : write(nu_diag_out+n-1,802) (pSin(k), k = 1,nilyr)
112 3700 : write(nu_diag_out+n-1,*) ' '
113 : endif
114 : enddo ! nx
115 :
116 : 802 format (f24.17,2x)
117 : 803 format (a25,2x)
118 : 899 format (43x,a24)
119 : 900 format (a25,2x,f24.17)
120 :
121 925 : end subroutine hbrine_diags
122 :
123 : !=======================================================================
124 : !
125 : ! Writes diagnostic info (max, min, global sums, etc) to standard out
126 : !
127 : ! authors: Nicole Jeffery, LANL
128 :
129 925 : subroutine bgc_diags ()
130 :
131 : use icedrv_arrays_column, only: ocean_bio, zfswin, fbio_atmice, fbio_snoice
132 : use icedrv_arrays_column, only: Zoo, grow_net, ice_bio_net, trcrn_sw
133 : use icedrv_domain_size, only: ncat, nblyr, n_algae, n_zaero
134 : use icedrv_domain_size, only: n_doc, n_don, n_fed, n_fep, nilyr, nslyr
135 : use icedrv_flux, only: flux_bio, flux_bio_atm
136 : use icedrv_state, only: vicen, vice, trcr
137 :
138 : ! local variables
139 :
140 : integer (kind=int_kind) :: &
141 : k, n, nn, kk, klev
142 :
143 : logical (kind=log_kind) :: &
144 : skl_bgc, z_tracers, dEdd_algae
145 :
146 : ! fields at diagnostic points
147 : real (kind=dbl_kind) :: &
148 1140 : pNit_sk, pAm_sk, pSil_sk, phum_sk, &
149 570 : pDMSPp_sk, pDMSPd_sk, pDMS_sk, &
150 1425 : pNit_ac, pAm_ac, pSil_ac, pDMSP_ac, pDMS_ac, &
151 855 : pflux_NO, pflux_Am, phum_ac, &
152 285 : pflux_snow_NO, pflux_snow_Am, &
153 570 : pflux_atm_NO, pflux_atm_Am, pgrow_net, &
154 285 : pflux_hum
155 :
156 : logical (kind=log_kind) :: &
157 : tr_bgc_DMS, tr_bgc_PON, &
158 : tr_bgc_N, tr_bgc_C, tr_bgc_DON, tr_zaero, tr_bgc_hum, &
159 : tr_bgc_Nit, tr_bgc_Am, tr_bgc_Sil, tr_bgc_Fe
160 :
161 : integer (kind=int_kind) :: &
162 : nt_fbri, nt_sice, nt_bgc_nit, nt_bgc_am, nt_bgc_sil, &
163 : nt_bgc_hum, nt_bgc_pon, nt_bgc_dmspp, nt_bgc_dmspd, nt_bgc_dms, &
164 : nlt_bgc_hum, nlt_bgc_Nit, nlt_bgc_Am, nlt_bgc_Sil, nlt_chl_sw, &
165 : nlt_bgc_DMSPp, nlt_bgc_DMS
166 : integer (kind=int_kind), dimension(icepack_max_algae) :: &
167 : nt_bgc_n, nlt_bgc_N
168 : integer (kind=int_kind), dimension(icepack_max_doc) :: &
169 : nt_bgc_doc, nlt_bgc_DOC
170 : integer (kind=int_kind), dimension(icepack_max_don) :: &
171 : nt_bgc_don, nlt_bgc_DON
172 : integer (kind=int_kind), dimension(icepack_max_aero) :: &
173 : nt_zaero, nlt_zaero, nlt_zaero_sw
174 : integer (kind=int_kind), dimension(icepack_max_fe) :: &
175 : nt_bgc_fed, nt_bgc_fep, nlt_bgc_Fed, nlt_bgc_Fep
176 :
177 : real (kind=dbl_kind), dimension(icepack_max_algae) :: &
178 3990 : pN_ac, pN_tot, pN_sk, pflux_N
179 : real (kind=dbl_kind), dimension(icepack_max_doc) :: &
180 1995 : pDOC_ac, pDOC_sk
181 : real (kind=dbl_kind), dimension(icepack_max_don) :: &
182 855 : pDON_ac, pDON_sk
183 : real (kind=dbl_kind), dimension(icepack_max_fe ) :: &
184 2850 : pFed_ac, pFed_sk, pFep_ac, pFep_sk
185 : real (kind=dbl_kind), dimension(icepack_max_aero) :: &
186 5700 : pflux_zaero, pflux_snow_zaero, pflux_atm_zaero, &
187 1995 : pflux_atm_zaero_s
188 :
189 : ! vertical fields of category 1 at diagnostic points for bgc layer model
190 : real (kind=dbl_kind), dimension(2) :: &
191 3420 : pNOs, pAms, pPONs, phums
192 : real (kind=dbl_kind), dimension(2,icepack_max_algae) :: &
193 2850 : pNs
194 : real (kind=dbl_kind), dimension(2,icepack_max_doc) :: &
195 2850 : pDOCs
196 : real (kind=dbl_kind), dimension(2,icepack_max_don) :: &
197 1140 : pDONs
198 : real (kind=dbl_kind), dimension(2,icepack_max_fe ) :: &
199 3990 : pFeds, pFeps
200 : real (kind=dbl_kind), dimension(2,icepack_max_aero) :: &
201 5415 : pzaeros
202 : real (kind=dbl_kind), dimension(nblyr+1) :: &
203 13557 : pNO, pAm, pPON, pzfswin, pZoo, phum
204 : real (kind=dbl_kind), dimension(nblyr+1,icepack_max_algae) :: &
205 7206 : pN
206 : real (kind=dbl_kind), dimension(nblyr+1,icepack_max_aero) :: &
207 14127 : pzaero
208 : real (kind=dbl_kind), dimension(nblyr+1,icepack_max_doc) :: &
209 7206 : pDOC
210 : real (kind=dbl_kind), dimension(nblyr+1,icepack_max_don) :: &
211 2592 : pDON
212 : real (kind=dbl_kind), dimension(nblyr+1,icepack_max_fe ) :: &
213 9798 : pFed, pFep
214 : real (kind=dbl_kind), dimension (nblyr+1) :: &
215 2307 : zspace
216 : real (kind=dbl_kind), dimension (nslyr+nilyr+2) :: &
217 3135 : pchlsw
218 : real (kind=dbl_kind), dimension(nslyr+nilyr+2,icepack_max_aero) :: &
219 19095 : pzaerosw
220 :
221 : character(len=*), parameter :: subname='(bgc_diags)'
222 :
223 : !-----------------------------------------------------------------
224 : ! query Icepack values
225 : !-----------------------------------------------------------------
226 :
227 : call icepack_query_parameters(skl_bgc_out=skl_bgc, &
228 925 : z_tracers_out=z_tracers, dEdd_algae_out=dEdd_algae)
229 : call icepack_query_tracer_flags( &
230 : tr_bgc_DMS_out=tr_bgc_DMS, tr_bgc_PON_out=tr_bgc_PON, &
231 : tr_bgc_N_out=tr_bgc_N, tr_bgc_C_out=tr_bgc_C, &
232 : tr_bgc_DON_out=tr_bgc_DON, tr_zaero_out=tr_zaero, &
233 : tr_bgc_hum_out=tr_bgc_hum,tr_bgc_Sil_out=tr_bgc_Sil, &
234 : tr_bgc_Nit_out=tr_bgc_Nit, tr_bgc_Am_out=tr_bgc_Am, &
235 925 : tr_bgc_Fe_out=tr_bgc_Fe)
236 : call icepack_query_tracer_indices( &
237 : nt_fbri_out=nt_fbri, nt_sice_out=nt_sice, nt_zaero_out=nt_zaero, &
238 : nt_bgc_n_out=nt_bgc_n, nt_bgc_doc_out=nt_bgc_doc, &
239 : nt_bgc_don_out=nt_bgc_don, nt_bgc_sil_out=nt_bgc_sil, &
240 : nt_bgc_fed_out=nt_bgc_fed, nt_bgc_fep_out=nt_bgc_fep, &
241 : nt_bgc_nit_out=nt_bgc_nit, nt_bgc_am_out=nt_bgc_am, &
242 : nt_bgc_hum_out=nt_bgc_hum, nt_bgc_pon_out=nt_bgc_pon, &
243 : nt_bgc_dmspp_out=nt_bgc_dmspp, nt_bgc_dmspd_out=nt_bgc_dmspd, &
244 : nt_bgc_dms_out=nt_bgc_dms, nlt_chl_sw_out=nlt_chl_sw, &
245 : nlt_bgc_N_out=nlt_bgc_N, nlt_zaero_out=nlt_zaero, &
246 : nlt_zaero_sw_out=nlt_zaero_sw, &
247 : nlt_bgc_Fed_out=nlt_bgc_Fed, nlt_bgc_Fep_out=nlt_bgc_Fep, &
248 : nlt_bgc_hum_out=nlt_bgc_hum, nlt_bgc_Nit_out=nlt_bgc_Nit, &
249 : nlt_bgc_Am_out=nlt_bgc_Am, nlt_bgc_Sil_out=nlt_bgc_Sil, &
250 : nlt_bgc_DOC_out=nlt_bgc_DOC, nlt_bgc_DON_out=nlt_bgc_DON, &
251 925 : nlt_bgc_DMSPp_out=nlt_bgc_DMSPp, nlt_bgc_DMS_out=nlt_bgc_DMS)
252 925 : call icepack_warnings_flush(nu_diag)
253 925 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
254 0 : file=__FILE__,line= __LINE__)
255 :
256 8067 : zspace(:) = c1/real(nblyr,kind=dbl_kind)
257 925 : zspace(1) = zspace(1)*p5
258 925 : zspace(nblyr+1) = zspace(nblyr+1)*p5
259 925 : klev = 1+nilyr+nslyr
260 :
261 : !-----------------------------------------------------------------
262 : ! biogeochemical state of the ice
263 : !-----------------------------------------------------------------
264 : ! NOTE these are computed for the last timestep only (not avg)
265 :
266 4625 : do n = 1, nx
267 3700 : pAm_ac = c0
268 3700 : pSil_ac = c0
269 3700 : phum_ac = c0
270 3700 : pDMSP_ac = c0
271 3700 : pDMS_ac = c0
272 3700 : pN_ac(:) = c0
273 3700 : pDOC_ac(:)= c0
274 3700 : pDON_ac(:)= c0
275 3700 : pFed_ac(:)= c0
276 3700 : pFep_ac(:)= c0
277 3700 : pNit_ac = c0
278 3700 : if (tr_bgc_N) then
279 14800 : do k = 1,n_algae
280 14800 : pN_ac(k) = ocean_bio(n,nlt_bgc_N(k))
281 : enddo !n_algae
282 : endif !tr_bgc_N
283 3700 : if (tr_bgc_C) then
284 11100 : do k = 1,n_doc
285 11100 : pDOC_ac(k) = ocean_bio(n,nlt_bgc_DOC(k))
286 : enddo !n_algae
287 : endif !tr_bgc_N
288 3700 : if (tr_bgc_DON) then
289 7400 : do k = 1,n_don
290 7400 : pDON_ac(k) = ocean_bio(n,nlt_bgc_DON(k))
291 : enddo
292 : endif
293 3700 : if (tr_bgc_Fe ) then
294 688 : do k = 1,n_fed
295 688 : pFed_ac (k) = ocean_bio(n,nlt_bgc_Fed (k))
296 : enddo
297 688 : do k = 1,n_fep
298 688 : pFep_ac (k) = ocean_bio(n,nlt_bgc_Fep (k))
299 : enddo
300 : endif
301 3700 : if (tr_bgc_Nit) pNit_ac = ocean_bio(n,nlt_bgc_Nit)
302 3700 : if (tr_bgc_Am) pAm_ac = ocean_bio(n,nlt_bgc_Am)
303 3700 : if (tr_bgc_Sil) pSil_ac = ocean_bio(n,nlt_bgc_Sil)
304 3700 : if (tr_bgc_hum) phum_ac = ocean_bio(n,nlt_bgc_hum)
305 3700 : if (tr_bgc_DMS) then
306 3700 : pDMSP_ac = ocean_bio(n,nlt_bgc_DMSPp)
307 3700 : pDMS_ac = ocean_bio(n,nlt_bgc_DMS)
308 : endif
309 :
310 : ! fluxes in mmol/m^2/d
311 : ! concentrations are bulk in mmol/m^3
312 : ! iron is in 10^-3 mmol/m^3 (nM)
313 :
314 3700 : if (skl_bgc) then
315 172 : pNit_sk = c0
316 172 : pAm_sk = c0
317 172 : pSil_sk = c0
318 172 : phum_sk = c0
319 172 : pDMSPp_sk = c0
320 172 : pDMSPd_sk = c0
321 172 : pDMS_sk = c0
322 172 : pN_sk (:) = c0
323 172 : pflux_N(:) = c0
324 172 : pDOC_sk(:) = c0
325 172 : pDON_sk(:) = c0
326 172 : pFed_sk(:) = c0
327 172 : pFep_sk(:) = c0
328 172 : pgrow_net = grow_net(n)
329 :
330 688 : do k = 1,n_algae
331 516 : pN_sk(k) = trcr (n, nt_bgc_N(k))
332 688 : pflux_N(k) = flux_bio(n,nlt_bgc_N(k))*mps_to_cmpdy/c100
333 : enddo
334 172 : if (tr_bgc_C) then
335 516 : do k = 1,n_doc
336 516 : pDOC_sk(k) = trcr (n,nt_bgc_DOC(k))
337 : enddo
338 : endif
339 172 : if (tr_bgc_DON) then
340 344 : do k = 1,n_don
341 344 : pDON_sk(k) = trcr (n,nt_bgc_DON(k))
342 : enddo
343 : endif
344 172 : if (tr_bgc_Fe ) then
345 344 : do k = 1,n_fed
346 344 : pFed_sk (k) = trcr (n,nt_bgc_Fed(k))
347 : enddo
348 344 : do k = 1,n_fep
349 344 : pFep_sk (k) = trcr (n,nt_bgc_Fep(k))
350 : enddo
351 : endif
352 172 : if (tr_bgc_Nit) then
353 172 : pNit_sk = trcr (n, nt_bgc_Nit)
354 172 : pflux_NO = flux_bio(n,nlt_bgc_Nit)*mps_to_cmpdy/c100
355 : endif
356 172 : if (tr_bgc_Am) then
357 172 : pAm_sk = trcr (n, nt_bgc_Am)
358 172 : pflux_Am = flux_bio(n,nlt_bgc_Am)*mps_to_cmpdy/c100
359 : endif
360 172 : if (tr_bgc_Sil) then
361 172 : pSil_sk = trcr (n, nt_bgc_Sil)
362 : endif
363 172 : if (tr_bgc_hum) then
364 172 : phum_sk = trcr (n, nt_bgc_hum)
365 172 : pflux_hum = flux_bio(n,nlt_bgc_hum)*mps_to_cmpdy/c100
366 : endif
367 172 : if (tr_bgc_DMS) then
368 172 : pDMSPp_sk = trcr (n,nt_bgc_DMSPp)
369 172 : pDMSPd_sk = trcr (n,nt_bgc_DMSPd)
370 172 : pDMS_sk = trcr (n,nt_bgc_DMS)
371 : endif
372 :
373 3528 : elseif (z_tracers) then ! zbgc
374 3528 : pflux_NO = c0
375 3528 : pN_tot (:) = c0
376 3528 : pflux_Am = c0
377 3528 : pflux_hum = c0
378 3528 : pflux_atm_Am = c0
379 3528 : pflux_snow_Am = c0
380 3528 : pflux_N (:) = c0
381 3528 : pflux_NO = c0
382 3528 : pflux_atm_NO = c0
383 3528 : pflux_snow_NO = c0
384 3528 : pflux_zaero (:) = c0
385 3528 : pflux_atm_zaero_s(:) = c0
386 3528 : pflux_atm_zaero (:) = c0
387 3528 : pflux_snow_zaero (:) = c0
388 :
389 3528 : if (tr_bgc_Nit) then
390 3528 : pflux_NO = flux_bio(n,nlt_bgc_Nit)*mps_to_cmpdy/c100
391 3528 : pflux_atm_NO = fbio_atmice(n,nlt_bgc_Nit)*mps_to_cmpdy/c100
392 3528 : pflux_snow_NO = fbio_snoice(n,nlt_bgc_Nit)*mps_to_cmpdy/c100
393 : endif
394 3528 : if (tr_bgc_Am) then
395 3528 : pflux_Am = flux_bio(n,nlt_bgc_Am)*mps_to_cmpdy/c100
396 3528 : pflux_atm_Am = fbio_atmice(n,nlt_bgc_Am)*mps_to_cmpdy/c100
397 3528 : pflux_snow_Am = fbio_snoice(n,nlt_bgc_Am)*mps_to_cmpdy/c100
398 : endif
399 3528 : if (tr_bgc_hum) then
400 3528 : pflux_hum = flux_bio(n,nlt_bgc_hum)*mps_to_cmpdy/c100
401 : endif
402 3528 : if (tr_bgc_N) then
403 14112 : do k = 1,n_algae
404 14112 : pflux_N(k) = flux_bio(n,nlt_bgc_N(k))*mps_to_cmpdy/c100
405 : enddo
406 : endif
407 3528 : if (tr_zaero) then
408 23492 : do k = 1,n_zaero
409 20136 : pflux_zaero(k) = flux_bio(n,nlt_zaero(k))*mps_to_cmpdy/c100
410 20136 : pflux_atm_zaero_s(k)= flux_bio_atm(n,nlt_zaero(k))*mps_to_cmpdy/c100 !*aice
411 20136 : pflux_atm_zaero(k) = fbio_atmice(n,nlt_zaero(k))*mps_to_cmpdy/c100
412 23492 : pflux_snow_zaero(k) = fbio_snoice(n,nlt_zaero(k))*mps_to_cmpdy/c100
413 : enddo
414 : endif
415 31752 : do k = 1, nblyr+1
416 28224 : pzfswin(k) = c0
417 28224 : pZoo (k) = c0
418 169344 : do nn = 1,ncat
419 141120 : pzfswin(k) = pzfswin(k) + zfswin(n,k,nn)*vicen(n,nn)
420 169344 : pZoo (k) = pZoo (k) + Zoo (n,k,nn)*vicen(n,nn)
421 : enddo !nn
422 28224 : if (vice(n) > c0) then
423 21048 : pzfswin(k) = pzfswin(k)/vice(n)
424 21048 : pZoo (k) = pZoo (k)/vice(n)
425 : endif !vice
426 28224 : pAm (k ) = c0
427 112896 : pN (k,:) = c0
428 112896 : pDOC (k,:) = c0
429 56448 : pDON (k,:) = c0
430 84672 : pFed (k,:) = c0
431 84672 : pFep (k,:) = c0
432 197568 : pzaero(k,:) = c0
433 28224 : pPON (k ) = c0
434 28224 : phum (k ) = c0
435 28224 : pNO (k ) = c0
436 28224 : if (tr_bgc_Nit) pNO(k) = trcr(n,nt_bgc_Nit+k-1)
437 28224 : if (tr_bgc_Am) pAm(k) = trcr(n,nt_bgc_Am+k-1)
438 28224 : if (tr_bgc_N) then
439 112896 : do nn = 1, n_algae
440 112896 : pN (k,nn) = trcr(n,nt_bgc_N(nn)+k-1)
441 : enddo
442 : endif
443 28224 : if (tr_bgc_C) then
444 84672 : do nn = 1, n_doc
445 84672 : pDOC (k,nn) = trcr(n,nt_bgc_DOC(nn)+k-1)
446 : enddo
447 : endif
448 28224 : if (tr_bgc_DON) then
449 56448 : do nn = 1, n_don
450 56448 : pDON (k,nn) = trcr(n,nt_bgc_DON(nn)+k-1)
451 : enddo
452 : endif
453 28224 : if (tr_bgc_Fe) then
454 2752 : do nn = 1, n_fed
455 2752 : pFed (k,nn) = trcr(n,nt_bgc_Fed(nn)+k-1)
456 : enddo
457 2752 : do nn = 1, n_fep
458 2752 : pFep (k,nn) = trcr(n,nt_bgc_Fep(nn)+k-1)
459 : enddo
460 : endif
461 28224 : if (tr_zaero) then
462 187936 : do nn = 1, n_zaero
463 187936 : pzaero(k,nn) = trcr(n,nt_zaero(nn)+k-1)
464 : enddo
465 : endif
466 28224 : if (tr_bgc_PON) pPON(k) = trcr(n,nt_bgc_PON+k-1)
467 31752 : if (tr_bgc_hum) phum(k) = trcr(n,nt_bgc_hum+k-1)
468 : enddo ! k
469 3528 : if (tr_bgc_N) then
470 14112 : do nn = 1,n_algae
471 14112 : pN_tot(nn) = ice_bio_net(n,nlt_bgc_N(nn))
472 : enddo
473 3528 : pgrow_net = grow_net(n)
474 : endif ! tr_bgc_N
475 10584 : do k = 1 , 2 ! snow concentration
476 7056 : pAms (k ) = c0
477 28224 : pNs (k,:) = c0
478 28224 : pDOCs (k,:) = c0
479 14112 : pDONs (k,:) = c0
480 21168 : pFeds (k,:)= c0
481 21168 : pFeps (k,:)= c0
482 49392 : pzaeros(k,:) = c0
483 7056 : pPONs (k ) = c0
484 7056 : phums (k ) = c0
485 7056 : pNOs (k ) = c0
486 7056 : if (tr_bgc_Nit) pNOs(k) = trcr(n,nt_bgc_Nit+nblyr+k)
487 7056 : if (tr_bgc_Am) pAms(k) = trcr(n,nt_bgc_Am +nblyr+k)
488 7056 : if (tr_bgc_N) then
489 28224 : do nn = 1, n_algae
490 28224 : pNs (k,nn) = trcr(n,nt_bgc_N(nn)+nblyr+k)
491 : enddo
492 : endif
493 7056 : if (tr_bgc_C) then
494 21168 : do nn = 1, n_doc
495 21168 : pDOCs (k,nn) = trcr(n,nt_bgc_DOC(nn)+nblyr+k)
496 : enddo
497 : endif
498 7056 : if (tr_bgc_DON) then
499 14112 : do nn = 1, n_don
500 14112 : pDONs (k,nn) = trcr(n,nt_bgc_DON(nn)+nblyr+k)
501 : enddo
502 : endif
503 7056 : if (tr_bgc_Fe ) then
504 688 : do nn = 1, n_fed
505 688 : pFeds (k,nn) = trcr(n,nt_bgc_Fed(nn)+nblyr+k)
506 : enddo
507 688 : do nn = 1, n_fep
508 688 : pFeps (k,nn) = trcr(n,nt_bgc_Fep(nn)+nblyr+k)
509 : enddo
510 : endif
511 7056 : if (tr_zaero) then
512 46984 : do nn = 1, n_zaero
513 46984 : pzaeros(k,nn) = trcr(n,nt_zaero(nn)+nblyr+k)
514 : enddo
515 : endif
516 7056 : if (tr_bgc_PON)pPONs(k) =trcr(n,nt_bgc_PON+nblyr+k)
517 10584 : if (tr_bgc_hum)phums(k) =trcr(n,nt_bgc_hum+nblyr+k)
518 : enddo ! k
519 : endif ! ztracers
520 3700 : pchlsw(:) = c0
521 3700 : pzaerosw(:,:) = c0
522 3700 : if (dEdd_algae) then
523 0 : do k = 0, klev
524 0 : if (tr_bgc_N) pchlsw(k+1) = trcrn_sw(n,nlt_chl_sw+k,1)
525 0 : if (tr_zaero) then
526 0 : do nn = 1,n_zaero
527 0 : pzaerosw(k+1,nn) = trcrn_sw(n,nlt_zaero_sw(nn) + k,1)
528 : enddo
529 : endif
530 : enddo
531 : endif ! dEdd_algae
532 :
533 : !-----------------------------------------------------------------
534 : ! start spewing
535 : !-----------------------------------------------------------------
536 :
537 3700 : if (z_tracers) then
538 3528 : write(nu_diag_out+n-1,803) 'zfswin PAR '
539 3528 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
540 31752 : write(nu_diag_out+n-1,802) (pzfswin(k), k = 1,nblyr+1)
541 3528 : write(nu_diag_out+n-1,*) ' '
542 3528 : write(nu_diag_out+n-1,803) 'Losses: Zoo(mmol/m^3) '
543 3528 : write(nu_diag_out+n-1,803) ' Brine Conc. ',' Brine Conc'
544 3528 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
545 31752 : write(nu_diag_out+n-1,802) (pZoo(k), k = 1,nblyr+1)
546 3528 : write(nu_diag_out+n-1,*) ' '
547 : endif
548 3700 : if (tr_bgc_Nit) then
549 3700 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
550 3700 : write(nu_diag_out+n-1,*) ' nitrate conc. (mmol/m^3) or flux (mmol/m^2/d)'
551 3700 : write(nu_diag_out+n-1,900) 'Ocean conc = ',pNit_ac
552 3700 : write(nu_diag_out+n-1,900) 'ice-ocean flux = ',pflux_NO
553 3700 : if (skl_bgc) then
554 172 : write(nu_diag_out+n-1,900) 'Bulk ice conc. = ',pNit_sk
555 3528 : elseif (z_tracers) then
556 3528 : write(nu_diag_out+n-1,900) 'atm-ice flux = ',pflux_atm_NO
557 3528 : write(nu_diag_out+n-1,900) 'snow-ice flux = ',pflux_snow_NO
558 3528 : write(nu_diag_out+n-1,*) ' snow + ice conc'
559 3528 : write(nu_diag_out+n-1,803) ' nitrate'
560 10584 : write(nu_diag_out+n-1,802) (pNOs(k), k = 1,2)
561 31752 : write(nu_diag_out+n-1,802) (pNO(k), k = 1,nblyr+1)
562 3528 : write(nu_diag_out+n-1,*) ' '
563 : endif
564 : endif
565 3700 : if (tr_bgc_PON .and. z_tracers) then
566 3528 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
567 3528 : write(nu_diag_out+n-1,*) ' PON snow + ice conc. (mmol/m^3)'
568 3528 : write(nu_diag_out+n-1,803) ' PON'
569 10584 : write(nu_diag_out+n-1,802) (pPONs(k), k = 1,2)
570 31752 : write(nu_diag_out+n-1,802) (pPON(k), k = 1,nblyr+1)
571 3528 : write(nu_diag_out+n-1,*) ' '
572 : endif
573 3700 : if (tr_bgc_hum) then
574 3700 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
575 3700 : write(nu_diag_out+n-1,*) ' hum snow + ice conc. (mmolC/m^3)'
576 3700 : write(nu_diag_out+n-1,900) 'Ocean conc = ',phum_ac
577 3700 : write(nu_diag_out+n-1,900) 'ice-ocean flux = ',pflux_hum
578 3700 : if (skl_bgc) then
579 172 : write(nu_diag_out+n-1,900) 'Bulk ice conc. = ',phum_sk
580 3528 : elseif (z_tracers) then
581 3528 : write(nu_diag_out+n-1,803) ' hum'
582 10584 : write(nu_diag_out+n-1,802) (phums(k), k = 1,2)
583 31752 : write(nu_diag_out+n-1,802) (phum(k), k = 1,nblyr+1)
584 3528 : write(nu_diag_out+n-1,*) ' '
585 : endif
586 : endif
587 3700 : if (tr_bgc_Am) then
588 3700 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
589 3700 : write(nu_diag_out+n-1,*) ' ammonium conc. (mmol/m^3) or flux (mmol/m^2/d)'
590 3700 : write(nu_diag_out+n-1,900) 'Ocean conc = ',pAm_ac
591 3700 : write(nu_diag_out+n-1,900) 'ice-ocean flux = ',pflux_Am
592 3700 : if (skl_bgc) then
593 172 : write(nu_diag_out+n-1,900) 'Bulk ice conc. = ',pAm_sk
594 3528 : elseif (z_tracers) then
595 3528 : write(nu_diag_out+n-1,900) 'atm-ice flux = ',pflux_atm_Am
596 3528 : write(nu_diag_out+n-1,900) 'snow-ice flux = ',pflux_snow_Am
597 3528 : write(nu_diag_out+n-1,*) ' snow + ice conc.'
598 3528 : write(nu_diag_out+n-1,803) ' ammonium'
599 10584 : write(nu_diag_out+n-1,802) (pAms(k), k = 1,2)
600 31752 : write(nu_diag_out+n-1,802) (pAm(k), k = 1,nblyr+1)
601 3528 : write(nu_diag_out+n-1,*) ' '
602 : endif
603 : endif
604 3700 : if (tr_bgc_N) then
605 3700 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
606 3700 : write(nu_diag_out+n-1,900) 'tot algal growth (1/d) = ',pgrow_net
607 14800 : do kk = 1,n_algae
608 11100 : write(nu_diag_out+n-1,*) ' algal conc. (mmol N/m^3) or flux (mmol N/m^2/d)'
609 11100 : write(nu_diag_out+n-1,1020) ' type:', kk
610 11100 : write(nu_diag_out+n-1,900) 'Ocean conc = ',pN_ac(kk)
611 11100 : write(nu_diag_out+n-1,900) 'ice-ocean flux = ',pflux_N(kk)
612 14800 : if (skl_bgc) then
613 516 : write(nu_diag_out+n-1,900) 'Bulk ice conc. = ',pN_sk(kk)
614 10584 : elseif (z_tracers) then
615 10584 : write(nu_diag_out+n-1,900) 'Tot ice (mmolN/m^2) = ',pN_tot(kk)
616 10584 : write(nu_diag_out+n-1,*) ' snow + ice conc.'
617 10584 : write(nu_diag_out+n-1,803) ' algal N'
618 31752 : write(nu_diag_out+n-1,802) (pNs(k,kk), k = 1,2)
619 95256 : write(nu_diag_out+n-1,802) (pN(k,kk), k = 1,nblyr+1)
620 10584 : write(nu_diag_out+n-1,*) ' '
621 : endif
622 : enddo
623 : endif
624 3700 : if (tr_bgc_C) then
625 7400 : do kk = 1,1 !n_doc
626 3700 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
627 3700 : write(nu_diag_out+n-1,*) ' DOC conc. (mmol C/m^3)'
628 3700 : write(nu_diag_out+n-1,1020) ' type:', kk
629 3700 : write(nu_diag_out+n-1,900) 'Ocean conc = ',pDOC_ac(kk)
630 7400 : if (skl_bgc) then
631 172 : write(nu_diag_out+n-1,900)'Bulk ice conc. = ',pDOC_sk(kk)
632 3528 : elseif (z_tracers) then
633 3528 : write(nu_diag_out+n-1,*) ' snow + ice conc.'
634 3528 : write(nu_diag_out+n-1,803) ' DOC'
635 10584 : write(nu_diag_out+n-1,802) (pDOCs(k,kk), k = 1,2)
636 31752 : write(nu_diag_out+n-1,802) (pDOC(k,kk), k = 1,nblyr+1)
637 3528 : write(nu_diag_out+n-1,*) ' '
638 : endif
639 : enddo
640 : endif
641 3700 : if (tr_bgc_DON) then
642 7400 : do kk = 1,n_don
643 3700 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
644 3700 : write(nu_diag_out+n-1,*) ' DON conc. (mmol N/m^3)'
645 3700 : write(nu_diag_out+n-1,1020) ' type:', kk
646 3700 : write(nu_diag_out+n-1,900) 'Ocean conc = ',pDON_ac(kk)
647 7400 : if (skl_bgc) then
648 172 : write(nu_diag_out+n-1,900)'Bulk ice conc. = ',pDON_sk(kk)
649 3528 : elseif (z_tracers) then
650 3528 : write(nu_diag_out+n-1,*) ' snow + ice conc.'
651 3528 : write(nu_diag_out+n-1,803) ' DON'
652 10584 : write(nu_diag_out+n-1,802) (pDONs(k,kk), k = 1,2)
653 31752 : write(nu_diag_out+n-1,802) (pDON(k,kk), k = 1,nblyr+1)
654 3528 : write(nu_diag_out+n-1,*) ' '
655 : endif
656 : enddo
657 : endif
658 3700 : if (tr_bgc_Fe ) then
659 688 : do kk = 1,n_fed
660 344 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
661 344 : write(nu_diag_out+n-1,*) ' dFe conc. (nM)'
662 344 : write(nu_diag_out+n-1,1020) ' type:', kk
663 344 : write(nu_diag_out+n-1,900) 'Ocean conc = ',pFed_ac (kk)
664 688 : if (skl_bgc) then
665 172 : write(nu_diag_out+n-1,900)'Bulk ice conc. = ',pFed_sk (kk)
666 172 : elseif (z_tracers) then
667 172 : write(nu_diag_out+n-1,*) ' snow + ice conc.'
668 172 : write(nu_diag_out+n-1,803) ' Fed'
669 516 : write(nu_diag_out+n-1,802) (pFeds (k,kk), k = 1,2)
670 1548 : write(nu_diag_out+n-1,802) (pFed (k,kk), k = 1,nblyr+1)
671 172 : write(nu_diag_out+n-1,*) ' '
672 : endif
673 : enddo
674 688 : do kk = 1,n_fep
675 344 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
676 344 : write(nu_diag_out+n-1,*) ' pFe conc. (nM)'
677 344 : write(nu_diag_out+n-1,1020) ' type:', kk
678 344 : write(nu_diag_out+n-1,900) 'Ocean conc = ',pFep_ac (kk)
679 688 : if (skl_bgc) then
680 172 : write(nu_diag_out+n-1,900)'Bulk ice conc. = ',pFep_sk (kk)
681 172 : elseif (z_tracers) then
682 172 : write(nu_diag_out+n-1,*) ' snow + ice conc.'
683 172 : write(nu_diag_out+n-1,803) ' Fep'
684 516 : write(nu_diag_out+n-1,802) (pFeps (k,kk), k = 1,2)
685 1548 : write(nu_diag_out+n-1,802) (pFep (k,kk), k = 1,nblyr+1)
686 172 : write(nu_diag_out+n-1,*) ' '
687 : endif
688 : enddo
689 : endif
690 3700 : if (tr_bgc_DMS) then
691 3700 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
692 3700 : write(nu_diag_out+n-1,*) ' DMS (mmol/m^3) '
693 3700 : write(nu_diag_out+n-1,900) 'Ocean DMSP = ',pDMSP_ac
694 3700 : write(nu_diag_out+n-1,900) 'Ocean DMS = ',pDMS_ac
695 3700 : if (skl_bgc) then
696 172 : write(nu_diag_out+n-1,900) 'Ice DMSPp = ',pDMSPp_sk
697 172 : write(nu_diag_out+n-1,900) 'Ice DMSPd = ',pDMSPd_sk
698 172 : write(nu_diag_out+n-1,900) 'Ice DMS = ',pDMS_sk
699 : endif
700 : endif
701 3700 : if (tr_zaero .and. z_tracers) then
702 23492 : do kk = 1,n_zaero
703 20136 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
704 20136 : write(nu_diag_out+n-1,*) ' aerosol conc. (kg/m^3) or flux (kg/m^2/d)'
705 20136 : write(nu_diag_out+n-1,1020) ' type: ',kk
706 20136 : write(nu_diag_out+n-1,900) 'Atm source flux = ',pflux_atm_zaero_s(kk)
707 20136 : write(nu_diag_out+n-1,900) 'ice-ocean flux*aice = ',pflux_zaero(kk)
708 20136 : write(nu_diag_out+n-1,900) 'atm-ice flux*aice = ',pflux_atm_zaero(kk)
709 20136 : write(nu_diag_out+n-1,900) 'snow-ice flux*aice = ',pflux_snow_zaero(kk)
710 20136 : write(nu_diag_out+n-1,*) ' snow + ice conc.'
711 20136 : write(nu_diag_out+n-1,803) ' aerosol'
712 60408 : write(nu_diag_out+n-1,802) (pzaeros(k,kk), k = 1,2)
713 181224 : write(nu_diag_out+n-1,802) (pzaero(k,kk), k = 1,nblyr+1)
714 23492 : write(nu_diag_out+n-1,*) ' '
715 : enddo
716 : endif
717 4625 : if (dEdd_algae) then
718 0 : if (tr_zaero) then
719 0 : do kk = 1,n_zaero
720 0 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
721 0 : write(nu_diag_out+n-1,*) ' Cat 1 aerosol conc. (kg/m^3) on delta-Eddington grid '
722 0 : write(nu_diag_out+n-1,802) (pzaerosw(k,kk), k = 1,klev +1)
723 : enddo
724 : endif
725 0 : if (tr_bgc_N) then
726 0 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
727 0 : write(nu_diag_out+n-1,*) ' Cat 1 chl (mg/m^3) on delta-Eddington grid '
728 0 : write(nu_diag_out+n-1,802) (pchlsw(k), k = 1,klev +1)
729 : endif
730 : endif
731 : enddo ! nx
732 :
733 : 802 format (f24.17,2x)
734 : 803 format (a25,2x)
735 : 900 format (a25,2x,f24.17)
736 : 1020 format (a30,2x,i6) ! integer
737 :
738 925 : end subroutine bgc_diags
739 :
740 : !=======================================================================
741 : !
742 : ! Writes diagnostic info (max, min, global sums, etc) to standard out
743 : !
744 : ! authors: Nicole Jeffery, LANL
745 :
746 0 : subroutine zsal_diags ()
747 :
748 : use icedrv_arrays_column, only: fzsal, fzsal_g, sice_rho, bTiz
749 : use icedrv_arrays_column, only: iDi, bphi, dhbr_top, dhbr_bot, darcy_V
750 : use icedrv_domain_size, only: nblyr, ncat, nilyr
751 : use icedrv_state, only: aicen, aice, vice, trcr, trcrn, vicen, vsno
752 :
753 : ! local variables
754 :
755 : integer (kind=int_kind) :: &
756 : k, n, nn
757 :
758 : ! fields at diagnostic points
759 : real (kind=dbl_kind), dimension(nx) :: &
760 0 : phinS, phinS1, &
761 0 : phbrn,pdh_top1,pdh_bot1, psice_rho, pfzsal, &
762 0 : pfzsal_g, pdarcy_V1
763 :
764 : ! vertical fields of category 1 at diagnostic points for bgc layer model
765 : real (kind=dbl_kind), dimension(nblyr+2) :: &
766 0 : pphin, pphin1
767 : real (kind=dbl_kind), dimension(nblyr) :: &
768 0 : pSin, pSice, pSin1
769 : real (kind=dbl_kind), dimension(nblyr+1) :: &
770 0 : pbTiz, piDin
771 :
772 : real (kind=dbl_kind) :: &
773 0 : rhosi, rhow, rhos
774 :
775 : logical (kind=log_kind) :: tr_brine
776 : integer (kind=int_kind) :: nt_fbri, nt_bgc_S, nt_sice
777 :
778 : character(len=*), parameter :: subname='(zsal_diags)'
779 :
780 : !-----------------------------------------------------------------
781 : ! query Icepack values
782 : !-----------------------------------------------------------------
783 :
784 0 : call icepack_query_parameters(rhosi_out=rhosi, rhow_out=rhow, rhos_out=rhos)
785 0 : call icepack_query_tracer_flags(tr_brine_out=tr_brine)
786 : call icepack_query_tracer_indices(nt_fbri_out=nt_fbri, &
787 0 : nt_bgc_S_out=nt_bgc_S, nt_sice_out=nt_sice)
788 0 : call icepack_warnings_flush(nu_diag)
789 0 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
790 0 : file=__FILE__,line= __LINE__)
791 :
792 : !-----------------------------------------------------------------
793 : ! salinity and microstructure of the ice
794 : !-----------------------------------------------------------------
795 : ! NOTE these are computed for the last timestep only (not avg)
796 :
797 0 : do n = 1, nx
798 0 : pfzsal = fzsal(n)
799 0 : pfzsal_g = fzsal_g(n)
800 0 : phinS = c0
801 0 : phinS1 = c0
802 0 : phbrn = c0
803 0 : psice_rho = c0
804 0 : pdh_top1 = c0
805 0 : pdh_bot1 = c0
806 0 : pdarcy_V1 = c0
807 0 : do nn = 1, ncat
808 0 : psice_rho = psice_rho(n) + sice_rho(n,nn)*aicen(n,nn)
809 : enddo
810 0 : if (aice(n) > c0) psice_rho = psice_rho/aice(n)
811 0 : if (tr_brine .and. aice(n) > c0) &
812 0 : phinS = trcr(n,nt_fbri)*vice(n)/aice(n)
813 0 : if (aicen(n,1) > c0) then
814 0 : if (tr_brine) phinS1 = trcrn(n,nt_fbri,1) &
815 0 : * vicen(n,1)/aicen(n,1)
816 0 : pdh_top1 = dhbr_top(n,1)
817 0 : pdh_bot1 = dhbr_bot(n,1)
818 0 : pdarcy_V1 = darcy_V(n,1)
819 : endif
820 0 : if (tr_brine .AND. aice(n) > c0) &
821 0 : phbrn = (c1 - rhosi/rhow)*vice(n)/aice(n) &
822 0 : - rhos /rhow *vsno(n)/aice(n)
823 0 : do k = 1, nblyr+1
824 0 : pbTiz(k) = c0
825 0 : piDin(k) = c0
826 0 : do nn = 1,ncat
827 0 : pbTiz(k) = pbTiz(k) + bTiz(n,k,nn)*vicen(n,nn)
828 0 : piDin(k) = piDin(k) + iDi(n,k,nn)*vicen(n,nn)
829 : enddo
830 0 : if (vice(n) > c0) then
831 0 : pbTiz(k) = pbTiz(k)/vice(n)
832 0 : piDin(k) = piDin(k)/vice(n)
833 : endif
834 : enddo ! k
835 0 : do k = 1, nblyr+2
836 0 : pphin (k) = c0
837 0 : pphin1(k) = c0
838 0 : if (aicen(n,1) > c0) pphin1(k) = bphi(n,k,1)
839 0 : do nn = 1,ncat
840 0 : pphin(k) = pphin(k) + bphi(n,k,nn)*vicen(n,nn)
841 : enddo
842 0 : if (vice(n) > c0) pphin(k) = pphin(k)/vice(n)
843 : enddo
844 0 : do k = 1,nblyr
845 0 : pSin (k) = c0
846 0 : pSin1(k) = c0
847 0 : pSin (k) = trcr(n,nt_bgc_S+k-1)
848 0 : if (aicen(n,1) > c0) pSin1(k) = trcrn(n,nt_bgc_S+k-1,1)
849 : enddo
850 0 : do k = 1,nilyr
851 0 : pSice(k) = trcr(n,nt_sice+k-1)
852 : enddo
853 : enddo ! nx
854 :
855 : !-----------------------------------------------------------------
856 : ! start spewing
857 : !-----------------------------------------------------------------
858 :
859 0 : write(nu_diag_out+n-1,*) ' '
860 0 : write(nu_diag_out+n-1,*) ' Brine height '
861 0 : write(nu_diag_out+n-1,900) 'hbrin = ',phinS
862 0 : write(nu_diag_out+n-1,900) 'hbrin cat 1 = ',phinS1
863 0 : write(nu_diag_out+n-1,900) 'Freeboard = ',phbrn
864 0 : write(nu_diag_out+n-1,900) 'dhbrin cat 1 top = ',pdh_top1
865 0 : write(nu_diag_out+n-1,900) 'dhbrin cat 1 bottom = ',pdh_bot1
866 0 : write(nu_diag_out+n-1,*) ' '
867 0 : write(nu_diag_out+n-1,*) ' zSalinity '
868 0 : write(nu_diag_out+n-1,900) 'Avg density (kg/m^3) = ',psice_rho
869 0 : write(nu_diag_out+n-1,900) 'Salt flux (kg/m^2/s) = ',pfzsal
870 0 : write(nu_diag_out+n-1,900) 'Grav. Drain. Salt flux = ',pfzsal_g
871 0 : write(nu_diag_out+n-1,900) 'Darcy V cat 1 (m/s) = ',pdarcy_V1
872 0 : write(nu_diag_out+n-1,*) ' '
873 0 : write(nu_diag_out+n-1,*) ' Top down bgc Layer Model'
874 0 : write(nu_diag_out+n-1,*) ' '
875 0 : write(nu_diag_out+n-1,803) 'bTiz ice temp'
876 0 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
877 0 : write(nu_diag_out+n-1,802) (pbTiz(k), k = 1,nblyr+1)
878 0 : write(nu_diag_out+n-1,*) ' '
879 0 : write(nu_diag_out+n-1,803) 'iDi diffusivity'
880 0 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
881 0 : write(nu_diag_out+n-1,802) (piDin(k), k = 1,nblyr+1)
882 0 : write(nu_diag_out+n-1,*) ' '
883 0 : write(nu_diag_out+n-1,803) 'bphi porosity'
884 0 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
885 0 : write(nu_diag_out+n-1,802) (pphin(k), k = 1,nblyr+1)
886 0 : write(nu_diag_out+n-1,*) ' '
887 0 : write(nu_diag_out+n-1,803) 'phi1 porosity'
888 0 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
889 0 : write(nu_diag_out+n-1,802) (pphin1(k), k = 1,nblyr+1)
890 0 : write(nu_diag_out+n-1,*) ' '
891 0 : write(nu_diag_out+n-1,803) 'zsal cat 1'
892 0 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
893 0 : write(nu_diag_out+n-1,802) (pSin1(k), k = 1,nblyr)
894 0 : write(nu_diag_out+n-1,*) ' '
895 0 : write(nu_diag_out+n-1,803) 'zsal Avg S'
896 0 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
897 0 : write(nu_diag_out+n-1,802) (pSin(k), k = 1,nblyr)
898 0 : write(nu_diag_out+n-1,*) ' '
899 0 : write(nu_diag_out+n-1,803) 'Sice Ice S'
900 0 : write(nu_diag_out+n-1,*) '---------------------------------------------------'
901 0 : write(nu_diag_out+n-1,802) (pSice(k), k = 1,nilyr)
902 0 : write(nu_diag_out+n-1,*) ' '
903 :
904 : 802 format (f24.17,2x)
905 : 803 format (a25,2x)
906 : 900 format (a25,2x,f24.17)
907 :
908 0 : end subroutine zsal_diags
909 :
910 : !=======================================================================
911 :
912 : end module icedrv_diagnostics_bgc
913 :
914 : !=======================================================================
|