Line data Source code
1 : !=======================================================================
2 : !
3 : ! Contains CICE component driver routines common to all drivers.
4 : !
5 : ! authors Elizabeth C. Hunke, LANL
6 : ! Philip W. Jones, LANL
7 : ! William H. Lipscomb, LANL
8 : !
9 : ! 2008 ECH: created module by moving subroutines from drivers/cice4/
10 : ! 2014 ECH: created column package
11 :
12 : module ice_step_mod
13 :
14 : use ice_kinds_mod
15 : use ice_blocks, only: block, get_block
16 : use ice_blocks, only: nx_block, ny_block
17 : use ice_constants, only: c0, c1, c1000, c4, p25
18 : use ice_constants, only: field_loc_center, field_loc_NEcorner, &
19 : field_loc_Nface, field_loc_Eface, & ! LCOV_EXCL_LINE
20 : field_type_scalar, field_type_vector
21 : use ice_domain, only: halo_info, nblocks, blocks_ice
22 : use ice_domain_size, only: max_blocks
23 : use ice_exit, only: abort_ice
24 : use ice_fileunits, only: nu_diag
25 : use icepack_intfc, only: icepack_warnings_flush, icepack_warnings_aborted
26 : use icepack_intfc, only: icepack_prep_radiation
27 : use icepack_intfc, only: icepack_step_therm1
28 : use icepack_intfc, only: icepack_step_therm2
29 : use icepack_intfc, only: icepack_aggregate
30 : use icepack_intfc, only: icepack_step_ridge
31 : use icepack_intfc, only: icepack_step_wavefracture
32 : use icepack_intfc, only: icepack_step_radiation
33 : use icepack_intfc, only: icepack_ocn_mixed_layer, icepack_atm_boundary
34 : use icepack_intfc, only: icepack_biogeochemistry, icepack_load_ocean_bio_array
35 : use icepack_intfc, only: icepack_max_algae, icepack_max_nbtrcr, icepack_max_don
36 : use icepack_intfc, only: icepack_max_doc, icepack_max_dic, icepack_max_aero
37 : use icepack_intfc, only: icepack_max_fe, icepack_max_iso
38 : use icepack_intfc, only: icepack_query_parameters
39 : use icepack_intfc, only: icepack_query_tracer_flags, icepack_query_tracer_sizes
40 : use icepack_intfc, only: icepack_query_tracer_indices
41 :
42 : implicit none
43 : private
44 :
45 : public :: step_therm1, step_therm2, step_dyn_horiz, step_dyn_ridge, &
46 : step_snow, prep_radiation, step_radiation, ocean_mixed_layer, & ! LCOV_EXCL_LINE
47 : update_state, biogeochemistry, step_dyn_wave, step_prep
48 :
49 : real (kind=dbl_kind), dimension (:,:,:), allocatable :: &
50 : uvelT_icep, & ! uvel for wind stress computation in icepack ! LCOV_EXCL_LINE
51 : vvelT_icep ! vvel for wind stress computation in icepack
52 :
53 : !=======================================================================
54 :
55 : contains
56 :
57 : !=======================================================================
58 :
59 5784 : subroutine save_init
60 : ! saves initial values for aice, aicen, vicen, vsnon
61 :
62 : use ice_state, only: aice, aicen, aice_init, aicen_init, &
63 : vicen, vicen_init, vsnon, vsnon_init
64 :
65 : character(len=*), parameter :: subname = '(save_init)'
66 :
67 : !-----------------------------------------------------------------
68 : ! Save the ice area passed to the coupler (so that history fields
69 : ! can be made consistent with coupler fields).
70 : ! Save the initial ice area and volume in each category.
71 : !-----------------------------------------------------------------
72 :
73 16227768 : aice_init = aice
74 81115632 : aicen_init = aicen
75 81115632 : vicen_init = vicen
76 81115632 : vsnon_init = vsnon
77 :
78 5784 : end subroutine save_init
79 :
80 : !=======================================================================
81 :
82 5784 : subroutine step_prep
83 : ! prep for step, called outside nblock loop
84 :
85 : use ice_flux, only: uatm, vatm, uatmT, vatmT
86 : use ice_grid, only: grid_atm_dynu, grid_atm_dynv, grid_average_X2Y
87 : use ice_state, only: uvel, vvel
88 :
89 : logical (kind=log_kind) :: &
90 : highfreq ! highfreq flag
91 :
92 : logical (kind=log_kind), save :: &
93 : first_call = .true. ! first call flag
94 :
95 : character(len=*), parameter :: subname = '(step_prep)'
96 :
97 : ! Save initial state
98 :
99 5784 : call save_init
100 :
101 : ! Compute uatmT, vatmT
102 :
103 5784 : call grid_average_X2Y('S',uatm,grid_atm_dynu,uatmT,'T')
104 5784 : call grid_average_X2Y('S',vatm,grid_atm_dynv,vatmT,'T')
105 :
106 : !-----------------------------------------------------------------
107 : ! Compute uvelT_icep, vvelT_icep
108 : !-----------------------------------------------------------------
109 :
110 5784 : if (first_call) then
111 37 : allocate(uvelT_icep(nx_block,ny_block,max_blocks))
112 37 : allocate(vvelT_icep(nx_block,ny_block,max_blocks))
113 111936 : uvelT_icep = c0
114 111936 : vvelT_icep = c0
115 : endif
116 :
117 5784 : call icepack_query_parameters(highfreq_out=highfreq)
118 :
119 5784 : if (highfreq) then
120 0 : call grid_average_X2Y('A', uvel, 'U', uvelT_icep, 'T')
121 0 : call grid_average_X2Y('A', vvel, 'U', vvelT_icep, 'T')
122 : endif
123 :
124 5784 : first_call = .false.
125 :
126 5784 : end subroutine step_prep
127 :
128 : !=======================================================================
129 : !
130 : ! Scales radiation fields computed on the previous time step.
131 : !
132 : ! authors: Elizabeth Hunke, LANL
133 :
134 22944 : subroutine prep_radiation (iblk)
135 :
136 : use ice_domain_size, only: ncat, nilyr, nslyr
137 : use ice_flux, only: scale_factor, swvdr, swvdf, swidr, swidf, &
138 : alvdr_ai, alvdf_ai, alidr_ai, alidf_ai, & ! LCOV_EXCL_LINE
139 : alvdr_init, alvdf_init, alidr_init, alidf_init
140 : use ice_arrays_column, only: fswsfcn, fswintn, &
141 : fswthrun, fswthrun_vdr, fswthrun_vdf, fswthrun_idr, fswthrun_idf, & ! LCOV_EXCL_LINE
142 : fswpenln, Sswabsn, Iswabsn
143 : use ice_state, only: aice, aicen
144 : use ice_timers, only: ice_timer_start, ice_timer_stop, timer_sw
145 :
146 : integer (kind=int_kind), intent(in) :: &
147 : iblk ! block index
148 :
149 : ! local variables
150 :
151 : integer (kind=int_kind) :: &
152 : ilo,ihi,jlo,jhi, & ! beginning and end of physical domain ! LCOV_EXCL_LINE
153 : i, j ! horizontal indices
154 :
155 : type (block) :: &
156 : this_block ! block information for current block
157 :
158 : character(len=*), parameter :: subname = '(prep_radiation)'
159 :
160 22944 : call ice_timer_start(timer_sw,iblk) ! shortwave
161 :
162 16186320 : alvdr_init(:,:,iblk) = c0
163 16186320 : alvdf_init(:,:,iblk) = c0
164 16186320 : alidr_init(:,:,iblk) = c0
165 16186320 : alidf_init(:,:,iblk) = c0
166 :
167 22944 : this_block = get_block(blocks_ice(iblk),iblk)
168 22944 : ilo = this_block%ilo
169 22944 : ihi = this_block%ihi
170 22944 : jlo = this_block%jlo
171 22944 : jhi = this_block%jhi
172 :
173 : !-----------------------------------------------------------------
174 : ! Compute netsw scaling factor (new netsw / old netsw)
175 : !-----------------------------------------------------------------
176 :
177 707688 : do j = jlo, jhi
178 13826928 : do i = ilo, ihi
179 :
180 13119240 : alvdr_init(i,j,iblk) = alvdr_ai(i,j,iblk)
181 13119240 : alvdf_init(i,j,iblk) = alvdf_ai(i,j,iblk)
182 13119240 : alidr_init(i,j,iblk) = alidr_ai(i,j,iblk)
183 13119240 : alidf_init(i,j,iblk) = alidf_ai(i,j,iblk)
184 :
185 0 : call icepack_prep_radiation (scale_factor=scale_factor(i,j,iblk), &
186 : aice = aice (i,j, iblk), aicen = aicen (i,j, :,iblk), & ! LCOV_EXCL_LINE
187 : swvdr = swvdr (i,j, iblk), swvdf = swvdf (i,j, iblk), & ! LCOV_EXCL_LINE
188 : swidr = swidr (i,j, iblk), swidf = swidf (i,j, iblk), & ! LCOV_EXCL_LINE
189 : alvdr_ai = alvdr_ai(i,j, iblk), alvdf_ai = alvdf_ai(i,j, iblk), & ! LCOV_EXCL_LINE
190 : alidr_ai = alidr_ai(i,j, iblk), alidf_ai = alidf_ai(i,j, iblk), & ! LCOV_EXCL_LINE
191 : fswsfcn = fswsfcn (i,j, :,iblk), fswintn = fswintn (i,j, :,iblk), & ! LCOV_EXCL_LINE
192 : fswthrun = fswthrun(i,j, :,iblk), & ! LCOV_EXCL_LINE
193 : fswthrun_vdr = fswthrun_vdr(i,j, :,iblk), & ! LCOV_EXCL_LINE
194 : fswthrun_vdf = fswthrun_vdf(i,j, :,iblk), & ! LCOV_EXCL_LINE
195 : fswthrun_idr = fswthrun_idr(i,j, :,iblk), & ! LCOV_EXCL_LINE
196 : fswthrun_idf = fswthrun_idf(i,j, :,iblk), & ! LCOV_EXCL_LINE
197 : fswpenln = fswpenln(i,j,:,:,iblk), & ! LCOV_EXCL_LINE
198 13803984 : Sswabsn = Sswabsn (i,j,:,:,iblk), Iswabsn = Iswabsn (i,j,:,:,iblk))
199 :
200 : enddo ! i
201 : enddo ! j
202 :
203 22944 : call icepack_warnings_flush(nu_diag)
204 22944 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
205 0 : file=__FILE__, line=__LINE__)
206 :
207 22944 : call ice_timer_stop(timer_sw,iblk) ! shortwave
208 :
209 22944 : end subroutine prep_radiation
210 :
211 : !=======================================================================
212 : !
213 : ! Driver for updating ice and snow internal temperatures and
214 : ! computing thermodynamic growth rates and coupler fluxes.
215 : !
216 : ! authors: William H. Lipscomb, LANL
217 :
218 22944 : subroutine step_therm1 (dt, iblk)
219 :
220 : use ice_arrays_column, only: ffracn, dhsn, &
221 : Cdn_ocn, Cdn_ocn_skin, Cdn_ocn_floe, Cdn_ocn_keel, Cdn_atm_ratio, & ! LCOV_EXCL_LINE
222 : Cdn_atm, Cdn_atm_skin, Cdn_atm_floe, Cdn_atm_rdg, Cdn_atm_pond, & ! LCOV_EXCL_LINE
223 : hfreebd, hdraft, hridge, distrdg, hkeel, dkeel, lfloe, dfloe, & ! LCOV_EXCL_LINE
224 : fswsfcn, fswintn, Sswabsn, Iswabsn, meltsliqn, meltsliq, & ! LCOV_EXCL_LINE
225 : fswthrun, fswthrun_vdr, fswthrun_vdf, fswthrun_idr, fswthrun_idf
226 : use ice_calendar, only: yday
227 : use ice_domain_size, only: ncat, nilyr, nslyr, n_iso, n_aero
228 : use ice_flux, only: frzmlt, sst, Tf, strocnxT_iavg, strocnyT_iavg, rside, fbot, Tbot, Tsnice, &
229 : meltsn, melttn, meltbn, congeln, snoicen, uatmT, vatmT, fside, wlat, & ! LCOV_EXCL_LINE
230 : wind, rhoa, potT, Qa, zlvl, zlvs, strax, stray, flatn, fsensn, fsurfn, fcondtopn, & ! LCOV_EXCL_LINE
231 : flw, fsnow, fpond, sss, mlt_onset, frz_onset, fcondbotn, fcondbot, fsloss, & ! LCOV_EXCL_LINE
232 : frain, Tair, strairxT, strairyT, fsurf, fcondtop, fsens, & ! LCOV_EXCL_LINE
233 : flat, fswabs, flwout, evap, evaps, evapi, Tref, Qref, Uref, fresh, fsalt, fhocn, & ! LCOV_EXCL_LINE
234 : fswthru, fswthru_vdr, fswthru_vdf, fswthru_idr, fswthru_idf, & ! LCOV_EXCL_LINE
235 : meltt, melts, meltb, congel, snoice, & ! LCOV_EXCL_LINE
236 : flatn_f, fsensn_f, fsurfn_f, fcondtopn_f, & ! LCOV_EXCL_LINE
237 : send_i2x_per_cat, fswthrun_ai, dsnow
238 : use ice_flux_bgc, only: dsnown, faero_atm, faero_ocn, fiso_atm, fiso_ocn, &
239 : Qa_iso, Qref_iso, fiso_evap, HDO_ocn, H2_16O_ocn, H2_18O_ocn
240 : use ice_grid, only: lmask_n, lmask_s, tmask
241 : use ice_state, only: aice, aicen, aicen_init, vicen_init, &
242 : vice, vicen, vsno, vsnon, trcrn, vsnon_init
243 : #ifdef CICE_IN_NEMO
244 : use ice_state, only: aice_init
245 : #endif
246 :
247 : #ifdef CESMCOUPLED
248 : use ice_prescribed_mod, only: prescribed_ice
249 : #else
250 : logical (kind=log_kind) :: &
251 : prescribed_ice ! if .true., use prescribed ice instead of computed
252 : #endif
253 : real (kind=dbl_kind), intent(in) :: &
254 : dt ! time step (s)
255 :
256 : integer (kind=int_kind), intent(in) :: &
257 : iblk ! block index
258 :
259 : ! local variables
260 : #ifdef CICE_IN_NEMO
261 : real (kind=dbl_kind) :: &
262 : raice ! reciprocal of ice concentration
263 : #endif
264 : integer (kind=int_kind) :: &
265 : ilo,ihi,jlo,jhi, & ! beginning and end of physical domain ! LCOV_EXCL_LINE
266 : i, j , & ! horizontal indices ! LCOV_EXCL_LINE
267 : n , & ! thickness category index ! LCOV_EXCL_LINE
268 : k, kk ! indices for aerosols
269 :
270 : integer (kind=int_kind) :: &
271 : ntrcr, nt_apnd, nt_hpnd, nt_ipnd, nt_alvl, nt_vlvl, nt_Tsfc, & ! LCOV_EXCL_LINE
272 : nt_iage, nt_FY, nt_qice, nt_sice, nt_aero, nt_qsno, & ! LCOV_EXCL_LINE
273 : nt_isosno, nt_isoice, nt_rsnw, nt_smice, nt_smliq
274 :
275 : logical (kind=log_kind) :: &
276 : tr_iage, tr_FY, tr_iso, tr_aero, tr_pond, & ! LCOV_EXCL_LINE
277 : tr_pond_lvl, tr_pond_topo, calc_Tsfc, snwgrain
278 :
279 : real (kind=dbl_kind) :: &
280 5760 : puny ! a very small number
281 :
282 : real (kind=dbl_kind), dimension(n_aero,2,ncat) :: &
283 212928 : aerosno, aeroice ! kg/m^2
284 :
285 : real (kind=dbl_kind), dimension(n_iso,ncat) :: &
286 97728 : isosno, isoice ! kg/m^2
287 :
288 : real (kind=dbl_kind), dimension(nslyr,ncat) :: &
289 212928 : rsnwn, smicen, smliqn
290 :
291 : type (block) :: &
292 : this_block ! block information for current block
293 :
294 : character(len=*), parameter :: subname = '(step_therm1)'
295 :
296 22944 : call icepack_query_parameters(puny_out=puny)
297 22944 : call icepack_query_parameters(calc_Tsfc_out=calc_Tsfc)
298 22944 : call icepack_query_parameters(snwgrain_out=snwgrain)
299 22944 : call icepack_query_tracer_sizes(ntrcr_out=ntrcr)
300 : call icepack_query_tracer_flags( &
301 : tr_iage_out=tr_iage, tr_FY_out=tr_FY, tr_iso_out=tr_iso, & ! LCOV_EXCL_LINE
302 : tr_aero_out=tr_aero, tr_pond_out=tr_pond, & ! LCOV_EXCL_LINE
303 22944 : tr_pond_lvl_out=tr_pond_lvl, tr_pond_topo_out=tr_pond_topo)
304 : call icepack_query_tracer_indices( &
305 : nt_apnd_out=nt_apnd, nt_hpnd_out=nt_hpnd, nt_ipnd_out=nt_ipnd, & ! LCOV_EXCL_LINE
306 : nt_alvl_out=nt_alvl, nt_vlvl_out=nt_vlvl, nt_Tsfc_out=nt_Tsfc, & ! LCOV_EXCL_LINE
307 : nt_iage_out=nt_iage, nt_FY_out=nt_FY, & ! LCOV_EXCL_LINE
308 : nt_qice_out=nt_qice, nt_sice_out=nt_sice, & ! LCOV_EXCL_LINE
309 : nt_aero_out=nt_aero, nt_qsno_out=nt_qsno, & ! LCOV_EXCL_LINE
310 : nt_rsnw_out=nt_rsnw, nt_smice_out=nt_smice, nt_smliq_out=nt_smliq, & ! LCOV_EXCL_LINE
311 22944 : nt_isosno_out=nt_isosno, nt_isoice_out=nt_isoice)
312 22944 : call icepack_warnings_flush(nu_diag)
313 22944 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
314 0 : file=__FILE__, line=__LINE__)
315 :
316 : #ifndef CESMCOUPLED
317 22944 : prescribed_ice = .false.
318 : #endif
319 :
320 252384 : rsnwn (:,:) = c0
321 252384 : smicen (:,:) = c0
322 252384 : smliqn (:,:) = c0
323 137664 : isoice (:,:) = c0
324 367104 : aerosno(:,:,:) = c0
325 367104 : aeroice(:,:,:) = c0
326 :
327 : #ifdef CICE_IN_NEMO
328 : do j = 1, ny_block
329 : do i = 1, nx_block
330 :
331 : !---------------------------------------------------------------
332 : ! Scale frain and fsnow by ice concentration as these fields
333 : ! are supplied by NEMO multiplied by ice concentration
334 : !---------------------------------------------------------------
335 :
336 : if (aice_init(i,j,iblk) > puny) then
337 : raice = c1 / aice_init(i,j,iblk)
338 : frain(i,j,iblk) = frain(i,j,iblk)*raice
339 : fsnow(i,j,iblk) = fsnow(i,j,iblk)*raice
340 : else
341 : frain(i,j,iblk) = c0
342 : fsnow(i,j,iblk) = c0
343 : endif
344 :
345 : enddo ! i
346 : enddo ! j
347 : #endif
348 :
349 22944 : this_block = get_block(blocks_ice(iblk),iblk)
350 22944 : ilo = this_block%ilo
351 22944 : ihi = this_block%ihi
352 22944 : jlo = this_block%jlo
353 22944 : jhi = this_block%jhi
354 :
355 707688 : do j = jlo, jhi
356 13826928 : do i = ilo, ihi
357 :
358 13119240 : if (snwgrain) then
359 0 : do n = 1, ncat
360 0 : do k = 1, nslyr
361 0 : rsnwn (k,n) = trcrn(i,j,nt_rsnw +k-1,n,iblk)
362 0 : smicen(k,n) = trcrn(i,j,nt_smice+k-1,n,iblk)
363 0 : smliqn(k,n) = trcrn(i,j,nt_smliq+k-1,n,iblk)
364 : enddo
365 : enddo
366 : endif ! snwgrain
367 :
368 13119240 : if (tr_iso) then ! trcrn(nt_iso*) has units kg/m^3
369 0 : do n=1,ncat
370 0 : do k=1,n_iso
371 0 : isosno(k,n) = trcrn(i,j,nt_isosno+k-1,n,iblk) * vsnon_init(i,j,n,iblk)
372 0 : isoice(k,n) = trcrn(i,j,nt_isoice+k-1,n,iblk) * vicen_init(i,j,n,iblk)
373 : enddo
374 : enddo
375 : endif ! tr_iso
376 :
377 13119240 : if (tr_aero) then ! trcrn(nt_aero) has units kg/m^3
378 0 : do n=1,ncat
379 0 : do k=1,n_aero
380 0 : aerosno (k,:,n) = &
381 : trcrn(i,j,nt_aero+(k-1)*4 :nt_aero+(k-1)*4+1,n,iblk) & ! LCOV_EXCL_LINE
382 0 : * vsnon_init(i,j,n,iblk)
383 0 : aeroice (k,:,n) = &
384 : trcrn(i,j,nt_aero+(k-1)*4+2:nt_aero+(k-1)*4+3,n,iblk) & ! LCOV_EXCL_LINE
385 0 : * vicen_init(i,j,n,iblk)
386 : enddo
387 : enddo
388 : endif ! tr_aero
389 :
390 13119240 : if (tmask(i,j,iblk)) then
391 :
392 : call icepack_step_therm1(dt=dt, ncat=ncat, &
393 : nilyr=nilyr, nslyr=nslyr, & ! LCOV_EXCL_LINE
394 : aicen_init = aicen_init (i,j,:,iblk), & ! LCOV_EXCL_LINE
395 : vicen_init = vicen_init (i,j,:,iblk), & ! LCOV_EXCL_LINE
396 : vsnon_init = vsnon_init (i,j,:,iblk), & ! LCOV_EXCL_LINE
397 : aice = aice (i,j, iblk), & ! LCOV_EXCL_LINE
398 : aicen = aicen (i,j,:,iblk), & ! LCOV_EXCL_LINE
399 : vice = vice (i,j, iblk), & ! LCOV_EXCL_LINE
400 : vicen = vicen (i,j,:,iblk), & ! LCOV_EXCL_LINE
401 : vsno = vsno (i,j, iblk), & ! LCOV_EXCL_LINE
402 : vsnon = vsnon (i,j,:,iblk), & ! LCOV_EXCL_LINE
403 : uvel = uvelT_icep (i,j, iblk), & ! LCOV_EXCL_LINE
404 : vvel = vvelT_icep (i,j, iblk), & ! LCOV_EXCL_LINE
405 : Tsfc = trcrn (i,j,nt_Tsfc,:,iblk), & ! LCOV_EXCL_LINE
406 : zqsn = trcrn (i,j,nt_qsno:nt_qsno+nslyr-1,:,iblk), & ! LCOV_EXCL_LINE
407 : zqin = trcrn (i,j,nt_qice:nt_qice+nilyr-1,:,iblk), & ! LCOV_EXCL_LINE
408 : zSin = trcrn (i,j,nt_sice:nt_sice+nilyr-1,:,iblk), & ! LCOV_EXCL_LINE
409 : alvl = trcrn (i,j,nt_alvl,:,iblk), & ! LCOV_EXCL_LINE
410 : vlvl = trcrn (i,j,nt_vlvl,:,iblk), & ! LCOV_EXCL_LINE
411 : apnd = trcrn (i,j,nt_apnd,:,iblk), & ! LCOV_EXCL_LINE
412 : hpnd = trcrn (i,j,nt_hpnd,:,iblk), & ! LCOV_EXCL_LINE
413 : ipnd = trcrn (i,j,nt_ipnd,:,iblk), & ! LCOV_EXCL_LINE
414 : iage = trcrn (i,j,nt_iage,:,iblk), & ! LCOV_EXCL_LINE
415 : FY = trcrn (i,j,nt_FY ,:,iblk), & ! LCOV_EXCL_LINE
416 : rsnwn = rsnwn (:,:), & ! LCOV_EXCL_LINE
417 : smicen = smicen (:,:), & ! LCOV_EXCL_LINE
418 : smliqn = smliqn (:,:), & ! LCOV_EXCL_LINE
419 : aerosno = aerosno (:,:,:), & ! LCOV_EXCL_LINE
420 : aeroice = aeroice (:,:,:), & ! LCOV_EXCL_LINE
421 : isosno = isosno (:,:), & ! LCOV_EXCL_LINE
422 : isoice = isoice (:,:), & ! LCOV_EXCL_LINE
423 : uatm = uatmT (i,j, iblk), & ! LCOV_EXCL_LINE
424 : vatm = vatmT (i,j, iblk), & ! LCOV_EXCL_LINE
425 : wind = wind (i,j, iblk), & ! LCOV_EXCL_LINE
426 : zlvl = zlvl (i,j, iblk), & ! LCOV_EXCL_LINE
427 : zlvs = zlvs (i,j, iblk), & ! LCOV_EXCL_LINE
428 : Qa = Qa (i,j, iblk), & ! LCOV_EXCL_LINE
429 : Qa_iso = Qa_iso (i,j,:,iblk), & ! LCOV_EXCL_LINE
430 : rhoa = rhoa (i,j, iblk), & ! LCOV_EXCL_LINE
431 : Tair = Tair (i,j, iblk), & ! LCOV_EXCL_LINE
432 : Tref = Tref (i,j, iblk), & ! LCOV_EXCL_LINE
433 : Qref = Qref (i,j, iblk), & ! LCOV_EXCL_LINE
434 : Qref_iso = Qref_iso (i,j,:,iblk), & ! LCOV_EXCL_LINE
435 : Uref = Uref (i,j, iblk), & ! LCOV_EXCL_LINE
436 : Cdn_atm_ratio= Cdn_atm_ratio(i,j, iblk), & ! LCOV_EXCL_LINE
437 : Cdn_ocn = Cdn_ocn (i,j, iblk), & ! LCOV_EXCL_LINE
438 : Cdn_ocn_skin = Cdn_ocn_skin(i,j, iblk), & ! LCOV_EXCL_LINE
439 : Cdn_ocn_floe = Cdn_ocn_floe(i,j, iblk), & ! LCOV_EXCL_LINE
440 : Cdn_ocn_keel = Cdn_ocn_keel(i,j, iblk), & ! LCOV_EXCL_LINE
441 : Cdn_atm = Cdn_atm (i,j, iblk), & ! LCOV_EXCL_LINE
442 : Cdn_atm_skin = Cdn_atm_skin(i,j, iblk), & ! LCOV_EXCL_LINE
443 : Cdn_atm_floe = Cdn_atm_floe(i,j, iblk), & ! LCOV_EXCL_LINE
444 : Cdn_atm_pond = Cdn_atm_pond(i,j, iblk), & ! LCOV_EXCL_LINE
445 : Cdn_atm_rdg = Cdn_atm_rdg (i,j, iblk), & ! LCOV_EXCL_LINE
446 : hfreebd = hfreebd (i,j, iblk), & ! LCOV_EXCL_LINE
447 : hdraft = hdraft (i,j, iblk), & ! LCOV_EXCL_LINE
448 : hridge = hridge (i,j, iblk), & ! LCOV_EXCL_LINE
449 : distrdg = distrdg (i,j, iblk), & ! LCOV_EXCL_LINE
450 : hkeel = hkeel (i,j, iblk), & ! LCOV_EXCL_LINE
451 : dkeel = dkeel (i,j, iblk), & ! LCOV_EXCL_LINE
452 : lfloe = lfloe (i,j, iblk), & ! LCOV_EXCL_LINE
453 : dfloe = dfloe (i,j, iblk), & ! LCOV_EXCL_LINE
454 : strax = strax (i,j, iblk), & ! LCOV_EXCL_LINE
455 : stray = stray (i,j, iblk), & ! LCOV_EXCL_LINE
456 : strairxT = strairxT (i,j, iblk), & ! LCOV_EXCL_LINE
457 : strairyT = strairyT (i,j, iblk), & ! LCOV_EXCL_LINE
458 : potT = potT (i,j, iblk), & ! LCOV_EXCL_LINE
459 : sst = sst (i,j, iblk), & ! LCOV_EXCL_LINE
460 : sss = sss (i,j, iblk), & ! LCOV_EXCL_LINE
461 : Tf = Tf (i,j, iblk), & ! LCOV_EXCL_LINE
462 : strocnxT = strocnxT_iavg(i,j, iblk), & ! LCOV_EXCL_LINE
463 : strocnyT = strocnyT_iavg(i,j, iblk), & ! LCOV_EXCL_LINE
464 : fbot = fbot (i,j, iblk), & ! LCOV_EXCL_LINE
465 : Tbot = Tbot (i,j, iblk), & ! LCOV_EXCL_LINE
466 : Tsnice = Tsnice (i,j, iblk), & ! LCOV_EXCL_LINE
467 : frzmlt = frzmlt (i,j, iblk), & ! LCOV_EXCL_LINE
468 : rside = rside (i,j, iblk), & ! LCOV_EXCL_LINE
469 : fside = fside (i,j, iblk), & ! LCOV_EXCL_LINE
470 : wlat = wlat (i,j, iblk), & ! LCOV_EXCL_LINE
471 : fsnow = fsnow (i,j, iblk), & ! LCOV_EXCL_LINE
472 : frain = frain (i,j, iblk), & ! LCOV_EXCL_LINE
473 : fpond = fpond (i,j, iblk), & ! LCOV_EXCL_LINE
474 : fsloss = fsloss (i,j, iblk), & ! LCOV_EXCL_LINE
475 : fsurf = fsurf (i,j, iblk), & ! LCOV_EXCL_LINE
476 : fsurfn = fsurfn (i,j,:,iblk), & ! LCOV_EXCL_LINE
477 : fcondtop = fcondtop (i,j, iblk), & ! LCOV_EXCL_LINE
478 : fcondtopn = fcondtopn (i,j,:,iblk), & ! LCOV_EXCL_LINE
479 : fcondbot = fcondbot (i,j, iblk), & ! LCOV_EXCL_LINE
480 : fcondbotn = fcondbotn (i,j,:,iblk), & ! LCOV_EXCL_LINE
481 : fswsfcn = fswsfcn (i,j,:,iblk), & ! LCOV_EXCL_LINE
482 : fswintn = fswintn (i,j,:,iblk), & ! LCOV_EXCL_LINE
483 : fswthrun = fswthrun (i,j,:,iblk), & ! LCOV_EXCL_LINE
484 : fswthrun_vdr = fswthrun_vdr (i,j,:,iblk),& ! LCOV_EXCL_LINE
485 : fswthrun_vdf = fswthrun_vdf (i,j,:,iblk),& ! LCOV_EXCL_LINE
486 : fswthrun_idr = fswthrun_idr (i,j,:,iblk),& ! LCOV_EXCL_LINE
487 : fswthrun_idf = fswthrun_idf (i,j,:,iblk),& ! LCOV_EXCL_LINE
488 : fswabs = fswabs (i,j, iblk), & ! LCOV_EXCL_LINE
489 : flwout = flwout (i,j, iblk), & ! LCOV_EXCL_LINE
490 : Sswabsn = Sswabsn (i,j,:,:,iblk), & ! LCOV_EXCL_LINE
491 : Iswabsn = Iswabsn (i,j,:,:,iblk), & ! LCOV_EXCL_LINE
492 : flw = flw (i,j, iblk), & ! LCOV_EXCL_LINE
493 : fsens = fsens (i,j, iblk), & ! LCOV_EXCL_LINE
494 : fsensn = fsensn (i,j,:,iblk), & ! LCOV_EXCL_LINE
495 : flat = flat (i,j, iblk), & ! LCOV_EXCL_LINE
496 : flatn = flatn (i,j,:,iblk), & ! LCOV_EXCL_LINE
497 : evap = evap (i,j, iblk), & ! LCOV_EXCL_LINE
498 : evaps = evaps (i,j, iblk), & ! LCOV_EXCL_LINE
499 : evapi = evapi (i,j, iblk), & ! LCOV_EXCL_LINE
500 : fresh = fresh (i,j, iblk), & ! LCOV_EXCL_LINE
501 : fsalt = fsalt (i,j, iblk), & ! LCOV_EXCL_LINE
502 : fhocn = fhocn (i,j, iblk), & ! LCOV_EXCL_LINE
503 : fswthru = fswthru (i,j, iblk), & ! LCOV_EXCL_LINE
504 : fswthru_vdr = fswthru_vdr (i,j, iblk), & ! LCOV_EXCL_LINE
505 : fswthru_vdf = fswthru_vdf (i,j, iblk), & ! LCOV_EXCL_LINE
506 : fswthru_idr = fswthru_idr (i,j, iblk), & ! LCOV_EXCL_LINE
507 : fswthru_idf = fswthru_idf (i,j, iblk), & ! LCOV_EXCL_LINE
508 : flatn_f = flatn_f (i,j,:,iblk), & ! LCOV_EXCL_LINE
509 : fsensn_f = fsensn_f (i,j,:,iblk), & ! LCOV_EXCL_LINE
510 : fsurfn_f = fsurfn_f (i,j,:,iblk), & ! LCOV_EXCL_LINE
511 : fcondtopn_f = fcondtopn_f (i,j,:,iblk), & ! LCOV_EXCL_LINE
512 : faero_atm = faero_atm (i,j,1:n_aero,iblk), & ! LCOV_EXCL_LINE
513 : faero_ocn = faero_ocn (i,j,1:n_aero,iblk), & ! LCOV_EXCL_LINE
514 : fiso_atm = fiso_atm (i,j,:,iblk), & ! LCOV_EXCL_LINE
515 : fiso_ocn = fiso_ocn (i,j,:,iblk), & ! LCOV_EXCL_LINE
516 : fiso_evap = fiso_evap (i,j,:,iblk), & ! LCOV_EXCL_LINE
517 : HDO_ocn = HDO_ocn (i,j, iblk), & ! LCOV_EXCL_LINE
518 : H2_16O_ocn = H2_16O_ocn (i,j, iblk), & ! LCOV_EXCL_LINE
519 : H2_18O_ocn = H2_18O_ocn (i,j, iblk), & ! LCOV_EXCL_LINE
520 : dhsn = dhsn (i,j,:,iblk), & ! LCOV_EXCL_LINE
521 : ffracn = ffracn (i,j,:,iblk), & ! LCOV_EXCL_LINE
522 : meltt = meltt (i,j, iblk), & ! LCOV_EXCL_LINE
523 : melttn = melttn (i,j,:,iblk), & ! LCOV_EXCL_LINE
524 : meltb = meltb (i,j, iblk), & ! LCOV_EXCL_LINE
525 : meltbn = meltbn (i,j,:,iblk), & ! LCOV_EXCL_LINE
526 : melts = melts (i,j, iblk), & ! LCOV_EXCL_LINE
527 : meltsn = meltsn (i,j,:,iblk), & ! LCOV_EXCL_LINE
528 : congel = congel (i,j, iblk), & ! LCOV_EXCL_LINE
529 : congeln = congeln (i,j,:,iblk), & ! LCOV_EXCL_LINE
530 : snoice = snoice (i,j, iblk), & ! LCOV_EXCL_LINE
531 : snoicen = snoicen (i,j,:,iblk), & ! LCOV_EXCL_LINE
532 : dsnow = dsnow (i,j, iblk), & ! LCOV_EXCL_LINE
533 : dsnown = dsnown (i,j,:,iblk), & ! LCOV_EXCL_LINE
534 : meltsliq = meltsliq (i,j, iblk), & ! LCOV_EXCL_LINE
535 : meltsliqn = meltsliqn (i,j,:,iblk), & ! LCOV_EXCL_LINE
536 : lmask_n = lmask_n (i,j, iblk), & ! LCOV_EXCL_LINE
537 : lmask_s = lmask_s (i,j, iblk), & ! LCOV_EXCL_LINE
538 : mlt_onset = mlt_onset (i,j, iblk), & ! LCOV_EXCL_LINE
539 : frz_onset = frz_onset (i,j, iblk), & ! LCOV_EXCL_LINE
540 10719744 : yday=yday, prescribed_ice=prescribed_ice)
541 :
542 : !-----------------------------------------------------------------
543 : ! handle per-category i2x fields, no merging
544 : !-----------------------------------------------------------------
545 :
546 10719744 : if (send_i2x_per_cat) then
547 0 : do n = 1, ncat
548 : ! TODO (mvertens, 2018-12-22): do we need to add the band separated quantities
549 : ! for MOM6 here also?
550 :
551 0 : fswthrun_ai(i,j,n,iblk) = fswthrun(i,j,n,iblk)*aicen_init(i,j,n,iblk)
552 : enddo ! ncat
553 : endif
554 :
555 : endif
556 :
557 13119240 : if (snwgrain) then
558 0 : do n = 1, ncat
559 0 : do k = 1, nslyr
560 0 : trcrn(i,j,nt_rsnw +k-1,n,iblk) = rsnwn (k,n)
561 0 : trcrn(i,j,nt_smice+k-1,n,iblk) = smicen(k,n)
562 0 : trcrn(i,j,nt_smliq+k-1,n,iblk) = smliqn(k,n)
563 : enddo
564 : enddo
565 : endif ! snwgrain
566 :
567 13119240 : if (tr_iso) then
568 0 : do n = 1, ncat
569 0 : if (vicen(i,j,n,iblk) > puny) &
570 0 : isoice(:,n) = isoice(:,n)/vicen(i,j,n,iblk)
571 0 : if (vsnon(i,j,n,iblk) > puny) &
572 0 : isosno(:,n) = isosno(:,n)/vsnon(i,j,n,iblk)
573 0 : do k = 1, n_iso
574 0 : trcrn(i,j,nt_isosno+k-1,n,iblk) = isosno(k,n)
575 0 : trcrn(i,j,nt_isoice+k-1,n,iblk) = isoice(k,n)
576 : enddo
577 : enddo
578 : endif ! tr_iso
579 :
580 13803984 : if (tr_aero) then
581 0 : do n = 1, ncat
582 0 : if (vicen(i,j,n,iblk) > puny) &
583 0 : aeroice(:,:,n) = aeroice(:,:,n)/vicen(i,j,n,iblk)
584 0 : if (vsnon(i,j,n,iblk) > puny) &
585 0 : aerosno(:,:,n) = aerosno(:,:,n)/vsnon(i,j,n,iblk)
586 0 : do k = 1, n_aero
587 0 : do kk = 1, 2
588 0 : trcrn(i,j,nt_aero+(k-1)*4+kk-1,n,iblk)=aerosno(k,kk,n)
589 0 : trcrn(i,j,nt_aero+(k-1)*4+kk+1,n,iblk)=aeroice(k,kk,n)
590 : enddo
591 : enddo
592 : enddo
593 : endif ! tr_aero
594 :
595 : enddo ! i
596 : enddo ! j
597 :
598 22944 : call icepack_warnings_flush(nu_diag)
599 22944 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
600 0 : file=__FILE__, line=__LINE__)
601 :
602 22944 : end subroutine step_therm1
603 :
604 : !=======================================================================
605 : ! Driver for thermodynamic changes not needed for coupling:
606 : ! transport in thickness space, lateral growth and melting.
607 : !
608 : ! authors: William H. Lipscomb, LANL
609 : ! Elizabeth C. Hunke, LANL
610 :
611 22944 : subroutine step_therm2 (dt, iblk)
612 :
613 : use ice_arrays_column, only: hin_max, ocean_bio, wave_sig_ht, &
614 : wave_spectrum, wavefreq, dwavefreq, & ! LCOV_EXCL_LINE
615 : first_ice, bgrid, cgrid, igrid, floe_rad_c, floe_binwidth, & ! LCOV_EXCL_LINE
616 : d_afsd_latg, d_afsd_newi, d_afsd_latm, d_afsd_weld
617 : use ice_calendar, only: yday
618 : use ice_domain_size, only: ncat, nilyr, nslyr, nblyr, nfsd
619 : use ice_flux, only: fresh, frain, fpond, frzmlt, frazil, frz_onset, &
620 : fsalt, Tf, sss, salinz, fhocn, rside, fside, wlat, & ! LCOV_EXCL_LINE
621 : meltl, frazil_diag
622 : use ice_flux_bgc, only: flux_bio, faero_ocn, &
623 : fiso_ocn, HDO_ocn, H2_16O_ocn, H2_18O_ocn
624 : use ice_grid, only: tmask
625 : use ice_state, only: aice, aicen, aice0, trcr_depend, &
626 : aicen_init, vicen_init, trcrn, vicen, vsnon, & ! LCOV_EXCL_LINE
627 : trcr_base, n_trcr_strata, nt_strata
628 :
629 : real (kind=dbl_kind), intent(in) :: &
630 : dt ! time step
631 :
632 : integer (kind=int_kind), intent(in) :: &
633 : iblk ! block index
634 :
635 : ! local variables
636 :
637 : integer (kind=int_kind) :: &
638 : ilo,ihi,jlo,jhi, & ! beginning and end of physical domain ! LCOV_EXCL_LINE
639 : i, j ! horizontal indices
640 :
641 : integer (kind=int_kind) :: &
642 : ntrcr, nbtrcr, nltrcr
643 :
644 : logical (kind=log_kind) :: &
645 : tr_fsd, & ! floe size distribution tracers ! LCOV_EXCL_LINE
646 : z_tracers ! vertical biogeochemistry
647 :
648 : type (block) :: &
649 : this_block ! block information for current block
650 :
651 : character(len=*), parameter :: subname = '(step_therm2)'
652 :
653 22944 : call icepack_query_parameters(z_tracers_out=z_tracers)
654 22944 : call icepack_query_tracer_sizes(ntrcr_out=ntrcr, nbtrcr_out=nbtrcr)
655 22944 : call icepack_query_tracer_flags(tr_fsd_out=tr_fsd)
656 22944 : call icepack_warnings_flush(nu_diag)
657 22944 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
658 0 : file=__FILE__, line=__LINE__)
659 :
660 : ! nltrcr is only used as a zbgc flag in icepack (number of zbgc tracers > 0)
661 22944 : if (z_tracers) then
662 0 : nltrcr = 1
663 : else
664 22944 : nltrcr = 0
665 : endif
666 :
667 22944 : this_block = get_block(blocks_ice(iblk),iblk)
668 22944 : ilo = this_block%ilo
669 22944 : ihi = this_block%ihi
670 22944 : jlo = this_block%jlo
671 22944 : jhi = this_block%jhi
672 :
673 707688 : do j = jlo, jhi
674 13826928 : do i = ilo, ihi
675 :
676 13803984 : if (tmask(i,j,iblk)) then
677 :
678 : ! significant wave height for FSD
679 10719744 : if (tr_fsd) &
680 0 : wave_sig_ht(i,j,iblk) = c4*SQRT(SUM(wave_spectrum(i,j,:,iblk)*dwavefreq(:)))
681 :
682 : call icepack_step_therm2(dt=dt, ncat=ncat, &
683 : nltrcr=nltrcr, nilyr=nilyr, nslyr=nslyr, nblyr=nblyr, & ! LCOV_EXCL_LINE
684 : hin_max = hin_max (:), & ! LCOV_EXCL_LINE
685 : aicen = aicen (i,j,:,iblk), & ! LCOV_EXCL_LINE
686 : vicen = vicen (i,j,:,iblk), & ! LCOV_EXCL_LINE
687 : vsnon = vsnon (i,j,:,iblk), & ! LCOV_EXCL_LINE
688 : aicen_init = aicen_init(i,j,:,iblk), & ! LCOV_EXCL_LINE
689 : vicen_init = vicen_init(i,j,:,iblk), & ! LCOV_EXCL_LINE
690 : trcrn = trcrn (i,j,:,:,iblk), & ! LCOV_EXCL_LINE
691 : aice0 = aice0 (i,j, iblk), & ! LCOV_EXCL_LINE
692 : aice = aice (i,j, iblk), & ! LCOV_EXCL_LINE
693 : trcr_depend= trcr_depend(:), & ! LCOV_EXCL_LINE
694 : trcr_base = trcr_base(:,:), & ! LCOV_EXCL_LINE
695 : n_trcr_strata = n_trcr_strata(:), & ! LCOV_EXCL_LINE
696 : nt_strata = nt_strata(:,:), & ! LCOV_EXCL_LINE
697 : Tf = Tf (i,j, iblk), & ! LCOV_EXCL_LINE
698 : sss = sss (i,j, iblk), & ! LCOV_EXCL_LINE
699 : salinz = salinz (i,j,:,iblk), & ! LCOV_EXCL_LINE
700 : rside = rside (i,j, iblk), & ! LCOV_EXCL_LINE
701 : meltl = meltl (i,j, iblk), & ! LCOV_EXCL_LINE
702 : fside = fside (i,j, iblk), & ! LCOV_EXCL_LINE
703 : wlat = wlat (i,j, iblk), & ! LCOV_EXCL_LINE
704 : frzmlt = frzmlt (i,j, iblk), & ! LCOV_EXCL_LINE
705 : frazil = frazil (i,j, iblk), & ! LCOV_EXCL_LINE
706 : frain = frain (i,j, iblk), & ! LCOV_EXCL_LINE
707 : fpond = fpond (i,j, iblk), & ! LCOV_EXCL_LINE
708 : fresh = fresh (i,j, iblk), & ! LCOV_EXCL_LINE
709 : fsalt = fsalt (i,j, iblk), & ! LCOV_EXCL_LINE
710 : fhocn = fhocn (i,j, iblk), & ! LCOV_EXCL_LINE
711 : bgrid = bgrid, & ! LCOV_EXCL_LINE
712 : cgrid = cgrid, & ! LCOV_EXCL_LINE
713 : igrid = igrid, & ! LCOV_EXCL_LINE
714 : faero_ocn = faero_ocn (i,j,:,iblk), & ! LCOV_EXCL_LINE
715 : first_ice = first_ice (i,j,:,iblk), & ! LCOV_EXCL_LINE
716 : flux_bio = flux_bio (i,j,1:nbtrcr,iblk), & ! LCOV_EXCL_LINE
717 : ocean_bio = ocean_bio (i,j,1:nbtrcr,iblk), & ! LCOV_EXCL_LINE
718 : frazil_diag= frazil_diag(i,j,iblk), & ! LCOV_EXCL_LINE
719 : frz_onset = frz_onset (i,j, iblk), & ! LCOV_EXCL_LINE
720 : yday = yday, & ! LCOV_EXCL_LINE
721 : fiso_ocn = fiso_ocn (i,j,:,iblk), & ! LCOV_EXCL_LINE
722 : HDO_ocn = HDO_ocn (i,j, iblk), & ! LCOV_EXCL_LINE
723 : H2_16O_ocn = H2_16O_ocn(i,j, iblk), & ! LCOV_EXCL_LINE
724 : H2_18O_ocn = H2_18O_ocn(i,j, iblk), & ! LCOV_EXCL_LINE
725 : nfsd = nfsd, & ! LCOV_EXCL_LINE
726 : wave_sig_ht= wave_sig_ht(i,j,iblk), & ! LCOV_EXCL_LINE
727 : wave_spectrum = wave_spectrum(i,j,:,iblk), & ! LCOV_EXCL_LINE
728 : wavefreq = wavefreq(:), & ! LCOV_EXCL_LINE
729 : dwavefreq = dwavefreq(:), & ! LCOV_EXCL_LINE
730 : d_afsd_latg= d_afsd_latg(i,j,:,iblk),& ! LCOV_EXCL_LINE
731 : d_afsd_newi= d_afsd_newi(i,j,:,iblk),& ! LCOV_EXCL_LINE
732 : d_afsd_latm= d_afsd_latm(i,j,:,iblk),& ! LCOV_EXCL_LINE
733 : d_afsd_weld= d_afsd_weld(i,j,:,iblk),& ! LCOV_EXCL_LINE
734 : floe_rad_c = floe_rad_c(:), & ! LCOV_EXCL_LINE
735 10719744 : floe_binwidth = floe_binwidth(:))
736 : endif ! tmask
737 :
738 : enddo ! i
739 : enddo ! j
740 :
741 22944 : call icepack_warnings_flush(nu_diag)
742 22944 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
743 0 : file=__FILE__, line=__LINE__)
744 :
745 22944 : end subroutine step_therm2
746 :
747 : !=======================================================================
748 : !
749 : ! finalize thermo updates
750 : !
751 : ! authors: Elizabeth Hunke, LANL
752 :
753 11568 : subroutine update_state (dt, daidt, dvidt, dagedt, offset)
754 :
755 : use ice_domain_size, only: ncat
756 : ! use ice_grid, only: tmask
757 : use ice_state, only: aicen, trcrn, vicen, vsnon, &
758 : aice, trcr, vice, vsno, aice0, trcr_depend, & ! LCOV_EXCL_LINE
759 : bound_state, trcr_base, nt_strata, n_trcr_strata
760 : use ice_flux, only: Tf
761 : use ice_timers, only: ice_timer_start, ice_timer_stop, timer_bound, timer_updstate
762 :
763 : real (kind=dbl_kind), intent(in) :: &
764 : dt ! time step
765 :
766 : real (kind=dbl_kind), dimension(:,:,:), intent(inout), optional :: &
767 : daidt, & ! change in ice area per time step ! LCOV_EXCL_LINE
768 : dvidt, & ! change in ice volume per time step ! LCOV_EXCL_LINE
769 : dagedt ! change in ice age per time step
770 :
771 : real (kind=dbl_kind), intent(in), optional :: &
772 : offset ! d(age)/dt time offset = dt for thermo, 0 for dyn
773 :
774 : integer (kind=int_kind) :: &
775 : iblk, & ! block index ! LCOV_EXCL_LINE
776 : i,j, & ! horizontal indices ! LCOV_EXCL_LINE
777 : ntrcr, & ! ! LCOV_EXCL_LINE
778 : nt_iage !
779 :
780 : logical (kind=log_kind) :: &
781 : tr_iage !
782 :
783 : character(len=*), parameter :: subname='(update_state)'
784 :
785 11568 : call ice_timer_start(timer_updstate)
786 11568 : call icepack_query_tracer_flags(tr_iage_out=tr_iage)
787 11568 : call icepack_query_tracer_sizes(ntrcr_out=ntrcr)
788 11568 : call icepack_query_tracer_indices(nt_iage_out=nt_iage)
789 11568 : call icepack_warnings_flush(nu_diag)
790 11568 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
791 0 : file=__FILE__, line=__LINE__)
792 :
793 : !-------------------------------------------------------------------
794 : ! Ghost cell updates for state variables.
795 : !-------------------------------------------------------------------
796 :
797 11568 : call ice_timer_start(timer_bound)
798 : call bound_state (aicen, &
799 : vicen, vsnon, & ! LCOV_EXCL_LINE
800 11568 : ntrcr, trcrn)
801 11568 : call ice_timer_stop(timer_bound)
802 :
803 5760 : !$OMP PARALLEL DO PRIVATE(iblk,i,j) SCHEDULE(runtime)
804 17376 : do iblk = 1, nblocks
805 420480 : do j = 1, ny_block
806 14303760 : do i = 1, nx_block
807 :
808 : !-----------------------------------------------------------------
809 : ! Aggregate the updated state variables (includes ghost cells).
810 : !-----------------------------------------------------------------
811 :
812 : ! if (tmask(i,j,iblk)) &
813 : call icepack_aggregate(ncat = ncat, & ! LCOV_EXCL_LINE
814 : aicen = aicen(i,j,:,iblk), & ! LCOV_EXCL_LINE
815 : trcrn = trcrn(i,j,:,:,iblk), & ! LCOV_EXCL_LINE
816 : vicen = vicen(i,j,:,iblk), & ! LCOV_EXCL_LINE
817 : vsnon = vsnon(i,j,:,iblk), & ! LCOV_EXCL_LINE
818 : aice = aice (i,j, iblk), & ! LCOV_EXCL_LINE
819 : trcr = trcr (i,j,:,iblk), & ! LCOV_EXCL_LINE
820 : vice = vice (i,j, iblk), & ! LCOV_EXCL_LINE
821 : vsno = vsno (i,j, iblk), & ! LCOV_EXCL_LINE
822 : aice0 = aice0(i,j, iblk), & ! LCOV_EXCL_LINE
823 : ntrcr = ntrcr, & ! LCOV_EXCL_LINE
824 : trcr_depend = trcr_depend(:), & ! LCOV_EXCL_LINE
825 : trcr_base = trcr_base(:,:), & ! LCOV_EXCL_LINE
826 : n_trcr_strata = n_trcr_strata(:), & ! LCOV_EXCL_LINE
827 : nt_strata = nt_strata(:,:), & ! LCOV_EXCL_LINE
828 13894848 : Tf = Tf(i,j,iblk))
829 :
830 14292192 : if (present(offset)) then
831 :
832 : !-----------------------------------------------------------------
833 : ! Compute thermodynamic area and volume tendencies.
834 : !-----------------------------------------------------------------
835 :
836 13894848 : daidt(i,j,iblk) = (aice(i,j,iblk) - daidt(i,j,iblk)) / dt
837 13894848 : dvidt(i,j,iblk) = (vice(i,j,iblk) - dvidt(i,j,iblk)) / dt
838 13894848 : if (tr_iage) then
839 13894848 : if (offset > c0) then ! thermo
840 6947424 : if (trcr(i,j,nt_iage,iblk) > c0) &
841 : dagedt(i,j,iblk) = (trcr(i,j,nt_iage,iblk) & ! LCOV_EXCL_LINE
842 6377816 : - dagedt(i,j,iblk) - offset) / dt
843 : else ! dynamics
844 : dagedt(i,j,iblk) = (trcr(i,j,nt_iage,iblk) &
845 6947424 : - dagedt(i,j,iblk)) / dt
846 : endif
847 : endif ! tr_iage
848 : endif ! present(offset)
849 :
850 : enddo ! i
851 : enddo ! j
852 : enddo ! iblk
853 : !$OMP END PARALLEL DO
854 :
855 11568 : call icepack_warnings_flush(nu_diag)
856 11568 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
857 0 : file=__FILE__, line=__LINE__)
858 11568 : call ice_timer_stop(timer_updstate)
859 :
860 11568 : end subroutine update_state
861 :
862 : !=======================================================================
863 : !
864 : ! Run one time step of wave-fracturing the floe size distribution
865 : !
866 : ! authors: Lettie Roach, NIWA
867 : ! Elizabeth C. Hunke, LANL
868 :
869 0 : subroutine step_dyn_wave (dt)
870 :
871 : use ice_arrays_column, only: wave_spectrum, &
872 : d_afsd_wave, floe_rad_l, floe_rad_c, wavefreq, dwavefreq
873 : use ice_domain_size, only: ncat, nfsd, nfreq
874 : use ice_state, only: trcrn, aicen, aice, vice
875 : use ice_timers, only: ice_timer_start, ice_timer_stop, timer_column, &
876 : timer_fsd
877 :
878 : real (kind=dbl_kind), intent(in) :: &
879 : dt ! time step
880 :
881 : ! local variables
882 :
883 : type (block) :: &
884 : this_block ! block information for current block
885 :
886 : integer (kind=int_kind) :: &
887 : ilo,ihi,jlo,jhi, & ! beginning and end of physical domain ! LCOV_EXCL_LINE
888 : iblk, & ! block index ! LCOV_EXCL_LINE
889 : i, j ! horizontal indices
890 :
891 : character (len=char_len) :: wave_spec_type
892 :
893 : character(len=*), parameter :: subname = '(step_dyn_wave)'
894 :
895 0 : call ice_timer_start(timer_column)
896 0 : call ice_timer_start(timer_fsd)
897 :
898 0 : call icepack_query_parameters(wave_spec_type_out=wave_spec_type)
899 0 : call icepack_warnings_flush(nu_diag)
900 0 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
901 0 : file=__FILE__, line=__LINE__)
902 :
903 0 : !$OMP PARALLEL DO PRIVATE(iblk,i,j,ilo,ihi,jlo,jhi,this_block)
904 0 : do iblk = 1, nblocks
905 :
906 0 : this_block = get_block(blocks_ice(iblk),iblk)
907 0 : ilo = this_block%ilo
908 0 : ihi = this_block%ihi
909 0 : jlo = this_block%jlo
910 0 : jhi = this_block%jhi
911 :
912 0 : do j = jlo, jhi
913 0 : do i = ilo, ihi
914 0 : d_afsd_wave(i,j,:,iblk) = c0
915 : call icepack_step_wavefracture (wave_spec_type, &
916 : dt, ncat, nfsd, nfreq, & ! LCOV_EXCL_LINE
917 : aice (i,j, iblk), & ! LCOV_EXCL_LINE
918 : vice (i,j, iblk), & ! LCOV_EXCL_LINE
919 : aicen (i,j,:, iblk), & ! LCOV_EXCL_LINE
920 : floe_rad_l(:), floe_rad_c(:), & ! LCOV_EXCL_LINE
921 : wave_spectrum (i,j,:, iblk), & ! LCOV_EXCL_LINE
922 : wavefreq(:), dwavefreq(:), & ! LCOV_EXCL_LINE
923 : trcrn (i,j,:,:,iblk), & ! LCOV_EXCL_LINE
924 0 : d_afsd_wave (i,j,:, iblk))
925 : end do ! i
926 : end do ! j
927 : end do ! iblk
928 : !$OMP END PARALLEL DO
929 :
930 0 : call ice_timer_stop(timer_fsd)
931 0 : call ice_timer_stop(timer_column)
932 :
933 0 : end subroutine step_dyn_wave
934 :
935 : !=======================================================================
936 : !
937 : ! Run one time step of dynamics and horizontal transport.
938 : ! NOTE: The evp and transport modules include boundary updates, so
939 : ! they cannot be done inside a single block loop.
940 : !
941 : ! authors: William H. Lipscomb, LANL
942 : ! Elizabeth C. Hunke, LANL
943 :
944 5784 : subroutine step_dyn_horiz (dt)
945 :
946 : use ice_boundary, only: ice_HaloUpdate
947 : use ice_dyn_evp, only: evp
948 : use ice_dyn_eap, only: eap
949 : use ice_dyn_vp, only: implicit_solver
950 : use ice_dyn_shared, only: kdyn
951 : use ice_flux, only: strocnxU, strocnyU, strocnxT_iavg, strocnyT_iavg
952 : use ice_flux, only: init_history_dyn
953 : use ice_grid, only: grid_average_X2Y
954 : use ice_state, only: aiU
955 : use ice_transport_driver, only: advection, transport_upwind, transport_remap
956 :
957 : real (kind=dbl_kind), intent(in) :: &
958 : dt ! dynamics time step
959 :
960 : ! local variables
961 :
962 : type (block) :: &
963 : this_block ! block information for current block
964 :
965 : integer (kind=int_kind) :: &
966 : ilo,ihi,jlo,jhi, & ! beginning and end of physical domain ! LCOV_EXCL_LINE
967 : iblk, & ! block index ! LCOV_EXCL_LINE
968 : i, j ! horizontal indices
969 :
970 : real (kind=dbl_kind), dimension (nx_block,ny_block,max_blocks) :: &
971 : work1, & ! temporary ! LCOV_EXCL_LINE
972 5015568 : work2 ! temporary
973 :
974 : character(len=*), parameter :: subname = '(step_dyn_horiz)'
975 :
976 5784 : call init_history_dyn ! initialize dynamic history variables
977 :
978 : !-----------------------------------------------------------------
979 : ! Ice dynamics (momentum equation)
980 : !-----------------------------------------------------------------
981 :
982 5784 : if (kdyn == 1) call evp (dt)
983 5784 : if (kdyn == 2) call eap (dt)
984 5784 : if (kdyn == 3) call implicit_solver (dt)
985 :
986 : !-----------------------------------------------------------------
987 : ! Compute strocnxT_iavg, strocnyT_iavg for thermo and coupling
988 : !-----------------------------------------------------------------
989 :
990 : ! strocn computed on U, N, E as needed. Map strocn U divided by aiU to T
991 : ! conservation requires aiU be divided before averaging
992 16221984 : work1 = c0
993 16221984 : work2 = c0
994 2880 : !$OMP PARALLEL DO PRIVATE(iblk,i,j,ilo,ihi,jlo,jhi,this_block)
995 8688 : do iblk = 1, nblocks
996 5784 : this_block = get_block(blocks_ice(iblk), iblk)
997 5784 : ilo = this_block%ilo
998 5784 : ihi = this_block%ihi
999 5784 : jlo = this_block%jlo
1000 5784 : jhi = this_block%jhi
1001 198672 : do j = jlo, jhi
1002 6369528 : do i = ilo, ihi
1003 6363744 : if (aiU(i,j,iblk) /= c0) then
1004 5840054 : work1(i,j,iblk) = strocnxU(i,j,iblk)/aiU(i,j,iblk)
1005 5840054 : work2(i,j,iblk) = strocnyU(i,j,iblk)/aiU(i,j,iblk)
1006 : endif
1007 : enddo
1008 : enddo
1009 : enddo
1010 : !$OMP END PARALLEL DO
1011 0 : call ice_HaloUpdate (work1, halo_info, &
1012 5784 : field_loc_NEcorner, field_type_vector)
1013 0 : call ice_HaloUpdate (work2, halo_info, &
1014 5784 : field_loc_NEcorner, field_type_vector)
1015 5784 : call grid_average_X2Y('F', work1, 'U', strocnxT_iavg, 'T') ! shift
1016 5784 : call grid_average_X2Y('F', work2, 'U', strocnyT_iavg, 'T')
1017 :
1018 : !-----------------------------------------------------------------
1019 : ! Horizontal ice transport
1020 : !-----------------------------------------------------------------
1021 :
1022 5784 : if (advection == 'upwind') then
1023 0 : call transport_upwind (dt) ! upwind
1024 5784 : elseif (advection == 'remap') then
1025 5784 : call transport_remap (dt) ! incremental remapping
1026 : endif
1027 :
1028 5784 : end subroutine step_dyn_horiz
1029 :
1030 : !=======================================================================
1031 : !
1032 : ! Run one time step of ridging.
1033 : !
1034 : ! authors: William H. Lipscomb, LANL
1035 : ! Elizabeth C. Hunke, LANL
1036 :
1037 22944 : subroutine step_dyn_ridge (dt, ndtd, iblk)
1038 :
1039 : use ice_arrays_column, only: hin_max, first_ice
1040 : use ice_domain_size, only: ncat, nilyr, nslyr, n_aero, nblyr
1041 : use ice_flux, only: &
1042 : rdg_conv, rdg_shear, dardg1dt, dardg2dt, & ! LCOV_EXCL_LINE
1043 : dvirdgdt, opening, fpond, fresh, fhocn, & ! LCOV_EXCL_LINE
1044 : aparticn, krdgn, aredistn, vredistn, dardg1ndt, dardg2ndt, & ! LCOV_EXCL_LINE
1045 : dvirdgndt, araftn, vraftn, fsalt, Tf
1046 : use ice_flux_bgc, only: flux_bio, faero_ocn, fiso_ocn
1047 : use ice_grid, only: tmask
1048 : use ice_state, only: trcrn, vsnon, aicen, vicen, &
1049 : aice, aice0, trcr_depend, n_trcr_strata, & ! LCOV_EXCL_LINE
1050 : trcr_base, nt_strata
1051 : use ice_timers, only: ice_timer_start, ice_timer_stop, timer_column, &
1052 : timer_ridge
1053 :
1054 : real (kind=dbl_kind), intent(in) :: &
1055 : dt ! time step
1056 :
1057 : integer (kind=int_kind), intent(in) :: &
1058 : ndtd, & ! number of dynamics subcycles ! LCOV_EXCL_LINE
1059 : iblk ! block index
1060 :
1061 : ! local variables
1062 :
1063 : type (block) :: &
1064 : this_block ! block information for current block
1065 :
1066 : integer (kind=int_kind) :: &
1067 : ilo,ihi,jlo,jhi, & ! beginning and end of physical domain ! LCOV_EXCL_LINE
1068 : i, j, & ! horizontal indices ! LCOV_EXCL_LINE
1069 : ntrcr, & ! ! LCOV_EXCL_LINE
1070 : nbtrcr !
1071 :
1072 : character(len=*), parameter :: subname = '(step_dyn_ridge)'
1073 :
1074 : !-----------------------------------------------------------------
1075 : ! Ridging
1076 : !-----------------------------------------------------------------
1077 :
1078 22944 : call ice_timer_start(timer_column,iblk)
1079 22944 : call ice_timer_start(timer_ridge,iblk)
1080 :
1081 22944 : call icepack_query_tracer_sizes(ntrcr_out=ntrcr, nbtrcr_out=nbtrcr)
1082 22944 : call icepack_warnings_flush(nu_diag)
1083 22944 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
1084 0 : file=__FILE__, line=__LINE__)
1085 :
1086 22944 : this_block = get_block(blocks_ice(iblk), iblk)
1087 22944 : ilo = this_block%ilo
1088 22944 : ihi = this_block%ihi
1089 22944 : jlo = this_block%jlo
1090 22944 : jhi = this_block%jhi
1091 :
1092 707688 : do j = jlo, jhi
1093 13826928 : do i = ilo, ihi
1094 :
1095 : !echmod: this changes the answers, continue using tmask for now
1096 : ! call aggregate_area (ncat, aicen(i,j,:,iblk), atmp, atmp0)
1097 : ! if (atmp > c0) then
1098 :
1099 13803984 : if (tmask(i,j,iblk)) then
1100 :
1101 : call icepack_step_ridge (dt=dt, ndtd=ndtd, &
1102 : nilyr=nilyr, nslyr=nslyr, nblyr=nblyr, & ! LCOV_EXCL_LINE
1103 : ncat=ncat, n_aero=n_aero, hin_max=hin_max(:), & ! LCOV_EXCL_LINE
1104 : trcr_depend = trcr_depend (:), & ! LCOV_EXCL_LINE
1105 : trcr_base = trcr_base (:,:), & ! LCOV_EXCL_LINE
1106 : n_trcr_strata = n_trcr_strata(:), & ! LCOV_EXCL_LINE
1107 : nt_strata = nt_strata (:,:), & ! LCOV_EXCL_LINE
1108 : trcrn = trcrn (i,j,:,:,iblk), & ! LCOV_EXCL_LINE
1109 : rdg_conv = rdg_conv (i,j, iblk), & ! LCOV_EXCL_LINE
1110 : rdg_shear = rdg_shear(i,j, iblk), & ! LCOV_EXCL_LINE
1111 : aicen = aicen (i,j,:,iblk), & ! LCOV_EXCL_LINE
1112 : vicen = vicen (i,j,:,iblk), & ! LCOV_EXCL_LINE
1113 : vsnon = vsnon (i,j,:,iblk), & ! LCOV_EXCL_LINE
1114 : aice0 = aice0 (i,j, iblk), & ! LCOV_EXCL_LINE
1115 : dardg1dt = dardg1dt (i,j, iblk), & ! LCOV_EXCL_LINE
1116 : dardg2dt = dardg2dt (i,j, iblk), & ! LCOV_EXCL_LINE
1117 : dvirdgdt = dvirdgdt (i,j, iblk), & ! LCOV_EXCL_LINE
1118 : opening = opening (i,j, iblk), & ! LCOV_EXCL_LINE
1119 : fpond = fpond (i,j, iblk), & ! LCOV_EXCL_LINE
1120 : fresh = fresh (i,j, iblk), & ! LCOV_EXCL_LINE
1121 : fhocn = fhocn (i,j, iblk), & ! LCOV_EXCL_LINE
1122 : faero_ocn = faero_ocn(i,j,:,iblk), & ! LCOV_EXCL_LINE
1123 : fiso_ocn = fiso_ocn (i,j,:,iblk), & ! LCOV_EXCL_LINE
1124 : aparticn = aparticn (i,j,:,iblk), & ! LCOV_EXCL_LINE
1125 : krdgn = krdgn (i,j,:,iblk), & ! LCOV_EXCL_LINE
1126 : aredistn = aredistn (i,j,:,iblk), & ! LCOV_EXCL_LINE
1127 : vredistn = vredistn (i,j,:,iblk), & ! LCOV_EXCL_LINE
1128 : dardg1ndt = dardg1ndt(i,j,:,iblk), & ! LCOV_EXCL_LINE
1129 : dardg2ndt = dardg2ndt(i,j,:,iblk), & ! LCOV_EXCL_LINE
1130 : dvirdgndt = dvirdgndt(i,j,:,iblk), & ! LCOV_EXCL_LINE
1131 : araftn = araftn (i,j,:,iblk), & ! LCOV_EXCL_LINE
1132 : vraftn = vraftn (i,j,:,iblk), & ! LCOV_EXCL_LINE
1133 : aice = aice (i,j, iblk), & ! LCOV_EXCL_LINE
1134 : fsalt = fsalt (i,j, iblk), & ! LCOV_EXCL_LINE
1135 : first_ice = first_ice(i,j,:,iblk), & ! LCOV_EXCL_LINE
1136 : flux_bio = flux_bio (i,j,1:nbtrcr,iblk), & ! LCOV_EXCL_LINE
1137 10719744 : Tf = Tf(i,j,iblk))
1138 :
1139 : endif ! tmask
1140 :
1141 : enddo ! i
1142 : enddo ! j
1143 :
1144 22944 : call icepack_warnings_flush(nu_diag)
1145 22944 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
1146 0 : file=__FILE__, line=__LINE__)
1147 :
1148 22944 : call ice_timer_stop(timer_ridge,iblk)
1149 22944 : call ice_timer_stop(timer_column,iblk)
1150 :
1151 22944 : end subroutine step_dyn_ridge
1152 :
1153 : !=======================================================================
1154 : !
1155 : ! Updates snow tracers
1156 : !
1157 : ! authors: Elizabeth C. Hunke, LANL
1158 : ! Nicole Jeffery, LANL
1159 :
1160 0 : subroutine step_snow (dt, iblk)
1161 :
1162 : use ice_calendar, only: nstreams
1163 : use ice_domain_size, only: ncat, nslyr, nilyr
1164 : use ice_flux, only: snwcnt, wind, fresh, fhocn, fsloss, fsnow
1165 : use ice_state, only: trcrn, vsno, vsnon, vicen, aicen, aice
1166 : use icepack_intfc, only: icepack_step_snow
1167 :
1168 : real (kind=dbl_kind), intent(in) :: &
1169 : dt ! time step
1170 :
1171 : integer (kind=int_kind), intent(in) :: &
1172 : iblk ! block index
1173 :
1174 : ! local variables
1175 :
1176 : integer (kind=int_kind) :: &
1177 : nt_smice, nt_smliq, nt_rsnw, & ! LCOV_EXCL_LINE
1178 : nt_Tsfc, nt_qice, nt_sice, nt_qsno, & ! LCOV_EXCL_LINE
1179 : nt_alvl, nt_vlvl, nt_rhos
1180 :
1181 : integer (kind=int_kind) :: &
1182 : ilo,ihi,jlo,jhi, & ! beginning and end of physical domain ! LCOV_EXCL_LINE
1183 : i, j, & ! horizontal indices ! LCOV_EXCL_LINE
1184 : ns ! history streams index
1185 :
1186 : real (kind=dbl_kind) :: &
1187 0 : puny
1188 :
1189 : real (kind=dbl_kind) :: &
1190 0 : fhs ! flag for presence of snow
1191 :
1192 : character(len=*), parameter :: subname = '(step_snow)'
1193 :
1194 : type (block) :: &
1195 : this_block ! block information for current block
1196 :
1197 0 : this_block = get_block(blocks_ice(iblk),iblk)
1198 0 : ilo = this_block%ilo
1199 0 : ihi = this_block%ihi
1200 0 : jlo = this_block%jlo
1201 0 : jhi = this_block%jhi
1202 :
1203 : !-----------------------------------------------------------------
1204 : ! query icepack values
1205 : !-----------------------------------------------------------------
1206 :
1207 0 : call icepack_query_parameters(puny_out=puny)
1208 : call icepack_query_tracer_indices( &
1209 : nt_smice_out=nt_smice, nt_smliq_out=nt_smliq, & ! LCOV_EXCL_LINE
1210 : nt_rsnw_out=nt_rsnw, nt_Tsfc_out=nt_Tsfc, & ! LCOV_EXCL_LINE
1211 : nt_qice_out=nt_qice, nt_sice_out=nt_sice, nt_qsno_out=nt_qsno, & ! LCOV_EXCL_LINE
1212 0 : nt_alvl_out=nt_alvl, nt_vlvl_out=nt_vlvl, nt_rhos_out=nt_rhos)
1213 0 : call icepack_warnings_flush(nu_diag)
1214 0 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
1215 0 : file=__FILE__, line=__LINE__)
1216 :
1217 : !-----------------------------------------------------------------
1218 : ! Snow redistribution and metamorphosis
1219 : !-----------------------------------------------------------------
1220 :
1221 0 : do j = jlo, jhi
1222 0 : do i = ilo, ihi
1223 :
1224 : call icepack_step_snow (dt, nilyr, &
1225 : nslyr, ncat, & ! LCOV_EXCL_LINE
1226 : wind (i,j, iblk), & ! LCOV_EXCL_LINE
1227 : aice (i,j, iblk), & ! LCOV_EXCL_LINE
1228 : aicen(i,j,:,iblk), & ! LCOV_EXCL_LINE
1229 : vicen(i,j,:,iblk), & ! LCOV_EXCL_LINE
1230 : vsnon(i,j,:,iblk), & ! LCOV_EXCL_LINE
1231 : trcrn(i,j,nt_Tsfc,:,iblk), & ! LCOV_EXCL_LINE
1232 : trcrn(i,j,nt_qice,:,iblk), & ! top layer only ! LCOV_EXCL_LINE
1233 : trcrn(i,j,nt_sice,:,iblk), & ! top layer only ! LCOV_EXCL_LINE
1234 : trcrn(i,j,nt_qsno:nt_qsno+nslyr-1,:,iblk), & ! LCOV_EXCL_LINE
1235 : trcrn(i,j,nt_alvl,:,iblk), & ! LCOV_EXCL_LINE
1236 : trcrn(i,j,nt_vlvl,:,iblk), & ! LCOV_EXCL_LINE
1237 : trcrn(i,j,nt_smice:nt_smice+nslyr-1,:,iblk), & ! LCOV_EXCL_LINE
1238 : trcrn(i,j,nt_smliq:nt_smliq+nslyr-1,:,iblk), & ! LCOV_EXCL_LINE
1239 : trcrn(i,j,nt_rsnw:nt_rsnw+nslyr-1,:,iblk), & ! LCOV_EXCL_LINE
1240 : trcrn(i,j,nt_rhos:nt_rhos+nslyr-1,:,iblk), & ! LCOV_EXCL_LINE
1241 : fresh (i,j,iblk), & ! LCOV_EXCL_LINE
1242 : fhocn (i,j,iblk), & ! LCOV_EXCL_LINE
1243 : fsloss (i,j,iblk), & ! LCOV_EXCL_LINE
1244 0 : fsnow (i,j,iblk))
1245 : enddo
1246 : enddo
1247 :
1248 : ! increment counter for history averaging
1249 0 : do j = jlo, jhi
1250 0 : do i = ilo, ihi
1251 0 : fhs = c0
1252 0 : if (vsno(i,j,iblk) > puny) fhs = c1
1253 0 : do ns = 1, nstreams
1254 0 : snwcnt(i,j,iblk,ns) = snwcnt(i,j,iblk,ns) + fhs
1255 : enddo
1256 : enddo
1257 : enddo
1258 :
1259 0 : end subroutine step_snow
1260 :
1261 : !=======================================================================
1262 : !
1263 : ! Computes radiation fields
1264 : !
1265 : ! authors: William H. Lipscomb, LANL
1266 : ! David Bailey, NCAR
1267 : ! Elizabeth C. Hunke, LANL
1268 :
1269 22944 : subroutine step_radiation (dt, iblk)
1270 :
1271 : use ice_arrays_column, only: ffracn, dhsn, &
1272 : fswsfcn, fswintn, fswpenln, Sswabsn, Iswabsn, & ! LCOV_EXCL_LINE
1273 : fswthrun, fswthrun_vdr, fswthrun_vdf, fswthrun_idr, fswthrun_idf, & ! LCOV_EXCL_LINE
1274 : albicen, albsnon, albpndn, & ! LCOV_EXCL_LINE
1275 : alvdrn, alidrn, alvdfn, alidfn, apeffn, trcrn_sw, snowfracn, & ! LCOV_EXCL_LINE
1276 : swgrid, igrid
1277 : use ice_calendar, only: calendar_type, days_per_year, nextsw_cday, yday, msec
1278 : use ice_domain_size, only: ncat, n_aero, nilyr, nslyr, n_zaero, n_algae, nblyr
1279 : use ice_flux, only: swvdr, swvdf, swidr, swidf, coszen, fsnow
1280 : use ice_grid, only: TLAT, TLON, tmask
1281 : use ice_state, only: aicen, vicen, vsnon, trcrn
1282 : use ice_timers, only: ice_timer_start, ice_timer_stop, timer_sw
1283 : use ice_communicate, only: my_task
1284 : use ice_diagnostics, only: npnt, print_points, pmloc, piloc, pjloc
1285 :
1286 : real (kind=dbl_kind), intent(in) :: &
1287 : dt ! time step
1288 :
1289 : integer (kind=int_kind), intent(in) :: &
1290 : iblk ! block index
1291 :
1292 : ! local variables
1293 :
1294 : integer (kind=int_kind) :: &
1295 : ilo,ihi,jlo,jhi, & ! beginning and end of physical domain ! LCOV_EXCL_LINE
1296 : i, j, n, k, & ! horizontal indices ! LCOV_EXCL_LINE
1297 : ipoint ! index for print diagnostic
1298 :
1299 : type (block) :: &
1300 : this_block ! block information for current block
1301 :
1302 : integer (kind=int_kind) :: &
1303 : nt_Tsfc, nt_alvl, nt_rsnw, & ! LCOV_EXCL_LINE
1304 : nt_apnd, nt_hpnd, nt_ipnd, nt_aero, nlt_chl_sw, & ! LCOV_EXCL_LINE
1305 : ntrcr, nbtrcr, nbtrcr_sw, nt_fbri
1306 :
1307 : integer (kind=int_kind), dimension(icepack_max_algae) :: &
1308 : nt_bgc_N
1309 :
1310 : integer (kind=int_kind), dimension(icepack_max_aero) :: &
1311 : nlt_zaero_sw, nt_zaero
1312 :
1313 : logical (kind=log_kind) :: &
1314 : tr_bgc_N, tr_zaero, tr_brine, dEdd_algae, modal_aero, snwgrain
1315 :
1316 : real (kind=dbl_kind), dimension(ncat) :: &
1317 68928 : fbri ! brine height to ice thickness
1318 :
1319 : real(kind= dbl_kind), dimension(:,:), allocatable :: &
1320 : ztrcr_sw, & ! zaerosols (kg/m^3) and chla (mg/m^3) ! LCOV_EXCL_LINE
1321 22944 : rsnow ! snow grain radius tracer (10^-6 m)
1322 :
1323 : logical (kind=log_kind) :: &
1324 : debug, & ! flag for printing debugging information ! LCOV_EXCL_LINE
1325 : l_print_point ! flag for printing debugging information
1326 :
1327 : character(len=*), parameter :: subname = '(step_radiation)'
1328 :
1329 22944 : call ice_timer_start(timer_sw,iblk) ! shortwave
1330 :
1331 : call icepack_query_tracer_sizes(ntrcr_out=ntrcr, &
1332 22944 : nbtrcr_out=nbtrcr, nbtrcr_sw_out=nbtrcr_sw)
1333 : call icepack_query_tracer_flags( &
1334 22944 : tr_brine_out=tr_brine, tr_bgc_N_out=tr_bgc_N, tr_zaero_out=tr_zaero)
1335 : call icepack_query_tracer_indices( &
1336 : nt_Tsfc_out=nt_Tsfc, nt_alvl_out=nt_alvl, nt_rsnw_out=nt_rsnw, & ! LCOV_EXCL_LINE
1337 : nt_apnd_out=nt_apnd, nt_hpnd_out=nt_hpnd, nt_ipnd_out=nt_ipnd, nt_aero_out=nt_aero, & ! LCOV_EXCL_LINE
1338 : nlt_chl_sw_out=nlt_chl_sw, nlt_zaero_sw_out=nlt_zaero_sw, & ! LCOV_EXCL_LINE
1339 22944 : nt_fbri_out=nt_fbri, nt_zaero_out=nt_zaero, nt_bgc_N_out=nt_bgc_N)
1340 : call icepack_query_parameters(dEdd_algae_out=dEdd_algae, modal_aero_out=modal_aero, &
1341 22944 : snwgrain_out=snwgrain)
1342 22944 : call icepack_warnings_flush(nu_diag)
1343 22944 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
1344 0 : file=__FILE__, line=__LINE__)
1345 :
1346 22944 : allocate(ztrcr_sw(nbtrcr_sw,ncat))
1347 22944 : allocate(rsnow(nslyr,ncat))
1348 :
1349 22944 : this_block = get_block(blocks_ice(iblk),iblk)
1350 22944 : ilo = this_block%ilo
1351 22944 : ihi = this_block%ihi
1352 22944 : jlo = this_block%jlo
1353 22944 : jhi = this_block%jhi
1354 :
1355 707688 : do j = jlo, jhi
1356 13826928 : do i = ilo, ihi
1357 :
1358 13119240 : l_print_point = .false.
1359 13119240 : debug = .false.
1360 13119240 : if (debug .and. print_points) then
1361 0 : do ipoint = 1, npnt
1362 0 : if (my_task == pmloc(ipoint) .and. &
1363 : i == piloc(ipoint) .and. & ! LCOV_EXCL_LINE
1364 : j == pjloc(ipoint)) & ! LCOV_EXCL_LINE
1365 0 : l_print_point = .true.
1366 0 : write (nu_diag, *) 'my_task = ',my_task
1367 : enddo ! ipoint
1368 : endif
1369 78715440 : fbri (:) = c0
1370 144311640 : ztrcr_sw(:,:) = c0
1371 144311640 : rsnow (:,:) = c0
1372 78715440 : do n = 1, ncat
1373 65596200 : if (tr_brine) fbri(n) = trcrn(i,j,nt_fbri,n,iblk)
1374 78715440 : if (snwgrain) then
1375 0 : do k = 1, nslyr
1376 0 : rsnow(k,n) = trcrn(i,j,nt_rsnw+k-1,n,iblk)
1377 : enddo
1378 : endif
1379 : enddo
1380 :
1381 13119240 : if (tmask(i,j,iblk)) then
1382 :
1383 : call icepack_step_radiation (dt=dt, &
1384 : swgrid=swgrid(:), igrid=igrid(:), & ! LCOV_EXCL_LINE
1385 : fbri=fbri(:), & ! LCOV_EXCL_LINE
1386 : aicen=aicen(i,j, :,iblk), & ! LCOV_EXCL_LINE
1387 : vicen=vicen(i,j, :,iblk), & ! LCOV_EXCL_LINE
1388 : vsnon=vsnon(i,j, :,iblk), & ! LCOV_EXCL_LINE
1389 : Tsfcn=trcrn(i,j,nt_Tsfc,:,iblk), & ! LCOV_EXCL_LINE
1390 : alvln=trcrn(i,j,nt_alvl,:,iblk), & ! LCOV_EXCL_LINE
1391 : apndn=trcrn(i,j,nt_apnd,:,iblk), & ! LCOV_EXCL_LINE
1392 : hpndn=trcrn(i,j,nt_hpnd,:,iblk), & ! LCOV_EXCL_LINE
1393 : ipndn=trcrn(i,j,nt_ipnd,:,iblk), & ! LCOV_EXCL_LINE
1394 : aeron=trcrn(i,j,nt_aero:nt_aero+4*n_aero-1,:,iblk), & ! LCOV_EXCL_LINE
1395 : bgcNn=trcrn(i,j,nt_bgc_N(1):nt_bgc_N(1)+n_algae*(nblyr+3)-1,:,iblk), & ! LCOV_EXCL_LINE
1396 : zaeron=trcrn(i,j,nt_zaero(1):nt_zaero(1)+n_zaero*(nblyr+3)-1,:,iblk), & ! LCOV_EXCL_LINE
1397 : trcrn_bgcsw=ztrcr_sw, & ! LCOV_EXCL_LINE
1398 : TLAT=TLAT(i,j,iblk), TLON=TLON(i,j,iblk), & ! LCOV_EXCL_LINE
1399 : calendar_type=calendar_type, & ! LCOV_EXCL_LINE
1400 : days_per_year=days_per_year, & ! LCOV_EXCL_LINE
1401 : nextsw_cday=nextsw_cday, yday=yday, & ! LCOV_EXCL_LINE
1402 : sec=msec, & ! LCOV_EXCL_LINE
1403 : swvdr =swvdr (i,j ,iblk), swvdf =swvdf (i,j ,iblk), & ! LCOV_EXCL_LINE
1404 : swidr =swidr (i,j ,iblk), swidf =swidf (i,j ,iblk), & ! LCOV_EXCL_LINE
1405 : coszen =coszen (i,j ,iblk), fsnow =fsnow (i,j ,iblk), & ! LCOV_EXCL_LINE
1406 : alvdrn =alvdrn (i,j,: ,iblk), alvdfn =alvdfn (i,j,: ,iblk), & ! LCOV_EXCL_LINE
1407 : alidrn =alidrn (i,j,: ,iblk), alidfn =alidfn (i,j,: ,iblk), & ! LCOV_EXCL_LINE
1408 : fswsfcn =fswsfcn (i,j,: ,iblk), fswintn =fswintn (i,j,: ,iblk), & ! LCOV_EXCL_LINE
1409 : fswthrun =fswthrun (i,j,: ,iblk), & ! LCOV_EXCL_LINE
1410 : fswthrun_vdr =fswthrun_vdr (i,j,: ,iblk), & ! LCOV_EXCL_LINE
1411 : fswthrun_vdf =fswthrun_vdf (i,j,: ,iblk), & ! LCOV_EXCL_LINE
1412 : fswthrun_idr =fswthrun_idr (i,j,: ,iblk), & ! LCOV_EXCL_LINE
1413 : fswthrun_idf =fswthrun_idf (i,j,: ,iblk), & ! LCOV_EXCL_LINE
1414 : fswpenln=fswpenln(i,j,:,:,iblk), & ! LCOV_EXCL_LINE
1415 : Sswabsn =Sswabsn (i,j,:,:,iblk), Iswabsn =Iswabsn (i,j,:,:,iblk), & ! LCOV_EXCL_LINE
1416 : albicen =albicen (i,j,: ,iblk), albsnon =albsnon (i,j,: ,iblk), & ! LCOV_EXCL_LINE
1417 : albpndn =albpndn (i,j,: ,iblk), apeffn =apeffn (i,j,: ,iblk), & ! LCOV_EXCL_LINE
1418 : snowfracn=snowfracn(i,j,: ,iblk), & ! LCOV_EXCL_LINE
1419 : dhsn =dhsn (i,j,: ,iblk), ffracn =ffracn(i,j,:,iblk), & ! LCOV_EXCL_LINE
1420 10719744 : rsnow =rsnow (:,:), l_print_point=l_print_point)
1421 : endif
1422 :
1423 13803984 : if (dEdd_algae .and. (tr_zaero .or. tr_bgc_N)) then
1424 0 : do n = 1, ncat
1425 0 : do k = 1, nbtrcr_sw
1426 0 : trcrn_sw(i,j,k,n,iblk) = ztrcr_sw(k,n)
1427 : enddo
1428 : enddo
1429 : endif
1430 :
1431 : enddo ! i
1432 : enddo ! j
1433 :
1434 22944 : call icepack_warnings_flush(nu_diag)
1435 22944 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
1436 0 : file=__FILE__, line=__LINE__)
1437 :
1438 22944 : deallocate(ztrcr_sw)
1439 22944 : deallocate(rsnow)
1440 :
1441 22944 : call ice_timer_stop(timer_sw,iblk) ! shortwave
1442 :
1443 68832 : end subroutine step_radiation
1444 :
1445 : !=======================================================================
1446 : ! Ocean mixed layer calculation (internal to sea ice model).
1447 : ! Allows heat storage in ocean for uncoupled runs.
1448 : !
1449 : ! authors: John Weatherly, CRREL
1450 : ! C.M. Bitz, UW
1451 : ! Elizabeth C. Hunke, LANL
1452 : ! Bruce P. Briegleb, NCAR
1453 : ! William H. Lipscomb, LANL
1454 :
1455 22944 : subroutine ocean_mixed_layer (dt, iblk)
1456 :
1457 : use ice_arrays_column, only: Cdn_atm, Cdn_atm_ratio
1458 : use ice_flux, only: sst, Tf, Qa, uatmT, vatmT, wind, potT, rhoa, zlvl, &
1459 : frzmlt, fhocn, fswthru, flw, flwout_ocn, fsens_ocn, flat_ocn, evap_ocn, & ! LCOV_EXCL_LINE
1460 : alvdr_ocn, alidr_ocn, alvdf_ocn, alidf_ocn, swidf, swvdf, swidr, swvdr, & ! LCOV_EXCL_LINE
1461 : qdp, hmix, strairx_ocn, strairy_ocn, Tref_ocn, Qref_ocn
1462 : use ice_grid, only: tmask
1463 : use ice_state, only: aice
1464 :
1465 : real (kind=dbl_kind), intent(in) :: &
1466 : dt ! time step
1467 :
1468 : integer (kind=int_kind), intent(in) :: &
1469 : iblk ! block index
1470 :
1471 : ! local variables
1472 :
1473 5760 : real (kind=dbl_kind) :: albocn
1474 :
1475 : real (kind=dbl_kind), parameter :: &
1476 : frzmlt_max = c1000 ! max magnitude of frzmlt (W/m^2)
1477 :
1478 : integer (kind=int_kind) :: &
1479 : i, j , & ! horizontal indices ! LCOV_EXCL_LINE
1480 : ij ! combined ij index
1481 :
1482 : real (kind=dbl_kind), dimension(nx_block,ny_block) :: &
1483 : delt , & ! potential temperature difference (K) ! LCOV_EXCL_LINE
1484 : delq , & ! specific humidity difference (kg/kg) ! LCOV_EXCL_LINE
1485 : shcoef, & ! transfer coefficient for sensible heat ! LCOV_EXCL_LINE
1486 5039808 : lhcoef ! transfer coefficient for latent heat
1487 :
1488 : integer (kind=int_kind) :: &
1489 : icells ! number of ocean cells
1490 :
1491 : integer (kind=int_kind), dimension(nx_block*ny_block) :: &
1492 40128 : indxi, indxj ! compressed indices for ocean cells
1493 :
1494 : character(len=*), parameter :: subname = '(ocn_mixed_layer)'
1495 :
1496 : !-----------------------------------------------------------------
1497 :
1498 22944 : call icepack_query_parameters(albocn_out=albocn)
1499 22944 : call icepack_warnings_flush(nu_diag)
1500 22944 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
1501 0 : file=__FILE__, line=__LINE__)
1502 :
1503 : !-----------------------------------------------------------------
1504 : ! Identify ocean cells.
1505 : ! Set fluxes to zero in land cells.
1506 : !-----------------------------------------------------------------
1507 :
1508 22944 : icells = 0
1509 15455688 : indxi(:) = 0
1510 15455688 : indxj(:) = 0
1511 :
1512 753576 : do j = 1, ny_block
1513 16186320 : do i = 1, nx_block
1514 :
1515 16163376 : if (tmask(i,j,iblk)) then
1516 12324840 : icells = icells + 1
1517 12324840 : indxi(icells) = i
1518 12324840 : indxj(icells) = j
1519 : else
1520 3107904 : sst (i,j,iblk) = c0
1521 3107904 : frzmlt (i,j,iblk) = c0
1522 3107904 : flwout_ocn(i,j,iblk) = c0
1523 3107904 : fsens_ocn (i,j,iblk) = c0
1524 3107904 : flat_ocn (i,j,iblk) = c0
1525 3107904 : evap_ocn (i,j,iblk) = c0
1526 : endif
1527 : enddo ! i
1528 : enddo ! j
1529 :
1530 : !-----------------------------------------------------------------
1531 : ! Compute boundary layer quantities
1532 : !-----------------------------------------------------------------
1533 :
1534 12347784 : do ij = 1, icells
1535 12324840 : i = indxi(ij)
1536 12324840 : j = indxj(ij)
1537 :
1538 : call icepack_atm_boundary(sfctype = 'ocn', &
1539 : Tsf = sst (i,j,iblk), & ! LCOV_EXCL_LINE
1540 : potT = potT (i,j,iblk), & ! LCOV_EXCL_LINE
1541 : uatm = uatmT (i,j,iblk), & ! LCOV_EXCL_LINE
1542 : vatm = vatmT (i,j,iblk), & ! LCOV_EXCL_LINE
1543 : wind = wind (i,j,iblk), & ! LCOV_EXCL_LINE
1544 : zlvl = zlvl (i,j,iblk), & ! LCOV_EXCL_LINE
1545 : Qa = Qa (i,j,iblk), & ! LCOV_EXCL_LINE
1546 : rhoa = rhoa (i,j,iblk), & ! LCOV_EXCL_LINE
1547 : strx = strairx_ocn(i,j,iblk), & ! LCOV_EXCL_LINE
1548 : stry = strairy_ocn(i,j,iblk), & ! LCOV_EXCL_LINE
1549 : Tref = Tref_ocn (i,j,iblk), & ! LCOV_EXCL_LINE
1550 : Qref = Qref_ocn (i,j,iblk), & ! LCOV_EXCL_LINE
1551 : delt = delt (i,j), & ! LCOV_EXCL_LINE
1552 : delq = delq (i,j), & ! LCOV_EXCL_LINE
1553 : lhcoef = lhcoef (i,j), & ! LCOV_EXCL_LINE
1554 : shcoef = shcoef (i,j), & ! LCOV_EXCL_LINE
1555 : Cdn_atm = Cdn_atm (i,j,iblk), & ! LCOV_EXCL_LINE
1556 15633504 : Cdn_atm_ratio_n = Cdn_atm_ratio(i,j,iblk))
1557 : enddo ! ij
1558 :
1559 22944 : call icepack_warnings_flush(nu_diag)
1560 22944 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
1561 0 : file=__FILE__, line=__LINE__)
1562 :
1563 : !-----------------------------------------------------------------
1564 : ! Ocean albedo
1565 : ! For now, assume albedo = albocn in each spectral band.
1566 : !-----------------------------------------------------------------
1567 :
1568 16186320 : alvdr_ocn(:,:,iblk) = albocn
1569 16186320 : alidr_ocn(:,:,iblk) = albocn
1570 16186320 : alvdf_ocn(:,:,iblk) = albocn
1571 16186320 : alidf_ocn(:,:,iblk) = albocn
1572 :
1573 : !-----------------------------------------------------------------
1574 : ! Compute ocean fluxes and update SST
1575 : !-----------------------------------------------------------------
1576 :
1577 12347784 : do ij = 1, icells
1578 12324840 : i = indxi(ij)
1579 12324840 : j = indxj(ij)
1580 :
1581 0 : call icepack_ocn_mixed_layer(alvdr_ocn=alvdr_ocn(i,j,iblk), swvdr =swvdr (i,j,iblk), &
1582 : alidr_ocn=alidr_ocn(i,j,iblk), swidr =swidr (i,j,iblk), & ! LCOV_EXCL_LINE
1583 : alvdf_ocn=alvdf_ocn(i,j,iblk), swvdf =swvdf (i,j,iblk), & ! LCOV_EXCL_LINE
1584 : alidf_ocn=alidf_ocn(i,j,iblk), swidf =swidf (i,j,iblk), & ! LCOV_EXCL_LINE
1585 : sst =sst (i,j,iblk), flwout_ocn=flwout_ocn(i,j,iblk), & ! LCOV_EXCL_LINE
1586 : fsens_ocn=fsens_ocn(i,j,iblk), shcoef=shcoef(i,j), & ! LCOV_EXCL_LINE
1587 : flat_ocn =flat_ocn (i,j,iblk), lhcoef=lhcoef(i,j), & ! LCOV_EXCL_LINE
1588 : evap_ocn =evap_ocn (i,j,iblk), flw =flw (i,j,iblk), & ! LCOV_EXCL_LINE
1589 : delt =delt (i,j), delq =delq (i,j), & ! LCOV_EXCL_LINE
1590 : aice =aice (i,j,iblk), fhocn =fhocn (i,j,iblk), & ! LCOV_EXCL_LINE
1591 : fswthru =fswthru (i,j,iblk), hmix =hmix (i,j,iblk), & ! LCOV_EXCL_LINE
1592 : Tf =Tf (i,j,iblk), qdp =qdp (i,j,iblk), & ! LCOV_EXCL_LINE
1593 22204944 : frzmlt =frzmlt (i,j,iblk), dt =dt)
1594 : enddo ! ij
1595 :
1596 22944 : call icepack_warnings_flush(nu_diag)
1597 22944 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
1598 0 : file=__FILE__, line=__LINE__)
1599 :
1600 22944 : end subroutine ocean_mixed_layer
1601 :
1602 : !=======================================================================
1603 :
1604 22944 : subroutine biogeochemistry (dt, iblk)
1605 :
1606 : use ice_arrays_column, only: upNO, upNH, iDi, iki, zfswin, &
1607 : darcy_V, grow_net, & ! LCOV_EXCL_LINE
1608 : PP_net, hbri,dhbr_bot, dhbr_top, Zoo,& ! LCOV_EXCL_LINE
1609 : fbio_snoice, fbio_atmice, ocean_bio, & ! LCOV_EXCL_LINE
1610 : first_ice, fswpenln, bphi, bTiz, ice_bio_net, & ! LCOV_EXCL_LINE
1611 : snow_bio_net, fswthrun, & ! LCOV_EXCL_LINE
1612 : ocean_bio_all, sice_rho, & ! LCOV_EXCL_LINE
1613 : bgrid, igrid, icgrid, cgrid
1614 : use ice_domain_size, only: nblyr, nilyr, nslyr, n_algae, n_zaero, ncat, &
1615 : n_doc, n_dic, n_don, n_fed, n_fep
1616 : use ice_flux, only: meltbn, melttn, congeln, snoicen, &
1617 : sst, sss, fsnow, meltsn
1618 : use ice_flux_bgc, only: hin_old, flux_bio, flux_bio_atm, faero_atm, &
1619 : nit, amm, sil, dmsp, dms, algalN, doc, don, dic, fed, fep, zaeros, hum
1620 : use ice_state, only: aicen_init, vicen_init, aicen, vicen, vsnon, &
1621 : trcrn, vsnon_init, aice0
1622 : use ice_timers, only: timer_bgc, ice_timer_start, ice_timer_stop
1623 :
1624 : real (kind=dbl_kind), intent(in) :: &
1625 : dt ! time step
1626 :
1627 : integer (kind=int_kind), intent(in) :: &
1628 : iblk ! block index
1629 :
1630 : ! local variables
1631 :
1632 : integer (kind=int_kind) :: &
1633 : i, j , & ! horizontal indices ! LCOV_EXCL_LINE
1634 : ilo,ihi,jlo,jhi, & ! beginning and end of physical domain ! LCOV_EXCL_LINE
1635 : mm ! tracer index
1636 :
1637 : type (block) :: &
1638 : this_block ! block information for current block
1639 :
1640 : integer (kind=int_kind) :: &
1641 : nbtrcr, ntrcr
1642 :
1643 : integer (kind=int_kind), dimension(icepack_max_aero) :: &
1644 : nlt_zaero
1645 :
1646 : integer (kind=int_kind), dimension(icepack_max_nbtrcr) :: &
1647 : bio_index_o
1648 :
1649 : logical (kind=log_kind) :: &
1650 : skl_bgc, tr_brine, tr_zaero
1651 :
1652 : character(len=*), parameter :: subname='(biogeochemistry)'
1653 :
1654 22944 : call icepack_query_tracer_flags(tr_brine_out=tr_brine)
1655 22944 : call icepack_query_parameters(skl_bgc_out=skl_bgc)
1656 22944 : call icepack_query_tracer_sizes(ntrcr_out=ntrcr, nbtrcr_out=nbtrcr)
1657 22944 : call icepack_query_tracer_flags(tr_zaero_out=tr_zaero)
1658 22944 : call icepack_query_tracer_indices(nlt_zaero_out=nlt_zaero)
1659 22944 : call icepack_query_tracer_indices(bio_index_o_out=bio_index_o)
1660 22944 : call icepack_warnings_flush(nu_diag)
1661 22944 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
1662 0 : file=__FILE__, line=__LINE__)
1663 :
1664 22944 : if (tr_brine .or. skl_bgc) then
1665 :
1666 0 : call ice_timer_start(timer_bgc,iblk) ! biogeochemistry
1667 :
1668 0 : this_block = get_block(blocks_ice(iblk),iblk)
1669 0 : ilo = this_block%ilo
1670 0 : ihi = this_block%ihi
1671 0 : jlo = this_block%jlo
1672 0 : jhi = this_block%jhi
1673 :
1674 : ! Define ocean concentrations for tracers used in simulation
1675 0 : do j = jlo, jhi
1676 0 : do i = ilo, ihi
1677 :
1678 : call icepack_load_ocean_bio_array(max_nbtrcr = icepack_max_nbtrcr, &
1679 : max_algae = icepack_max_algae, max_don = icepack_max_don, & ! LCOV_EXCL_LINE
1680 : max_doc = icepack_max_doc, max_dic = icepack_max_dic, & ! LCOV_EXCL_LINE
1681 : max_aero = icepack_max_aero, max_fe = icepack_max_fe, & ! LCOV_EXCL_LINE
1682 : nit = nit(i,j, iblk), amm = amm (i,j, iblk), & ! LCOV_EXCL_LINE
1683 : sil = sil(i,j, iblk), dmsp = dmsp (i,j, iblk), & ! LCOV_EXCL_LINE
1684 : dms = dms(i,j, iblk), algalN = algalN(i,j,:,iblk), & ! LCOV_EXCL_LINE
1685 : doc = doc(i,j,:,iblk), don = don (i,j,:,iblk), & ! LCOV_EXCL_LINE
1686 : dic = dic(i,j,:,iblk), fed = fed (i,j,:,iblk), & ! LCOV_EXCL_LINE
1687 : fep = fep(i,j,:,iblk), zaeros = zaeros(i,j,:,iblk), & ! LCOV_EXCL_LINE
1688 : hum = hum(i,j, iblk), & ! LCOV_EXCL_LINE
1689 0 : ocean_bio_all = ocean_bio_all(i,j,:,iblk))
1690 :
1691 0 : do mm = 1,nbtrcr
1692 0 : ocean_bio(i,j,mm,iblk) = ocean_bio_all(i,j,bio_index_o(mm),iblk)
1693 : enddo ! mm
1694 0 : if (tr_zaero) then
1695 0 : do mm = 1, n_zaero ! update aerosols
1696 0 : flux_bio_atm(i,j,nlt_zaero(mm),iblk) = faero_atm(i,j,mm,iblk)
1697 : enddo ! mm
1698 : endif
1699 :
1700 : call icepack_biogeochemistry(dt=dt, ntrcr=ntrcr, nbtrcr=nbtrcr,&
1701 : bgrid=bgrid, igrid=igrid, icgrid=icgrid, cgrid=cgrid, & ! LCOV_EXCL_LINE
1702 : nblyr=nblyr, nilyr=nilyr, nslyr=nslyr, n_algae=n_algae, n_zaero=n_zaero, & ! LCOV_EXCL_LINE
1703 : ncat=ncat, n_doc=n_doc, n_dic=n_dic, n_don=n_don, n_fed=n_fed, n_fep=n_fep, & ! LCOV_EXCL_LINE
1704 : upNO = upNO (i,j, iblk), & ! LCOV_EXCL_LINE
1705 : upNH = upNH (i,j, iblk), & ! LCOV_EXCL_LINE
1706 : iDi = iDi (i,j,:,:, iblk), & ! LCOV_EXCL_LINE
1707 : iki = iki (i,j,:,:, iblk), & ! LCOV_EXCL_LINE
1708 : zfswin = zfswin (i,j,:,:, iblk), & ! LCOV_EXCL_LINE
1709 : darcy_V = darcy_V (i,j,:, iblk), & ! LCOV_EXCL_LINE
1710 : grow_net = grow_net (i,j, iblk), & ! LCOV_EXCL_LINE
1711 : PP_net = PP_net (i,j, iblk), & ! LCOV_EXCL_LINE
1712 : hbri = hbri (i,j, iblk), & ! LCOV_EXCL_LINE
1713 : dhbr_bot = dhbr_bot (i,j,:, iblk), & ! LCOV_EXCL_LINE
1714 : dhbr_top = dhbr_top (i,j,:, iblk), & ! LCOV_EXCL_LINE
1715 : Zoo = Zoo (i,j,:,:, iblk), & ! LCOV_EXCL_LINE
1716 : fbio_snoice = fbio_snoice (i,j,:, iblk), & ! LCOV_EXCL_LINE
1717 : fbio_atmice = fbio_atmice (i,j,:, iblk), & ! LCOV_EXCL_LINE
1718 : ocean_bio = ocean_bio (i,j,1:nbtrcr, iblk), & ! LCOV_EXCL_LINE
1719 : first_ice = first_ice (i,j,:, iblk), & ! LCOV_EXCL_LINE
1720 : fswpenln = fswpenln (i,j,:,:, iblk), & ! LCOV_EXCL_LINE
1721 : bphi = bphi (i,j,:,:, iblk), & ! LCOV_EXCL_LINE
1722 : bTiz = bTiz (i,j,:,:, iblk), & ! LCOV_EXCL_LINE
1723 : ice_bio_net = ice_bio_net (i,j,1:nbtrcr, iblk), & ! LCOV_EXCL_LINE
1724 : snow_bio_net = snow_bio_net(i,j,1:nbtrcr, iblk), & ! LCOV_EXCL_LINE
1725 : fswthrun = fswthrun (i,j,:, iblk), & ! LCOV_EXCL_LINE
1726 : sice_rho = sice_rho (i,j,:, iblk), & ! LCOV_EXCL_LINE
1727 : meltbn = meltbn (i,j,:, iblk), & ! LCOV_EXCL_LINE
1728 : melttn = melttn (i,j,:, iblk), & ! LCOV_EXCL_LINE
1729 : congeln = congeln (i,j,:, iblk), & ! LCOV_EXCL_LINE
1730 : snoicen = snoicen (i,j,:, iblk), & ! LCOV_EXCL_LINE
1731 : sst = sst (i,j, iblk), & ! LCOV_EXCL_LINE
1732 : sss = sss (i,j, iblk), & ! LCOV_EXCL_LINE
1733 : fsnow = fsnow (i,j, iblk), & ! LCOV_EXCL_LINE
1734 : meltsn = meltsn (i,j,:, iblk), & ! LCOV_EXCL_LINE
1735 : hin_old = hin_old (i,j,:, iblk), & ! LCOV_EXCL_LINE
1736 : flux_bio = flux_bio (i,j,1:nbtrcr, iblk), & ! LCOV_EXCL_LINE
1737 : flux_bio_atm = flux_bio_atm(i,j,1:nbtrcr, iblk), & ! LCOV_EXCL_LINE
1738 : aicen_init = aicen_init (i,j,:, iblk), & ! LCOV_EXCL_LINE
1739 : vicen_init = vicen_init (i,j,:, iblk), & ! LCOV_EXCL_LINE
1740 : aicen = aicen (i,j,:, iblk), & ! LCOV_EXCL_LINE
1741 : vicen = vicen (i,j,:, iblk), & ! LCOV_EXCL_LINE
1742 : vsnon = vsnon (i,j,:, iblk), & ! LCOV_EXCL_LINE
1743 : aice0 = aice0 (i,j, iblk), & ! LCOV_EXCL_LINE
1744 : trcrn = trcrn (i,j,:,:, iblk), & ! LCOV_EXCL_LINE
1745 : vsnon_init = vsnon_init (i,j,:, iblk), & ! LCOV_EXCL_LINE
1746 0 : skl_bgc = skl_bgc)
1747 :
1748 : enddo ! i
1749 : enddo ! j
1750 :
1751 0 : call icepack_warnings_flush(nu_diag)
1752 0 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
1753 0 : file=__FILE__, line=__LINE__)
1754 :
1755 0 : call ice_timer_stop(timer_bgc,iblk) ! biogeochemistry
1756 :
1757 : endif ! tr_brine .or. skl_bgc
1758 :
1759 22944 : end subroutine biogeochemistry
1760 :
1761 : !=======================================================================
1762 :
1763 : end module ice_step_mod
1764 :
1765 : !=======================================================================
|