MOM6
user_shelf_init Module Reference

Data Types

type  user_ice_shelf_cs
 

Functions/Subroutines

subroutine, public user_initialize_shelf_mass (mass_shelf, area_shelf_h, h_shelf, hmask, G, CS, param_file, new_sim)
 
subroutine, public user_init_ice_thickness (h_shelf, area_shelf_h, hmask, G, param_file)
 
subroutine, public user_update_shelf_mass (mass_shelf, area_shelf_h, h_shelf, hmask, G, CS, Time, new_sim)
 
subroutine write_user_log (param_file)
 

Variables

logical first_call = .true.
 

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(param_file_type), intent(in)  param_file 
)
Parameters
[in]gThe ocean's grid structure
[in]param_fileA structure to parse for run-time parameters

Definition at line 173 of file user_shelf_init.F90.

References user_initialize_shelf_mass().

Referenced by mom_ice_shelf_initialize::initialize_ice_thickness().

173  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
174  real, dimension(SZI_(G),SZJ_(G)), intent(out) :: area_shelf_h, hmask, h_shelf
175  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
176 
177  ! This subroutine initializes the ice shelf thickness. Currently it does so
178  ! calling USER_initialize_shelf_mass, but this can be revised as needed.
179  real, dimension(SZI_(G),SZJ_(G)) :: mass_shelf
180  type(user_ice_shelf_cs), pointer :: cs => null()
181 
182  call user_initialize_shelf_mass(mass_shelf, area_shelf_h, h_shelf, hmask, g, cs, param_file, .true.)
183 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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(user_ice_shelf_cs), pointer  CS,
type(param_file_type), intent(in)  param_file,
logical  new_sim 
)
Parameters
[in]gThe ocean's grid structure
[in]param_fileA structure to parse for run-time parameters

Definition at line 115 of file user_shelf_init.F90.

References first_call, user_update_shelf_mass(), and write_user_log().

Referenced by user_init_ice_thickness().

115 
116  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
117  real, dimension(SZI_(G),SZJ_(G)), intent(out) :: mass_shelf, area_shelf_h, hmask, h_shelf
118  type(user_ice_shelf_cs), pointer :: cs
119  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
120  logical :: new_sim
121 
122 ! Arguments: mass_shelf - The mass per unit area averaged over the full ocean
123 ! cell, in kg m-2. (Intent out)
124 ! (out) area_shelf_h - The area of the ocean cell that is covered by the
125 ! rigid ice shelf, in m2.
126 ! (in) G - The ocean's grid structure.
127 ! (in) param_file - A structure indicating the open file to parse for
128 ! model parameter values.
129 
130 
131 ! just check for cvs
132 ! This subroutine sets up the initial mass and area covered by the ice shelf.
133  real :: rho_ocean ! The ocean's typical density, in kg m-3.
134  real :: max_draft ! The maximum ocean draft of the ice shelf, in m.
135  real :: min_draft ! The minimum ocean draft of the ice shelf, in m.
136  real :: flat_shelf_width ! The range over which the shelf is min_draft thick.
137  real :: c1 ! The maximum depths in m.
138  character(len=40) :: mdl = "USER_initialize_shelf_mass" ! This subroutine's name.
139  integer :: i, j
140 
141  ! call MOM_error(FATAL, "USER_shelf_init.F90, USER_set_shelf_mass: " // &
142  ! "Unmodified user routine called - you must edit the routine to use it")
143 
144  if (.not.associated(cs)) allocate(cs)
145 
146  ! Read all relevant parameters and write them to the model log.
147  if (first_call) call write_user_log(param_file)
148  call get_param(param_file, mdl, "RHO_0", cs%Rho_ocean, &
149  "The mean ocean density used with BOUSSINESQ true to \n"//&
150  "calculate accelerations and the mass for conservation \n"//&
151  "properties, or with BOUSSINSEQ false to convert some \n"//&
152  "parameters from vertical units of m to kg m-2.", &
153  units="kg m-3", default=1035.0)
154  call get_param(param_file, mdl, "SHELF_MAX_DRAFT", cs%max_draft, &
155  units="m", default=1.0)
156  call get_param(param_file, mdl, "SHELF_MIN_DRAFT", cs%min_draft, &
157  units="m", default=1.0)
158  call get_param(param_file, mdl, "FLAT_SHELF_WIDTH", cs%flat_shelf_width, &
159  units="axis_units", default=0.0)
160  call get_param(param_file, mdl, "SHELF_SLOPE_SCALE", cs%shelf_slope_scale, &
161  units="axis_units", default=0.0)
162  call get_param(param_file, mdl, "SHELF_EDGE_POS_0", cs%pos_shelf_edge_0, &
163  units="axis_units", default=0.0)
164  call get_param(param_file, mdl, "SHELF_SPEED", cs%shelf_speed, &
165  units="axis_units day-1", default=0.0)
166 
167  call user_update_shelf_mass(mass_shelf, area_shelf_h, h_shelf, hmask, g, cs, set_time(0,0), new_sim)
168 
169 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the call graph for this function:
Here is the caller graph for this function:

◆ user_update_shelf_mass()

subroutine, public user_shelf_init::user_update_shelf_mass ( real, dimension(szi_(g),szj_(g)), intent(inout)  mass_shelf,
real, dimension(szi_(g),szj_(g)), intent(inout)  area_shelf_h,
real, dimension(szi_(g),szj_(g)), intent(inout)  h_shelf,
real, dimension(szi_(g),szj_(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 
)
Parameters
[in]gThe ocean's grid structure

Definition at line 187 of file user_shelf_init.F90.

Referenced by user_initialize_shelf_mass().

187  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
188  real, dimension(SZI_(G),SZJ_(G)), intent(inout) :: mass_shelf, area_shelf_h, hmask, h_shelf
189  type(user_ice_shelf_cs), pointer :: cs
190  type(time_type), intent(in) :: time
191  logical, intent(in) :: new_sim
192 
193 ! Arguments: mass_shelf - The mass per unit area averaged over the full ocean
194 ! cell, in kg m-2. (Intent out)
195 ! (out) area_shelf_h - The area of the ocean cell that is covered by the
196 ! rigid ice shelf, in m2.
197 ! (in) G - The ocean's grid structure.
198 ! (in) param_file - A structure indicating the open file to parse for
199 ! model parameter values.
200 
201  real :: c1, edge_pos, slope_pos
202  integer :: i, j
203 
204  edge_pos = cs%pos_shelf_edge_0 + cs%shelf_speed*(time_type_to_real(time) / 86400.0)
205 
206  slope_pos = edge_pos - cs%flat_shelf_width
207  c1 = 0.0 ; if (cs%shelf_slope_scale > 0.0) c1 = 1.0 / cs%shelf_slope_scale
208 
209 
210  do j=g%jsd,g%jed ;
211 
212  if (((j+g%jdg_offset) .le. g%domain%njglobal+g%domain%njhalo) .AND. &
213  ((j+g%jdg_offset) .ge. g%domain%njhalo+1)) then
214 
215  do i=g%isc,g%iec
216 
217 ! if (((i+G%idg_offset) <= G%domain%niglobal+G%domain%nihalo) .AND. &
218 ! ((i+G%idg_offset) >= G%domain%nihalo+1)) then
219 
220  if ((j.ge.g%jsc) .and. (j.le.g%jec)) then
221 
222  if (new_sim) then ; if (g%geoLonCu(i-1,j) >= edge_pos) then
223  ! Everything past the edge is open ocean.
224  mass_shelf(i,j) = 0.0
225  area_shelf_h(i,j) = 0.0
226  hmask(i,j) = 0.0
227  h_shelf(i,j) = 0.0
228  else
229  if (g%geoLonCu(i,j) > edge_pos) then
230  area_shelf_h(i,j) = g%areaT(i,j) * (edge_pos - g%geoLonCu(i-1,j)) / &
231  (g%geoLonCu(i,j) - g%geoLonCu(i-1,j))
232  hmask(i,j) = 2.0
233  else
234  area_shelf_h(i,j) = g%areaT(i,j)
235  hmask(i,j) = 1.0
236  endif
237 
238  if (g%geoLonT(i,j) > slope_pos) then
239  h_shelf(i,j) = cs%min_draft
240  mass_shelf(i,j) = cs%Rho_ocean * cs%min_draft
241  else
242  mass_shelf(i,j) = cs%Rho_ocean * (cs%min_draft + &
243  (cs%max_draft - cs%min_draft) * &
244  min(1.0, (c1*(slope_pos - g%geoLonT(i,j)))**2) )
245  h_shelf(i,j) = (cs%min_draft + &
246  (cs%max_draft - cs%min_draft) * &
247  min(1.0, (c1*(slope_pos - g%geoLonT(i,j)))**2) )
248  endif
249 
250  endif ; endif ; endif
251 
252  if ((i+g%idg_offset) .eq. g%domain%nihalo+1) then
253  hmask(i-1,j) = 3.0
254  endif
255 
256  enddo ; endif ; enddo
257 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the caller graph for this function:

◆ write_user_log()

subroutine user_shelf_init::write_user_log ( type(param_file_type), intent(in)  param_file)
private
Parameters
[in]param_fileA structure to parse for run-time parameters

Definition at line 261 of file user_shelf_init.F90.

References first_call.

Referenced by user_initialize_shelf_mass().

261  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
262 
263  character(len=128) :: version = '$Id: user_shelf_init.F90,v 1.1.2.7 2012/06/19 22:15:52 Robert.Hallberg Exp $'
264  character(len=128) :: tagname = '$Name: MOM_ogrp $'
265  character(len=40) :: mdl = "user_shelf_init" ! This module's name.
266 
267  call log_version(param_file, mdl, version, tagname)
268  first_call = .false.
269 
Here is the caller graph for this function:

Variable Documentation

◆ first_call

logical user_shelf_init::first_call = .true.
private

Definition at line 100 of file user_shelf_init.F90.

Referenced by user_initialize_shelf_mass(), and write_user_log().

100 logical :: first_call = .true.