Line data Source code
1 : !=======================================================================
2 :
3 : ! Floe size distribution history output
4 : !
5 : ! 2019 Elizabeth Hunke
6 :
7 : module ice_history_fsd
8 :
9 : use ice_kinds_mod
10 : use ice_domain_size, only: max_nstrm, nfsd
11 : use ice_constants, only: c0, c1, c8, c100, mps_to_cmpdy
12 : use ice_fileunits, only: nu_nml, nml_filename, &
13 : get_fileunit, release_fileunit
14 : use ice_fileunits, only: nu_diag
15 : use ice_exit, only: abort_ice
16 : use icepack_intfc, only: icepack_warnings_flush, icepack_warnings_aborted
17 : use icepack_intfc, only: icepack_query_parameters, &
18 : icepack_query_tracer_flags, icepack_query_tracer_indices
19 :
20 : implicit none
21 : private
22 : public :: accum_hist_fsd, init_hist_fsd_2D, init_hist_fsd_3Df, &
23 : init_hist_fsd_4Df
24 :
25 : !---------------------------------------------------------------
26 : ! flags: write to output file if true or histfreq value
27 : !---------------------------------------------------------------
28 :
29 : character (len=max_nstrm), public :: &
30 : f_afsd = 'm', f_afsdn = 'm', & ! LCOV_EXCL_LINE
31 : f_dafsd_newi = 'm', f_dafsd_latg = 'm', & ! LCOV_EXCL_LINE
32 : f_dafsd_latm = 'm', f_dafsd_wave = 'm', & ! LCOV_EXCL_LINE
33 : f_dafsd_weld = 'm', f_wave_sig_ht = 'm', & ! LCOV_EXCL_LINE
34 : f_aice_ww = 'x', f_diam_ww = 'x', & ! LCOV_EXCL_LINE
35 : f_hice_ww = 'x', f_fsdrad = 'x', & ! LCOV_EXCL_LINE
36 : f_fsdperim = 'x'
37 :
38 : !---------------------------------------------------------------
39 : ! namelist variables
40 : !---------------------------------------------------------------
41 :
42 : namelist / icefields_fsd_nml / &
43 : f_afsd , f_afsdn , & ! LCOV_EXCL_LINE
44 : f_dafsd_newi, f_dafsd_latg , & ! LCOV_EXCL_LINE
45 : f_dafsd_latm, f_dafsd_wave , & ! LCOV_EXCL_LINE
46 : f_dafsd_weld, f_wave_sig_ht, & ! LCOV_EXCL_LINE
47 : f_aice_ww , f_diam_ww , & ! LCOV_EXCL_LINE
48 : f_hice_ww , f_fsdrad , & ! LCOV_EXCL_LINE
49 : f_fsdperim
50 :
51 : !---------------------------------------------------------------
52 : ! field indices
53 : !---------------------------------------------------------------
54 :
55 : integer (kind=int_kind), dimension(max_nstrm) :: &
56 : n_afsd , n_afsdn , & ! LCOV_EXCL_LINE
57 : n_dafsd_newi, n_dafsd_latg , & ! LCOV_EXCL_LINE
58 : n_dafsd_latm, n_dafsd_wave , & ! LCOV_EXCL_LINE
59 : n_dafsd_weld, n_wave_sig_ht, & ! LCOV_EXCL_LINE
60 : n_aice_ww , n_diam_ww , & ! LCOV_EXCL_LINE
61 : n_hice_ww , n_fsdrad , & ! LCOV_EXCL_LINE
62 : n_fsdperim
63 :
64 : !=======================================================================
65 :
66 : contains
67 :
68 : !=======================================================================
69 :
70 : ! Initialize history files
71 : ! authors Elizabeth C. Hunke, LANL
72 :
73 37 : subroutine init_hist_fsd_2D
74 :
75 : use ice_broadcast, only: broadcast_scalar
76 : use ice_calendar, only: nstreams
77 : use ice_communicate, only: my_task, master_task
78 : use ice_history_shared, only: tstr2D, tcstr, define_hist_field
79 : use ice_fileunits, only: goto_nml
80 :
81 : integer (kind=int_kind) :: ns
82 : integer (kind=int_kind) :: nml_error ! namelist i/o error flag
83 : logical (kind=log_kind) :: tr_fsd, wave_spec
84 : character (len=char_len_long) :: tmpstr2 ! test namelist
85 : character(len=char_len) :: nml_name ! text namelist name
86 :
87 : character(len=*), parameter :: subname = '(init_hist_fsd_2D)'
88 :
89 37 : call icepack_query_tracer_flags(tr_fsd_out=tr_fsd)
90 37 : call icepack_query_parameters(wave_spec_out=wave_spec)
91 37 : call icepack_warnings_flush(nu_diag)
92 37 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
93 0 : file=__FILE__, line=__LINE__)
94 :
95 37 : if (tr_fsd) then
96 :
97 : !-----------------------------------------------------------------
98 : ! read namelist
99 : !-----------------------------------------------------------------
100 :
101 0 : if (my_task == master_task) then
102 0 : nml_name = 'icefields_fsd_nml'
103 0 : write(nu_diag,*) subname,' Reading ', trim(nml_name)
104 :
105 : ! open namelist file
106 0 : call get_fileunit(nu_nml)
107 0 : open (nu_nml, file=trim(nml_filename), status='old',iostat=nml_error)
108 0 : if (nml_error /= 0) then
109 : call abort_ice(subname//'ERROR: '//trim(nml_name)//' open file '// &
110 : trim(nml_filename), & ! LCOV_EXCL_LINE
111 0 : file=__FILE__, line=__LINE__)
112 : endif
113 :
114 : ! goto this namelist
115 0 : call goto_nml(nu_nml,trim(nml_name),nml_error)
116 0 : if (nml_error /= 0) then
117 : call abort_ice(subname//'ERROR: searching for '// trim(nml_name), &
118 0 : file=__FILE__, line=__LINE__)
119 : endif
120 :
121 : ! read namelist
122 0 : nml_error = 1
123 0 : do while (nml_error > 0)
124 0 : read(nu_nml, nml=icefields_fsd_nml,iostat=nml_error)
125 : ! check if error
126 0 : if (nml_error /= 0) then
127 : ! backspace and re-read erroneous line
128 0 : backspace(nu_nml)
129 0 : read(nu_nml,fmt='(A)') tmpstr2
130 : call abort_ice(subname//'ERROR: ' // trim(nml_name) // ' reading ' // &
131 0 : trim(tmpstr2), file=__FILE__, line=__LINE__)
132 : endif
133 : end do
134 :
135 0 : close(nu_nml)
136 0 : call release_fileunit(nu_nml)
137 : endif
138 :
139 0 : call broadcast_scalar (f_afsd, master_task)
140 0 : call broadcast_scalar (f_afsdn, master_task)
141 0 : call broadcast_scalar (f_dafsd_newi, master_task)
142 0 : call broadcast_scalar (f_dafsd_latg, master_task)
143 0 : call broadcast_scalar (f_dafsd_latm, master_task)
144 0 : call broadcast_scalar (f_dafsd_wave, master_task)
145 0 : call broadcast_scalar (f_dafsd_weld, master_task)
146 0 : call broadcast_scalar (f_wave_sig_ht, master_task)
147 0 : call broadcast_scalar (f_aice_ww, master_task)
148 0 : call broadcast_scalar (f_diam_ww, master_task)
149 0 : call broadcast_scalar (f_hice_ww, master_task)
150 0 : call broadcast_scalar (f_fsdrad, master_task)
151 0 : call broadcast_scalar (f_fsdperim, master_task)
152 :
153 : ! 2D variables
154 :
155 0 : do ns = 1, nstreams
156 :
157 0 : if (f_wave_sig_ht(1:1) /= 'x') &
158 : call define_hist_field(n_wave_sig_ht,"wave_sig_ht","m",tstr2D, tcstr, & ! LCOV_EXCL_LINE
159 : "significant height of wind and swell waves", & ! LCOV_EXCL_LINE
160 : "from attenuated spectrum in ice", c1, c0, & ! LCOV_EXCL_LINE
161 0 : ns, f_wave_sig_ht)
162 0 : if (f_aice_ww(1:1) /= 'x') &
163 : call define_hist_field(n_aice_ww,"aice_ww","1",tstr2D, tcstr, & ! LCOV_EXCL_LINE
164 : "Concentration of floes > Dmin", & ! LCOV_EXCL_LINE
165 : "for waves", c1, c0, & ! LCOV_EXCL_LINE
166 0 : ns, f_aice_ww)
167 0 : if (f_diam_ww(1:1) /= 'x') &
168 : call define_hist_field(n_diam_ww,"diam_ww","m",tstr2D, tcstr, & ! LCOV_EXCL_LINE
169 : "Average (number) diameter of floes > Dmin", & ! LCOV_EXCL_LINE
170 : "for waves", c1, c0, & ! LCOV_EXCL_LINE
171 0 : ns, f_diam_ww)
172 0 : if (f_hice_ww(1:1) /= 'x') &
173 : call define_hist_field(n_hice_ww,"hice_ww","m",tstr2D, tcstr, & ! LCOV_EXCL_LINE
174 : "Thickness of floes > Dmin", & ! LCOV_EXCL_LINE
175 : "for waves", c1, c0, & ! LCOV_EXCL_LINE
176 0 : ns, f_hice_ww)
177 0 : if (f_fsdrad(1:1) /= 'x') &
178 : call define_hist_field(n_fsdrad,"fsdrad","m",tstr2D, tcstr, & ! LCOV_EXCL_LINE
179 : "floe size distribution, representative radius", & ! LCOV_EXCL_LINE
180 0 : " ", c1, c0, ns, f_fsdrad)
181 0 : if (f_fsdperim(1:1) /= 'x') &
182 : call define_hist_field(n_fsdperim,"fsdperim","1/m",tstr2D, tcstr, & ! LCOV_EXCL_LINE
183 : "floe size distribution, perimeter", & ! LCOV_EXCL_LINE
184 0 : "per unit ice area", c1, c0, ns, f_fsdperim)
185 :
186 :
187 : enddo ! nstreams
188 :
189 : else ! tr_fsd
190 :
191 37 : f_afsd = 'x'
192 37 : f_afsdn = 'x'
193 37 : f_dafsd_newi = 'x'
194 37 : f_dafsd_latg = 'x'
195 37 : f_dafsd_latm = 'x'
196 37 : f_dafsd_wave = 'x'
197 37 : f_dafsd_weld = 'x'
198 37 : f_wave_sig_ht = 'x'
199 37 : f_fsdrad = 'x'
200 37 : f_fsdperim = 'x'
201 37 : if (.not. wave_spec) then
202 37 : f_aice_ww = 'x'
203 37 : f_diam_ww = 'x'
204 37 : f_hice_ww = 'x'
205 : endif
206 :
207 : endif ! tr_fsd
208 :
209 37 : end subroutine init_hist_fsd_2D
210 :
211 : !=======================================================================
212 :
213 37 : subroutine init_hist_fsd_3Df
214 :
215 : use ice_calendar, only: nstreams, histfreq
216 : use ice_history_shared, only: tstr3Df, tcstr, define_hist_field
217 :
218 : logical (kind=log_kind) :: tr_fsd
219 : integer (kind=int_kind) :: ns
220 : character(len=*), parameter :: subname = '(init_hist_fsd_3Df)'
221 :
222 37 : call icepack_query_tracer_flags(tr_fsd_out=tr_fsd)
223 37 : call icepack_warnings_flush(nu_diag)
224 37 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
225 0 : file=__FILE__, line=__LINE__)
226 :
227 : !-----------------------------------------------------------------
228 : ! 3D (category) variables must be looped separately
229 : !-----------------------------------------------------------------
230 :
231 37 : if (tr_fsd) then
232 :
233 0 : do ns = 1, nstreams
234 0 : if (histfreq(ns) /= 'x') then
235 :
236 0 : if (f_afsd(1:1) /= 'x') &
237 : call define_hist_field(n_afsd,"afsd", "1/m", tstr3Df, tcstr, & ! LCOV_EXCL_LINE
238 : "areal floe size distribution", & ! LCOV_EXCL_LINE
239 0 : "per unit bin width ", c1, c0, ns, f_afsd)
240 0 : if (f_dafsd_newi(1:1) /= 'x') &
241 : call define_hist_field(n_dafsd_newi,"dafsd_newi","1/s",tstr3Df, tcstr, & ! LCOV_EXCL_LINE
242 : "Change in fsd: new ice", & ! LCOV_EXCL_LINE
243 0 : "Avg over freq period", c1, c0, ns, f_dafsd_newi)
244 0 : if (f_dafsd_latg(1:1) /= 'x') &
245 : call define_hist_field(n_dafsd_latg,"dafsd_latg","1/s",tstr3Df, tcstr, & ! LCOV_EXCL_LINE
246 : "Change in fsd: lateral growth", & ! LCOV_EXCL_LINE
247 0 : "Avg over freq period", c1, c0, ns, f_dafsd_latg)
248 0 : if (f_dafsd_latm(1:1) /= 'x') &
249 : call define_hist_field(n_dafsd_latm,"dafsd_latm","1/s",tstr3Df, tcstr, & ! LCOV_EXCL_LINE
250 : "Change in fsd: lateral melt", & ! LCOV_EXCL_LINE
251 0 : "Avg over freq period", c1, c0, ns, f_dafsd_latm)
252 0 : if (f_dafsd_wave(1:1) /= 'x') &
253 : call define_hist_field(n_dafsd_wave,"dafsd_wave","1/s",tstr3Df, tcstr, & ! LCOV_EXCL_LINE
254 : "Change in fsd: waves", & ! LCOV_EXCL_LINE
255 0 : "Avg over freq period", c1, c0, ns, f_dafsd_wave)
256 0 : if (f_dafsd_weld(1:1) /= 'x') &
257 : call define_hist_field(n_dafsd_weld,"dafsd_weld","1/s",tstr3Df, tcstr, & ! LCOV_EXCL_LINE
258 : "Change in fsd: welding", & ! LCOV_EXCL_LINE
259 0 : "Avg over freq period", c1, c0, ns, f_dafsd_weld)
260 : endif ! if (histfreq(ns) /= 'x')
261 : enddo ! ns
262 :
263 : endif ! tr_fsd
264 :
265 37 : end subroutine init_hist_fsd_3Df
266 :
267 : !=======================================================================
268 :
269 37 : subroutine init_hist_fsd_4Df
270 :
271 : use ice_calendar, only: nstreams, histfreq
272 : use ice_history_shared, only: tstr4Df, tcstr, define_hist_field
273 :
274 : logical (kind=log_kind) :: tr_fsd
275 : integer (kind=int_kind) :: ns
276 : character(len=*), parameter :: subname = '(init_hist_fsd_4Df)'
277 :
278 37 : call icepack_query_tracer_flags(tr_fsd_out=tr_fsd)
279 37 : call icepack_warnings_flush(nu_diag)
280 37 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
281 0 : file=__FILE__, line=__LINE__)
282 :
283 : !-----------------------------------------------------------------
284 : ! 4D (floe size, thickness category) variables
285 : !-----------------------------------------------------------------
286 :
287 37 : if (tr_fsd) then
288 :
289 0 : do ns = 1, nstreams
290 0 : if (histfreq(ns) /= 'x') then
291 :
292 0 : if (f_afsdn(1:1) /= 'x') &
293 : call define_hist_field(n_afsdn,"afsdn","1/m",tstr4Df, tcstr, & ! LCOV_EXCL_LINE
294 : "areal floe size and thickness distribution", & ! LCOV_EXCL_LINE
295 0 : "per unit bin width", c1, c0, ns, f_afsdn)
296 :
297 : endif ! if (histfreq(ns) /= 'x') then
298 : enddo ! ns
299 :
300 : endif ! tr_fsd
301 :
302 37 : end subroutine init_hist_fsd_4Df
303 :
304 : !=======================================================================
305 :
306 : ! accumulate average ice quantities or snapshots
307 : ! author: Elizabeth C. Hunke, LANL
308 :
309 23104 : subroutine accum_hist_fsd (dt, iblk)
310 :
311 : use ice_blocks, only: nx_block, ny_block
312 : use ice_constants, only: c0, c1, c2, c4
313 : use ice_history_shared, only: a2D, a3Df, a4Df, nfsd_hist, &
314 : ncat_hist, accum_hist_field, n3Dacum, n4Dscum
315 : use ice_state, only: trcrn, aicen, vicen, aice
316 : use ice_arrays_column, only: wave_sig_ht, floe_rad_c, floe_binwidth, &
317 : d_afsd_newi, d_afsd_latg, d_afsd_latm, d_afsd_wave, d_afsd_weld
318 :
319 : real (kind=dbl_kind), intent(in) :: &
320 : dt ! time step
321 :
322 : integer (kind=int_kind), intent(in) :: &
323 : iblk ! block index
324 :
325 : ! local variables
326 :
327 : integer (kind=int_kind) :: &
328 : i, j, n, k, & ! loop indices ! LCOV_EXCL_LINE
329 : nt_fsd ! fsd tracer index
330 : logical (kind=log_kind) :: tr_fsd
331 5792 : real (kind=dbl_kind) :: floeshape, puny
332 :
333 5792 : real (kind=dbl_kind) :: workb, workc
334 5067872 : real (kind=dbl_kind), dimension(nx_block,ny_block) :: worka
335 5073664 : real (kind=dbl_kind), dimension(nx_block,ny_block,nfsd_hist) :: worke
336 25241408 : real (kind=dbl_kind), dimension(nx_block,ny_block,nfsd_hist,ncat_hist) :: workd
337 :
338 : character(len=*), parameter :: subname = '(accum_hist_fsd)'
339 :
340 23104 : call icepack_query_parameters(floeshape_out=floeshape, puny_out=puny)
341 23104 : call icepack_query_tracer_flags(tr_fsd_out=tr_fsd)
342 23104 : call icepack_query_tracer_indices(nt_fsd_out=nt_fsd)
343 23104 : call icepack_warnings_flush(nu_diag)
344 23104 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
345 0 : file=__FILE__, line=__LINE__)
346 :
347 23104 : if (tr_fsd) then
348 :
349 : !---------------------------------------------------------------
350 : ! increment field
351 : !---------------------------------------------------------------
352 :
353 : ! 2D fields
354 0 : if (allocated(a2D)) then
355 :
356 0 : if (f_wave_sig_ht(1:1)/= 'x') &
357 : call accum_hist_field(n_wave_sig_ht, iblk, & ! LCOV_EXCL_LINE
358 0 : wave_sig_ht(:,:,iblk), a2D)
359 :
360 0 : if (f_aice_ww(1:1)/= 'x') then
361 0 : do j = 1, ny_block
362 0 : do i = 1, nx_block
363 0 : worka(i,j) = c0
364 0 : do n = 1, ncat_hist
365 0 : do k = 1, nfsd_hist
366 0 : worka(i,j) = worka(i,j) + aicen(i,j,n,iblk)*trcrn(i,j,nt_fsd+k-1,n,iblk)
367 : end do
368 : end do
369 : end do
370 : end do
371 0 : call accum_hist_field(n_aice_ww, iblk, worka(:,:), a2D)
372 : endif
373 :
374 0 : if (f_diam_ww (1:1) /= 'x') then
375 0 : do j = 1, ny_block
376 0 : do i = 1, nx_block
377 0 : worka(i,j) = c0
378 0 : workb = c0
379 0 : do n = 1, ncat_hist
380 0 : do k = 1, nfsd_hist
381 0 : workc = aicen(i,j,n,iblk)*trcrn(i,j,nt_fsd+k-1,n,iblk) &
382 0 : / (c4*floeshape*floe_rad_c(k)**2)
383 : ! number-mean radius
384 0 : worka(i,j) = worka(i,j) + workc * floe_rad_c(k)
385 : ! normalization factor
386 0 : workb = workb + workc
387 : end do
388 : end do
389 : ! number-mean diameter, following Eqn. 5 in HT2017
390 0 : if (workb > c0) worka(i,j) = c2*worka(i,j) / workb
391 0 : worka(i,j) = MAX(c2*floe_rad_c(1),worka(i,j)) ! >= smallest resolved diameter
392 : end do
393 : end do
394 0 : call accum_hist_field(n_diam_ww, iblk, worka(:,:), a2D)
395 : endif
396 :
397 0 : if (f_hice_ww (1:1) /= 'x') then
398 0 : do j = 1, ny_block
399 0 : do i = 1, nx_block
400 0 : worka(i,j) = c0
401 0 : workb = c0
402 0 : do n = 1, ncat_hist
403 0 : do k = 1, nfsd_hist
404 0 : workb = workb + aicen(i,j,n,iblk)*trcrn(i,j,nt_fsd+k-1,n,iblk)
405 0 : worka(i,j) = worka(i,j) + vicen(i,j,n,iblk)*trcrn(i,j,nt_fsd+k-1,n,iblk)
406 : end do
407 : end do
408 0 : if (workb > puny) then
409 0 : worka(i,j) = worka(i,j) / workb
410 : else
411 0 : worka(i,j) = c1
412 : endif
413 : end do
414 : end do
415 0 : call accum_hist_field(n_hice_ww, iblk, worka(:,:), a2D)
416 : endif
417 :
418 0 : if (f_fsdrad(1:1) /= 'x') then
419 0 : do j = 1, ny_block
420 0 : do i = 1, nx_block
421 0 : worka(i,j) = c0
422 0 : if (aice(i,j,iblk) > puny) then
423 0 : do k = 1, nfsd_hist
424 0 : do n = 1, ncat_hist
425 0 : worka(i,j) = worka(i,j) &
426 : + (trcrn(i,j,nt_fsd+k-1,n,iblk) * floe_rad_c(k) & ! LCOV_EXCL_LINE
427 0 : * aicen(i,j,n,iblk)/aice(i,j,iblk))
428 : end do
429 : end do
430 : endif
431 : end do
432 : end do
433 0 : call accum_hist_field(n_fsdrad, iblk, worka(:,:), a2D)
434 : endif
435 :
436 0 : if (f_fsdperim(1:1) /= 'x') then
437 0 : do j = 1, ny_block
438 0 : do i = 1, nx_block
439 0 : worka(i,j) = c0
440 0 : if (aice(i,j,iblk) > puny) then
441 0 : do k = 1, nfsd_hist
442 0 : do n = 1, ncat_hist
443 0 : worka(i,j) = worka(i,j) &
444 : + (c8*floeshape*trcrn(i,j,nt_fsd+k-1,n,iblk)*floe_rad_c(k) & ! LCOV_EXCL_LINE
445 0 : *aicen(i,j,n,iblk)/(c4*floeshape*floe_rad_c(k)**2 *aice(i,j,iblk)))
446 : end do
447 : end do
448 : endif
449 : end do
450 : end do
451 0 : call accum_hist_field(n_fsdperim, iblk, worka, a2D)
452 : endif
453 :
454 : endif ! a2D allocated
455 :
456 : ! 3D category fields
457 0 : if (allocated(a3Df)) then
458 :
459 0 : if (f_afsd(1:1) /= 'x') then
460 0 : do j = 1, ny_block
461 0 : do i = 1, nx_block
462 0 : do k = 1, nfsd_hist
463 0 : worke(i,j,k)=c0
464 0 : do n = 1, ncat_hist
465 0 : worke(i,j,k) = worke(i,j,k) + (trcrn(i,j,nt_fsd+k-1,n,iblk) &
466 0 : * aicen(i,j,n,iblk)/floe_binwidth(k))
467 : end do
468 : end do
469 : end do
470 : end do
471 0 : call accum_hist_field(n_afsd-n3Dacum, iblk, nfsd_hist, worke, a3Df)
472 : endif
473 :
474 0 : if (f_dafsd_newi(1:1)/= 'x') &
475 : call accum_hist_field(n_dafsd_newi-n3Dacum, iblk, nfsd_hist, & ! LCOV_EXCL_LINE
476 0 : d_afsd_newi(:,:,1:nfsd_hist,iblk)/dt, a3Df)
477 0 : if (f_dafsd_latg(1:1)/= 'x') &
478 : call accum_hist_field(n_dafsd_latg-n3Dacum, iblk, nfsd_hist, & ! LCOV_EXCL_LINE
479 0 : d_afsd_latg(:,:,1:nfsd_hist,iblk)/dt, a3Df)
480 0 : if (f_dafsd_latm(1:1)/= 'x') &
481 : call accum_hist_field(n_dafsd_latm-n3Dacum, iblk, nfsd_hist, & ! LCOV_EXCL_LINE
482 0 : d_afsd_latm(:,:,1:nfsd_hist,iblk)/dt, a3Df)
483 0 : if (f_dafsd_wave(1:1)/= 'x') &
484 : call accum_hist_field(n_dafsd_wave-n3Dacum, iblk, nfsd_hist, & ! LCOV_EXCL_LINE
485 0 : d_afsd_wave(:,:,1:nfsd_hist,iblk)/dt, a3Df)
486 0 : if (f_dafsd_weld(1:1)/= 'x') &
487 : call accum_hist_field(n_dafsd_weld-n3Dacum, iblk, nfsd_hist, & ! LCOV_EXCL_LINE
488 0 : d_afsd_weld(:,:,1:nfsd_hist,iblk)/dt, a3Df)
489 : endif ! a3Df allocated
490 :
491 : ! 4D floe size, thickness category fields
492 0 : if (allocated(a4Df)) then
493 :
494 0 : if (f_afsdn(1:1) /= 'x') then
495 0 : do n = 1, ncat_hist
496 0 : do k = 1, nfsd_hist
497 0 : do j = 1, ny_block
498 0 : do i = 1, nx_block
499 0 : workd(i,j,k,n) = trcrn(i,j,nt_fsd+k-1,n,iblk) &
500 0 : * aicen(i,j,n,iblk)/floe_binwidth(k)
501 : end do
502 : end do
503 : end do
504 : end do
505 : call accum_hist_field(n_afsdn-n4Dscum, iblk, nfsd_hist, ncat_hist, &
506 0 : workd(:,:,1:nfsd_hist,1:ncat_hist), a4Df)
507 : endif
508 :
509 : endif ! a4Df allocated
510 :
511 : endif ! tr_fsd
512 :
513 23104 : end subroutine accum_hist_fsd
514 :
515 : !=======================================================================
516 :
517 : end module ice_history_fsd
518 :
519 : !=======================================================================
|