seasonal_glacier_boundary Derived Type

type, public, extends(glacier_boundary) :: seasonal_glacier_boundary

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.


Inherits

type~~seasonal_glacier_boundary~~InheritsGraph type~seasonal_glacier_boundary seasonal_glacier_boundary type~glacier_boundary glacier_boundary type~seasonal_glacier_boundary->type~glacier_boundary

Contents


Components

TypeVisibility AttributesNameInitial
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


Constructor

public interface seasonal_glacier_boundary

  • private pure function constructor(thickness, frequency, amplitude, mean, chi, square) result(this)

    Author
    Chris MacMackin
    Date
    October 2017

    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.

    Arguments

    Type IntentOptional AttributesName
    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.

    Return Value type(seasonal_glacier_boundary)


Type-Bound Procedures

procedure, public :: thickness_upper_bound => bound_array

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.

  • private pure function bound_array(this)

    Author
    Chris MacMackin
    Date
    September 2016

    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).

    Arguments

    Type IntentOptional AttributesName
    class(glacier_boundary), intent(in) :: this

    Return Value integer, dimension(2)

procedure, public :: thickness_upper_type => bound_type

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].

  • private pure function bound_type(this)

    Author
    Chris MacMackin
    Date
    January 2017

    Default implementation of the methods getting the boundary types for a glacier. It returns an array which indicates free boundaries.

    Arguments

    Type IntentOptional AttributesName
    class(glacier_boundary), intent(in) :: this

    Return Value integer, dimension(:), allocatable

procedure, public :: 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.

  • private pure function seasonal_lower_bound(this) result(bound_array)

    Author
    Chris MacMackin
    Date
    October 2017

    Indicates that one layer of cells at the lower boundary in the first dimension should be omitted. This is appropriate for Dirichlet boundary conditions.

    Arguments

    Type IntentOptional AttributesName
    class(seasonal_glacier_boundary), intent(in) :: this

    Return Value integer, dimension(2)

procedure, public :: 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.

  • private pure function seasonal_lower_bound(this) result(bound_array)

    Author
    Chris MacMackin
    Date
    October 2017

    Indicates that one layer of cells at the lower boundary in the first dimension should be omitted. This is appropriate for Dirichlet boundary conditions.

    Arguments

    Type IntentOptional AttributesName
    class(seasonal_glacier_boundary), intent(in) :: this

    Return Value integer, dimension(2)

procedure, public :: 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.

  • private pure function seasonal_upper_bound(this) result(bound_array)

    Author
    Chris MacMackin
    Date
    October 2017

    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.

    Arguments

    Type IntentOptional AttributesName
    class(seasonal_glacier_boundary), intent(in) :: this

    Return Value integer, dimension(2)

procedure, public :: 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].

  • private pure function seasonal_lower_type(this) result(bound_type)

    Author
    Chris MacMackin
    Date
    January 2017

    Specifies that the lower boundary in the first dimension has Dirichlet boundary conditions.

    Arguments

    Type IntentOptional AttributesName
    class(seasonal_glacier_boundary), intent(in) :: this

    Return Value integer, dimension(:), allocatable

procedure, public :: 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].

  • private pure function seasonal_lower_type(this) result(bound_type)

    Author
    Chris MacMackin
    Date
    January 2017

    Specifies that the lower boundary in the first dimension has Dirichlet boundary conditions.

    Arguments

    Type IntentOptional AttributesName
    class(seasonal_glacier_boundary), intent(in) :: this

    Return Value integer, dimension(:), allocatable

procedure, public :: 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].

  • private pure function seasonal_upper_type(this) result(bound_type)

    Author
    Chris MacMackin
    Date
    January 2017

    Specifies that the upper boundary in the first dimension has Neumann boundary conditions.

    Arguments

    Type IntentOptional AttributesName
    class(seasonal_glacier_boundary), intent(in) :: this

    Return Value integer, dimension(:), allocatable

procedure, public :: 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.

  • private function seasonal_residuals(this, thickness, velocity, viscosity, t) result(residuals)

    Author
    Chris MacMackin
    Date
    October 2017

    Returns the difference between the glacier conditions of the plume and the Dirichlet conditions prescribed in the model of Dallaston et al. (2015)

    Arguments

    Type IntentOptional AttributesName
    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.

    Return Value real(kind=r8), allocatable, dimension(:)

    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.

Source Code

  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