MOM6
MOM_marine_ice.F90
1 !> Routines incorporating the effects of marine ice (sea-ice and icebergs) into
2 !! the ocean model dynamics and thermodynamics.
4 
5 ! This file is part of MOM6. See LICENSE.md for the license.
6 
7 use mom_constants, only : hlf
8 use mom_diag_mediator, only : post_data, query_averaging_enabled, diag_ctrl
9 use mom_domains, only : pass_var, pass_vector, agrid, bgrid_ne, cgrid_ne
10 use mom_domains, only : to_all, omit_corners
11 use mom_error_handler, only : mom_error, fatal, warning
12 use mom_error_handler, only : mom_error, mom_mesg, fatal, warning
16 use mom_grid, only : ocean_grid_type
17 use mom_time_manager, only : time_type
19 use mom_variables, only : surface
20 
21 implicit none ; private
22 
23 #include <MOM_memory.h>
24 
25 public iceberg_forces, iceberg_fluxes, marine_ice_init
26 
27 !> Control structure for MOM_marine_ice
28 type, public :: marine_ice_cs ; private
29  real :: kv_iceberg !< The viscosity of the icebergs [L4 Z-2 T-1 ~> m2 s-1] (for ice rigidity)
30  real :: berg_area_threshold !< Fraction of grid cell which iceberg must occupy
31  !! so that fluxes below are set to zero. (0.5 is a
32  !! good value to use.) Not applied for negative values.
33  real :: latent_heat_fusion !< Latent heat of fusion [Q ~> J kg-1]
34  real :: density_iceberg !< A typical density of icebergs [R ~> kg m-3] (for ice rigidity)
35 
36  type(time_type), pointer :: time !< A pointer to the ocean model's clock.
37  type(diag_ctrl), pointer :: diag !< A structure that is used to regulate the timing of diagnostic output.
38 end type marine_ice_cs
39 
40 contains
41 
42 !> add_berg_flux_to_shelf adds rigidity and ice-area coverage due to icebergs
43 !! to the forces type fields, and adds ice-areal coverage and modifies various
44 !! thermodynamic fluxes due to the presence of icebergs.
45 subroutine iceberg_forces(G, forces, use_ice_shelf, sfc_state, time_step, CS)
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 
97 end subroutine iceberg_forces
98 
99 !> iceberg_fluxes adds ice-area-coverage and modifies various
100 !! thermodynamic fluxes due to the presence of icebergs.
101 subroutine iceberg_fluxes(G, US, fluxes, use_ice_shelf, sfc_state, time_step, CS)
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 
170 end subroutine iceberg_fluxes
171 
172 !> Initialize control structure for MOM_marine_ice
173 subroutine marine_ice_init(Time, G, param_file, diag, CS)
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 
203 end subroutine marine_ice_init
204 
205 end module mom_marine_ice
mom_time_manager
Wraps the FMS time manager functions.
Definition: MOM_time_manager.F90:2
mom_forcing_type::mech_forcing
Structure that contains pointers to the mechanical forcing at the surface used to drive the liquid oc...
Definition: MOM_forcing_type.F90:204
mom_variables::surface
Pointers to various fields which may be used describe the surface state of MOM, and which will be ret...
Definition: MOM_variables.F90:41
mom_marine_ice
Routines incorporating the effects of marine ice (sea-ice and icebergs) into the ocean model dynamics...
Definition: MOM_marine_ice.F90:3
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_constants
Provides a few physical constants.
Definition: MOM_constants.F90:2
mom_diag_mediator
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Definition: MOM_diag_mediator.F90:3
mom_domains::pass_var
Do a halo update on an array.
Definition: MOM_domains.F90:54
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
mom_diag_mediator::post_data
Make a diagnostic available for averaging or output.
Definition: MOM_diag_mediator.F90:70
mom_domains::pass_vector
Do a halo update on a pair of arrays representing the two components of a vector.
Definition: MOM_domains.F90:59
mom_forcing_type
This module implements boundary forcing for MOM6.
Definition: MOM_forcing_type.F90:2
mom_forcing_type::allocate_forcing_type
Allocate the fields of a (flux) forcing type, based on either a set of input flags for each group of ...
Definition: MOM_forcing_type.F90:43
mom_domains
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
mom_variables
Provides transparent structures with groups of MOM6 variables and supporting routines.
Definition: MOM_variables.F90:2
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_forcing_type::forcing
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
Definition: MOM_forcing_type.F90:66
mom_marine_ice::marine_ice_cs
Control structure for MOM_marine_ice.
Definition: MOM_marine_ice.F90:28
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26
mom_diag_mediator::diag_ctrl
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Definition: MOM_diag_mediator.F90:240