Line data Source code
1 : !=======================================================================
2 : !
3 : ! Main driver for time stepping of Icepack
4 : !
5 : ! authors Elizabeth C. Hunke, LANL
6 :
7 : module icedrv_RunMod
8 :
9 : use icedrv_kinds
10 : use icedrv_constants, only: c0, c1, nu_diag
11 : use icepack_intfc, only: icepack_warnings_flush
12 : use icepack_intfc, only: icepack_warnings_aborted
13 : use icepack_intfc, only: icepack_query_parameters
14 : use icepack_intfc, only: icepack_query_tracer_flags
15 : use icepack_intfc, only: icepack_query_tracer_sizes
16 : use icedrv_system, only: icedrv_system_abort
17 :
18 : implicit none
19 : private
20 : public :: icedrv_run, ice_step
21 :
22 : !=======================================================================
23 :
24 : contains
25 :
26 : !=======================================================================
27 : !
28 : ! This is the main driver routine for advancing CICE forward in time.
29 : !
30 : ! author Elizabeth C. Hunke, LANL
31 :
32 45 : subroutine icedrv_run
33 :
34 : use icedrv_calendar, only: istep, istep1, time, dt, stop_now, calendar
35 : use icedrv_forcing, only: get_forcing, get_wave_spec
36 : use icedrv_forcing_bgc, only: faero_default, fiso_default, get_forcing_bgc
37 : use icedrv_flux, only: init_flux_atm_ocn
38 :
39 : logical (kind=log_kind) :: skl_bgc, z_tracers, tr_aero, tr_zaero, &
40 : wave_spec, tr_fsd, tr_iso
41 :
42 : character(len=*), parameter :: subname='(icedrv_run)'
43 :
44 : !--------------------------------------------------------------------
45 : ! timestep loop
46 : !--------------------------------------------------------------------
47 :
48 : call icepack_query_tracer_flags(tr_aero_out=tr_aero, tr_zaero_out=tr_zaero, &
49 45 : tr_fsd_out=tr_fsd, tr_iso_out=tr_iso)
50 45 : call icepack_warnings_flush(nu_diag)
51 45 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
52 0 : file=__FILE__,line= __LINE__)
53 :
54 342495 : timeLoop: do
55 :
56 342540 : call ice_step
57 :
58 342540 : istep = istep + 1 ! update time step counters
59 342540 : istep1 = istep1 + 1
60 342540 : time = time + dt ! determine the time and date
61 :
62 342540 : call calendar(time) ! at the end of the timestep
63 :
64 342540 : if (stop_now >= 1) exit timeLoop
65 :
66 : call icepack_query_parameters(skl_bgc_out=skl_bgc, z_tracers_out=z_tracers,&
67 342495 : wave_spec_out=wave_spec)
68 342495 : call icepack_warnings_flush(nu_diag)
69 342495 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
70 0 : file=__FILE__,line= __LINE__)
71 :
72 342495 : if (tr_fsd .and. wave_spec) call get_wave_spec ! wave spectrum in ice
73 342495 : call get_forcing(istep1) ! get forcing from data arrays
74 :
75 : ! biogeochemistry forcing
76 342495 : if (tr_iso) call fiso_default ! default values
77 342495 : if (tr_aero .or. tr_zaero) call faero_default ! default values
78 342495 : if (skl_bgc .or. z_tracers) call get_forcing_bgc ! biogeochemistry
79 :
80 342540 : call init_flux_atm_ocn ! initialize atmosphere, ocean fluxes
81 :
82 : enddo timeLoop
83 :
84 45 : end subroutine icedrv_run
85 :
86 : !=======================================================================
87 : !
88 : ! Calls drivers for physics components, some initialization, and output
89 : !
90 : ! author Elizabeth C. Hunke, LANL
91 :
92 342540 : subroutine ice_step
93 :
94 : use icedrv_calendar, only: dt, dt_dyn, ndtd, diagfreq, write_restart, istep
95 : use icedrv_diagnostics, only: runtime_diags, init_mass_diags
96 : ! use icedrv_diagnostics, only: icedrv_diagnostics_debug
97 : use icedrv_diagnostics_bgc, only: hbrine_diags, zsal_diags, bgc_diags
98 : use icedrv_flux, only: init_history_therm, init_history_bgc, &
99 : daidtt, daidtd, dvidtt, dvidtd, dagedtt, dagedtd, init_history_dyn
100 : use icedrv_restart, only: dumpfile, final_restart
101 : use icedrv_restart_bgc, only: write_restart_bgc
102 : use icedrv_step, only: prep_radiation, step_therm1, step_therm2, &
103 : update_state, step_dyn_ridge, step_radiation, &
104 : biogeochemistry, step_dyn_wave
105 :
106 : integer (kind=int_kind) :: &
107 : k ! dynamics supercycling index
108 :
109 : logical (kind=log_kind) :: &
110 : calc_Tsfc, skl_bgc, solve_zsal, z_tracers, tr_brine, & ! from icepack
111 : tr_fsd, wave_spec
112 :
113 : real (kind=dbl_kind) :: &
114 118572 : offset ! d(age)/dt time offset
115 :
116 : character(len=*), parameter :: subname='(ice_step)'
117 :
118 : ! call icedrv_diagnostics_debug ('beginning time step')
119 :
120 : !-----------------------------------------------------------------
121 : ! query Icepack values
122 : !-----------------------------------------------------------------
123 :
124 342540 : call icepack_query_parameters(skl_bgc_out=skl_bgc, z_tracers_out=z_tracers)
125 : call icepack_query_parameters(solve_zsal_out=solve_zsal, &
126 : calc_Tsfc_out=calc_Tsfc, &
127 342540 : wave_spec_out=wave_spec)
128 342540 : call icepack_query_tracer_flags(tr_brine_out=tr_brine,tr_fsd_out=tr_fsd)
129 342540 : call icepack_warnings_flush(nu_diag)
130 342540 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
131 0 : file=__FILE__,line= __LINE__)
132 :
133 : !-----------------------------------------------------------------
134 : ! initialize diagnostics
135 : !-----------------------------------------------------------------
136 :
137 342540 : call init_mass_diags ! diagnostics per timestep
138 342540 : call init_history_therm
139 342540 : call init_history_bgc
140 :
141 : !-----------------------------------------------------------------
142 : ! Scale radiation fields
143 : !-----------------------------------------------------------------
144 :
145 342540 : if (calc_Tsfc) call prep_radiation ()
146 :
147 : ! call icedrv_diagnostics_debug ('post prep_radiation')
148 :
149 : !-----------------------------------------------------------------
150 : ! thermodynamics and biogeochemistry
151 : !-----------------------------------------------------------------
152 :
153 342540 : call step_therm1 (dt) ! vertical thermodynamics
154 342540 : call biogeochemistry (dt) ! biogeochemistry
155 342540 : call step_therm2 (dt) ! ice thickness distribution thermo
156 :
157 : ! clean up, update tendency diagnostics
158 342540 : offset = dt
159 342540 : call update_state (dt, daidtt, dvidtt, dagedtt, offset)
160 :
161 : ! call icedrv_diagnostics_debug ('post thermo')
162 :
163 : !-----------------------------------------------------------------
164 : ! dynamics, transport, ridging
165 : !-----------------------------------------------------------------
166 :
167 342540 : call init_history_dyn
168 :
169 : ! wave fracture of the floe size distribution
170 : ! note this is called outside of the dynamics subcycling loop
171 342540 : if (tr_fsd .and. wave_spec) call step_dyn_wave(dt)
172 :
173 709212 : do k = 1, ndtd
174 :
175 : ! ridging
176 366672 : call step_dyn_ridge (dt_dyn, ndtd)
177 :
178 : ! clean up, update tendency diagnostics
179 366672 : offset = c0
180 709212 : call update_state (dt_dyn, daidtd, dvidtd, dagedtd, offset)
181 :
182 : enddo
183 :
184 : ! call icedrv_diagnostics_debug ('post dynamics')
185 :
186 : !-----------------------------------------------------------------
187 : ! albedo, shortwave radiation
188 : !-----------------------------------------------------------------
189 :
190 342540 : call step_radiation (dt)
191 :
192 : !-----------------------------------------------------------------
193 : ! get ready for coupling and the next time step
194 : !-----------------------------------------------------------------
195 :
196 342540 : call coupling_prep
197 :
198 : ! call icedrv_diagnostics_debug ('post step_rad, cpl')
199 :
200 : !-----------------------------------------------------------------
201 : ! write data
202 : !-----------------------------------------------------------------
203 :
204 342540 : if (mod(istep,diagfreq) == 0) then
205 37392 : call runtime_diags(dt) ! log file
206 37392 : if (solve_zsal) call zsal_diags
207 37392 : if (skl_bgc .or. z_tracers) call bgc_diags
208 37392 : if (tr_brine) call hbrine_diags
209 : endif
210 :
211 342540 : if (write_restart == 1) then
212 811 : call dumpfile ! core variables for restarting
213 811 : if (solve_zsal .or. skl_bgc .or. z_tracers) &
214 25 : call write_restart_bgc ! biogeochemistry
215 811 : call final_restart
216 : endif
217 :
218 342540 : end subroutine ice_step
219 :
220 : !=======================================================================
221 : !
222 : ! Prepare for coupling
223 : !
224 : ! authors: Elizabeth C. Hunke, LANL
225 :
226 342540 : subroutine coupling_prep
227 :
228 : use icedrv_arrays_column, only: alvdfn, alidfn, alvdrn, alidrn, &
229 : albicen, albsnon, albpndn, apeffn, fzsal_g, fzsal, snowfracn
230 : use icedrv_calendar, only: dt
231 : use icedrv_domain_size, only: ncat, nx
232 : use icedrv_flux, only: alvdf, alidf, alvdr, alidr, albice, albsno, &
233 : albpnd, apeff_ai, fpond, fresh, l_mpond_fresh, &
234 : alvdf_ai, alidf_ai, alvdr_ai, alidr_ai, fhocn_ai, &
235 : fresh_ai, fsalt_ai, fsalt, &
236 : fswthru_ai, fhocn, fswthru, scale_factor, snowfrac, &
237 : swvdr, swidr, swvdf, swidf, &
238 : frzmlt_init, frzmlt, &
239 : fzsal_ai, fzsal_g_ai, flux_bio, flux_bio_ai
240 : use icedrv_forcing, only: oceanmixed_ice
241 : use icedrv_state, only: aicen
242 : use icedrv_step, only: ocean_mixed_layer
243 :
244 : ! local variables
245 :
246 : integer (kind=int_kind) :: &
247 : n , & ! thickness category index
248 : i , & ! horizontal index
249 : k , & ! tracer index
250 : nbtrcr
251 :
252 : real (kind=dbl_kind) :: &
253 118572 : netsw, & ! flag for shortwave radiation presence
254 118572 : rhofresh, & !
255 118572 : puny !
256 :
257 : character(len=*), parameter :: subname='(coupling_prep)'
258 :
259 : !-----------------------------------------------------------------
260 : ! Save current value of frzmlt for diagnostics.
261 : ! Update mixed layer with heat and radiation from ice.
262 : !-----------------------------------------------------------------
263 :
264 342540 : call icepack_query_parameters(puny_out=puny, rhofresh_out=rhofresh)
265 342540 : call icepack_query_tracer_sizes(nbtrcr_out=nbtrcr)
266 342540 : call icepack_warnings_flush(nu_diag)
267 342540 : if (icepack_warnings_aborted()) call icedrv_system_abort(string=subname, &
268 0 : file=__FILE__,line= __LINE__)
269 :
270 1712700 : do i = 1, nx
271 1712700 : frzmlt_init (i) = frzmlt(i)
272 : enddo
273 :
274 342540 : if (oceanmixed_ice) &
275 342540 : call ocean_mixed_layer (dt) ! ocean surface fluxes and sst
276 :
277 : !-----------------------------------------------------------------
278 : ! Aggregate albedos
279 : !-----------------------------------------------------------------
280 :
281 1712700 : do i = 1, nx
282 1370160 : alvdf(i) = c0
283 1370160 : alidf(i) = c0
284 1370160 : alvdr(i) = c0
285 1370160 : alidr(i) = c0
286 :
287 1370160 : albice(i) = c0
288 1370160 : albsno(i) = c0
289 1370160 : albpnd(i) = c0
290 1370160 : apeff_ai(i) = c0
291 1712700 : snowfrac(i) = c0
292 : enddo
293 1958712 : do n = 1, ncat
294 8423400 : do i = 1, nx
295 8080860 : if (aicen(i,n) > puny) then
296 :
297 4361008 : alvdf(i) = alvdf(i) + alvdfn(i,n)*aicen(i,n)
298 4361008 : alidf(i) = alidf(i) + alidfn(i,n)*aicen(i,n)
299 4361008 : alvdr(i) = alvdr(i) + alvdrn(i,n)*aicen(i,n)
300 4361008 : alidr(i) = alidr(i) + alidrn(i,n)*aicen(i,n)
301 :
302 4361008 : netsw = swvdr(i) + swidr(i) + swvdf(i) + swidf(i)
303 4361008 : if (netsw > puny) then ! sun above horizon
304 3430369 : albice(i) = albice(i) + albicen(i,n)*aicen(i,n)
305 3430369 : albsno(i) = albsno(i) + albsnon(i,n)*aicen(i,n)
306 3430369 : albpnd(i) = albpnd(i) + albpndn(i,n)*aicen(i,n)
307 : endif
308 :
309 4361008 : apeff_ai(i) = apeff_ai(i) + apeffn(i,n)*aicen(i,n) ! for history
310 4361008 : snowfrac(i) = snowfrac(i) + snowfracn(i,n)*aicen(i,n) ! for history
311 :
312 : endif ! aicen > puny
313 : enddo
314 : enddo
315 :
316 1712700 : do i = 1, nx
317 :
318 : !-----------------------------------------------------------------
319 : ! reduce fresh by fpond for coupling
320 : !-----------------------------------------------------------------
321 :
322 1370160 : if (l_mpond_fresh) then
323 96528 : fpond(i) = fpond(i) * rhofresh/dt
324 96528 : fresh(i) = fresh(i) - fpond(i)
325 : endif
326 :
327 : !----------------------------------------------------------------
328 : ! Store grid box mean albedos and fluxes before scaling by aice
329 : !----------------------------------------------------------------
330 :
331 1370160 : alvdf_ai (i) = alvdf (i)
332 1370160 : alidf_ai (i) = alidf (i)
333 1370160 : alvdr_ai (i) = alvdr (i)
334 1370160 : alidr_ai (i) = alidr (i)
335 1370160 : fresh_ai (i) = fresh (i)
336 1370160 : fsalt_ai (i) = fsalt (i)
337 1370160 : fhocn_ai (i) = fhocn (i)
338 1370160 : fswthru_ai(i) = fswthru(i)
339 1370160 : fzsal_ai (i) = fzsal (i)
340 1370160 : fzsal_g_ai(i) = fzsal_g(i)
341 :
342 1370160 : if (nbtrcr > 0) then
343 1832784 : do k = 1, nbtrcr
344 1832784 : flux_bio_ai (i,k) = flux_bio (i,k)
345 : enddo
346 : endif
347 :
348 : !-----------------------------------------------------------------
349 : ! Save net shortwave for scaling factor in scale_factor
350 : !-----------------------------------------------------------------
351 474288 : scale_factor(i) = &
352 474288 : swvdr(i)*(c1 - alvdr_ai(i)) &
353 474288 : + swvdf(i)*(c1 - alvdf_ai(i)) &
354 474288 : + swidr(i)*(c1 - alidr_ai(i)) &
355 1712700 : + swidf(i)*(c1 - alidf_ai(i))
356 : enddo
357 :
358 342540 : end subroutine coupling_prep
359 :
360 : !=======================================================================
361 :
362 : end module icedrv_RunMod
363 :
364 : !=======================================================================
|