MOM6
mom_marine_ice Module Reference

Detailed Description

Routines incorporating the effects of marine ice (sea-ice and icebergs) into the ocean model dynamics and thermodynamics.

Data Types

type  marine_ice_cs
 Control structure for MOM_marine_ice. More...
 

Functions/Subroutines

subroutine, public iceberg_forces (G, forces, use_ice_shelf, sfc_state, time_step, CS)
 add_berg_flux_to_shelf adds rigidity and ice-area coverage due to icebergs to the forces type fields, and adds ice-areal coverage and modifies various thermodynamic fluxes due to the presence of icebergs. More...
 
subroutine, public iceberg_fluxes (G, US, fluxes, use_ice_shelf, sfc_state, time_step, CS)
 iceberg_fluxes adds ice-area-coverage and modifies various thermodynamic fluxes due to the presence of icebergs. More...
 
subroutine, public marine_ice_init (Time, G, param_file, diag, CS)
 Initialize control structure for MOM_marine_ice. More...
 

Function/Subroutine Documentation

◆ iceberg_fluxes()

subroutine, public mom_marine_ice::iceberg_fluxes ( type(ocean_grid_type), intent(inout)  G,
type(unit_scale_type), intent(in)  US,
type(forcing), intent(inout)  fluxes,
logical, intent(in)  use_ice_shelf,
type(surface), intent(inout)  sfc_state,
real, intent(in)  time_step,
type(marine_ice_cs), pointer  CS 
)

iceberg_fluxes adds ice-area-coverage and modifies various thermodynamic fluxes due to the presence of icebergs.

Parameters
[in,out]gThe ocean's grid structure
[in]usA dimensional unit scaling type
[in,out]fluxesA structure with pointers to themodynamic, tracer and mass exchange forcing fields
[in,out]sfc_stateA structure containing fields that describe the surface state of the ocean.
[in]use_ice_shelfIf true, this configuration uses ice shelves.
[in]time_stepThe coupling time step [s].
csPointer to the control structure for MOM_marine_ice

Definition at line 101 of file MOM_marine_ice.F90.

102  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
103  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
104  type(forcing), intent(inout) :: fluxes !< A structure with pointers to themodynamic,
105  !! tracer and mass exchange forcing fields
106  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
107  !! describe the surface state of the ocean.
108  logical, intent(in) :: use_ice_shelf !< If true, this configuration uses ice shelves.
109  real, intent(in) :: time_step !< The coupling time step [s].
110  type(marine_ice_CS), pointer :: CS !< Pointer to the control structure for MOM_marine_ice
111 
112  real :: fraz ! refreezing rate [R Z T-1 ~> kg m-2 s-1]
113  real :: I_dt_LHF ! The inverse of the timestep times the latent heat of fusion times [Q-1 T-1 ~> kg J-1 s-1].
114  integer :: i, j, is, ie, js, je, isd, ied, jsd, jed
115  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
116  isd = g%isd ; jsd = g%jsd ; ied = g%ied ; jed = g%jed
117  ! This routine adds iceberg data to the ice shelf data (if ice shelf is used)
118  ! which can then be used to change the top of ocean boundary condition used in
119  ! the ocean model. This routine is taken from the add_shelf_flux subroutine
120  ! within the ice shelf model.
121 
122  if (.not.associated(cs)) return
123  if (.not.(associated(fluxes%area_berg) .and. associated(fluxes%ustar_berg) .and. &
124  associated(fluxes%mass_berg) ) ) return
125  if (.not.(associated(fluxes%frac_shelf_h) .and. associated(fluxes%ustar_shelf)) ) return
126 
127 
128  if (.not.(associated(fluxes%area_berg) .and. associated(fluxes%ustar_berg) .and. &
129  associated(fluxes%mass_berg) ) ) return
130  if (.not. use_ice_shelf) then
131  fluxes%frac_shelf_h(:,:) = 0.
132  fluxes%ustar_shelf(:,:) = 0.
133  endif
134  do j=jsd,jed ; do i=isd,ied ; if (g%areaT(i,j) > 0.0) then
135  fluxes%frac_shelf_h(i,j) = fluxes%frac_shelf_h(i,j) + fluxes%area_berg(i,j)
136  fluxes%ustar_shelf(i,j) = fluxes%ustar_shelf(i,j) + fluxes%ustar_berg(i,j)
137  endif ; enddo ; enddo
138 
139  !Zero'ing out other fluxes under the tabular icebergs
140  if (cs%berg_area_threshold >= 0.) then
141  i_dt_lhf = 1.0 / (us%s_to_T*time_step * cs%latent_heat_fusion)
142  do j=jsd,jed ; do i=isd,ied
143  if (fluxes%frac_shelf_h(i,j) > cs%berg_area_threshold) then
144  ! Only applying for ice shelf covering most of cell.
145 
146  if (associated(fluxes%sw)) fluxes%sw(i,j) = 0.0
147  if (associated(fluxes%lw)) fluxes%lw(i,j) = 0.0
148  if (associated(fluxes%latent)) fluxes%latent(i,j) = 0.0
149  if (associated(fluxes%evap)) fluxes%evap(i,j) = 0.0
150 
151  ! Add frazil formation diagnosed by the ocean model [Q R Z ~> J m-2] in the
152  ! form of surface layer evaporation [R Z T-1 ~> kg m-2 s-1]. Update lprec in the
153  ! control structure for diagnostic purposes.
154 
155  if (allocated(sfc_state%frazil)) then
156  fraz = sfc_state%frazil(i,j) * i_dt_lhf
157  if (associated(fluxes%evap)) fluxes%evap(i,j) = fluxes%evap(i,j) - fraz
158  ! if (associated(fluxes%lprec)) fluxes%lprec(i,j) = fluxes%lprec(i,j) - fraz
159  sfc_state%frazil(i,j) = 0.0
160  endif
161 
162  !Alon: Should these be set to zero too?
163  if (associated(fluxes%sens)) fluxes%sens(i,j) = 0.0
164  if (associated(fluxes%salt_flux)) fluxes%salt_flux(i,j) = 0.0
165  if (associated(fluxes%lprec)) fluxes%lprec(i,j) = 0.0
166  endif
167  enddo ; enddo
168  endif
169 

◆ iceberg_forces()

subroutine, public mom_marine_ice::iceberg_forces ( type(ocean_grid_type), intent(inout)  G,
type(mech_forcing), intent(inout)  forces,
logical, intent(in)  use_ice_shelf,
type(surface), intent(inout)  sfc_state,
real, intent(in)  time_step,
type(marine_ice_cs), pointer  CS 
)

add_berg_flux_to_shelf adds rigidity and ice-area coverage due to icebergs to the forces type fields, and adds ice-areal coverage and modifies various thermodynamic fluxes due to the presence of icebergs.

Parameters
[in,out]gThe ocean's grid structure
[in,out]forcesA structure with the driving mechanical forces
[in,out]sfc_stateA structure containing fields that describe the surface state of the ocean.
[in]use_ice_shelfIf true, this configuration uses ice shelves.
[in]time_stepThe coupling time step [s].
csPointer to the control structure for MOM_marine_ice

Definition at line 45 of file MOM_marine_ice.F90.

46  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
47  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
48  type(surface), intent(inout) :: sfc_state !< A structure containing fields that
49  !! describe the surface state of the ocean.
50  logical, intent(in) :: use_ice_shelf !< If true, this configuration uses ice shelves.
51  real, intent(in) :: time_step !< The coupling time step [s].
52  type(marine_ice_CS), pointer :: CS !< Pointer to the control structure for MOM_marine_ice
53 
54  real :: kv_rho_ice ! The viscosity of ice divided by its density [L4 Z-2 T-1 R-1 ~> m5 kg-1 s-1].
55  integer :: i, j, is, ie, js, je
56  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
57  !This routine adds iceberg data to the ice shelf data (if ice shelf is used)
58  !which can then be used to change the top of ocean boundary condition used in
59  !the ocean model. This routine is taken from the add_shelf_flux subroutine
60  !within the ice shelf model.
61 
62  if (.not.associated(cs)) return
63 
64  if (.not.(associated(forces%area_berg) .and. associated(forces%mass_berg) ) ) return
65 
66  if (.not.(associated(forces%frac_shelf_u) .and. associated(forces%frac_shelf_v) .and. &
67  associated(forces%rigidity_ice_u) .and. associated(forces%rigidity_ice_v)) ) return
68 
69  ! This section sets or augments the values of fields in forces.
70  if (.not. use_ice_shelf) then
71  forces%frac_shelf_u(:,:) = 0.0 ; forces%frac_shelf_v(:,:) = 0.0
72  endif
73  if (.not. forces%accumulate_rigidity) then
74  forces%rigidity_ice_u(:,:) = 0.0 ; forces%rigidity_ice_v(:,:) = 0.0
75  endif
76 
77  call pass_var(forces%area_berg, g%domain, to_all+omit_corners, halo=1, complete=.false.)
78  call pass_var(forces%mass_berg, g%domain, to_all+omit_corners, halo=1, complete=.true.)
79  kv_rho_ice = cs%kv_iceberg / cs%density_iceberg
80  do j=js,je ; do i=is-1,ie
81  if ((g%areaT(i,j) + g%areaT(i+1,j) > 0.0)) & ! .and. (G%dxdy_u(I,j) > 0.0)) &
82  forces%frac_shelf_u(i,j) = forces%frac_shelf_u(i,j) + &
83  (forces%area_berg(i,j)*g%areaT(i,j) + forces%area_berg(i+1,j)*g%areaT(i+1,j)) / &
84  (g%areaT(i,j) + g%areaT(i+1,j))
85  forces%rigidity_ice_u(i,j) = forces%rigidity_ice_u(i,j) + kv_rho_ice * &
86  min(forces%mass_berg(i,j), forces%mass_berg(i+1,j))
87  enddo ; enddo
88  do j=js-1,je ; do i=is,ie
89  if ((g%areaT(i,j) + g%areaT(i,j+1) > 0.0)) & ! .and. (G%dxdy_v(i,J) > 0.0)) &
90  forces%frac_shelf_v(i,j) = forces%frac_shelf_v(i,j) + &
91  (forces%area_berg(i,j)*g%areaT(i,j) + forces%area_berg(i,j+1)*g%areaT(i,j+1)) / &
92  (g%areaT(i,j) + g%areaT(i,j+1))
93  forces%rigidity_ice_v(i,j) = forces%rigidity_ice_v(i,j) + kv_rho_ice * &
94  min(forces%mass_berg(i,j), forces%mass_berg(i,j+1))
95  enddo ; enddo
96 

◆ marine_ice_init()

subroutine, public mom_marine_ice::marine_ice_init ( type(time_type), intent(in), target  Time,
type(ocean_grid_type), intent(in)  G,
type(param_file_type), intent(in)  param_file,
type(diag_ctrl), intent(inout), target  diag,
type(marine_ice_cs), pointer  CS 
)

Initialize control structure for MOM_marine_ice.

Parameters
[in]timeCurrent model time
[in]gOcean grid structure
[in]param_fileRuntime parameter handles
[in,out]diagDiagnostics control structure
csPointer to the control structure for MOM_marine_ice

Definition at line 173 of file MOM_marine_ice.F90.

174  type(time_type), target, intent(in) :: Time !< Current model time
175  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
176  type(param_file_type), intent(in) :: param_file !< Runtime parameter handles
177  type(diag_ctrl), target, intent(inout) :: diag !< Diagnostics control structure
178  type(marine_ice_CS), pointer :: CS !< Pointer to the control structure for MOM_marine_ice
179 ! This include declares and sets the variable "version".
180 #include "version_variable.h"
181  character(len=40) :: mdl = "MOM_marine_ice" ! This module's name.
182 
183  if (associated(cs)) then
184  call mom_error(warning, "marine_ice_init called with an associated control structure.")
185  return
186  else ; allocate(cs) ; endif
187 
188  ! Write all relevant parameters to the model log.
189  call log_version(mdl, version)
190 
191  call get_param(param_file, mdl, "KV_ICEBERG", cs%kv_iceberg, &
192  "The viscosity of the icebergs", &
193  units="m2 s-1", default=1.0e10, scale=g%US%Z_to_L**2*g%US%m_to_L**2*g%US%T_to_s)
194  call get_param(param_file, mdl, "DENSITY_ICEBERGS", cs%density_iceberg, &
195  "A typical density of icebergs.", units="kg m-3", default=917.0, scale=g%US%kg_m3_to_R)
196  call get_param(param_file, mdl, "LATENT_HEAT_FUSION", cs%latent_heat_fusion, &
197  "The latent heat of fusion.", units="J/kg", default=hlf, scale=g%US%J_kg_to_Q)
198  call get_param(param_file, mdl, "BERG_AREA_THRESHOLD", cs%berg_area_threshold, &
199  "Fraction of grid cell which iceberg must occupy, so that fluxes "//&
200  "below berg are set to zero. Not applied for negative "//&
201  "values.", units="non-dim", default=-1.0)
202