MOM6
mom_ice_shelf_initialize Module Reference

Detailed Description

Initialize ice shelf variables.

Functions/Subroutines

subroutine, public initialize_ice_thickness (h_shelf, area_shelf_h, hmask, G, US, PF)
 Initialize ice shelf thickness. More...
 
subroutine initialize_ice_thickness_from_file (h_shelf, area_shelf_h, hmask, G, US, PF)
 Initialize ice shelf thickness from file. More...
 
subroutine initialize_ice_thickness_channel (h_shelf, area_shelf_h, hmask, G, US, PF)
 Initialize ice shelf thickness for a channel configuration. More...
 

Function/Subroutine Documentation

◆ initialize_ice_thickness()

subroutine, public mom_ice_shelf_initialize::initialize_ice_thickness ( real, dimension(szdi_(g),szdj_(g)), intent(inout)  h_shelf,
real, dimension(szdi_(g),szdj_(g)), intent(inout)  area_shelf_h,
real, dimension(szdi_(g),szdj_(g)), intent(inout)  hmask,
type(ocean_grid_type), intent(in)  G,
type(unit_scale_type), intent(in)  US,
type(param_file_type), intent(in)  PF 
)

Initialize ice shelf thickness.

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

Definition at line 29 of file MOM_ice_shelf_initialize.F90.

29  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
30  real, dimension(SZDI_(G),SZDJ_(G)), &
31  intent(inout) :: h_shelf !< The ice shelf thickness [Z ~> m].
32  real, dimension(SZDI_(G),SZDJ_(G)), &
33  intent(inout) :: area_shelf_h !< The area per cell covered by the ice shelf [m2].
34  real, dimension(SZDI_(G),SZDJ_(G)), &
35  intent(inout) :: hmask !< A mask indicating which tracer points are
36  !! partly or fully covered by an ice-shelf
37  type(unit_scale_type), intent(in) :: US !< A structure containing unit conversion factors
38  type(param_file_type), intent(in) :: PF !< A structure to parse for run-time parameters
39 
40  integer :: i, j
41  character(len=40) :: mdl = "initialize_ice_thickness" ! This subroutine's name.
42  character(len=200) :: config
43 
44  call get_param(pf, mdl, "ICE_PROFILE_CONFIG", config, &
45  "This specifies how the initial ice profile is specified. "//&
46  "Valid values are: CHANNEL, FILE, and USER.", &
47  fail_if_missing=.true.)
48 
49  select case ( trim(config) )
50  case ("CHANNEL"); call initialize_ice_thickness_channel (h_shelf, area_shelf_h, hmask, g, us, pf)
51  case ("FILE"); call initialize_ice_thickness_from_file (h_shelf, area_shelf_h, hmask, g, us, pf)
52  case ("USER"); call user_init_ice_thickness (h_shelf, area_shelf_h, hmask, g, us, pf)
53  case default ; call mom_error(fatal,"MOM_initialize: "// &
54  "Unrecognized ice profile setup "//trim(config))
55  end select
56 

References initialize_ice_thickness_channel(), initialize_ice_thickness_from_file(), mom_error_handler::mom_error(), and user_shelf_init::user_init_ice_thickness().

Referenced by mom_ice_shelf::initialize_ice_shelf().

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

◆ initialize_ice_thickness_channel()

subroutine mom_ice_shelf_initialize::initialize_ice_thickness_channel ( real, dimension(szdi_(g),szdj_(g)), intent(inout)  h_shelf,
real, dimension(szdi_(g),szdj_(g)), intent(inout)  area_shelf_h,
real, dimension(szdi_(g),szdj_(g)), intent(inout)  hmask,
type(ocean_grid_type), intent(in)  G,
type(unit_scale_type), intent(in)  US,
type(param_file_type), intent(in)  PF 
)
private

Initialize ice shelf thickness for a channel configuration.

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

Definition at line 148 of file MOM_ice_shelf_initialize.F90.

148  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
149  real, dimension(SZDI_(G),SZDJ_(G)), &
150  intent(inout) :: h_shelf !< The ice shelf thickness [Z ~> m].
151  real, dimension(SZDI_(G),SZDJ_(G)), &
152  intent(inout) :: area_shelf_h !< The area per cell covered by the ice shelf [m2].
153  real, dimension(SZDI_(G),SZDJ_(G)), &
154  intent(inout) :: hmask !< A mask indicating which tracer points are
155  !! partly or fully covered by an ice-shelf
156  type(unit_scale_type), intent(in) :: US !< A structure containing unit conversion factors
157  type(param_file_type), intent(in) :: PF !< A structure to parse for run-time parameters
158 
159  character(len=40) :: mdl = "initialize_ice_shelf_thickness_channel" ! This subroutine's name.
160  real :: max_draft, min_draft, flat_shelf_width, c1, slope_pos
161  real :: edge_pos, shelf_slope_scale, Rho_ocean
162  integer :: i, j, jsc, jec, jsd, jed, jedg, nyh, isc, iec, isd, ied
163  integer :: j_off
164 
165  jsc = g%jsc ; jec = g%jec ; isc = g%isc ; iec = g%iec
166  jsd = g%jsd ; jed = g%jed ; isd = g%isd ; ied = g%ied
167  nyh = g%domain%njhalo ; jedg = g%domain%njglobal+nyh
168  j_off = g%jdg_offset
169 
170  call mom_mesg(mdl//": setting thickness")
171 
172  call get_param(pf, mdl, "SHELF_MAX_DRAFT", max_draft, &
173  units="m", default=1.0, scale=us%m_to_Z)
174  call get_param(pf, mdl, "SHELF_MIN_DRAFT", min_draft, &
175  units="m", default=1.0, scale=us%m_to_Z)
176  call get_param(pf, mdl, "FLAT_SHELF_WIDTH", flat_shelf_width, &
177  units="axis_units", default=0.0)
178  call get_param(pf, mdl, "SHELF_SLOPE_SCALE", shelf_slope_scale, &
179  units="axis_units", default=0.0)
180  call get_param(pf, mdl, "SHELF_EDGE_POS_0", edge_pos, &
181  units="axis_units", default=0.0)
182 ! call get_param(param_file, mdl, "RHO_0", Rho_ocean, &
183 ! "The mean ocean density used with BOUSSINESQ true to "//&
184 ! "calculate accelerations and the mass for conservation "//&
185 ! "properties, or with BOUSSINSEQ false to convert some "//&
186 ! "parameters from vertical units of m to kg m-2.", &
187 ! units="kg m-3", default=1035.0, scale=US%Z_to_m)
188 
189  slope_pos = edge_pos - flat_shelf_width
190  c1 = 0.0 ; if (shelf_slope_scale > 0.0) c1 = 1.0 / shelf_slope_scale
191 
192 
193  do j=g%jsd,g%jed
194 
195  if (((j+j_off) <= jedg) .AND. ((j+j_off) >= nyh+1)) then
196 
197  do i=g%isc,g%iec
198 
199  if ((j >= jsc) .and. (j <= jec)) then
200 
201  if (g%geoLonCu(i-1,j) >= edge_pos) then
202  ! Everything past the edge is open ocean.
203 ! mass_shelf(i,j) = 0.0
204  area_shelf_h(i,j) = 0.0
205  hmask(i,j) = 0.0
206  h_shelf(i,j) = 0.0
207  else
208  if (g%geoLonCu(i,j) > edge_pos) then
209  area_shelf_h(i,j) = us%L_to_m**2*g%areaT(i,j) * (edge_pos - g%geoLonCu(i-1,j)) / &
210  (g%geoLonCu(i,j) - g%geoLonCu(i-1,j))
211  hmask(i,j) = 2.0
212  else
213  area_shelf_h(i,j) = us%L_to_m**2*g%areaT(i,j)
214  hmask(i,j) = 1.0
215  endif
216 
217  if (g%geoLonT(i,j) > slope_pos) then
218  h_shelf(i,j) = min_draft
219 ! mass_shelf(i,j) = Rho_ocean * min_draft
220  else
221 ! mass_shelf(i,j) = Rho_ocean * (min_draft + &
222 ! (CS%max_draft - CS%min_draft) * &
223 ! min(1.0, (c1*(slope_pos - G%geoLonT(i,j)))**2) )
224  h_shelf(i,j) = (min_draft + &
225  (max_draft - min_draft) * &
226  min(1.0, (c1*(slope_pos - g%geoLonT(i,j)))**2) )
227  endif
228 
229  endif
230  endif
231 
232  if ((i+g%idg_offset) == g%domain%nihalo+1) then
233  hmask(i-1,j) = 3.0
234  endif
235 
236  enddo
237  endif ; enddo
238 

References mom_error_handler::mom_mesg().

Referenced by initialize_ice_thickness().

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

◆ initialize_ice_thickness_from_file()

subroutine mom_ice_shelf_initialize::initialize_ice_thickness_from_file ( real, dimension(szdi_(g),szdj_(g)), intent(inout)  h_shelf,
real, dimension(szdi_(g),szdj_(g)), intent(inout)  area_shelf_h,
real, dimension(szdi_(g),szdj_(g)), intent(inout)  hmask,
type(ocean_grid_type), intent(in)  G,
type(unit_scale_type), intent(in)  US,
type(param_file_type), intent(in)  PF 
)
private

Initialize ice shelf thickness from file.

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

Definition at line 61 of file MOM_ice_shelf_initialize.F90.

61  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
62  real, dimension(SZDI_(G),SZDJ_(G)), &
63  intent(inout) :: h_shelf !< The ice shelf thickness [m].
64  real, dimension(SZDI_(G),SZDJ_(G)), &
65  intent(inout) :: area_shelf_h !< The area per cell covered by the ice shelf [m2].
66  real, dimension(SZDI_(G),SZDJ_(G)), &
67  intent(inout) :: hmask !< A mask indicating which tracer points are
68  !! partly or fully covered by an ice-shelf
69  type(unit_scale_type), intent(in) :: US !< A structure containing unit conversion factors
70  type(param_file_type), intent(in) :: PF !< A structure to parse for run-time parameters
71 
72  ! This subroutine reads ice thickness and area from a file and puts it into
73  ! h_shelf and area_shelf_h in m (and dimensionless) and updates hmask
74  character(len=200) :: filename,thickness_file,inputdir ! Strings for file/path
75  character(len=200) :: thickness_varname, area_varname ! Variable name in file
76  character(len=40) :: mdl = "initialize_ice_thickness_from_file" ! This subroutine's name.
77  integer :: i, j, isc, jsc, iec, jec
78  real :: len_sidestress, mask, udh
79 
80  call mom_mesg(" MOM_ice_shelf_init_profile.F90, initialize_thickness_from_file: reading thickness")
81 
82  call get_param(pf, mdl, "INPUTDIR", inputdir, default=".")
83  inputdir = slasher(inputdir)
84  call get_param(pf, mdl, "ICE_THICKNESS_FILE", thickness_file, &
85  "The file from which the bathymetry is read.", &
86  default="ice_shelf_h.nc")
87  call get_param(pf, mdl, "LEN_SIDE_STRESS", len_sidestress, &
88  "position past which shelf sides are stress free.", &
89  default=0.0, units="axis_units")
90 
91  filename = trim(inputdir)//trim(thickness_file)
92  call log_param(pf, mdl, "INPUTDIR/THICKNESS_FILE", filename)
93  call get_param(pf, mdl, "ICE_THICKNESS_VARNAME", thickness_varname, &
94  "The name of the thickness variable in ICE_THICKNESS_FILE.", &
95  default="h_shelf")
96  call get_param(pf, mdl, "ICE_AREA_VARNAME", area_varname, &
97  "The name of the area variable in ICE_THICKNESS_FILE.", &
98  default="area_shelf_h")
99 
100  if (.not.file_exists(filename, g%Domain)) call mom_error(fatal, &
101  " initialize_topography_from_file: Unable to open "//trim(filename))
102 
103  call mom_read_data(filename, trim(thickness_varname), h_shelf, g%Domain, scale=us%m_to_Z)
104  call mom_read_data(filename,trim(area_varname),area_shelf_h,g%Domain)
105 
106 ! call get_param(PF, mdl, "ICE_BOUNDARY_CONFIG", config, &
107 ! "This specifies how the ice domain boundary is specified", &
108 ! fail_if_missing=.true.)
109 
110  isc = g%isc ; jsc = g%jsc ; iec = g%iec ; jec = g%jec
111 
112  do j=jsc,jec
113  do i=isc,iec
114 
115  ! taper ice shelf in area where there is no sidestress -
116  ! but do not interfere with hmask
117 
118  if ((g%geoLonCv(i,j) > len_sidestress).and. &
119  (len_sidestress > 0.)) then
120  udh = exp(-(g%geoLonCv(i,j)-len_sidestress)/5.0) * h_shelf(i,j)
121  if (udh <= 25.0) then
122  h_shelf(i,j) = 0.0
123  area_shelf_h(i,j) = 0.0
124  else
125  h_shelf(i,j) = udh
126  endif
127  endif
128 
129  ! update thickness mask
130 
131  if (area_shelf_h(i,j) >= us%L_to_m**2*g%areaT(i,j)) then
132  hmask(i,j) = 1.
133  elseif (area_shelf_h(i,j) == 0.0) then
134  hmask(i,j) = 0.
135  elseif ((area_shelf_h(i,j) > 0) .and. (area_shelf_h(i,j) <= us%L_to_m**2*g%areaT(i,j))) then
136  hmask(i,j) = 2.
137  else
138  call mom_error(fatal,mdl// " AREA IN CELL OUT OF RANGE")
139  endif
140  enddo
141  enddo
142 
143 

References mom_error_handler::mom_error(), and mom_error_handler::mom_mesg().

Referenced by initialize_ice_thickness().

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