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 108 of file MOM_marine_ice.F90.

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

Referenced by mom_ocean_model_mct::update_ocean_model(), and mom_ocean_model_nuopc::update_ocean_model().

Here is the caller graph for this function:

◆ 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 47 of file MOM_marine_ice.F90.

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

References mom_domains::to_all.

Referenced by mom_ocean_model_mct::update_ocean_model(), and mom_ocean_model_nuopc::update_ocean_model().

Here is the caller graph for this function:

◆ 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 181 of file MOM_marine_ice.F90.

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

References mom_error_handler::mom_error().

Referenced by mom_ocean_model_mct::ocean_model_init(), and mom_ocean_model_nuopc::ocean_model_init().

Here is the call graph for this function:
Here is the caller graph for this function: