A type with procedures for getting the boundary conditions of the ice shelf model used by Dallaston et al. (2015), but with the ice flux at the grounding line varying sinusoidally in time. There are Dirichlet conditions at the lower bound of the first condition as well as, for thickness, the upper bound.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
real(kind=r8), | private | :: | thickness | = | 1.0_r8 | The thickness of the glacier at the inflowing boundary |
|
real(kind=r8), | private | :: | frequency | = | 1.0_r8 | The angular frequency of the oscillations in ice flux |
|
real(kind=r8), | private | :: | amplitude | = | 0.5_r8 | The amplitude of the oscillations in ice flux |
|
real(kind=r8), | private | :: | mean | = | 1.0_r8 | The time-average of the ice flux, about which it oscillates |
|
real(kind=r8), | private | :: | chi | = | 1.0_r8 | The dimensionless ratio $\chi \equiv \frac{\rho_igh_0x_x}{2\eta_0u_0}$ |
|
logical, | private | :: | square | = | .false. | If true, produce a square wave, otherwise produce a sinusoid |
Constructs a boundary condition object for an ice shelf based on the conditions used in Dallaston et al. (2015), but with the ice flux at the grounding line varying in time.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=r8), | intent(in), | optional | :: | thickness | The ice thickness at the inflowing ice shelf boundary, defaults to 1.0 |
|
real(kind=r8), | intent(in), | optional | :: | frequency | The angular frequency of the oscillations in ice flux, defaults to 0.5 |
|
real(kind=r8), | intent(in), | optional | :: | amplitude | The amplitude of the oscillations in ice flux, defaults to 1.0 |
|
real(kind=r8), | intent(in), | optional | :: | mean | The time-average of the discharge, about which it oscillates, defaulting to 1.0 |
|
real(kind=r8), | intent(in), | optional | :: | chi | The dimensionless ratio , defaults to 1.0 |
|
logical, | intent(in), | optional | :: | square | If present and true, produce a square wave. Otherwise produce a sinusoid. |
Returns a 1D array which should be passed as the
exclude_upper_bound
/provide_upper_bound
argument when
getting or setting the raw representation of the thickness
field.
Default implementation of the method getting lower and upper boundary information, which is then passed to the methods for getting and setting raw representations of fields. It returns a 1D array of length 2, indicating free boundaries (the raw data should represent all cells contained in the field, not excluding any near the boundaries).
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(glacier_boundary), | intent(in) | :: | this |
Returns an array indicating what type of boundary conditions apply for thickness at the upper boundary of each dimension. The types are specified using the parameters in [boundary_types_mod].
Default implementation of the methods getting the boundary types for a glacier. It returns an array which indicates free boundaries.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(glacier_boundary), | intent(in) | :: | this |
Returns a 1D array which should be passed as the
exclude_lower_bound
/provide_lower_bound
argument when
getting or setting the raw representation of the thickness
field.
Indicates that one layer of cells at the lower boundary in the first dimension should be omitted. This is appropriate for Dirichlet boundary conditions.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(seasonal_glacier_boundary), | intent(in) | :: | this |
Returns a 1D array which should be passed as the
exclude_lower_bound
/provide_lower_bound
argument when
getting or setting the raw representation of the velocity
field.
Indicates that one layer of cells at the lower boundary in the first dimension should be omitted. This is appropriate for Dirichlet boundary conditions.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(seasonal_glacier_boundary), | intent(in) | :: | this |
Returns a 1D array which should be passed as the
exclude_upper_bound
/provide_upper_bound
argument when
getting or setting the raw representation of the velocity
field.
Indicates that one layer of cells at the upper boundary in the first dimension should be omitted for thickness. This is appropriate for Dirichlet boundary conditions.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(seasonal_glacier_boundary), | intent(in) | :: | this |
Returns an array indicating what type of boundary conditions apply for thickness at the lower boundary of each dimension. The types are specified using the parameters in [boundary_types_mod].
Specifies that the lower boundary in the first dimension has Dirichlet boundary conditions.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(seasonal_glacier_boundary), | intent(in) | :: | this |
Returns an array indicating what type of boundary conditions apply for velocity at the lower boundary of each dimension. The types are specified using the parameters in [boundary_types_mod].
Specifies that the lower boundary in the first dimension has Dirichlet boundary conditions.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(seasonal_glacier_boundary), | intent(in) | :: | this |
Returns an array indicating what type of boundary conditions apply for velocity at the upper boundary of each dimension. The types are specified using the parameters in [boundary_types_mod].
Specifies that the upper boundary in the first dimension has Neumann boundary conditions.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(seasonal_glacier_boundary), | intent(in) | :: | this |
Returns an array consisting of the difference between the required boundary values and those which actually exist. This can then be appended to a glacier's state vector. The order in which these are listed is as follows: lower thickness boundary, upper thickness boundary, lower velocity boundary, and upper velocity boundary.
Returns the difference between the glacier conditions of the plume and the Dirichlet conditions prescribed in the model of Dallaston et al. (2015)
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(seasonal_glacier_boundary), | intent(in) | :: | this | |||
class(scalar_field), | intent(in) | :: | thickness | A field containing the thickness of the glacier |
||
class(vector_field), | intent(in) | :: | velocity | A field containing the flow velocity of the glacier |
||
class(scalar_field), | intent(in) | :: | viscosity | A field containing the viscosity of the ice in the glacier. |
||
real(kind=r8), | intent(in) | :: | t | The time at which the boundary conditions are to be calculated. |
An array containing the difference between the required boundary values and those which are actually present. They are stored in the order: lower thickness boundary, upper thickness boundary, lower velocity boundary, and upper velocity boundary.
type, extends(glacier_boundary), public :: seasonal_glacier_boundary
!* Author: Chris MacMackin
! Date: October 2017
!
! A type with procedures for getting the boundary conditions of
! the ice shelf model used by Dallaston et al. (2015), but with
! the ice flux at the grounding line varying sinusoidally in
! time. There are Dirichlet conditions at the lower bound of the
! first condition as well as, for thickness, the upper bound.
!
private
real(r8) :: thickness = 1.0_r8
!! The thickness of the glacier at the inflowing boundary
real(r8) :: frequency = 1.0_r8
!! The angular frequency of the oscillations in ice flux
real(r8) :: amplitude = 0.5_r8
!! The amplitude of the oscillations in ice flux
real(r8) :: mean = 1.0_r8
!! The time-average of the ice flux, about which it oscillates
real(r8) :: chi = 1.0_r8
!! The dimensionless ratio
!! $\chi \equiv \frac{\rho_igh_0x_x}{2\eta_0u_0}$
logical :: square = .false.
!! If true, produce a square wave, otherwise produce a sinusoid
contains
procedure :: thickness_lower_bound => seasonal_lower_bound
!! Returns a 1D array which should be passed as the
!! `exclude_lower_bound`/`provide_lower_bound` argument when
!! getting or setting the raw representation of the thickness
!! field.
procedure :: velocity_lower_bound => seasonal_lower_bound
!! Returns a 1D array which should be passed as the
!! `exclude_lower_bound`/`provide_lower_bound` argument when
!! getting or setting the raw representation of the velocity
!! field.
procedure :: velocity_upper_bound => seasonal_upper_bound
!! Returns a 1D array which should be passed as the
!! `exclude_upper_bound`/`provide_upper_bound` argument when
!! getting or setting the raw representation of the velocity
!! field.
procedure :: thickness_lower_type => seasonal_lower_type
!! Returns an array indicating what type of boundary conditions
!! apply for thickness at the lower boundary of each
!! dimension. The types are specified using the parameters in
!! [boundary_types_mod].
procedure :: velocity_lower_type => seasonal_lower_type
!! Returns an array indicating what type of boundary conditions
!! apply for velocity at the lower boundary of each
!! dimension. The types are specified using the parameters in
!! [boundary_types_mod].
procedure :: velocity_upper_type => seasonal_upper_type
!! Returns an array indicating what type of boundary conditions
!! apply for velocity at the upper boundary of each
!! dimension. The types are specified using the parameters in
!! [boundary_types_mod].
procedure :: boundary_residuals => seasonal_residuals
!! Returns an array consisting of the difference between the
!! required boundary values and those which actually exist. This
!! can then be appended to a glacier's state vector. The order
!! in which these are listed is as follows: lower thickness
!! boundary, upper thickness boundary, lower velocity boundary,
!! and upper velocity boundary.
end type seasonal_glacier_boundary