Line data Source code
1 : !=======================================================================
2 : !
3 : ! Main driver for time stepping of CICE.
4 : !
5 : ! authors Elizabeth C. Hunke, LANL
6 : ! Philip W. Jones, LANL
7 : ! William H. Lipscomb, LANL
8 : !
9 : ! 2006 ECH: moved exit timeLoop to prevent execution of unnecessary timestep
10 : ! 2006 ECH: Streamlined for efficiency
11 : ! 2006 ECH: Converted to free source form (F90)
12 : ! 2007 BPB: Modified Delta-Eddington shortwave interface
13 : ! 2008 ECH: moved ESMF code to its own driver
14 :
15 : module CICE_RunMod
16 :
17 : use ice_kinds_mod
18 : use ice_communicate, only: my_task, master_task
19 : use ice_fileunits, only: nu_diag
20 : use ice_arrays_column, only: oceanmixed_ice
21 : use ice_constants, only: c0, c1
22 : use ice_constants, only: field_loc_center, field_type_scalar
23 : use ice_exit, only: abort_ice
24 : use ice_memusage, only: ice_memusage_print
25 : use icepack_intfc, only: icepack_warnings_flush, icepack_warnings_aborted
26 : use icepack_intfc, only: icepack_max_iso, icepack_max_aero
27 : use icepack_intfc, only: icepack_query_parameters
28 : use icepack_intfc, only: icepack_query_tracer_flags, icepack_query_tracer_sizes
29 :
30 : implicit none
31 : private
32 : public :: CICE_Run, ice_step
33 :
34 : !=======================================================================
35 :
36 : contains
37 :
38 : !=======================================================================
39 : !
40 : ! This is the main driver routine for advancing CICE forward in time.
41 : !
42 : ! author Elizabeth C. Hunke, LANL
43 : ! Philip W. Jones, LANL
44 : ! William H. Lipscomb, LANL
45 :
46 37 : subroutine CICE_Run
47 :
48 : use ice_calendar, only: dt, stop_now, advance_timestep
49 : use ice_forcing, only: get_forcing_atmo, get_forcing_ocn, &
50 : get_wave_spec
51 : use ice_forcing_bgc, only: get_forcing_bgc, get_atm_bgc, &
52 : fiso_default, faero_default
53 : use ice_flux, only: init_flux_atm, init_flux_ocn
54 : use ice_timers, only: ice_timer_start, ice_timer_stop, &
55 : timer_couple, timer_step
56 : logical (kind=log_kind) :: &
57 : tr_iso, tr_aero, tr_zaero, skl_bgc, z_tracers, wave_spec, tr_fsd
58 : character(len=*), parameter :: subname = '(CICE_Run)'
59 :
60 : !--------------------------------------------------------------------
61 : ! initialize error code and step timer
62 : !--------------------------------------------------------------------
63 :
64 37 : call ice_timer_start(timer_step) ! start timing entire run
65 :
66 : call icepack_query_parameters(skl_bgc_out=skl_bgc, &
67 : z_tracers_out=z_tracers, & ! LCOV_EXCL_LINE
68 37 : wave_spec_out=wave_spec)
69 : call icepack_query_tracer_flags(tr_iso_out=tr_iso, &
70 : tr_aero_out=tr_aero, & ! LCOV_EXCL_LINE
71 : tr_zaero_out=tr_zaero, & ! LCOV_EXCL_LINE
72 37 : tr_fsd_out=tr_fsd)
73 37 : call icepack_warnings_flush(nu_diag)
74 37 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
75 0 : file=__FILE__, line=__LINE__)
76 :
77 : #ifndef CICE_IN_NEMO
78 : !--------------------------------------------------------------------
79 : ! timestep loop
80 : !--------------------------------------------------------------------
81 :
82 5747 : timeLoop: do
83 : #endif
84 :
85 5784 : call ice_step
86 :
87 : ! tcraig, use advance_timestep now
88 : ! istep = istep + 1 ! update time step counters
89 : ! istep1 = istep1 + 1
90 : ! time = time + dt ! determine the time and date
91 : ! call calendar(time) ! at the end of the timestep
92 5784 : call advance_timestep() ! advance time
93 :
94 : #ifndef CICE_IN_NEMO
95 5784 : if (stop_now >= 1) exit timeLoop
96 : #endif
97 :
98 5747 : call ice_timer_start(timer_couple) ! atm/ocn coupling
99 :
100 : ! for now, wave_spectrum is constant in time
101 : ! if (tr_fsd .and. wave_spec) call get_wave_spec ! wave spectrum in ice
102 5747 : call get_forcing_atmo ! atmospheric forcing from data
103 5747 : call get_forcing_ocn(dt) ! ocean forcing from data
104 :
105 : ! isotopes
106 5747 : if (tr_iso) call fiso_default ! default values
107 : ! aerosols
108 : ! if (tr_aero) call faero_data ! data file
109 : ! if (tr_zaero) call fzaero_data ! data file (gx1)
110 5747 : if (tr_aero .or. tr_zaero) call faero_default ! default values
111 :
112 5747 : if (skl_bgc .or. z_tracers) call get_forcing_bgc ! biogeochemistry
113 5747 : if (z_tracers) call get_atm_bgc ! biogeochemistry
114 :
115 5747 : call init_flux_atm ! Initialize atmosphere fluxes sent to coupler
116 5747 : call init_flux_ocn ! initialize ocean fluxes sent to coupler
117 :
118 5747 : call ice_timer_stop(timer_couple) ! atm/ocn coupling
119 :
120 : #ifndef CICE_IN_NEMO
121 : enddo timeLoop
122 : #endif
123 :
124 : !--------------------------------------------------------------------
125 : ! end of timestep loop
126 : !--------------------------------------------------------------------
127 :
128 37 : call ice_timer_stop(timer_step) ! end timestepping loop timer
129 :
130 37 : end subroutine CICE_Run
131 :
132 : !=======================================================================
133 : !
134 : ! Calls drivers for physics components, some initialization, and output
135 : !
136 : ! author Elizabeth C. Hunke, LANL
137 : ! William H. Lipscomb, LANL
138 :
139 5784 : subroutine ice_step
140 :
141 : use ice_boundary, only: ice_HaloUpdate
142 : use ice_calendar, only: dt, dt_dyn, ndtd, diagfreq, write_restart, istep
143 : use ice_diagnostics, only: init_mass_diags, runtime_diags, debug_model, debug_ice
144 : use ice_diagnostics_bgc, only: hbrine_diags, bgc_diags
145 : use ice_domain, only: halo_info, nblocks
146 : use ice_dyn_eap, only: write_restart_eap
147 : use ice_dyn_shared, only: kdyn, kridge
148 : use ice_flux, only: scale_factor, init_history_therm, &
149 : daidtt, daidtd, dvidtt, dvidtd, dagedtt, dagedtd
150 : use ice_history, only: accum_hist
151 : use ice_history_bgc, only: init_history_bgc
152 : use ice_restart, only: final_restart
153 : use ice_restart_column, only: write_restart_age, write_restart_FY, &
154 : write_restart_lvl, write_restart_pond_lvl, & ! LCOV_EXCL_LINE
155 : write_restart_pond_topo, write_restart_aero, write_restart_fsd, & ! LCOV_EXCL_LINE
156 : write_restart_iso, write_restart_bgc, write_restart_hbrine, & ! LCOV_EXCL_LINE
157 : write_restart_snow
158 : use ice_restart_driver, only: dumpfile
159 : use ice_restoring, only: restore_ice, ice_HaloRestore
160 : use ice_step_mod, only: prep_radiation, step_therm1, step_therm2, &
161 : update_state, step_dyn_horiz, step_dyn_ridge, step_radiation, & ! LCOV_EXCL_LINE
162 : biogeochemistry, step_prep, step_dyn_wave, step_snow
163 : use ice_timers, only: ice_timer_start, ice_timer_stop, &
164 : timer_diags, timer_column, timer_thermo, timer_bound, & ! LCOV_EXCL_LINE
165 : timer_hist, timer_readwrite
166 :
167 : integer (kind=int_kind) :: &
168 : iblk , & ! block index ! LCOV_EXCL_LINE
169 : k , & ! dynamics supercycling index ! LCOV_EXCL_LINE
170 : ktherm ! thermodynamics is off when ktherm = -1
171 :
172 : real (kind=dbl_kind) :: &
173 1440 : offset ! d(age)/dt time offset
174 :
175 : logical (kind=log_kind) :: &
176 : tr_iage, tr_FY, tr_lvl, tr_fsd, tr_snow, & ! LCOV_EXCL_LINE
177 : tr_pond_lvl, tr_pond_topo, tr_brine, tr_iso, tr_aero, & ! LCOV_EXCL_LINE
178 : calc_Tsfc, skl_bgc, z_tracers, wave_spec
179 :
180 : character(len=*), parameter :: subname = '(ice_step)'
181 :
182 : character (len=char_len) :: plabeld
183 :
184 5784 : if (debug_model) then
185 0 : plabeld = 'beginning time step'
186 0 : do iblk = 1, nblocks
187 0 : call debug_ice (iblk, plabeld)
188 : enddo
189 : endif
190 :
191 : call icepack_query_parameters(calc_Tsfc_out=calc_Tsfc, skl_bgc_out=skl_bgc, &
192 5784 : z_tracers_out=z_tracers, ktherm_out=ktherm, wave_spec_out=wave_spec)
193 : call icepack_query_tracer_flags(tr_iage_out=tr_iage, tr_FY_out=tr_FY, &
194 : tr_lvl_out=tr_lvl, tr_pond_lvl_out=tr_pond_lvl, & ! LCOV_EXCL_LINE
195 : tr_pond_topo_out=tr_pond_topo, tr_brine_out=tr_brine, tr_aero_out=tr_aero, & ! LCOV_EXCL_LINE
196 5784 : tr_iso_out=tr_iso, tr_fsd_out=tr_fsd, tr_snow_out=tr_snow)
197 5784 : call icepack_warnings_flush(nu_diag)
198 5784 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
199 0 : file=__FILE__, line=__LINE__)
200 :
201 : !-----------------------------------------------------------------
202 : ! restoring on grid boundaries
203 : !-----------------------------------------------------------------
204 :
205 5784 : if (restore_ice) call ice_HaloRestore
206 :
207 : !-----------------------------------------------------------------
208 : ! initialize diagnostics and save initial state values
209 : !-----------------------------------------------------------------
210 :
211 5784 : call ice_timer_start(timer_diags) ! diagnostics/history
212 5784 : call init_mass_diags ! diagnostics per timestep
213 5784 : call init_history_therm
214 5784 : call init_history_bgc
215 5784 : call ice_timer_stop(timer_diags) ! diagnostics/history
216 :
217 5784 : call ice_timer_start(timer_column) ! column physics
218 5784 : call ice_timer_start(timer_thermo) ! thermodynamics
219 :
220 5784 : call step_prep
221 :
222 5784 : if (ktherm >= 0) then
223 2880 : !$OMP PARALLEL DO PRIVATE(iblk) SCHEDULE(runtime)
224 8688 : do iblk = 1, nblocks
225 :
226 : !-----------------------------------------------------------------
227 : ! scale radiation fields
228 : !-----------------------------------------------------------------
229 :
230 5784 : if (calc_Tsfc) call prep_radiation (iblk)
231 :
232 5784 : if (debug_model) then
233 0 : plabeld = 'post prep_radiation'
234 0 : call debug_ice (iblk, plabeld)
235 : endif
236 :
237 : !-----------------------------------------------------------------
238 : ! thermodynamics and biogeochemistry
239 : !-----------------------------------------------------------------
240 :
241 5784 : call step_therm1 (dt, iblk) ! vertical thermodynamics
242 :
243 5784 : if (debug_model) then
244 0 : plabeld = 'post step_therm1'
245 0 : call debug_ice (iblk, plabeld)
246 : endif
247 :
248 5784 : call biogeochemistry (dt, iblk) ! biogeochemistry
249 :
250 5784 : if (debug_model) then
251 0 : plabeld = 'post biogeochemistry'
252 0 : call debug_ice (iblk, plabeld)
253 : endif
254 :
255 5784 : call step_therm2 (dt, iblk) ! ice thickness distribution thermo
256 :
257 11568 : if (debug_model) then
258 0 : plabeld = 'post step_therm2'
259 0 : call debug_ice (iblk, plabeld)
260 : endif
261 :
262 : enddo
263 : !$OMP END PARALLEL DO
264 : endif ! ktherm > 0
265 :
266 : ! clean up, update tendency diagnostics
267 5784 : offset = dt
268 5784 : call update_state (dt, daidtt, dvidtt, dagedtt, offset)
269 :
270 5784 : call ice_timer_stop(timer_thermo) ! thermodynamics
271 5784 : call ice_timer_stop(timer_column) ! column physics
272 :
273 : !-----------------------------------------------------------------
274 : ! dynamics, transport, ridging
275 : !-----------------------------------------------------------------
276 :
277 : ! wave fracture of the floe size distribution
278 : ! note this is called outside of the dynamics subcycling loop
279 5784 : if (tr_fsd .and. wave_spec) call step_dyn_wave(dt)
280 :
281 11568 : do k = 1, ndtd
282 :
283 : ! momentum, stress, transport
284 5784 : call step_dyn_horiz (dt_dyn)
285 :
286 5784 : if (debug_model) then
287 0 : plabeld = 'post step_dyn_horiz'
288 0 : do iblk = 1, nblocks
289 0 : call debug_ice (iblk, plabeld)
290 : enddo ! iblk
291 : endif
292 :
293 : ! ridging
294 2880 : !$OMP PARALLEL DO PRIVATE(iblk) SCHEDULE(runtime)
295 8688 : do iblk = 1, nblocks
296 11568 : if (kridge > 0) call step_dyn_ridge (dt_dyn, ndtd, iblk)
297 : enddo
298 : !$OMP END PARALLEL DO
299 :
300 5784 : if (debug_model) then
301 0 : plabeld = 'post step_dyn_ridge'
302 0 : do iblk = 1, nblocks
303 0 : call debug_ice (iblk, plabeld)
304 : enddo ! iblk
305 : endif
306 :
307 : ! clean up, update tendency diagnostics
308 5784 : offset = c0
309 11568 : call update_state (dt_dyn, daidtd, dvidtd, dagedtd, offset)
310 :
311 : enddo
312 :
313 5784 : if (debug_model) then
314 0 : plabeld = 'post dynamics'
315 0 : do iblk = 1, nblocks
316 0 : call debug_ice (iblk, plabeld)
317 : enddo
318 : endif
319 :
320 5784 : call ice_timer_start(timer_column) ! column physics
321 5784 : call ice_timer_start(timer_thermo) ! thermodynamics
322 :
323 : !-----------------------------------------------------------------
324 : ! snow redistribution and metamorphosis
325 : !-----------------------------------------------------------------
326 :
327 5784 : if (tr_snow) then ! advanced snow physics
328 0 : !$OMP PARALLEL DO PRIVATE(iblk) SCHEDULE(runtime)
329 0 : do iblk = 1, nblocks
330 0 : call step_snow (dt, iblk)
331 : enddo
332 : !$OMP END PARALLEL DO
333 0 : call update_state (dt) ! clean up
334 : endif
335 :
336 2880 : !$OMP PARALLEL DO PRIVATE(iblk) SCHEDULE(runtime)
337 8688 : do iblk = 1, nblocks
338 :
339 : !-----------------------------------------------------------------
340 : ! albedo, shortwave radiation
341 : !-----------------------------------------------------------------
342 :
343 5784 : if (ktherm >= 0) call step_radiation (dt, iblk)
344 :
345 5784 : if (debug_model) then
346 0 : plabeld = 'post step_radiation'
347 0 : call debug_ice (iblk, plabeld)
348 : endif
349 :
350 : !-----------------------------------------------------------------
351 : ! get ready for coupling and the next time step
352 : !-----------------------------------------------------------------
353 :
354 5784 : call coupling_prep (iblk)
355 :
356 11568 : if (debug_model) then
357 0 : plabeld = 'post coupling_prep'
358 0 : call debug_ice (iblk, plabeld)
359 : endif
360 :
361 : enddo ! iblk
362 : !$OMP END PARALLEL DO
363 :
364 5784 : call ice_timer_start(timer_bound)
365 : call ice_HaloUpdate (scale_factor, halo_info, &
366 5784 : field_loc_center, field_type_scalar)
367 5784 : call ice_timer_stop(timer_bound)
368 :
369 5784 : call ice_timer_stop(timer_thermo) ! thermodynamics
370 5784 : call ice_timer_stop(timer_column) ! column physics
371 :
372 : !-----------------------------------------------------------------
373 : ! write data
374 : !-----------------------------------------------------------------
375 :
376 5784 : call ice_timer_start(timer_diags) ! diagnostics
377 5784 : if (mod(istep,diagfreq) == 0) then
378 5784 : call runtime_diags(dt) ! log file
379 5784 : if (skl_bgc .or. z_tracers) call bgc_diags
380 5784 : if (tr_brine) call hbrine_diags
381 5784 : if (my_task == master_task) then
382 984 : call ice_memusage_print(nu_diag,subname)
383 : endif
384 : endif
385 5784 : call ice_timer_stop(timer_diags) ! diagnostics
386 :
387 5784 : call ice_timer_start(timer_hist) ! history
388 5784 : call accum_hist (dt) ! history file
389 5784 : call ice_timer_stop(timer_hist) ! history
390 :
391 5784 : call ice_timer_start(timer_readwrite) ! reading/writing
392 5784 : if (write_restart == 1) then
393 49 : call dumpfile ! core variables for restarting
394 49 : if (tr_iage) call write_restart_age
395 49 : if (tr_FY) call write_restart_FY
396 49 : if (tr_lvl) call write_restart_lvl
397 49 : if (tr_pond_lvl) call write_restart_pond_lvl
398 49 : if (tr_pond_topo) call write_restart_pond_topo
399 49 : if (tr_snow) call write_restart_snow
400 49 : if (tr_fsd) call write_restart_fsd
401 49 : if (tr_iso) call write_restart_iso
402 49 : if (tr_aero) call write_restart_aero
403 49 : if (skl_bgc .or. z_tracers) &
404 0 : call write_restart_bgc
405 49 : if (tr_brine) call write_restart_hbrine
406 49 : if (kdyn == 2) call write_restart_eap
407 49 : call final_restart
408 : endif
409 5784 : call ice_timer_stop(timer_readwrite) ! reading/writing
410 :
411 5784 : end subroutine ice_step
412 :
413 : !=======================================================================
414 : !
415 : ! Prepare for coupling
416 : !
417 : ! authors: Elizabeth C. Hunke, LANL
418 :
419 22944 : subroutine coupling_prep (iblk)
420 :
421 : use ice_arrays_column, only: alvdfn, alidfn, alvdrn, alidrn, &
422 : albicen, albsnon, albpndn, apeffn, snowfracn
423 : use ice_blocks, only: nx_block, ny_block, get_block, block
424 : use ice_domain, only: blocks_ice
425 : use ice_calendar, only: dt, nstreams
426 : use ice_domain_size, only: ncat
427 : use ice_flux, only: alvdf, alidf, alvdr, alidr, albice, albsno, &
428 : albpnd, albcnt, apeff_ai, fpond, fresh, l_mpond_fresh, & ! LCOV_EXCL_LINE
429 : alvdf_ai, alidf_ai, alvdr_ai, alidr_ai, fhocn_ai, & ! LCOV_EXCL_LINE
430 : fresh_ai, fsalt_ai, fsalt, & ! LCOV_EXCL_LINE
431 : fswthru_ai, fhocn, scale_factor, snowfrac, & ! LCOV_EXCL_LINE
432 : fswthru, fswthru_vdr, fswthru_vdf, fswthru_idr, fswthru_idf, & ! LCOV_EXCL_LINE
433 : swvdr, swidr, swvdf, swidf, Tf, Tair, Qa, strairxT, strairyT, & ! LCOV_EXCL_LINE
434 : fsens, flat, fswabs, flwout, evap, Tref, Qref, & ! LCOV_EXCL_LINE
435 : scale_fluxes, frzmlt_init, frzmlt
436 : use ice_flux_bgc, only: faero_ocn, fiso_ocn, Qref_iso, fiso_evap, &
437 : flux_bio, flux_bio_ai
438 : use ice_grid, only: tmask
439 : use ice_state, only: aicen, aice
440 : #ifdef CICE_IN_NEMO
441 : use ice_state, only: aice_init
442 : use ice_flux, only: flatn_f, fsurfn_f
443 : #endif
444 : use ice_step_mod, only: ocean_mixed_layer
445 : use ice_timers, only: timer_couple, ice_timer_start, ice_timer_stop
446 :
447 : integer (kind=int_kind), intent(in) :: &
448 : iblk ! block index
449 :
450 : ! local variables
451 :
452 : integer (kind=int_kind) :: &
453 : ilo,ihi,jlo,jhi, & ! beginning and end of physical domain ! LCOV_EXCL_LINE
454 : n , & ! thickness category index ! LCOV_EXCL_LINE
455 : i,j , & ! horizontal indices ! LCOV_EXCL_LINE
456 : k , & ! tracer index ! LCOV_EXCL_LINE
457 : nbtrcr !
458 :
459 : type (block) :: &
460 : this_block ! block information for current block
461 :
462 : logical (kind=log_kind) :: &
463 : calc_Tsfc !
464 :
465 : real (kind=dbl_kind) :: &
466 : cszn , & ! counter for history averaging ! LCOV_EXCL_LINE
467 : puny , & ! ! LCOV_EXCL_LINE
468 : rhofresh , & ! ! LCOV_EXCL_LINE
469 5760 : netsw ! flag for shortwave radiation presence
470 :
471 : character(len=*), parameter :: subname = '(coupling_prep)'
472 :
473 22944 : call icepack_query_parameters(puny_out=puny, rhofresh_out=rhofresh)
474 22944 : call icepack_query_tracer_sizes(nbtrcr_out=nbtrcr)
475 22944 : call icepack_query_parameters(calc_Tsfc_out=calc_Tsfc)
476 22944 : call icepack_warnings_flush(nu_diag)
477 22944 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
478 0 : file=__FILE__, line=__LINE__)
479 :
480 : !-----------------------------------------------------------------
481 : ! Save current value of frzmlt for diagnostics.
482 : ! Update mixed layer with heat and radiation from ice.
483 : !-----------------------------------------------------------------
484 :
485 753576 : do j = 1, ny_block
486 16186320 : do i = 1, nx_block
487 16163376 : frzmlt_init (i,j,iblk) = frzmlt(i,j,iblk)
488 : enddo
489 : enddo
490 :
491 22944 : call ice_timer_start(timer_couple,iblk) ! atm/ocn coupling
492 :
493 22944 : if (oceanmixed_ice) &
494 22944 : call ocean_mixed_layer (dt,iblk) ! ocean surface fluxes and sst
495 :
496 : !-----------------------------------------------------------------
497 : ! Aggregate albedos
498 : !-----------------------------------------------------------------
499 :
500 753576 : do j = 1, ny_block
501 16186320 : do i = 1, nx_block
502 15432744 : alvdf(i,j,iblk) = c0
503 15432744 : alidf(i,j,iblk) = c0
504 15432744 : alvdr(i,j,iblk) = c0
505 15432744 : alidr(i,j,iblk) = c0
506 :
507 15432744 : albice(i,j,iblk) = c0
508 15432744 : albsno(i,j,iblk) = c0
509 15432744 : albpnd(i,j,iblk) = c0
510 15432744 : apeff_ai(i,j,iblk) = c0
511 15432744 : snowfrac(i,j,iblk) = c0
512 :
513 : ! for history averaging
514 15432744 : cszn = c0
515 15432744 : netsw = swvdr(i,j,iblk)+swidr(i,j,iblk)+swvdf(i,j,iblk)+swidf(i,j,iblk)
516 15432744 : if (netsw > puny) cszn = c1
517 31596120 : do n = 1, nstreams
518 30865488 : albcnt(i,j,iblk,n) = albcnt(i,j,iblk,n) + cszn
519 : enddo
520 : enddo
521 : enddo
522 :
523 22944 : this_block = get_block(blocks_ice(iblk),iblk)
524 22944 : ilo = this_block%ilo
525 22944 : ihi = this_block%ihi
526 22944 : jlo = this_block%jlo
527 22944 : jhi = this_block%jhi
528 :
529 137664 : do n = 1, ncat
530 3561384 : do j = jlo, jhi
531 69134640 : do i = ilo, ihi
532 69019920 : if (aicen(i,j,n,iblk) > puny) then
533 :
534 0 : alvdf(i,j,iblk) = alvdf(i,j,iblk) &
535 33915691 : + alvdfn(i,j,n,iblk)*aicen(i,j,n,iblk)
536 0 : alidf(i,j,iblk) = alidf(i,j,iblk) &
537 33915691 : + alidfn(i,j,n,iblk)*aicen(i,j,n,iblk)
538 0 : alvdr(i,j,iblk) = alvdr(i,j,iblk) &
539 33915691 : + alvdrn(i,j,n,iblk)*aicen(i,j,n,iblk)
540 0 : alidr(i,j,iblk) = alidr(i,j,iblk) &
541 33915691 : + alidrn(i,j,n,iblk)*aicen(i,j,n,iblk)
542 :
543 0 : netsw = swvdr(i,j,iblk) + swidr(i,j,iblk) &
544 33915691 : + swvdf(i,j,iblk) + swidf(i,j,iblk)
545 33915691 : if (netsw > puny) then ! sun above horizon
546 0 : albice(i,j,iblk) = albice(i,j,iblk) &
547 2702976 : + albicen(i,j,n,iblk)*aicen(i,j,n,iblk)
548 0 : albsno(i,j,iblk) = albsno(i,j,iblk) &
549 2702976 : + albsnon(i,j,n,iblk)*aicen(i,j,n,iblk)
550 0 : albpnd(i,j,iblk) = albpnd(i,j,iblk) &
551 2702976 : + albpndn(i,j,n,iblk)*aicen(i,j,n,iblk)
552 : endif
553 :
554 0 : apeff_ai(i,j,iblk) = apeff_ai(i,j,iblk) & ! for history
555 33915691 : + apeffn(i,j,n,iblk)*aicen(i,j,n,iblk)
556 0 : snowfrac(i,j,iblk) = snowfrac(i,j,iblk) & ! for history
557 33915691 : + snowfracn(i,j,n,iblk)*aicen(i,j,n,iblk)
558 :
559 : endif ! aicen > puny
560 : enddo
561 : enddo
562 : enddo
563 :
564 753576 : do j = 1, ny_block
565 16186320 : do i = 1, nx_block
566 :
567 : !-----------------------------------------------------------------
568 : ! reduce fresh by fpond for coupling
569 : !-----------------------------------------------------------------
570 :
571 15432744 : if (l_mpond_fresh) then
572 0 : fpond(i,j,iblk) = fpond(i,j,iblk) * rhofresh/dt
573 0 : fresh(i,j,iblk) = fresh(i,j,iblk) - fpond(i,j,iblk)
574 : endif
575 :
576 : !----------------------------------------------------------------
577 : ! Store grid box mean albedos and fluxes before scaling by aice
578 : !----------------------------------------------------------------
579 :
580 15432744 : alvdf_ai (i,j,iblk) = alvdf (i,j,iblk)
581 15432744 : alidf_ai (i,j,iblk) = alidf (i,j,iblk)
582 15432744 : alvdr_ai (i,j,iblk) = alvdr (i,j,iblk)
583 15432744 : alidr_ai (i,j,iblk) = alidr (i,j,iblk)
584 15432744 : fresh_ai (i,j,iblk) = fresh (i,j,iblk)
585 15432744 : fsalt_ai (i,j,iblk) = fsalt (i,j,iblk)
586 15432744 : fhocn_ai (i,j,iblk) = fhocn (i,j,iblk)
587 15432744 : fswthru_ai(i,j,iblk) = fswthru(i,j,iblk)
588 :
589 15432744 : if (nbtrcr > 0) then
590 0 : do k = 1, nbtrcr
591 0 : flux_bio_ai (i,j,k,iblk) = flux_bio (i,j,k,iblk)
592 : enddo
593 : endif
594 :
595 : !-----------------------------------------------------------------
596 : ! Save net shortwave for scaling factor in scale_factor
597 : !-----------------------------------------------------------------
598 0 : scale_factor(i,j,iblk) = &
599 : swvdr(i,j,iblk)*(c1 - alvdr_ai(i,j,iblk)) & ! LCOV_EXCL_LINE
600 : + swvdf(i,j,iblk)*(c1 - alvdf_ai(i,j,iblk)) & ! LCOV_EXCL_LINE
601 : + swidr(i,j,iblk)*(c1 - alidr_ai(i,j,iblk)) & ! LCOV_EXCL_LINE
602 16163376 : + swidf(i,j,iblk)*(c1 - alidf_ai(i,j,iblk))
603 :
604 : enddo
605 : enddo
606 :
607 : !-----------------------------------------------------------------
608 : ! Divide fluxes by ice area
609 : ! - the CESM coupler assumes fluxes are per unit ice area
610 : ! - also needed for global budget in diagnostics
611 : !-----------------------------------------------------------------
612 :
613 : call scale_fluxes (nx_block, ny_block, &
614 : tmask (:,:,iblk), nbtrcr, & ! LCOV_EXCL_LINE
615 : icepack_max_aero, & ! LCOV_EXCL_LINE
616 : aice (:,:,iblk), Tf (:,:,iblk), & ! LCOV_EXCL_LINE
617 : Tair (:,:,iblk), Qa (:,:,iblk), & ! LCOV_EXCL_LINE
618 : strairxT (:,:,iblk), strairyT(:,:,iblk), & ! LCOV_EXCL_LINE
619 : fsens (:,:,iblk), flat (:,:,iblk), & ! LCOV_EXCL_LINE
620 : fswabs (:,:,iblk), flwout (:,:,iblk), & ! LCOV_EXCL_LINE
621 : evap (:,:,iblk), & ! LCOV_EXCL_LINE
622 : Tref (:,:,iblk), Qref (:,:,iblk), & ! LCOV_EXCL_LINE
623 : fresh (:,:,iblk), fsalt (:,:,iblk), & ! LCOV_EXCL_LINE
624 : fhocn (:,:,iblk), & ! LCOV_EXCL_LINE
625 : fswthru (:,:,iblk), & ! LCOV_EXCL_LINE
626 : fswthru_vdr (:,:,iblk), & ! LCOV_EXCL_LINE
627 : fswthru_vdf (:,:,iblk), & ! LCOV_EXCL_LINE
628 : fswthru_idr (:,:,iblk), & ! LCOV_EXCL_LINE
629 : fswthru_idf (:,:,iblk), & ! LCOV_EXCL_LINE
630 : faero_ocn(:,:,:,iblk), & ! LCOV_EXCL_LINE
631 : alvdr (:,:,iblk), alidr (:,:,iblk), & ! LCOV_EXCL_LINE
632 : alvdf (:,:,iblk), alidf (:,:,iblk), & ! LCOV_EXCL_LINE
633 : flux_bio =flux_bio (:,:,1:nbtrcr,iblk), & ! LCOV_EXCL_LINE
634 : Qref_iso =Qref_iso (:,:,:,iblk), & ! LCOV_EXCL_LINE
635 : fiso_evap=fiso_evap(:,:,:,iblk), & ! LCOV_EXCL_LINE
636 22944 : fiso_ocn =fiso_ocn (:,:,:,iblk))
637 :
638 : #ifdef CICE_IN_NEMO
639 : !echmod - comment this out for efficiency, if .not. calc_Tsfc
640 : if (.not. calc_Tsfc) then
641 :
642 : !---------------------------------------------------------------
643 : ! If surface fluxes were provided, conserve these fluxes at ice
644 : ! free points by passing to ocean.
645 : !---------------------------------------------------------------
646 :
647 : call sfcflux_to_ocn &
648 : (nx_block, ny_block, & ! LCOV_EXCL_LINE
649 : tmask (:,:,iblk), aice_init(:,:,iblk), & ! LCOV_EXCL_LINE
650 : fsurfn_f (:,:,:,iblk), flatn_f(:,:,:,iblk), & ! LCOV_EXCL_LINE
651 : fresh (:,:,iblk), fhocn (:,:,iblk))
652 : endif
653 : !echmod
654 : #endif
655 22944 : call ice_timer_stop(timer_couple,iblk) ! atm/ocn coupling
656 :
657 22944 : end subroutine coupling_prep
658 :
659 : #ifdef CICE_IN_NEMO
660 :
661 : !=======================================================================
662 : !
663 : ! If surface heat fluxes are provided to CICE instead of CICE calculating
664 : ! them internally (i.e. .not. calc_Tsfc), then these heat fluxes can
665 : ! be provided at points which do not have ice. (This is could be due to
666 : ! the heat fluxes being calculated on a lower resolution grid or the
667 : ! heat fluxes not recalculated at every CICE timestep.) At ice free points,
668 : ! conserve energy and water by passing these fluxes to the ocean.
669 : !
670 : ! author: A. McLaren, Met Office
671 :
672 : subroutine sfcflux_to_ocn(nx_block, ny_block, &
673 : tmask, aice, & ! LCOV_EXCL_LINE
674 : fsurfn_f, flatn_f, & ! LCOV_EXCL_LINE
675 : fresh, fhocn)
676 :
677 : use ice_domain_size, only: ncat
678 :
679 : integer (kind=int_kind), intent(in) :: &
680 : nx_block, ny_block ! block dimensions
681 :
682 : logical (kind=log_kind), dimension (nx_block,ny_block), intent(in) :: &
683 : tmask ! land/boundary mask, thickness (T-cell)
684 :
685 : real (kind=dbl_kind), dimension(nx_block,ny_block), intent(in):: &
686 : aice ! initial ice concentration
687 :
688 : real (kind=dbl_kind), dimension(nx_block,ny_block,ncat), intent(in) :: &
689 : fsurfn_f, & ! net surface heat flux (provided as forcing) ! LCOV_EXCL_LINE
690 : flatn_f ! latent heat flux (provided as forcing)
691 :
692 : real (kind=dbl_kind), dimension(nx_block,ny_block), intent(inout):: &
693 : fresh , & ! fresh water flux to ocean (kg/m2/s) ! LCOV_EXCL_LINE
694 : fhocn ! actual ocn/ice heat flx (W/m**2)
695 :
696 :
697 : ! local variables
698 : integer (kind=int_kind) :: &
699 : i, j, n ! horizontal indices
700 :
701 : real (kind=dbl_kind) :: &
702 : puny, & ! ! LCOV_EXCL_LINE
703 : Lsub, & ! ! LCOV_EXCL_LINE
704 : rLsub ! 1/Lsub
705 :
706 : character(len=*), parameter :: subname = '(sfcflux_to_ocn)'
707 :
708 : call icepack_query_parameters(puny_out=puny, Lsub_out=Lsub)
709 : call icepack_warnings_flush(nu_diag)
710 : if (icepack_warnings_aborted()) call abort_ice(error_message=subname, &
711 : file=__FILE__, line=__LINE__)
712 : rLsub = c1 / Lsub
713 :
714 : do n = 1, ncat
715 : do j = 1, ny_block
716 : do i = 1, nx_block
717 : if (tmask(i,j) .and. aice(i,j) <= puny) then
718 : fhocn(i,j) = fhocn(i,j) &
719 : + fsurfn_f(i,j,n) + flatn_f(i,j,n)
720 : fresh(i,j) = fresh(i,j) &
721 : + flatn_f(i,j,n) * rLsub
722 : endif
723 : enddo ! i
724 : enddo ! j
725 : enddo ! n
726 :
727 :
728 : end subroutine sfcflux_to_ocn
729 :
730 : #endif
731 :
732 : !=======================================================================
733 :
734 : end module CICE_RunMod
735 :
736 : !=======================================================================
|