Line data Source code
1 : !=======================================================================
2 : !
3 : ! Contains Icepack component driver routines common to all drivers.
4 : !
5 : ! authors Elizabeth C. Hunke, LANL
6 :
7 : module icedrv_step
8 :
9 : use icedrv_constants, only: c0, nu_diag, c4
10 : use icedrv_kinds
11 : ! use icedrv_calendar, only: istep1
12 : use icedrv_forcing, only: ocn_data_type
13 : use icedrv_system, only: icedrv_system_abort
14 : use icepack_intfc, only: icepack_warnings_flush
15 : use icepack_intfc, only: icepack_warnings_aborted
16 : use icepack_intfc, only: icepack_query_tracer_flags
17 : use icepack_intfc, only: icepack_query_tracer_indices
18 : use icepack_intfc, only: icepack_query_tracer_sizes
19 : use icepack_intfc, only: icepack_query_parameters
20 :
21 : implicit none
22 : private
23 :
24 : public :: step_therm1, step_therm2, step_dyn_ridge, &
25 : prep_radiation, step_radiation, ocean_mixed_layer, &
26 : update_state, biogeochemistry, step_dyn_wave
27 :
28 : !=======================================================================
29 :
30 : contains
31 :
32 : !=======================================================================
33 : !
34 : ! Scales radiation fields computed on the previous time step.
35 : !
36 : ! authors: Elizabeth Hunke, LANL
37 :
38 346176 : subroutine prep_radiation ()
39 :
40 : use icedrv_domain_size, only: ncat, nilyr, nslyr, nx
41 : use icedrv_flux, only: scale_factor, swvdr, swvdf, swidr, swidf
42 : use icedrv_flux, only: alvdr_ai, alvdf_ai, alidr_ai, alidf_ai
43 : use icedrv_flux, only: alvdr_init, alvdf_init, alidr_init, alidf_init
44 : use icedrv_arrays_column, only: fswsfcn, fswintn
45 : use icedrv_arrays_column, only: fswthrun, fswthrun_vdr, fswthrun_vdf, fswthrun_idr, fswthrun_idf
46 : use icedrv_arrays_column, only: fswpenln, Sswabsn, Iswabsn
47 : use icedrv_state, only: aice, aicen
48 :
49 : ! column package includes
50 : use icepack_intfc, only: icepack_prep_radiation
51 :
52 : ! local variables
53 :
54 : integer (kind=int_kind) :: &
55 : i ! horizontal indices
56 :
57 : character(len=*), parameter :: subname='(prep_radiation)'
58 :
59 : !-----------------------------------------------------------------
60 : ! Compute netsw scaling factor (new netsw / old netsw)
61 : !-----------------------------------------------------------------
62 :
63 1730880 : do i = 1, nx
64 :
65 1384704 : alvdr_init(i) = alvdr_ai(i)
66 1384704 : alvdf_init(i) = alvdf_ai(i)
67 1384704 : alidr_init(i) = alidr_ai(i)
68 1384704 : alidf_init(i) = alidf_ai(i)
69 :
70 : call icepack_prep_radiation(ncat=ncat, nilyr=nilyr, nslyr=nslyr, &
71 509328 : aice=aice(i), aicen=aicen(i,:), &
72 509328 : swvdr=swvdr(i), swvdf=swvdf(i), &
73 509328 : swidr=swidr(i), swidf=swidf(i), &
74 509328 : alvdr_ai=alvdr_ai(i), alvdf_ai=alvdf_ai(i), &
75 509328 : alidr_ai=alidr_ai(i), alidf_ai=alidf_ai(i), &
76 509328 : scale_factor=scale_factor(i), &
77 : fswsfcn=fswsfcn(i,:), fswintn=fswintn(i,:), &
78 : fswthrun=fswthrun(i,:), &
79 : fswthrun_vdr=fswthrun_vdr(i,:), &
80 : fswthrun_vdf=fswthrun_vdf(i,:), &
81 : fswthrun_idr=fswthrun_idr(i,:), &
82 : fswthrun_idf=fswthrun_idf(i,:), &
83 : fswpenln=fswpenln(i,:,:), &
84 2240208 : Sswabsn=Sswabsn(i,:,:), Iswabsn=Iswabsn(i,:,:))
85 :
86 : enddo ! i
87 346176 : call icepack_warnings_flush(nu_diag)
88 346176 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
89 0 : file=__FILE__, line=__LINE__)
90 :
91 346176 : end subroutine prep_radiation
92 :
93 : !=======================================================================
94 : !
95 : ! Driver for updating ice and snow internal temperatures and
96 : ! computing thermodynamic growth rates and coupler fluxes.
97 : !
98 : ! authors: William H. Lipscomb, LANL
99 :
100 346176 : subroutine step_therm1 (dt)
101 :
102 : use icedrv_arrays_column, only: ffracn, dhsn
103 : use icedrv_arrays_column, only: Cdn_ocn, Cdn_ocn_skin, Cdn_ocn_floe
104 : use icedrv_arrays_column, only: Cdn_ocn_keel, Cdn_atm_ratio
105 : use icedrv_arrays_column, only: Cdn_atm, Cdn_atm_skin, Cdn_atm_floe
106 : use icedrv_arrays_column, only: Cdn_atm_rdg, Cdn_atm_pond
107 : use icedrv_arrays_column, only: hfreebd, hdraft, hridge, distrdg
108 : use icedrv_arrays_column, only: hkeel, dkeel, lfloe, dfloe
109 : use icedrv_arrays_column, only: fswsfcn, fswintn, Sswabsn, Iswabsn
110 : use icedrv_arrays_column, only: fswthrun, fswthrun_vdr, fswthrun_vdf, fswthrun_idr, fswthrun_idf
111 : use icedrv_calendar, only: yday
112 : use icedrv_domain_size, only: ncat, nilyr, nslyr, n_aero, n_iso, nx
113 : use icedrv_flux, only: frzmlt, sst, Tf, strocnxT, strocnyT, rside, fside, &
114 : fbot, Tbot, Tsnice
115 : use icedrv_flux, only: meltsn, melttn, meltbn, congeln, snoicen, uatm, vatm
116 : use icedrv_flux, only: wind, rhoa, potT, Qa, Qa_iso, zlvl, strax, stray, flatn
117 : use icedrv_flux, only: fsensn, fsurfn, fcondtopn, fcondbotn
118 : use icedrv_flux, only: flw, fsnow, fpond, sss, mlt_onset, frz_onset
119 : use icedrv_flux, only: frain, Tair, strairxT, strairyT, fsurf
120 : use icedrv_flux, only: fcondtop, fcondbot, fsens, fresh, fsalt, fhocn
121 : use icedrv_flux, only: flat, fswabs, flwout, evap, evaps, evapi
122 : use icedrv_flux, only: Tref, Qref, Qref_iso, Uref
123 : use icedrv_flux, only: meltt, melts, meltb, congel, snoice
124 : use icedrv_flux, only: fswthru, fswthru_vdr, fswthru_vdf, fswthru_idr, fswthru_idf
125 : use icedrv_flux, only: flatn_f, fsensn_f, fsurfn_f, fcondtopn_f
126 : use icedrv_flux, only: dsnown, faero_atm, faero_ocn
127 : use icedrv_flux, only: fiso_atm, fiso_ocn, fiso_evap
128 : use icedrv_flux, only: HDO_ocn, H2_16O_ocn, H2_18O_ocn
129 : use icedrv_init, only: lmask_n, lmask_s
130 : use icedrv_state, only: aice, aicen, aice_init, aicen_init, vicen_init
131 : use icedrv_state, only: vice, vicen, vsno, vsnon, trcrn, uvel, vvel, vsnon_init
132 :
133 : ! column packge includes
134 : use icepack_intfc, only: icepack_step_therm1
135 :
136 : logical (kind=log_kind) :: &
137 : prescribed_ice ! if .true., use prescribed ice instead of computed
138 :
139 : real (kind=dbl_kind), intent(in) :: &
140 : dt ! time step
141 :
142 : ! local variables
143 :
144 : integer (kind=int_kind) :: &
145 : i , & ! horizontal indices
146 : n , & ! thickness category index
147 : k, kk ! indices for aerosols
148 :
149 : integer (kind=int_kind) :: &
150 : ntrcr, nt_apnd, nt_hpnd, nt_ipnd, nt_alvl, nt_vlvl, nt_Tsfc, &
151 : nt_iage, nt_FY, nt_qice, nt_sice, nt_qsno, &
152 : nt_aero, nt_isosno, nt_isoice
153 :
154 : logical (kind=log_kind) :: &
155 : tr_iage, tr_FY, tr_aero, tr_iso, tr_pond, tr_pond_cesm, &
156 : tr_pond_lvl, tr_pond_topo, calc_Tsfc
157 :
158 : real (kind=dbl_kind), dimension(n_aero,2,ncat) :: &
159 6006732 : aerosno, aeroice ! kg/m^2
160 :
161 : real (kind=dbl_kind), dimension(n_iso,ncat) :: &
162 1330572 : isosno, isoice ! kg/m^2
163 :
164 : real (kind=dbl_kind) :: &
165 127332 : puny
166 :
167 : character(len=*), parameter :: subname='(step_therm1)'
168 :
169 : !-----------------------------------------------------------------
170 : ! query icepack values
171 : !-----------------------------------------------------------------
172 :
173 346176 : call icepack_query_parameters(puny_out=puny)
174 346176 : call icepack_query_parameters(calc_Tsfc_out=calc_Tsfc)
175 346176 : call icepack_warnings_flush(nu_diag)
176 346176 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
177 0 : file=__FILE__,line= __LINE__)
178 :
179 : call icepack_query_tracer_sizes( &
180 346176 : ntrcr_out=ntrcr)
181 346176 : call icepack_warnings_flush(nu_diag)
182 346176 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
183 0 : file=__FILE__,line= __LINE__)
184 :
185 : call icepack_query_tracer_flags( &
186 : tr_iage_out=tr_iage, tr_FY_out=tr_FY, &
187 : tr_aero_out=tr_aero, tr_iso_out=tr_iso, &
188 : tr_pond_out=tr_pond, tr_pond_cesm_out=tr_pond_cesm, &
189 346176 : tr_pond_lvl_out=tr_pond_lvl, tr_pond_topo_out=tr_pond_topo)
190 346176 : call icepack_warnings_flush(nu_diag)
191 346176 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
192 0 : file=__FILE__,line= __LINE__)
193 :
194 : call icepack_query_tracer_indices( &
195 : nt_apnd_out=nt_apnd, nt_hpnd_out=nt_hpnd, nt_ipnd_out=nt_ipnd, &
196 : nt_alvl_out=nt_alvl, nt_vlvl_out=nt_vlvl, nt_Tsfc_out=nt_Tsfc, &
197 : nt_iage_out=nt_iage, nt_FY_out=nt_FY, &
198 : nt_qice_out=nt_qice, nt_sice_out=nt_sice, &
199 : nt_aero_out=nt_aero, nt_qsno_out=nt_qsno, &
200 346176 : nt_isosno_out=nt_isosno, nt_isoice_out=nt_isoice)
201 346176 : call icepack_warnings_flush(nu_diag)
202 346176 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
203 0 : file=__FILE__,line= __LINE__)
204 :
205 : !-----------------------------------------------------------------
206 :
207 346176 : prescribed_ice = .false.
208 346176 : aerosno(:,:,:) = c0
209 346176 : aeroice(:,:,:) = c0
210 346176 : isosno (:,:) = c0
211 346176 : isoice (:,:) = c0
212 :
213 1730880 : do i = 1, nx
214 :
215 : !-----------------------------------------------------------------
216 : ! Save the ice area passed to the coupler (so that history fields
217 : ! can be made consistent with coupler fields).
218 : ! Save the initial ice area and volume in each category.
219 : !-----------------------------------------------------------------
220 :
221 1384704 : aice_init (i) = aice (i)
222 :
223 8268288 : do n = 1, ncat
224 6537408 : aicen_init(i,n) = aicen(i,n)
225 6537408 : vicen_init(i,n) = vicen(i,n)
226 7922112 : vsnon_init(i,n) = vsnon(i,n)
227 : enddo
228 :
229 : enddo ! i
230 :
231 1730880 : do i = 1, nx
232 1384704 : if (tr_aero) then
233 : ! trcrn(nt_aero) has units kg/m^3
234 368928 : do n=1,ncat
235 676368 : do k=1,n_aero
236 0 : aerosno (k,:,n) = &
237 0 : trcrn(i,nt_aero+(k-1)*4 :nt_aero+(k-1)*4+1,n) &
238 922320 : * vsnon_init(i,n)
239 0 : aeroice (k,:,n) = &
240 0 : trcrn(i,nt_aero+(k-1)*4+2:nt_aero+(k-1)*4+3,n) &
241 1229760 : * vicen_init(i,n)
242 : enddo
243 : enddo
244 : endif ! tr_aero
245 :
246 1384704 : if (tr_iso) then
247 : ! trcrn(nt_isosno/ice) has units kg/m^3
248 368928 : do n=1,ncat
249 1291248 : do k=1,n_iso
250 922320 : isosno(k,n) = trcrn(i,nt_isosno+k-1,n) * vsnon_init(i,n)
251 1229760 : isoice(k,n) = trcrn(i,nt_isoice+k-1,n) * vicen_init(i,n)
252 : enddo
253 : enddo
254 : endif ! tr_iso
255 :
256 : call icepack_step_therm1(dt=dt, ncat=ncat, nilyr=nilyr, nslyr=nslyr, &
257 : aicen_init = aicen_init(i,:), &
258 : vicen_init = vicen_init(i,:), &
259 : vsnon_init = vsnon_init(i,:), &
260 509328 : aice = aice(i), aicen = aicen(i,:), &
261 509328 : vice = vice(i), vicen = vicen(i,:), &
262 509328 : vsno = vsno(i), vsnon = vsnon(i,:), &
263 509328 : uvel = uvel(i), vvel = vvel(i), &
264 : Tsfc = trcrn(i,nt_Tsfc,:), &
265 0 : zqsn = trcrn(i,nt_qsno:nt_qsno+nslyr-1,:), &
266 0 : zqin = trcrn(i,nt_qice:nt_qice+nilyr-1,:), &
267 0 : zSin = trcrn(i,nt_sice:nt_sice+nilyr-1,:), &
268 : alvl = trcrn(i,nt_alvl,:), &
269 : vlvl = trcrn(i,nt_vlvl,:), &
270 : apnd = trcrn(i,nt_apnd,:), &
271 : hpnd = trcrn(i,nt_hpnd,:), &
272 : ipnd = trcrn(i,nt_ipnd,:), &
273 : iage = trcrn(i,nt_iage,:), &
274 : FY = trcrn(i,nt_FY,:), &
275 : aerosno = aerosno(:,:,:), &
276 : aeroice = aeroice(:,:,:), &
277 : isosno = isosno(:,:), &
278 : isoice = isoice(:,:), &
279 509328 : uatm = uatm(i), vatm = vatm(i), &
280 509328 : wind = wind(i), zlvl = zlvl(i), &
281 509328 : Qa = Qa(i), rhoa = rhoa(i), &
282 : Qa_iso = Qa_iso(i,:), &
283 509328 : Tair = Tair(i), Tref = Tref(i), &
284 509328 : Qref = Qref(i), Uref = Uref(i), &
285 : Qref_iso = Qref_iso(i,:), &
286 509328 : Cdn_atm_ratio = Cdn_atm_ratio(i),&
287 509328 : Cdn_ocn = Cdn_ocn(i), &
288 509328 : Cdn_ocn_skin = Cdn_ocn_skin(i), &
289 509328 : Cdn_ocn_floe = Cdn_ocn_floe(i), &
290 509328 : Cdn_ocn_keel = Cdn_ocn_keel(i), &
291 509328 : Cdn_atm = Cdn_atm(i), &
292 509328 : Cdn_atm_skin = Cdn_atm_skin(i), &
293 509328 : Cdn_atm_floe = Cdn_atm_floe(i), &
294 509328 : Cdn_atm_pond = Cdn_atm_pond(i), &
295 509328 : Cdn_atm_rdg = Cdn_atm_rdg(i), &
296 1018656 : hfreebd = hfreebd(i), hkeel = hkeel(i), &
297 509328 : hdraft = hdraft(i), hridge = hridge(i), &
298 1018656 : distrdg = distrdg(i), dkeel = dkeel(i), &
299 509328 : lfloe = lfloe(i), dfloe = dfloe(i), &
300 509328 : strax = strax(i), stray = stray(i), &
301 509328 : strairxT = strairxT(i), strairyT = strairyT(i), &
302 509328 : potT = potT(i), sst = sst(i), &
303 509328 : sss = sss(i), Tf = Tf(i), &
304 509328 : strocnxT = strocnxT(i), strocnyT = strocnyT(i), &
305 1018656 : fbot = fbot(i), frzmlt = frzmlt(i), &
306 509328 : Tbot = Tbot(i), Tsnice = Tsnice(i), &
307 509328 : rside = rside(i), fside = fside(i), &
308 509328 : fsnow = fsnow(i), frain = frain(i), &
309 509328 : fpond = fpond(i), &
310 509328 : fsurf = fsurf(i), fsurfn = fsurfn(i,:), &
311 509328 : fcondtop = fcondtop(i), fcondtopn = fcondtopn(i,:), &
312 509328 : fcondbot = fcondbot(i), fcondbotn = fcondbotn(i,:), &
313 : fswsfcn = fswsfcn(i,:), fswintn = fswintn(i,:), &
314 : fswthrun = fswthrun(i,:), &
315 : fswthrun_vdr = fswthrun_vdr(i,:), &
316 : fswthrun_vdf = fswthrun_vdf(i,:), &
317 : fswthrun_idr = fswthrun_idr(i,:), &
318 : fswthrun_idf = fswthrun_idf(i,:), &
319 509328 : fswabs = fswabs(i), &
320 1018656 : flwout = flwout(i), flw = flw(i), &
321 509328 : fsens = fsens(i), fsensn = fsensn(i,:), &
322 509328 : flat = flat(i), flatn = flatn(i,:), &
323 509328 : fresh = fresh(i), fsalt = fsalt(i), &
324 509328 : fhocn = fhocn(i), &
325 509328 : fswthru = fswthru(i), &
326 509328 : fswthru_vdr= fswthru_vdr(i), &
327 509328 : fswthru_vdf= fswthru_vdf(i), &
328 509328 : fswthru_idr= fswthru_idr(i), &
329 509328 : fswthru_idf= fswthru_idf(i), &
330 : flatn_f = flatn_f(i,:), fsensn_f = fsensn_f(i,:), &
331 : fsurfn_f = fsurfn_f(i,:), &
332 : fcondtopn_f = fcondtopn_f(i,:), &
333 : faero_atm = faero_atm(i,1:n_aero), &
334 : faero_ocn = faero_ocn(i,1:n_aero), &
335 : fiso_atm = fiso_atm (i,:), &
336 : fiso_ocn = fiso_ocn (i,:), &
337 : fiso_evap = fiso_evap (i,:), &
338 509328 : HDO_ocn = HDO_ocn (i), &
339 509328 : H2_16O_ocn = H2_16O_ocn (i), &
340 509328 : H2_18O_ocn = H2_18O_ocn (i), &
341 : Sswabsn = Sswabsn(i,:,:),Iswabsn = Iswabsn(i,:,:), &
342 509328 : evap = evap(i), evaps = evaps(i), evapi = evapi(i), &
343 : dhsn = dhsn(i,:), ffracn = ffracn(i,:), &
344 509328 : meltt = meltt(i), melttn = melttn(i,:), &
345 509328 : meltb = meltb(i), meltbn = meltbn(i,:), &
346 509328 : melts = melts(i), meltsn = meltsn(i,:), &
347 509328 : congel = congel(i), congeln = congeln(i,:), &
348 509328 : snoice = snoice(i), snoicen = snoicen(i,:), &
349 : dsnown = dsnown(i,:), &
350 509328 : lmask_n = lmask_n(i), lmask_s = lmask_s(i), &
351 509328 : mlt_onset=mlt_onset(i), frz_onset = frz_onset(i), &
352 11571264 : yday = yday, prescribed_ice = prescribed_ice)
353 :
354 1384704 : if (tr_aero) then
355 368928 : do n = 1, ncat
356 307440 : if (vicen(i,n) > puny) &
357 1148695 : aeroice(:,:,n) = aeroice(:,:,n)/vicen(i,n)
358 307440 : if (vsnon(i,n) > puny) &
359 1132865 : aerosno(:,:,n) = aerosno(:,:,n)/vsnon(i,n)
360 676368 : do k = 1, n_aero
361 1229760 : do kk = 1, 2
362 614880 : trcrn(i,nt_aero+(k-1)*4+kk-1,n)=aerosno(k,kk,n)
363 922320 : trcrn(i,nt_aero+(k-1)*4+kk+1,n)=aeroice(k,kk,n)
364 : enddo
365 : enddo
366 : enddo
367 : endif ! tr_aero
368 :
369 1730880 : if (tr_iso) then
370 368928 : do n = 1, ncat
371 996717 : if (vicen(i,n) > puny) isoice(:,n) = isoice(:,n)/vicen(i,n)
372 985359 : if (vsnon(i,n) > puny) isosno(:,n) = isosno(:,n)/vsnon(i,n)
373 1291248 : do k = 1, n_iso
374 922320 : trcrn(i,nt_isosno+k-1,n) = isosno(k,n)
375 1229760 : trcrn(i,nt_isoice+k-1,n) = isoice(k,n)
376 : enddo
377 : enddo
378 : endif ! tr_iso
379 :
380 : enddo ! i
381 346176 : call icepack_warnings_flush(nu_diag)
382 346176 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
383 0 : file=__FILE__, line=__LINE__)
384 :
385 346176 : end subroutine step_therm1
386 :
387 : !=======================================================================
388 : ! Driver for thermodynamic changes not needed for coupling:
389 : ! transport in thickness space, lateral growth and melting.
390 : !
391 : ! authors: William H. Lipscomb, LANL
392 : ! Elizabeth C. Hunke, LANL
393 :
394 346176 : subroutine step_therm2 (dt)
395 :
396 : use icedrv_arrays_column, only: hin_max, fzsal, ocean_bio, &
397 : wave_sig_ht, wave_spectrum, &
398 : wavefreq, dwavefreq, &
399 : floe_rad_c, floe_binwidth, &
400 : d_afsd_latg, d_afsd_newi, d_afsd_latm, d_afsd_weld
401 : use icedrv_arrays_column, only: first_ice, bgrid, cgrid, igrid
402 : use icedrv_calendar, only: yday
403 : use icedrv_domain_size, only: ncat, nilyr, nslyr, n_aero, nblyr, &
404 : nltrcr, nx, nfsd
405 : use icedrv_flux, only: fresh, frain, fpond, frzmlt, frazil, frz_onset
406 : use icedrv_flux, only: update_ocn_f, fsalt, Tf, sss, salinz, fhocn, rside, fside
407 : use icedrv_flux, only: meltl, frazil_diag, flux_bio, faero_ocn, fiso_ocn
408 : use icedrv_flux, only: HDO_ocn, H2_16O_ocn, H2_18O_ocn
409 : use icedrv_init, only: tmask
410 : use icedrv_state, only: aice, aicen, aice0, trcr_depend
411 : use icedrv_state, only: aicen_init, vicen_init, trcrn, vicen, vsnon
412 : use icedrv_state, only: trcr_base, n_trcr_strata, nt_strata
413 :
414 : ! column package_includes
415 : use icepack_intfc, only: icepack_step_therm2
416 :
417 : real (kind=dbl_kind), intent(in) :: &
418 : dt ! time step
419 :
420 : ! local variables
421 :
422 : integer (kind=int_kind) :: &
423 : i ! horizontal index
424 :
425 : integer (kind=int_kind) :: &
426 : ntrcr, nbtrcr
427 :
428 : logical (kind=log_kind) :: &
429 : tr_fsd ! floe size distribution tracers
430 :
431 : character(len=*), parameter :: subname='(step_therm2)'
432 :
433 : !-----------------------------------------------------------------
434 : ! query icepack values
435 : !-----------------------------------------------------------------
436 :
437 346176 : call icepack_query_tracer_sizes(ntrcr_out=ntrcr, nbtrcr_out=nbtrcr)
438 346176 : call icepack_query_tracer_flags(tr_fsd_out=tr_fsd)
439 346176 : call icepack_warnings_flush(nu_diag)
440 346176 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
441 0 : file=__FILE__,line= __LINE__)
442 :
443 : !-----------------------------------------------------------------
444 :
445 1730880 : do i = 1, nx
446 :
447 1730880 : if (tmask(i)) then
448 : ! wave_sig_ht - compute here to pass to add new ice
449 1038528 : if (tr_fsd) &
450 2165904 : wave_sig_ht(i) = c4*SQRT(SUM(wave_spectrum(i,:)*dwavefreq(:)))
451 :
452 : call icepack_step_therm2(dt=dt, ncat=ncat, &
453 : nltrcr=nltrcr, nilyr=nilyr, nslyr=nslyr, &
454 : hin_max=hin_max(:), nblyr=nblyr, &
455 : aicen=aicen(i,:), &
456 : vicen=vicen(i,:), &
457 : vsnon=vsnon(i,:), &
458 : aicen_init=aicen_init(i,:), &
459 : vicen_init=vicen_init(i,:), &
460 0 : trcrn=trcrn(i,1:ntrcr,:), &
461 381996 : aice0=aice0(i), &
462 381996 : aice =aice(i), &
463 0 : trcr_depend=trcr_depend(1:ntrcr), &
464 0 : trcr_base=trcr_base(1:ntrcr,:), &
465 0 : n_trcr_strata=n_trcr_strata(1:ntrcr), &
466 0 : nt_strata=nt_strata(1:ntrcr,:), &
467 381996 : Tf=Tf(i), sss=sss(i), &
468 381996 : salinz=salinz(i,:), fside=fside(i), &
469 381996 : rside=rside(i), meltl=meltl(i), &
470 381996 : frzmlt=frzmlt(i), frazil=frazil(i), &
471 381996 : frain=frain(i), fpond=fpond(i), &
472 381996 : fresh=fresh(i), fsalt=fsalt(i), &
473 381996 : fhocn=fhocn(i), update_ocn_f=update_ocn_f, &
474 : bgrid=bgrid, cgrid=cgrid, &
475 : igrid=igrid, faero_ocn=faero_ocn(i,:), &
476 : first_ice=first_ice(i,:), &
477 381996 : fzsal=fzsal(i), &
478 0 : flux_bio=flux_bio(i,1:nbtrcr), &
479 0 : ocean_bio=ocean_bio(i,1:nbtrcr), &
480 381996 : frazil_diag=frazil_diag(i), &
481 381996 : frz_onset=frz_onset(i), &
482 : yday=yday, &
483 : fiso_ocn=fiso_ocn(i,:), &
484 381996 : HDO_ocn=HDO_ocn(i), &
485 381996 : H2_16O_ocn=H2_16O_ocn(i), &
486 381996 : H2_18O_ocn=H2_18O_ocn(i), &
487 381996 : nfsd=nfsd, wave_sig_ht=wave_sig_ht(i), &
488 : wave_spectrum=wave_spectrum(i,:), &
489 : wavefreq=wavefreq(:), &
490 : dwavefreq=dwavefreq(:), &
491 : d_afsd_latg=d_afsd_latg(i,:), &
492 : d_afsd_newi=d_afsd_newi(i,:), &
493 : d_afsd_latm=d_afsd_latm(i,:), &
494 : d_afsd_weld=d_afsd_weld(i,:), &
495 : floe_rad_c=floe_rad_c(:), &
496 3330504 : floe_binwidth=floe_binwidth(:))
497 :
498 : endif ! tmask
499 :
500 : enddo ! i
501 346176 : call icepack_warnings_flush(nu_diag)
502 346176 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
503 0 : file=__FILE__, line=__LINE__)
504 :
505 346176 : end subroutine step_therm2
506 :
507 : !=======================================================================
508 : !
509 : ! finalize thermo updates
510 : !
511 : ! authors: Elizabeth Hunke, LANL
512 :
513 716484 : subroutine update_state (dt, daidt, dvidt, dagedt, offset)
514 :
515 : use icedrv_domain_size, only: ncat, nx
516 : use icedrv_init, only: tmask
517 : use icedrv_state, only: aicen, trcrn, vicen, vsnon
518 : use icedrv_state, only: aice, trcr, vice, vsno, aice0, trcr_depend
519 : use icedrv_state, only: trcr_base, nt_strata, n_trcr_strata
520 :
521 : ! column package includes
522 : use icepack_intfc, only: icepack_aggregate
523 :
524 : real (kind=dbl_kind), intent(in) :: &
525 : dt , & ! time step
526 : offset ! d(age)/dt time offset = dt for thermo, 0 for dyn
527 :
528 : real (kind=dbl_kind), dimension(:), intent(inout) :: &
529 : daidt, & ! change in ice area per time step
530 : dvidt, & ! change in ice volume per time step
531 : dagedt ! change in ice age per time step
532 :
533 : integer (kind=int_kind) :: &
534 : i, & ! horizontal indices
535 : ntrcr, & !
536 : nt_iage !
537 :
538 : logical (kind=log_kind) :: &
539 : tr_iage ! ice age tracer
540 :
541 : character(len=*), parameter :: subname='(update_state)'
542 :
543 : !-----------------------------------------------------------------
544 : ! query icepack values
545 : !-----------------------------------------------------------------
546 :
547 716484 : call icepack_query_tracer_sizes(ntrcr_out=ntrcr)
548 716484 : call icepack_warnings_flush(nu_diag)
549 716484 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
550 0 : file=__FILE__,line= __LINE__)
551 :
552 716484 : call icepack_query_tracer_indices(nt_iage_out=nt_iage)
553 716484 : call icepack_warnings_flush(nu_diag)
554 716484 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
555 0 : file=__FILE__,line= __LINE__)
556 :
557 716484 : call icepack_query_tracer_flags(tr_iage_out=tr_iage)
558 716484 : call icepack_warnings_flush(nu_diag)
559 716484 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
560 0 : file=__FILE__,line= __LINE__)
561 :
562 : !$OMP PARALLEL DO PRIVATE(i)
563 3582420 : do i = 1, nx
564 :
565 : !-----------------------------------------------------------------
566 : ! Aggregate the updated state variables (includes ghost cells).
567 : !-----------------------------------------------------------------
568 :
569 2865936 : if (tmask(i)) then
570 : call icepack_aggregate (ncat=ncat, &
571 0 : aicen=aicen(i,:), trcrn=trcrn(i,1:ntrcr,:), &
572 : vicen=vicen(i,:), vsnon=vsnon(i,:), &
573 790272 : aice =aice (i), trcr =trcr (i,1:ntrcr), &
574 790272 : vice =vice (i), vsno =vsno (i), &
575 790272 : aice0=aice0(i), &
576 : ntrcr=ntrcr, &
577 0 : trcr_depend=trcr_depend (1:ntrcr), &
578 0 : trcr_base=trcr_base (1:ntrcr,:), &
579 0 : n_trcr_strata=n_trcr_strata(1:ntrcr), &
580 3729996 : nt_strata=nt_strata (1:ntrcr,:))
581 : endif
582 :
583 : !-----------------------------------------------------------------
584 : ! Compute thermodynamic area and volume tendencies.
585 : !-----------------------------------------------------------------
586 :
587 2865936 : daidt(i) = (aice(i) - daidt(i)) / dt
588 2865936 : dvidt(i) = (vice(i) - dvidt(i)) / dt
589 3582420 : if (tr_iage) then
590 122976 : if (offset > c0) then ! thermo
591 61488 : if (trcr(i,nt_iage) > c0) &
592 0 : dagedt(i) = (trcr(i,nt_iage) &
593 46111 : - dagedt(i) - offset) / dt
594 : else ! dynamics
595 0 : dagedt(i) = (trcr(i,nt_iage) &
596 61488 : - dagedt(i)) / dt
597 : endif
598 : endif
599 :
600 : enddo ! i
601 : !$OMP END PARALLEL DO
602 716484 : call icepack_warnings_flush(nu_diag)
603 716484 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
604 0 : file=__FILE__, line=__LINE__)
605 :
606 716484 : end subroutine update_state
607 :
608 : !=======================================================================
609 : !
610 : ! Run one time step of wave-fracturing the floe size distribution
611 : !
612 : ! authors: Lettie Roach, NIWA
613 : ! Elizabeth C. Hunke, LANL
614 :
615 19008 : subroutine step_dyn_wave (dt)
616 :
617 : use icedrv_arrays_column, only: wave_spectrum, wave_sig_ht, &
618 : d_afsd_wave, floe_rad_l, floe_rad_c, wavefreq, dwavefreq
619 : use icedrv_domain_size, only: ncat, nfsd, nfreq, nx
620 : use icedrv_state, only: trcrn, aicen, aice, vice
621 : use icepack_intfc, only: icepack_step_wavefracture
622 :
623 : real (kind=dbl_kind), intent(in) :: &
624 : dt ! time step
625 :
626 : ! local variables
627 :
628 : integer (kind=int_kind) :: &
629 : i, j, & ! horizontal indices
630 : ntrcr, & !
631 : nbtrcr !
632 :
633 : character (len=char_len) :: wave_spec_type
634 :
635 : character(len=*), parameter :: subname = '(step_dyn_wave)'
636 :
637 19008 : call icepack_query_parameters(wave_spec_type_out=wave_spec_type)
638 19008 : call icepack_warnings_flush(nu_diag)
639 19008 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
640 0 : file=__FILE__,line= __LINE__)
641 :
642 95040 : do i = 1, nx
643 988416 : d_afsd_wave(i,:) = c0
644 : call icepack_step_wavefracture (wave_spec_type=wave_spec_type, &
645 : dt=dt, ncat=ncat, nfsd=nfsd, nfreq=nfreq, &
646 35040 : aice = aice (i), &
647 35040 : vice = vice (i), &
648 : aicen = aicen (i,:), &
649 : floe_rad_l = floe_rad_l (:), &
650 : floe_rad_c = floe_rad_c (:), &
651 : wave_spectrum = wave_spectrum(i,:), &
652 : wavefreq = wavefreq (:), &
653 : dwavefreq = dwavefreq (:), &
654 : trcrn = trcrn (i,:,:), &
655 95040 : d_afsd_wave = d_afsd_wave (i,:))
656 : end do ! i
657 :
658 19008 : call icepack_warnings_flush(nu_diag)
659 19008 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
660 0 : file=__FILE__,line= __LINE__)
661 :
662 19008 : end subroutine step_dyn_wave
663 :
664 : !=======================================================================
665 : !
666 : ! Run one time step of ridging.
667 : !
668 : ! authors: William H. Lipscomb, LANL
669 : ! Elizabeth C. Hunke, LANL
670 :
671 370308 : subroutine step_dyn_ridge (dt, ndtd)
672 :
673 : use icedrv_arrays_column, only: hin_max, fzsal, first_ice
674 : use icedrv_domain_size, only: ncat, nilyr, nslyr, n_aero, nblyr, nx
675 : use icedrv_flux, only: rdg_conv, rdg_shear, dardg1dt, dardg2dt
676 : use icedrv_flux, only: dvirdgdt, opening, closing, fpond, fresh, fhocn
677 : use icedrv_flux, only: aparticn, krdgn, aredistn, vredistn, dardg1ndt, dardg2ndt
678 : use icedrv_flux, only: dvirdgndt, araftn, vraftn, fsalt, flux_bio, faero_ocn, fiso_ocn
679 : use icedrv_init, only: tmask
680 : use icedrv_state, only: trcrn, vsnon, aicen, vicen
681 : use icedrv_state, only: aice, aice0, trcr_depend, n_trcr_strata
682 : use icedrv_state, only: trcr_base, nt_strata
683 :
684 : ! column package includes
685 : use icepack_intfc, only: icepack_step_ridge
686 :
687 : real (kind=dbl_kind), intent(in) :: &
688 : dt ! time step
689 :
690 : integer (kind=int_kind), intent(in) :: &
691 : ndtd ! number of dynamics subcycles
692 :
693 : ! local variables
694 :
695 : integer (kind=int_kind) :: &
696 : i, & ! horizontal indices
697 : ntrcr, & !
698 : nbtrcr !
699 :
700 : character(len=*), parameter :: subname='(step_dyn_ridge)'
701 :
702 : !-----------------------------------------------------------------
703 : ! query icepack values
704 : !-----------------------------------------------------------------
705 :
706 370308 : call icepack_query_tracer_sizes(ntrcr_out=ntrcr, nbtrcr_out=nbtrcr)
707 370308 : call icepack_warnings_flush(nu_diag)
708 370308 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
709 0 : file=__FILE__,line= __LINE__)
710 :
711 : !-----------------------------------------------------------------
712 : ! Ridging
713 : !-----------------------------------------------------------------
714 :
715 370308 : if (trim(ocn_data_type) == "SHEBA") then
716 :
717 1619820 : do i = 1, nx
718 :
719 : !echmod: this changes the answers, continue using tmask for now
720 : ! call aggregate_area (ncat, aicen(:), atmp, atmp0)
721 : ! if (atmp > c0) then
722 :
723 1619820 : if (tmask(i)) then
724 :
725 : call icepack_step_ridge(dt=dt, ndtd=ndtd, &
726 : nilyr=nilyr, nslyr=nslyr, &
727 : nblyr=nblyr, &
728 : ncat=ncat, hin_max=hin_max(:), &
729 361476 : rdg_conv=rdg_conv(i), rdg_shear=rdg_shear(i), &
730 : aicen=aicen(i,:), &
731 0 : trcrn=trcrn(i,1:ntrcr,:), &
732 : vicen=vicen(i,:), vsnon=vsnon(i,:), &
733 361476 : aice0=aice0(i), &
734 0 : trcr_depend=trcr_depend(1:ntrcr), &
735 0 : trcr_base=trcr_base(1:ntrcr,:), &
736 0 : n_trcr_strata=n_trcr_strata(1:ntrcr), &
737 0 : nt_strata=nt_strata(1:ntrcr,:), &
738 361476 : dardg1dt=dardg1dt(i), dardg2dt=dardg2dt(i), &
739 361476 : dvirdgdt=dvirdgdt(i), opening=opening(i), &
740 361476 : fpond=fpond(i), &
741 361476 : fresh=fresh(i), fhocn=fhocn(i), &
742 : n_aero=n_aero, &
743 : faero_ocn=faero_ocn(i,:), fiso_ocn=fiso_ocn(i,:), &
744 : aparticn=aparticn(i,:), krdgn=krdgn(i,:), &
745 : aredistn=aredistn(i,:), vredistn=vredistn(i,:), &
746 : dardg1ndt=dardg1ndt(i,:), dardg2ndt=dardg2ndt(i,:), &
747 : dvirdgndt=dvirdgndt(i,:), &
748 : araftn=araftn(i,:), vraftn=vraftn(i,:), &
749 361476 : aice=aice(i), fsalt=fsalt(i), &
750 361476 : first_ice=first_ice(i,:), fzsal=fzsal(i), &
751 0 : flux_bio=flux_bio(i,1:nbtrcr), &
752 2417796 : closing=closing(i) )
753 :
754 : endif ! tmask
755 :
756 : enddo ! i
757 :
758 : else ! closing not read in
759 :
760 231720 : do i = 1, nx
761 :
762 : !echmod: this changes the answers, continue using tmask for now
763 : ! call aggregate_area (ncat, aicen(:), atmp, atmp0)
764 : ! if (atmp > c0) then
765 :
766 231720 : if (tmask(i)) then
767 :
768 : call icepack_step_ridge (dt=dt, ndtd=ndtd, &
769 : nilyr=nilyr, nslyr=nslyr, &
770 : nblyr=nblyr, &
771 : ncat=ncat, hin_max=hin_max(:), &
772 46800 : rdg_conv=rdg_conv(i), rdg_shear=rdg_shear(i), &
773 : aicen=aicen(i,:), &
774 0 : trcrn=trcrn(i,1:ntrcr,:), &
775 : vicen=vicen(i,:), vsnon=vsnon(i,:), &
776 46800 : aice0=aice0(i), &
777 0 : trcr_depend=trcr_depend(1:ntrcr), &
778 0 : trcr_base=trcr_base(1:ntrcr,:), &
779 0 : n_trcr_strata=n_trcr_strata(1:ntrcr), &
780 0 : nt_strata=nt_strata(1:ntrcr,:), &
781 46800 : dardg1dt=dardg1dt(i), dardg2dt=dardg2dt(i), &
782 46800 : dvirdgdt=dvirdgdt(i), opening=opening(i), &
783 46800 : fpond=fpond(i), &
784 46800 : fresh=fresh(i), fhocn=fhocn(i), &
785 : n_aero=n_aero, &
786 : faero_ocn=faero_ocn(i,:), fiso_ocn=fiso_ocn(i,:), &
787 : aparticn=aparticn(i,:), krdgn=krdgn(i,:), &
788 : aredistn=aredistn(i,:), vredistn=vredistn(i,:), &
789 : dardg1ndt=dardg1ndt(i,:), dardg2ndt=dardg2ndt(i,:), &
790 : dvirdgndt=dvirdgndt(i,:), &
791 : araftn=araftn(i,:), vraftn=vraftn(i,:), &
792 46800 : aice=aice(i), fsalt=fsalt(i), &
793 46800 : first_ice=first_ice(i,:), fzsal=fzsal(i), &
794 326232 : flux_bio=flux_bio(i,1:nbtrcr))
795 :
796 : endif ! tmask
797 :
798 : enddo ! i
799 :
800 : endif
801 370308 : call icepack_warnings_flush(nu_diag)
802 370308 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
803 0 : file=__FILE__, line=__LINE__)
804 :
805 370308 : end subroutine step_dyn_ridge
806 :
807 : !=======================================================================
808 : !
809 : ! Computes radiation fields
810 : !
811 : ! authors: William H. Lipscomb, LANL
812 : ! David Bailey, NCAR
813 : ! Elizabeth C. Hunke, LANL
814 :
815 346176 : subroutine step_radiation (dt)
816 :
817 : use icedrv_arrays_column, only: ffracn, dhsn
818 : use icedrv_arrays_column, only: fswsfcn, fswintn, fswpenln, Sswabsn, Iswabsn
819 : use icedrv_arrays_column, only: fswthrun, fswthrun_vdr, fswthrun_vdf, fswthrun_idr, fswthrun_idf
820 : use icedrv_arrays_column, only: albicen, albsnon, albpndn
821 : use icedrv_arrays_column, only: alvdrn, alidrn, alvdfn, alidfn, apeffn, trcrn_sw, snowfracn
822 : use icedrv_arrays_column, only: kaer_tab, waer_tab, gaer_tab, kaer_bc_tab, waer_bc_tab
823 : use icedrv_arrays_column, only: gaer_bc_tab, bcenh, swgrid, igrid
824 : use icedrv_calendar, only: calendar_type, days_per_year, nextsw_cday, yday, sec
825 : use icedrv_domain_size, only: ncat, n_aero, nilyr, nslyr, n_zaero, n_algae, nblyr, nx
826 : use icedrv_flux, only: swvdr, swvdf, swidr, swidf, coszen, fsnow
827 : use icedrv_init, only: TLAT, TLON, tmask
828 : use icedrv_state, only: aicen, vicen, vsnon, trcrn
829 :
830 : ! column package includes
831 : use icepack_intfc, only: icepack_step_radiation
832 :
833 : real (kind=dbl_kind), intent(in) :: &
834 : dt ! time step
835 :
836 : ! local variables
837 :
838 : integer (kind=int_kind) :: &
839 : i, n, k ! horizontal indices
840 :
841 : integer (kind=int_kind) :: &
842 : max_aero, max_algae, nt_Tsfc, nt_alvl, &
843 : nt_apnd, nt_hpnd, nt_ipnd, nt_aero, nlt_chl_sw, &
844 : ntrcr, nbtrcr_sw, nt_fbri
845 :
846 : integer (kind=int_kind), dimension(:), allocatable :: &
847 346176 : nlt_zaero_sw, nt_zaero, nt_bgc_N
848 :
849 : logical (kind=log_kind) :: &
850 : tr_bgc_N, tr_zaero, tr_brine, dEdd_algae, modal_aero
851 :
852 : real (kind=dbl_kind), dimension(ncat) :: &
853 728952 : fbri ! brine height to ice thickness
854 :
855 : real(kind= dbl_kind), dimension(:,:), allocatable :: &
856 346176 : ztrcr_sw
857 :
858 : logical (kind=log_kind) :: &
859 : l_print_point ! flag for printing debugging information
860 :
861 : character(len=*), parameter :: subname='(step_radiation)'
862 :
863 : !-----------------------------------------------------------------
864 : ! query icepack values
865 : !-----------------------------------------------------------------
866 :
867 : call icepack_query_tracer_sizes( &
868 346176 : max_aero_out=max_aero, max_algae_out=max_algae)
869 346176 : call icepack_warnings_flush(nu_diag)
870 346176 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
871 0 : file=__FILE__,line= __LINE__)
872 346176 : allocate(nlt_zaero_sw(max_aero))
873 346176 : allocate(nt_zaero(max_aero))
874 346176 : allocate(nt_bgc_N(max_algae))
875 :
876 346176 : call icepack_query_tracer_sizes(ntrcr_out=ntrcr, nbtrcr_sw_out=nbtrcr_sw)
877 346176 : call icepack_warnings_flush(nu_diag)
878 346176 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
879 0 : file=__FILE__,line= __LINE__)
880 :
881 : call icepack_query_tracer_flags( &
882 346176 : tr_brine_out=tr_brine, tr_bgc_N_out=tr_bgc_N, tr_zaero_out=tr_zaero)
883 346176 : call icepack_warnings_flush(nu_diag)
884 346176 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
885 0 : file=__FILE__,line= __LINE__)
886 :
887 : call icepack_query_tracer_indices( &
888 : nt_Tsfc_out=nt_Tsfc, nt_alvl_out=nt_alvl, nt_apnd_out=nt_apnd, &
889 : nt_hpnd_out=nt_hpnd, nt_ipnd_out=nt_ipnd, nt_aero_out=nt_aero, &
890 : nlt_chl_sw_out=nlt_chl_sw, nlt_zaero_sw_out=nlt_zaero_sw, &
891 346176 : nt_fbri_out=nt_fbri, nt_zaero_out=nt_zaero, nt_bgc_N_out=nt_bgc_N)
892 346176 : call icepack_warnings_flush(nu_diag)
893 346176 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
894 0 : file=__FILE__,line= __LINE__)
895 :
896 346176 : call icepack_query_parameters(dEdd_algae_out=dEdd_algae, modal_aero_out=modal_aero)
897 346176 : call icepack_warnings_flush(nu_diag)
898 346176 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
899 0 : file=__FILE__,line= __LINE__)
900 :
901 : !-----------------------------------------------------------------
902 :
903 346176 : allocate(ztrcr_sw(nbtrcr_sw,ncat))
904 :
905 346176 : l_print_point = .false.
906 :
907 1730880 : do i = 1, nx
908 :
909 1384704 : fbri(:) = c0
910 7922112 : ztrcr_sw(:,:) = c0
911 7922112 : do n = 1, ncat
912 7922112 : if (tr_brine) fbri(n) = trcrn(i,nt_fbri,n)
913 : enddo
914 :
915 1384704 : if (tmask(i)) then
916 :
917 : call icepack_step_radiation(dt=dt, ncat=ncat, &
918 : nblyr=nblyr, nilyr=nilyr, &
919 : nslyr=nslyr, dEdd_algae=dEdd_algae, &
920 : swgrid=swgrid(:), igrid=igrid(:), &
921 : fbri=fbri(:), &
922 : aicen=aicen(i,:), vicen=vicen(i,:), &
923 : vsnon=vsnon(i,:), &
924 : Tsfcn=trcrn(i,nt_Tsfc,:), &
925 : alvln=trcrn(i,nt_alvl,:), &
926 : apndn=trcrn(i,nt_apnd,:), &
927 : hpndn=trcrn(i,nt_hpnd,:), &
928 : ipndn=trcrn(i,nt_ipnd,:), &
929 0 : aeron=trcrn(i,nt_aero:nt_aero+4*n_aero-1,:), &
930 0 : bgcNn=trcrn(i,nt_bgc_N(1):nt_bgc_N(1)+n_algae*(nblyr+3)-1,:), &
931 0 : zaeron=trcrn(i,nt_zaero(1):nt_zaero(1)+n_zaero*(nblyr+3)-1,:), &
932 : trcrn_bgcsw=ztrcr_sw, &
933 381996 : TLAT=TLAT(i), TLON=TLON(i), &
934 : calendar_type=calendar_type, &
935 : days_per_year=days_per_year, sec=sec, &
936 : nextsw_cday=nextsw_cday, yday=yday, &
937 : kaer_tab=kaer_tab, kaer_bc_tab=kaer_bc_tab(:,:), &
938 : waer_tab=waer_tab, waer_bc_tab=waer_bc_tab(:,:), &
939 : gaer_tab=gaer_tab, gaer_bc_tab=gaer_bc_tab(:,:), &
940 : bcenh=bcenh(:,:,:), modal_aero=modal_aero, &
941 381996 : swvdr=swvdr(i), swvdf=swvdf(i), &
942 381996 : swidr=swidr(i), swidf=swidf(i), &
943 381996 : coszen=coszen(i), fsnow=fsnow(i), &
944 : alvdrn=alvdrn(i,:), alvdfn=alvdfn(i,:), &
945 : alidrn=alidrn(i,:), alidfn=alidfn(i,:), &
946 : fswsfcn=fswsfcn(i,:), fswintn=fswintn(i,:), &
947 : fswthrun=fswthrun(i,:), &
948 : fswthrun_vdr=fswthrun_vdr(i,:), &
949 : fswthrun_vdf=fswthrun_vdf(i,:), &
950 : fswthrun_idr=fswthrun_idr(i,:), &
951 : fswthrun_idf=fswthrun_idf(i,:), &
952 : fswpenln=fswpenln(i,:,:), &
953 : Sswabsn=Sswabsn(i,:,:), Iswabsn=Iswabsn(i,:,:), &
954 : albicen=albicen(i,:), albsnon=albsnon(i,:), &
955 : albpndn=albpndn(i,:), apeffn=apeffn(i,:), &
956 : snowfracn=snowfracn(i,:), &
957 : dhsn=dhsn(i,:), ffracn=ffracn(i,:), &
958 1802520 : l_print_point=l_print_point)
959 :
960 : endif ! tmask
961 :
962 1730880 : if (dEdd_algae .and. (tr_zaero .or. tr_bgc_N)) then
963 0 : do n = 1, ncat
964 1384704 : do k = 1, nbtrcr_sw
965 0 : trcrn_sw(i,k,n) = ztrcr_sw(k,n)
966 : enddo
967 : enddo
968 : endif
969 :
970 : enddo ! i
971 346176 : call icepack_warnings_flush(nu_diag)
972 346176 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
973 0 : file=__FILE__, line=__LINE__)
974 :
975 346176 : deallocate(ztrcr_sw)
976 346176 : deallocate(nlt_zaero_sw)
977 346176 : deallocate(nt_zaero)
978 346176 : deallocate(nt_bgc_N)
979 :
980 346176 : end subroutine step_radiation
981 :
982 : !=======================================================================
983 : ! Ocean mixed layer calculation (internal to sea ice model).
984 : ! Allows heat storage in ocean for uncoupled runs.
985 : !
986 : ! authors: John Weatherly, CRREL
987 : ! C.M. Bitz, UW
988 : ! Elizabeth C. Hunke, LANL
989 : ! Bruce P. Briegleb, NCAR
990 : ! William H. Lipscomb, LANL
991 :
992 346176 : subroutine ocean_mixed_layer (dt)
993 :
994 : use icedrv_arrays_column, only: Cdn_atm, Cdn_atm_ratio
995 : use icepack_intfc, only: icepack_ocn_mixed_layer, icepack_atm_boundary
996 : use icedrv_init, only: tmask
997 : use icedrv_domain_size, only: nx
998 : use icedrv_flux, only: sst, Tf, Qa, uatm, vatm, wind, potT, rhoa, zlvl
999 : use icedrv_flux, only: frzmlt, fhocn, fswthru, flw, flwout_ocn, fsens_ocn, flat_ocn, evap_ocn
1000 : use icedrv_flux, only: alvdr_ocn, alidr_ocn, alvdf_ocn, alidf_ocn, swidf, swvdf, swidr, swvdr
1001 : use icedrv_flux, only: qdp, hmix, strairx_ocn, strairy_ocn, Tref_ocn, Qref_ocn
1002 : use icedrv_state, only: aice
1003 :
1004 : real (kind=dbl_kind), intent(in) :: &
1005 : dt ! time step
1006 :
1007 : ! local variables
1008 :
1009 : integer (kind=int_kind) :: &
1010 : i ! horizontal indices
1011 :
1012 : real (kind=dbl_kind) :: &
1013 127332 : albocn
1014 :
1015 : real (kind=dbl_kind), dimension(nx) :: &
1016 636660 : delt , & ! potential temperature difference (K)
1017 636660 : delq , & ! specific humidity difference (kg/kg)
1018 636660 : shcoef, & ! transfer coefficient for sensible heat
1019 636660 : lhcoef ! transfer coefficient for latent heat
1020 :
1021 : character(len=*), parameter :: subname='(ocean_mixed_layer)'
1022 :
1023 : !-----------------------------------------------------------------
1024 : ! query icepack values
1025 : !-----------------------------------------------------------------
1026 :
1027 346176 : call icepack_query_parameters(albocn_out=albocn)
1028 346176 : call icepack_warnings_flush(nu_diag)
1029 346176 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
1030 0 : file=__FILE__, line=__LINE__)
1031 :
1032 : !-----------------------------------------------------------------
1033 : ! Identify ocean cells.
1034 : ! Set fluxes to zero in land cells.
1035 : !-----------------------------------------------------------------
1036 :
1037 1730880 : do i = 1, nx
1038 1730880 : if (.not.tmask(i)) then
1039 346176 : sst (i) = c0
1040 346176 : frzmlt (i) = c0
1041 346176 : flwout_ocn(i) = c0
1042 346176 : fsens_ocn (i) = c0
1043 346176 : flat_ocn (i) = c0
1044 346176 : evap_ocn (i) = c0
1045 : endif
1046 : enddo ! i
1047 :
1048 : !-----------------------------------------------------------------
1049 : ! Compute boundary layer quantities
1050 : !-----------------------------------------------------------------
1051 :
1052 1730880 : do i = 1, nx
1053 1730880 : if (tmask(i)) then
1054 : call icepack_atm_boundary(sfctype = 'ocn', &
1055 381996 : Tsf = sst(i), &
1056 381996 : potT = potT(i), &
1057 381996 : uatm = uatm(i), &
1058 381996 : vatm = vatm(i), &
1059 381996 : wind = wind(i), &
1060 381996 : zlvl = zlvl(i), &
1061 381996 : Qa = Qa(i), &
1062 381996 : rhoa = rhoa(i), &
1063 381996 : strx = strairx_ocn(i), &
1064 381996 : stry = strairy_ocn(i), &
1065 381996 : Tref = Tref_ocn(i), &
1066 381996 : Qref = Qref_ocn(i), &
1067 381996 : delt = delt(i), &
1068 381996 : delq = delq(i), &
1069 381996 : lhcoef = lhcoef(i), &
1070 381996 : shcoef = shcoef(i), &
1071 381996 : Cdn_atm = Cdn_atm(i), &
1072 1038528 : Cdn_atm_ratio_n = Cdn_atm_ratio(i))
1073 : endif
1074 : enddo ! i
1075 346176 : call icepack_warnings_flush(nu_diag)
1076 346176 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
1077 0 : file=__FILE__, line=__LINE__)
1078 :
1079 : !-----------------------------------------------------------------
1080 : ! Ocean albedo
1081 : ! For now, assume albedo = albocn in each spectral band.
1082 : !-----------------------------------------------------------------
1083 :
1084 1730880 : alvdr_ocn(:) = albocn
1085 1730880 : alidr_ocn(:) = albocn
1086 1730880 : alvdf_ocn(:) = albocn
1087 1730880 : alidf_ocn(:) = albocn
1088 :
1089 : !-----------------------------------------------------------------
1090 : ! Compute ocean fluxes and update SST
1091 : !-----------------------------------------------------------------
1092 1730880 : do i = 1, nx
1093 1730880 : if (tmask(i)) then
1094 381996 : call icepack_ocn_mixed_layer(alvdr_ocn=alvdr_ocn(i), swvdr=swvdr(i), &
1095 381996 : alidr_ocn=alidr_ocn(i), swidr=swidr(i), &
1096 381996 : alvdf_ocn=alvdf_ocn(i), swvdf=swvdf(i), &
1097 381996 : alidf_ocn=alidf_ocn(i), swidf=swidf(i), &
1098 381996 : flwout_ocn=flwout_ocn(i),sst=sst(i), &
1099 381996 : fsens_ocn=fsens_ocn(i), shcoef=shcoef(i), &
1100 381996 : flat_ocn=flat_ocn(i), lhcoef=lhcoef(i), &
1101 381996 : evap_ocn=evap_ocn(i), flw=flw(i), &
1102 381996 : delt=delt(i), delq=delq(i), &
1103 381996 : aice=aice(i), fhocn=fhocn(i), &
1104 381996 : fswthru=fswthru(i), hmix=hmix(i), &
1105 381996 : Tf=Tf(i), qdp=qdp(i), &
1106 1038528 : frzmlt=frzmlt(i), dt=dt)
1107 : endif
1108 : enddo ! i
1109 346176 : call icepack_warnings_flush(nu_diag)
1110 346176 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
1111 0 : file=__FILE__, line=__LINE__)
1112 :
1113 346176 : end subroutine ocean_mixed_layer
1114 :
1115 : !=======================================================================
1116 :
1117 346176 : subroutine biogeochemistry (dt)
1118 :
1119 : use icedrv_arrays_column, only: upNO, upNH, iDi, iki, zfswin
1120 : use icedrv_arrays_column, only: zsal_tot, darcy_V, grow_net
1121 : use icedrv_arrays_column, only: PP_net, hbri,dhbr_bot, dhbr_top, Zoo
1122 : use icedrv_arrays_column, only: fbio_snoice, fbio_atmice, ocean_bio
1123 : use icedrv_arrays_column, only: first_ice, fswpenln, bphi, bTiz, ice_bio_net
1124 : use icedrv_arrays_column, only: snow_bio_net, fswthrun, Rayleigh_criteria
1125 : use icedrv_arrays_column, only: ocean_bio_all, sice_rho, fzsal, fzsal_g
1126 : use icedrv_arrays_column, only: bgrid, igrid, icgrid, cgrid
1127 : use icepack_intfc, only: icepack_biogeochemistry, icepack_load_ocean_bio_array
1128 : use icedrv_domain_size, only: nblyr, nilyr, nslyr, n_algae, n_zaero, ncat
1129 : use icedrv_domain_size, only: n_doc, n_dic, n_don, n_fed, n_fep, nx
1130 : use icedrv_flux, only: meltbn, melttn, congeln, snoicen
1131 : use icedrv_flux, only: sst, sss, fsnow, meltsn
1132 : use icedrv_flux, only: hin_old, flux_bio, flux_bio_atm, faero_atm
1133 : use icedrv_flux, only: nit, amm, sil, dmsp, dms, algalN, doc, don, dic, fed, fep, zaeros, hum
1134 : use icedrv_state, only: aicen_init, vicen_init, aicen, vicen, vsnon
1135 : use icedrv_state, only: trcrn, vsnon_init, aice0
1136 :
1137 : real (kind=dbl_kind), intent(in) :: &
1138 : dt ! time step
1139 :
1140 : ! local variables
1141 :
1142 : integer (kind=int_kind) :: &
1143 : i , & ! horizontal indices
1144 : mm ! tracer index
1145 :
1146 : integer (kind=int_kind) :: &
1147 : max_algae, max_nbtrcr, max_don, &
1148 : max_doc, max_dic, max_aero, max_fe, &
1149 : nbtrcr, ntrcr
1150 :
1151 : integer (kind=int_kind), dimension(:), allocatable :: &
1152 346176 : nlt_zaero
1153 :
1154 : integer (kind=int_kind), allocatable :: &
1155 346176 : bio_index_o(:)
1156 :
1157 : logical (kind=log_kind) :: &
1158 : skl_bgc, tr_brine, tr_zaero
1159 :
1160 : character(len=*), parameter :: subname='(biogeochemistry)'
1161 :
1162 : !-----------------------------------------------------------------
1163 : ! query icepack values
1164 : !-----------------------------------------------------------------
1165 :
1166 346176 : call icepack_query_tracer_flags(tr_brine_out=tr_brine)
1167 346176 : call icepack_warnings_flush(nu_diag)
1168 346176 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
1169 0 : file=__FILE__,line= __LINE__)
1170 :
1171 346176 : call icepack_query_parameters(skl_bgc_out=skl_bgc)
1172 346176 : call icepack_warnings_flush(nu_diag)
1173 346176 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
1174 0 : file=__FILE__,line= __LINE__)
1175 :
1176 : !-----------------------------------------------------------------
1177 :
1178 346176 : if (tr_brine .or. skl_bgc) then
1179 :
1180 : !-----------------------------------------------------------------
1181 :
1182 : call icepack_query_tracer_sizes( &
1183 : max_algae_out=max_algae, max_nbtrcr_out=max_nbtrcr, max_don_out=max_don, &
1184 22212 : max_doc_out=max_doc, max_dic_out=max_dic, max_aero_out=max_aero, max_fe_out=max_fe)
1185 22212 : call icepack_warnings_flush(nu_diag)
1186 22212 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
1187 0 : file=__FILE__,line= __LINE__)
1188 :
1189 : !-----------------------------------------------------------------
1190 :
1191 22212 : allocate(bio_index_o(max_nbtrcr))
1192 22212 : allocate(nlt_zaero(max_aero))
1193 :
1194 : !-----------------------------------------------------------------
1195 :
1196 22212 : call icepack_query_tracer_sizes(ntrcr_out=ntrcr, nbtrcr_out=nbtrcr)
1197 22212 : call icepack_query_tracer_flags(tr_zaero_out=tr_zaero)
1198 22212 : call icepack_query_tracer_indices(nlt_zaero_out=nlt_zaero)
1199 22212 : call icepack_query_tracer_indices(bio_index_o_out=bio_index_o)
1200 22212 : call icepack_warnings_flush(nu_diag)
1201 22212 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
1202 0 : file=__FILE__,line= __LINE__)
1203 :
1204 : !-----------------------------------------------------------------
1205 :
1206 : ! Define ocean concentrations for tracers used in simulation
1207 111060 : do i = 1, nx
1208 :
1209 : call icepack_load_ocean_bio_array(max_nbtrcr=max_nbtrcr,&
1210 : max_algae = max_algae, max_don = max_don, &
1211 : max_doc = max_doc, max_dic = max_dic, &
1212 : max_aero = max_aero, max_fe = max_fe, &
1213 27360 : nit = nit(i), amm = amm(i), &
1214 27360 : sil = sil(i), dmsp = dmsp(i), &
1215 27360 : dms = dms(i), algalN = algalN(i,:), &
1216 : doc = doc(i,:), don = don(i,:), &
1217 : dic = dic(i,:), fed = fed(i,:), &
1218 : fep = fep(i,:), zaeros = zaeros(i,:), &
1219 : ocean_bio_all=ocean_bio_all(i,:), &
1220 88848 : hum=hum(i))
1221 : ! call icepack_warnings_flush(nu_diag)
1222 : ! if (icepack_warnings_aborted()) call icedrv_system_abort(i, istep1, subname, &
1223 : ! file=__FILE__,line= __LINE__)
1224 :
1225 1832784 : do mm = 1,nbtrcr
1226 1832784 : ocean_bio(i,mm) = ocean_bio_all(i,bio_index_o(mm))
1227 : enddo ! mm
1228 88848 : if (tr_zaero) then
1229 564144 : do mm = 1, n_zaero ! update aerosols
1230 564144 : flux_bio_atm(i,nlt_zaero(mm)) = faero_atm(i,mm)
1231 : enddo ! mm
1232 : endif
1233 :
1234 : call icepack_biogeochemistry(dt=dt, ntrcr=ntrcr, nbtrcr=nbtrcr, &
1235 : ncat=ncat, nblyr=nblyr, nilyr=nilyr, nslyr=nslyr, &
1236 : n_algae=n_algae, n_zaero=n_zaero, &
1237 : n_doc=n_doc, n_dic=n_dic, n_don=n_don, &
1238 : n_fed=n_fed, n_fep=n_fep, &
1239 : bgrid=bgrid, igrid=igrid, icgrid=icgrid, cgrid=cgrid, &
1240 27360 : upNO = upNO(i), &
1241 27360 : upNH = upNH(i), &
1242 : iDi = iDi(i,:,:), &
1243 : iki = iki(i,:,:), &
1244 : zfswin = zfswin(i,:,:), &
1245 27360 : zsal_tot = zsal_tot(i), &
1246 : darcy_V = darcy_V(i,:), &
1247 27360 : grow_net = grow_net(i), &
1248 27360 : PP_net = PP_net(i), &
1249 27360 : hbri = hbri(i), &
1250 : dhbr_bot = dhbr_bot(i,:), &
1251 : dhbr_top = dhbr_top(i,:), &
1252 : Zoo = Zoo(i,:,:), &
1253 : fbio_snoice = fbio_snoice(i,:), &
1254 : fbio_atmice = fbio_atmice(i,:), &
1255 0 : ocean_bio = ocean_bio(i,1:nbtrcr), &
1256 : first_ice = first_ice(i,:), &
1257 : fswpenln = fswpenln(i,:,:), &
1258 : bphi = bphi(i,:,:), &
1259 : bTiz = bTiz(i,:,:), &
1260 0 : ice_bio_net = ice_bio_net(i,1:nbtrcr), &
1261 0 : snow_bio_net = snow_bio_net(i,1:nbtrcr), &
1262 : fswthrun = fswthrun(i,:), &
1263 27360 : Rayleigh_criteria = Rayleigh_criteria(i), &
1264 : sice_rho = sice_rho(i,:), &
1265 27360 : fzsal = fzsal(i), &
1266 27360 : fzsal_g = fzsal_g(i), &
1267 : meltbn = meltbn(i,:), &
1268 : melttn = melttn(i,:), &
1269 : congeln = congeln(i,:), &
1270 : snoicen = snoicen(i,:), &
1271 27360 : sst = sst(i), &
1272 27360 : sss = sss(i), &
1273 27360 : fsnow = fsnow(i), &
1274 : meltsn = meltsn(i,:), &
1275 : hin_old = hin_old(i,:), &
1276 0 : flux_bio = flux_bio(i,1:nbtrcr), &
1277 0 : flux_bio_atm = flux_bio_atm(i,1:nbtrcr), &
1278 : aicen_init = aicen_init(i,:), &
1279 : vicen_init = vicen_init(i,:), &
1280 : aicen = aicen(i,:), &
1281 : vicen = vicen(i,:), &
1282 : vsnon = vsnon(i,:), &
1283 27360 : aice0 = aice0(i), &
1284 0 : trcrn = trcrn(i,1:ntrcr,:), &
1285 : vsnon_init = vsnon_init(i,:), &
1286 275220 : skl_bgc = skl_bgc)
1287 :
1288 : ! call icepack_warnings_flush(nu_diag)
1289 : ! if (icepack_warnings_aborted()) call icedrv_system_abort(i, istep1, subname, &
1290 : ! __FILE__, __LINE__)
1291 :
1292 : enddo ! i
1293 22212 : call icepack_warnings_flush(nu_diag)
1294 22212 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
1295 0 : file=__FILE__, line=__LINE__)
1296 :
1297 22212 : deallocate(nlt_zaero)
1298 22212 : deallocate(bio_index_o)
1299 :
1300 : endif ! tr_brine .or. skl_bgc
1301 :
1302 346176 : end subroutine biogeochemistry
1303 :
1304 : !=======================================================================
1305 :
1306 : end module icedrv_step
1307 :
1308 : !=======================================================================
|