glacier Derived Type

type, public, abstract :: glacier

An abstract data type which represents large masses of ice, such as ice shelves and ice sheets.


Inherited by

type~~glacier~~InheritedByGraph type~glacier glacier type~ice_sheet ice_sheet type~ice_sheet->type~glacier type~cryosphere cryosphere type~cryosphere->type~glacier ice type~ice_shelf ice_shelf type~ice_shelf->type~glacier

Contents

Source Code


Type-Bound Procedures

procedure(get_scalar), public, deferred :: ice_thickness

Returns the thickness of the ice in the glacier across the domain.

  • function get_scalar(this) result(property) Prototype

    Arguments

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

    Return Value class(scalar_field), pointer

    The value of whatever property of the glacier is being returned.

procedure(get_r8), public, deferred :: ice_density

Returns the density of the ice, which is assumed to be uniform.

  • pure function get_r8(this) result(property) Prototype

    Arguments

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

    Return Value real(kind=r8)

    The value of whatever property of the glacier is being returned.

procedure(get_r8), public, deferred :: ice_temperature

Returns the temperature of the ice, which is assumed to be uniform.

  • pure function get_r8(this) result(property) Prototype

    Arguments

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

    Return Value real(kind=r8)

    The value of whatever property of the glacier is being returned.

procedure(get_residual), public, deferred :: residual

Computes the residual of the system of equations describing the glacier.

  • function get_residual(this, previous_states, melt_rate, basal_drag_parameter, water_density) result(residual) Prototype

    Arguments

    Type IntentOptional AttributesName
    class(glacier), intent(in) :: this
    class(glacier), intent(in), dimension(:):: previous_states

    The states of the glacier in the previous time steps. The first element of the array should be the most recent. The default implementation will only make use of the most recent state, but the fact that this is an array allows overriding methods to use older states for higher-order integration methods.

    class(scalar_field), intent(in) :: melt_rate

    Thickness of the ice above the glacier

    class(scalar_field), intent(in) :: basal_drag_parameter

    A paramter, e.g. coefficient of friction, needed to calculate the drag on basal surface of the glacier.

    real(kind=r8), intent(in) :: water_density

    The density of the water below the glacier

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

    The residual of the system of equations describing the glacier

procedure(precond), public, deferred :: precondition

Applies a preconditioner to the passed state vector.

  • function precond(this, previous_states, melt_rate, basal_drag_parameter, water_density, delta_state) result(preconditioned) Prototype

    Arguments

    Type IntentOptional AttributesName
    class(glacier), intent(inout) :: this
    class(glacier), intent(in), dimension(:):: previous_states

    The states of the glacier in the previous time steps. The first element of the array should be the most recent. The default implementation will only make use of the most recent state, but the fact that this is an array allows overriding methods to use older states for higher-order integration methods.

    class(scalar_field), intent(in) :: melt_rate

    Thickness of the ice above the glacier

    class(scalar_field), intent(in) :: basal_drag_parameter

    A paramter, e.g. coefficient of friction, needed to calculate the drag on basal surface of the glacier.

    real(kind=r8), intent(in) :: water_density

    The density of the water below the glacier

    real(kind=r8), intent(in), dimension(:):: delta_state

    The change to the state vector which is being preconditioned.

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

    The result of applying the preconditioner to delta_state.

procedure(solve_vel), public, deferred :: solve_velocity

Solves for the velocity field using the current thickness.

  • subroutine solve_vel(this, basal_drag, success) Prototype

    Arguments

    Type IntentOptional AttributesName
    class(glacier), intent(inout) :: this
    class(scalar_field), intent(in) :: basal_drag

    A paramter, e.g. coefficient of friction, needed to calculate the drag on basal surface of the glacier.

    logical, intent(out) :: success

    True if the integration is successful, false otherwise

procedure(setter), public, deferred :: update

Sets the state of the glacier.

  • subroutine setter(this, state_vector) Prototype

    Arguments

    Type IntentOptional AttributesName
    class(glacier), intent(inout) :: this
    real(kind=r8), intent(in), dimension(:):: state_vector

    A real array containing the data describing the state of the glacier.

procedure(time_setter), public, deferred :: set_time

Sets the time record for this glacier.

  • subroutine time_setter(this, time) Prototype

    Arguments

    Type IntentOptional AttributesName
    class(glacier), intent(inout) :: this
    real(kind=r8), intent(in) :: time

    The time at which the glacier is in the present state.

procedure(get_i), public, deferred :: data_size

Returns the number of elements in the glacier's state vector

  • pure function get_i(this) result(property) Prototype

    Arguments

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

    Return Value integer

    The value of whatever property of the glacier is being returned.

procedure(get_r81d), public, deferred :: state_vector

Returns the glacier's state vector, a 1D array with all necessary data to describe its state.

  • function get_r81d(this) result(state_vector) Prototype

    Arguments

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

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

    The state vector of the glacier

procedure(read_dat), public, deferred :: read_data

Read the glacier data from an HDF5 file on the disc.

  • subroutine read_dat(this, file_id, group_name, error) Prototype

    Arguments

    Type IntentOptional AttributesName
    class(glacier), intent(inout) :: this
    integer(kind=hid_t), intent(in) :: file_id

    The identifier for the HDF5 file/group from which the data will be read.

    character(len=*), intent(in) :: group_name

    The name of the group in the HDF5 file from which to read glacier's data.

    integer, intent(out) :: error

    Flag indicating whether routine ran without error. If no error occurs then has value 0.

procedure(write_dat), public, deferred :: write_data

Writes the data describing the glacier to the disc as an HDF5 file.

  • subroutine write_dat(this, file_id, group_name, error) Prototype

    Arguments

    Type IntentOptional AttributesName
    class(glacier), intent(in) :: this
    integer(kind=hid_t), intent(in) :: file_id

    The identifier for the HDF5 file/group in which this data is meant to be written.

    character(len=*), intent(in) :: group_name

    The name to give the group in the HDF5 file storing the glacier's data.

    integer, intent(out) :: error

    Flag indicating whether routine ran without error. If no error occurs then has value 0.

procedure(t_step), public, deferred :: time_step

Calculates the appropriate time step for integration.

  • function t_step(this) Prototype

    Arguments

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

    Return Value real(kind=r8)

    A time step which will allow integration of the ice shelf without causing numerical instability.

procedure(assign_ice), private, deferred :: assign

Copies the data from one glacier into another. This is only needed due to a bug in gfortran which means that the intrinsic assignment for glacier types is not using the appropriate defined assignment for the field components.

  • subroutine assign_ice(this, rhs) Prototype

    Arguments

    Type IntentOptional AttributesName
    class(glacier), intent(out) :: this
    class(glacier), intent(in) :: rhs

generic, public :: assignment(=) => assign

  • private subroutine shelf_assign(this, rhs)

    Author
    Chris MacMackin
    Date
    February 2017

    Copies the data from one ice shelf into another. This is only needed due to a bug in gfortran which means that the intrinsic assignment for glacier types is not using the appropriate defined assignment for the field components.

    Read more…

    Arguments

    Type IntentOptional AttributesName
    class(ice_shelf), intent(out) :: this
    class(glacier), intent(in) :: rhs

    The ice shelf to be assigned to this one.

procedure, public :: integrate => glacier_integrate

Performs a time-step of the integration, taking the state of the glacier to the specified time using the provided melt-rate data.

  • private subroutine glacier_integrate(this, old_states, basal_melt, basal_drag, water_density, time, success)

    Author
    Chris MacMackin
    Date
    November 2016

    Integrates the glacier's state to time. This is done using the NITSOL package of iterative Krylov solvers. If a different algorithm for the integration is desired, then this method may be overridden in the concrete implementations of the glacier type.

    Read more…

    Arguments

    Type IntentOptional AttributesName
    class(glacier), intent(inout) :: this
    class(glacier), intent(in), dimension(:):: old_states

    Previous states of the glacier, with the most recent one first.

    class(scalar_field), intent(in) :: basal_melt

    The melt rate that the bottom of the glacier experiences during this time step.

    class(scalar_field), intent(in) :: basal_drag

    A paramter, e.g. coefficient of friction, needed to calculate the drag on basal surface of the glacier.

    real(kind=r8), intent(in) :: water_density

    The density of the water below the glacier.

    real(kind=r8), intent(in) :: time

    The time to which the glacier should be integrated

    logical, intent(out) :: success

    True if the integration is successful, false otherwise

procedure, public :: integrate_layers => glacier_integrate_layers

Dummy routine which can be over-ridden to integrate internal layers of the glacier to the specified time.

  • private subroutine glacier_integrate_layers(this, old_states, time, success)

    Author
    Chris MacMackin
    Date
    September 2018

    Dummy routine which does nothing.

    Arguments

    Type IntentOptional AttributesName
    class(glacier), intent(inout) :: this
    class(glacier), intent(in), dimension(:):: old_states

    Previous states of the glacier, with the most recent one first.

    real(kind=r8), intent(in) :: time

    The time to which the glacier should be integrated

    logical, intent(out) :: success

    True if the integration is successful, false otherwise

Source Code

  type, abstract, public :: glacier
    !* Author: Christopehr MacMackin
    !  Date: April 2016
    !
    ! An abstract data type which represents large masses of ice, such
    ! as ice shelves and ice sheets.
    !
  contains
    procedure(get_scalar), deferred   :: ice_thickness
      !! Returns the thickness of the ice in the glacier across the domain.
    procedure(get_r8), deferred       :: ice_density
      !! Returns the density of the ice, which is assumed to be uniform.
    procedure(get_r8), deferred       :: ice_temperature
      !! Returns the temperature of the ice, which is assumed to be uniform.
    procedure(get_residual), deferred :: residual
      !! Computes the residual of the system of equations describing the
      !! glacier.
    procedure(precond), deferred      :: precondition
      !! Applies a preconditioner to the passed state vector.
    procedure(solve_vel), deferred    :: solve_velocity
      !! Solves for the velocity field using the current thickness.
    procedure(setter), deferred       :: update
      !! Sets the state of the glacier.
    procedure(time_setter), deferred  :: set_time
      !! Sets the time record for this glacier.
    procedure(get_i), deferred        :: data_size
      !! Returns the number of elements in the glacier's state vector
    procedure(get_r81d), deferred     :: state_vector
      !! Returns the glacier's state vector, a 1D array with all necessary 
      !! data to describe its state.
    procedure(read_dat), deferred     :: read_data
      !! Read the glacier data from an HDF5 file on the disc.
    procedure(write_dat), deferred    :: write_data
      !! Writes the data describing the glacier to the disc as an HDF5 file.
    procedure(t_step), deferred       :: time_step
      !! Calculates the appropriate time step for integration.
    procedure(assign_ice), private, deferred :: assign
      !! Copies the data from one glacier into another. This is only
      !! needed due to a bug in gfortran which means that the
      !! intrinsic assignment for glacier types is not using the
      !! appropriate defined assignment for the field components.
    generic                           :: assignment(=) => assign
    procedure                         :: integrate => glacier_integrate
      !! Performs a time-step of the integration, taking the state of
      !! the glacier to the specified time using the provided
      !! melt-rate data.
    procedure                         :: integrate_layers => glacier_integrate_layers
      !! Dummy routine which can be over-ridden to integrate internal
      !! layers of the glacier to the specified time.
  end type glacier