MOM6
user_shelf_init Module Reference

Detailed Description

This module specifies the initial values and evolving properties of the MOM6 ice shelf, using user-provided code.

Data Types

type  user_ice_shelf_cs
 The control structure for the user_ice_shelf module. More...
 

Functions/Subroutines

subroutine, public user_initialize_shelf_mass (mass_shelf, area_shelf_h, h_shelf, hmask, G, US, CS, param_file, new_sim)
 This subroutine sets up the initial mass and area covered by the ice shelf, based on user-provided code. More...
 
subroutine, public user_init_ice_thickness (h_shelf, area_shelf_h, hmask, G, US, param_file)
 This subroutine updates the ice shelf thickness, as specified by user-provided code. More...
 
subroutine, public user_update_shelf_mass (mass_shelf, area_shelf_h, h_shelf, hmask, G, CS, Time, new_sim)
 This subroutine updates the ice shelf mass, as specified by user-provided code. More...
 
subroutine write_user_log (param_file)
 This subroutine writes out the user ice shelf code version number to the model log. More...
 

Function/Subroutine Documentation

◆ user_init_ice_thickness()

subroutine, public user_shelf_init::user_init_ice_thickness ( real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(out)  h_shelf,
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(out)  area_shelf_h,
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(out)  hmask,
type(ocean_grid_type), intent(in)  G,
type(unit_scale_type), intent(in)  US,
type(param_file_type), intent(in)  param_file 
)

This subroutine updates the ice shelf thickness, as specified by user-provided code.

Parameters
[in]gThe ocean's grid structure
[out]h_shelfThe ice shelf thickness [m].
[out]area_shelf_hThe area per cell covered by the ice shelf [L2 ~> m2].
[out]hmaskA mask indicating which tracer points are
[in]usA structure containing unit conversion factors
[in]param_fileA structure to parse for run-time parameters

Definition at line 103 of file user_shelf_init.F90.

103  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
104  real, dimension(SZDI_(G),SZDJ_(G)), &
105  intent(out) :: h_shelf !< The ice shelf thickness [m].
106  real, dimension(SZDI_(G),SZDJ_(G)), &
107  intent(out) :: area_shelf_h !< The area per cell covered by the ice shelf [L2 ~> m2].
108  real, dimension(SZDI_(G),SZDJ_(G)), &
109  intent(out) :: hmask !< A mask indicating which tracer points are
110  !! partly or fully covered by an ice-shelf
111  type(unit_scale_type), intent(in) :: us !< A structure containing unit conversion factors
112  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
113 
114  ! This subroutine initializes the ice shelf thickness. Currently it does so
115  ! calling USER_initialize_shelf_mass, but this can be revised as needed.
116  real, dimension(SZI_(G),SZJ_(G)) :: mass_shelf
117  type(user_ice_shelf_cs), pointer :: cs => null()
118 
119  call user_initialize_shelf_mass(mass_shelf, area_shelf_h, h_shelf, hmask, g, us, cs, param_file, .true.)
120 

◆ user_initialize_shelf_mass()

subroutine, public user_shelf_init::user_initialize_shelf_mass ( real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(out)  mass_shelf,
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(out)  area_shelf_h,
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(out)  h_shelf,
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(out)  hmask,
type(ocean_grid_type), intent(in)  G,
type(unit_scale_type), intent(in)  US,
type(user_ice_shelf_cs), pointer  CS,
type(param_file_type), intent(in)  param_file,
logical, intent(in)  new_sim 
)

This subroutine sets up the initial mass and area covered by the ice shelf, based on user-provided code.

Parameters
[in]gThe ocean's grid structure
[out]mass_shelfThe ice shelf mass per unit area averaged
[out]h_shelfThe ice shelf thickness [Z ~> m].
[out]area_shelf_hThe area per cell covered by the ice shelf [L2 ~> m2].
[out]hmaskA mask indicating which tracer points are
[in]usA structure containing unit conversion factors
csA pointer to the user ice shelf control structure
[in]param_fileA structure to parse for run-time parameters
[in]new_simIf true, this is a new run; otherwise it is being started from a restart file.

Definition at line 44 of file user_shelf_init.F90.

44 
45  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
46  real, dimension(SZDI_(G),SZDJ_(G)), &
47  intent(out) :: mass_shelf !< The ice shelf mass per unit area averaged
48  !! over the full ocean cell [R Z ~> kg m-2].
49  real, dimension(SZDI_(G),SZDJ_(G)), &
50  intent(out) :: h_shelf !< The ice shelf thickness [Z ~> m].
51  real, dimension(SZDI_(G),SZDJ_(G)), &
52  intent(out) :: area_shelf_h !< The area per cell covered by the ice shelf [L2 ~> m2].
53  real, dimension(SZDI_(G),SZDJ_(G)), &
54  intent(out) :: hmask !< A mask indicating which tracer points are
55  !! partly or fully covered by an ice-shelf
56  type(unit_scale_type), intent(in) :: us !< A structure containing unit conversion factors
57  type(user_ice_shelf_cs), pointer :: cs !< A pointer to the user ice shelf control structure
58  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
59  logical, intent(in) :: new_sim !< If true, this is a new run; otherwise it is
60  !! being started from a restart file.
61 
62 ! This subroutine sets up the initial mass and area covered by the ice shelf.
63  real :: max_draft ! The maximum ocean draft of the ice shelf [Z ~> m].
64  real :: min_draft ! The minimum ocean draft of the ice shelf [Z ~> m].
65  real :: flat_shelf_width ! The range over which the shelf is min_draft thick.
66  real :: c1 ! The maximum depths in m.
67  character(len=40) :: mdl = "USER_initialize_shelf_mass" ! This subroutine's name.
68  integer :: i, j
69 
70  ! call MOM_error(FATAL, "USER_shelf_init.F90, USER_set_shelf_mass: " // &
71  ! "Unmodified user routine called - you must edit the routine to use it")
72 
73  if (.not.associated(cs)) allocate(cs)
74 
75  ! Read all relevant parameters and write them to the model log.
76  if (cs%first_call) call write_user_log(param_file)
77  cs%first_call = .false.
78  call get_param(param_file, mdl, "RHO_0", cs%Rho_ocean, &
79  "The mean ocean density used with BOUSSINESQ true to "//&
80  "calculate accelerations and the mass for conservation "//&
81  "properties, or with BOUSSINSEQ false to convert some "//&
82  "parameters from vertical units of m to kg m-2.", &
83  units="kg m-3", default=1035.0, scale=us%kg_m3_to_R)
84  call get_param(param_file, mdl, "SHELF_MAX_DRAFT", cs%max_draft, &
85  units="m", default=1.0, scale=us%m_to_Z)
86  call get_param(param_file, mdl, "SHELF_MIN_DRAFT", cs%min_draft, &
87  units="m", default=1.0, scale=us%m_to_Z)
88  call get_param(param_file, mdl, "FLAT_SHELF_WIDTH", cs%flat_shelf_width, &
89  units="axis_units", default=0.0)
90  call get_param(param_file, mdl, "SHELF_SLOPE_SCALE", cs%shelf_slope_scale, &
91  units="axis_units", default=0.0)
92  call get_param(param_file, mdl, "SHELF_EDGE_POS_0", cs%pos_shelf_edge_0, &
93  units="axis_units", default=0.0)
94  call get_param(param_file, mdl, "SHELF_SPEED", cs%shelf_speed, &
95  units="axis_units day-1", default=0.0)
96 
97  call user_update_shelf_mass(mass_shelf, area_shelf_h, h_shelf, hmask, g, cs, set_time(0,0), new_sim)
98 

◆ user_update_shelf_mass()

subroutine, public user_shelf_init::user_update_shelf_mass ( real, dimension(szdi_(g),szdj_(g)), intent(inout)  mass_shelf,
real, dimension(szdi_(g),szdj_(g)), intent(inout)  area_shelf_h,
real, dimension(szdi_(g),szdj_(g)), intent(inout)  h_shelf,
real, dimension(szdi_(g),szdj_(g)), intent(inout)  hmask,
type(ocean_grid_type), intent(in)  G,
type(user_ice_shelf_cs), pointer  CS,
type(time_type), intent(in)  Time,
logical, intent(in)  new_sim 
)

This subroutine updates the ice shelf mass, as specified by user-provided code.

Parameters
[in]gThe ocean's grid structure
[in,out]mass_shelfThe ice shelf mass per unit area averaged
[in,out]area_shelf_hThe area per cell covered by the ice shelf [L2 ~> m2].
[in,out]h_shelfThe ice shelf thickness [Z ~> m].
[in,out]hmaskA mask indicating which tracer points are
csA pointer to the user ice shelf control structure
[in]timeThe current model time
[in]new_simIf true, this the start of a new run.

Definition at line 125 of file user_shelf_init.F90.

125  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
126  real, dimension(SZDI_(G),SZDJ_(G)), &
127  intent(inout) :: mass_shelf !< The ice shelf mass per unit area averaged
128  !! over the full ocean cell [R Z ~> kg m-2].
129  real, dimension(SZDI_(G),SZDJ_(G)), &
130  intent(inout) :: area_shelf_h !< The area per cell covered by the ice shelf [L2 ~> m2].
131  real, dimension(SZDI_(G),SZDJ_(G)), &
132  intent(inout) :: h_shelf !< The ice shelf thickness [Z ~> m].
133  real, dimension(SZDI_(G),SZDJ_(G)), &
134  intent(inout) :: hmask !< A mask indicating which tracer points are
135  !! partly or fully covered by an ice-shelf
136  type(user_ice_shelf_cs), pointer :: cs !< A pointer to the user ice shelf control structure
137  type(time_type), intent(in) :: time !< The current model time
138  logical, intent(in) :: new_sim !< If true, this the start of a new run.
139 
140 
141  real :: c1, edge_pos, slope_pos
142  integer :: i, j
143 
144  edge_pos = cs%pos_shelf_edge_0 + cs%shelf_speed*(time_type_to_real(time) / 86400.0)
145 
146  slope_pos = edge_pos - cs%flat_shelf_width
147  c1 = 0.0 ; if (cs%shelf_slope_scale > 0.0) c1 = 1.0 / cs%shelf_slope_scale
148 
149 
150  do j=g%jsd,g%jed
151 
152  if (((j+g%jdg_offset) <= g%domain%njglobal+g%domain%njhalo) .AND. &
153  ((j+g%jdg_offset) >= g%domain%njhalo+1)) then
154 
155  do i=g%isc,g%iec
156 
157 ! if (((i+G%idg_offset) <= G%domain%niglobal+G%domain%nihalo) .AND. &
158 ! ((i+G%idg_offset) >= G%domain%nihalo+1)) then
159 
160  if ((j >= g%jsc) .and. (j <= g%jec)) then
161 
162  if (new_sim) then ; if (g%geoLonCu(i-1,j) >= edge_pos) then
163  ! Everything past the edge is open ocean.
164  mass_shelf(i,j) = 0.0
165  area_shelf_h(i,j) = 0.0
166  hmask(i,j) = 0.0
167  h_shelf(i,j) = 0.0
168  else
169  if (g%geoLonCu(i,j) > edge_pos) then
170  area_shelf_h(i,j) = g%areaT(i,j) * (edge_pos - g%geoLonCu(i-1,j)) / &
171  (g%geoLonCu(i,j) - g%geoLonCu(i-1,j))
172  hmask(i,j) = 2.0
173  else
174  area_shelf_h(i,j) = g%areaT(i,j)
175  hmask(i,j) = 1.0
176  endif
177 
178  if (g%geoLonT(i,j) > slope_pos) then
179  h_shelf(i,j) = cs%min_draft
180  mass_shelf(i,j) = cs%Rho_ocean * cs%min_draft
181  else
182  mass_shelf(i,j) = cs%Rho_ocean * (cs%min_draft + &
183  (cs%max_draft - cs%min_draft) * &
184  min(1.0, (c1*(slope_pos - g%geoLonT(i,j)))**2) )
185  h_shelf(i,j) = (cs%min_draft + &
186  (cs%max_draft - cs%min_draft) * &
187  min(1.0, (c1*(slope_pos - g%geoLonT(i,j)))**2) )
188  endif
189 
190  endif ; endif ; endif
191 
192  if ((i+g%idg_offset) == g%domain%nihalo+1) then
193  hmask(i-1,j) = 3.0
194  endif
195 
196  enddo ; endif ; enddo
197 

◆ write_user_log()

subroutine user_shelf_init::write_user_log ( type(param_file_type), intent(in)  param_file)
private

This subroutine writes out the user ice shelf code version number to the model log.

Parameters
[in]param_fileA structure to parse for run-time parameters

Definition at line 202 of file user_shelf_init.F90.

202  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
203 
204  character(len=128) :: version = '$Id: user_shelf_init.F90,v 1.1.2.7 2012/06/19 22:15:52 Robert.Hallberg Exp $'
205  character(len=128) :: tagname = '$Name: MOM_ogrp $'
206  character(len=40) :: mdl = "user_shelf_init" ! This module's name.
207 
208  call log_version(param_file, mdl, version, tagname)
209