!||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||| module ice !BOP ! !MODULE: ice ! ! !DESCRIPTION: ! This module currently contains routines for computing ice ! formation and the heat flux associated with ice formation. ! This heat flux is sent to the ice model via the flux coupler. ! In the future, this module could contain the driver for a ! subroutinized ice model. ! ! !REVISION HISTORY: ! SVN:$Id: ice.F90 17759 2009-08-12 20:20:36Z njn01 $ ! !USES: use TIMCOM_KindsMod use TIMCOM_IOUnitsMod use kinds_mod, only: int_kind, log_kind, char_len, r8 use blocks, only: nx_block, ny_block, block use domain_size use constants, only: cp_sw, latent_heat_fusion, delim_fmt, blank_fmt, & sea_ice_salinity, ppt_to_salt, ocn_ref_salinity, c0, p5, hflux_factor,& grav, eps2, ndelim_fmt use broadcast, only: broadcast_scalar use communicate, only: my_task, master_task use io_types, only: nml_in, nml_filename, stdout use time_management, only: freq_opt_never, freq_opt_nyear, freq_opt_nday, & freq_opt_nhour, freq_opt_nsecond, freq_opt_nstep, init_time_flag, & max_blocks_clinic, km, nt, avg_ts, back_to_back, dtt, check_time_flag,& partial_bottom_cells, KMT, DZT, DZ, freq_opt_nmonth, dt, tmix_matsuno,& tmix_iopt, ice_ts, access_time_flag use exit_mod, only: sigAbort, exit_timcom, flushm use prognostic ! use passive_tracers, only: tracer_ref_val use grid!, only: sfc_layer_varthick, sfc_layer_type implicit none private save ! !PUBLIC MEMBER FUNCTIONS: public :: init_ice, & increment_tlast_ice,& ! ice_formation, & ! tfreez, & ! ice_flx_to_coupler, & ! tmelt tavg_increment_sum_qflux ! !PUBLIC DATA MEMBERS: logical (log_kind), public :: & liceform ! flag to turn on/off ice formation logical (log_kind), public :: & lactive_ice ! T ==> ocn is coupled to an active ice model ! F ==> ocn is coupled to a dummy ice model integer (int_kind), public :: & ice_cpl_flag ! time flag id for coupled timestep real (r8), public :: & tlast_ice, &! time since last ice flux computed time_weight, & salice ! sea ice salinity in msu real (r8), dimension(:,:,:), allocatable, public :: & AQICE, &! sum of accumulated ice heat flux since tlast QFLUX, &! internal ocn heat flux due to ice formation QICE ! tot column cooling from ice form (in C*cm) real (r8), dimension(nx_block,ny_block,max_blocks_clinic), public :: & FW_FREEZE ! water flux at T points due to frazil ice formation real (r8), public :: & cp_over_lhfusion ! cp_sw/latent_heat_fusion real (r8), public :: & tavg_sum_qflux !EOP !BOC !----------------------------------------------------------------------- ! ! module variables ! !----------------------------------------------------------------------- integer (int_kind) :: & ice_flag, &! time flag id for ice formation kmxice ! lowest level from which to integrate ! ice formation real (r8) :: & salref ! ocean ref salinity in msu !EOC !*********************************************************************** contains !*********************************************************************** !BOP ! !IROUTINE: ! !INTERFACE: subroutine init_ice ! !DESCRIPTION: ! This routine initializes ice variables. It must be called ! before initializing restarts because this module add the accumulated ! ice heat flux to the restart file. ! ! !REVISION HISTORY: ! same as module !EOP !BOC !----------------------------------------------------------------------- ! ! local variables ! !----------------------------------------------------------------------- integer (int_kind) :: & nml_error ! namelist i/o error flag character (char_len) :: & ice_freq_opt ! option for frequency of computing ice integer (int_kind) :: & ice_freq_iopt, &! int option for freq units ice_freq ! freq for computing ice in above units namelist /ice_nml/ kmxice, ice_freq_opt, ice_freq, lactive_ice !----------------------------------------------------------------------- ! ! read input namelists ! !----------------------------------------------------------------------- cp_over_lhfusion = rho_sw*cp_sw/(latent_heat_fusion*rho_fw) kmxice = 1 ice_freq_opt = 'never' ice_freq = 100000 lactive_ice = .false. if (my_task == master_task) then open (nml_in, file=nml_filename, status='old',iostat=nml_error) if (nml_error /= 0) then nml_error = -1 else nml_error = 1 endif do while (nml_error > 0) read(nml_in, nml=ice_nml,iostat=nml_error) end do if (nml_error == 0) close(nml_in) endif call broadcast_scalar(nml_error, master_task) if (nml_error /= 0) then call exit_TIMCOM(sigAbort,'ERROR reading ice_nml') endif if (my_task == master_task) then write(stdout,blank_fmt) write(stdout,ndelim_fmt) write(stdout,blank_fmt) write(stdout,*) 'Ice:' write(stdout,blank_fmt) write(stdout,*) ' ice_nml namelist settings:' write(stdout,blank_fmt) write(stdout,ice_nml) write(stdout,delim_fmt) !*** !*** define salice and salref in msu !*** salice = sea_ice_salinity*ppt_to_salt salref = ocn_ref_salinity*ppt_to_salt liceform = .true. select case (ice_freq_opt) case ('never') write(stdout,'(a22)') 'Ice formation disabled' liceform = .false. ice_freq_iopt = freq_opt_never case ('coupled') write(stdout,'(a44)') & 'Ice formation computed on coupler time steps' ice_freq_iopt = freq_opt_never ! check coupler flag instead case ('nyear') write(stdout,'(a29,i4,a9)') 'Ice formation computed every ', & ice_freq,' years. ' ice_freq_iopt = freq_opt_nyear case ('nmonth') write(stdout,'(a29,i4,a9)') 'Ice formation computed every ', & ice_freq,' months. ' ice_freq_iopt = freq_opt_nmonth case ('nday') write(stdout,'(a29,i4,a9)') 'Ice formation computed every ', & ice_freq,' days. ' ice_freq_iopt = freq_opt_nday case ('nhour') write(stdout,'(a29,i4,a9)') 'Ice formation computed every ', & ice_freq,' hours. ' ice_freq_iopt = freq_opt_nhour case ('nsecond') write(stdout,'(a29,i4,a9)') 'Ice formation computed every ', & ice_freq,' seconds.' ice_freq_iopt = freq_opt_nsecond case ('nstep') write(stdout,'(a29,i4,a9)') 'Ice formation computed every ', & ice_freq,' steps. ' ice_freq_iopt = freq_opt_nstep case default ice_freq_iopt = -1000 end select if (liceform) then write(stdout,'(a20,1pe10.3)') 'Ice salinity(ppt) = ', & sea_ice_salinity write(stdout,'(a30,i3,a13)') 'Ice formation computed in top ', & kmxice, ' levels only.' endif call TIMCOM_IOUnitsFlush(TIMCOM_stdout) ; call TIMCOM_IOUnitsFlush(stdout) endif call broadcast_scalar(ice_freq_iopt, master_task) call broadcast_scalar(lactive_ice, master_task) if (ice_freq_iopt == -1000) then call exit_TIMCOM(sigAbort,'unknown restart frequency option') endif call broadcast_scalar(liceform, master_task) !----------------------------------------------------------------------- ! ! if ice turned on, broadcast remaining vars and allocate memory ! !----------------------------------------------------------------------- if (liceform) then call broadcast_scalar(ice_freq, master_task) call broadcast_scalar(kmxice, master_task) call broadcast_scalar(salice, master_task) call broadcast_scalar(salref, master_task) !*** !*** set up ice time flag and get coupled_ts id for local use !*** call init_time_flag('ice',ice_flag, default=.false., & owner = 'init_ice', & freq_opt = ice_freq_iopt, & freq = ice_freq) call access_time_flag('coupled_ts', ice_cpl_flag) !*** !*** must keep track of time since last ice flux computed !*** tlast_ice = c0 !*** !*** allocate and initialize ice flux arrays !*** allocate( QICE(nx_block,ny_block,max_blocks_clinic), & AQICE(nx_block,ny_block,max_blocks_clinic), & QFLUX(nx_block,ny_block,max_blocks_clinic)) QICE = c0 AQICE = c0 QFLUX = c0 endif FW_FREEZE = c0 !----------------------------------------------------------------------- !EOC call flushm (stdout) end subroutine init_ice !*********************************************************************** !BOP ! !IROUTINE: ! !INTERFACE: subroutine increment_tlast_ice use var, only: daodt ! !DESCRIPTION: ! This subroutine increments tlast_ice in a nonthreaded region. ! !REVISION HISTORY: ! same as module !EOP !BOC !----------------------------------------------------------------------- ! ! increment time since last evaluation ! !----------------------------------------------------------------------- ! if (avg_ts .or. back_to_back) then ! time_weight = p5 ! else time_weight = c1 ! endif ! tlast_ice = tlast_ice + dtt*time_weight tlast_ice = tlast_ice + (86400.d0/daodt)*time_weight !----------------------------------------------------------------------- !EOC end subroutine increment_tlast_ice !*********************************************************************** !*********************************************************************** !BOP ! !IROUTINE: tavg_increment_sum_qflux ! !INTERFACE: subroutine tavg_increment_sum_qflux(const) ! !DESCRIPTION: ! Increment the scalar tavg_sum_qflux with the weight for this timestep. ! ! !REVISION HISTORY: ! same as module ! !INPUT PARAMETERS: real (r8), intent(in) :: & const !EOP !BOC !----------------------------------------------------------------------- tavg_sum_qflux = tavg_sum_qflux + const ! const = tlast_ice !----------------------------------------------------------------------- !EOC end subroutine tavg_increment_sum_qflux !*********************************************************************** !*********************************************************************** end module ice !|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||