glacier_mod Module

Provides an abstract type to represent large masses of ice, such as ice sheets and ice shelves.


Uses

  • module~~glacier_mod~~UsesGraph module~glacier_mod glacier_mod logger_mod logger_mod module~glacier_mod->logger_mod factual_mod factual_mod module~glacier_mod->factual_mod hdf5 hdf5 module~glacier_mod->hdf5 iso_fortran_env iso_fortran_env module~glacier_mod->iso_fortran_env penf penf module~glacier_mod->penf module~nitsol_mod nitsol_mod module~glacier_mod->module~nitsol_mod module~nitsol_mod->iso_fortran_env

Used by

  • module~~glacier_mod~~UsedByGraph module~glacier_mod glacier_mod module~ice_shelf_mod ice_shelf_mod module~ice_shelf_mod->module~glacier_mod module~cryosphere_mod cryosphere_mod module~cryosphere_mod->module~glacier_mod module~ice_sheet_mod ice_sheet_mod module~ice_sheet_mod->module~glacier_mod

Contents


Variables

TypeVisibility AttributesNameInitial
character(len=12), public, parameter:: hdf_type_attr ='glacier_type'

Abstract Interfaces

abstract interface

  • private function get_scalar(this) result(property)

    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.

abstract interface

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

    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.

abstract interface

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

    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

abstract interface

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

    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.

abstract interface

  • private function get_r81d(this) result(state_vector)

    Arguments

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

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

    The state vector of the glacier

abstract interface

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

    Arguments

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

    Return Value integer

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

abstract interface

  • private function t_step(this)

    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.

abstract interface

  • private subroutine solve_vel(this, basal_drag, success)

    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

abstract interface

  • private subroutine setter(this, state_vector)

    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.

abstract interface

  • private subroutine time_setter(this, time)

    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.

abstract interface

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

    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.

abstract interface

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

    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.

abstract interface

  • private subroutine assign_ice(this, rhs)

    Arguments

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

abstract interface

  • public pure function thickness_func(location) result(thickness)

    Abstract interface for function providing the glacier thickness when a concrete object is being instantiated.

    Arguments

    Type IntentOptional AttributesName
    real(kind=r8), intent(in), dimension(:):: location

    The position $\vec{x}$ at which to compute the thickness

    Return Value real(kind=r8)

    The thickness of the glacier at location

abstract interface

  • public pure function velocity_func(location) result(velocity)

    Abstract interface for function providing the glacier velocity when a concrete object is being instantiated.

    Arguments

    Type IntentOptional AttributesName
    real(kind=r8), intent(in), dimension(:):: location

    The position $\vec{x}$ at which to compute the thickness

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

    The velocity vector of the ice in the glacier at location


Derived Types

type, public, abstract :: glacier

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

Type-Bound Procedures

procedure(get_scalar), public :: ice_thickness

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

procedure(get_r8), public :: ice_density

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

procedure(get_r8), public :: ice_temperature

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

procedure(get_residual), public :: residual

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

procedure(precond), public :: precondition

Applies a preconditioner to the passed state vector.

procedure(solve_vel), public :: solve_velocity

Solves for the velocity field using the current thickness.

procedure(setter), public :: update

Sets the state of the glacier.

procedure(time_setter), public :: set_time

Sets the time record for this glacier.

procedure(get_i), public :: data_size

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

procedure(get_r81d), public :: state_vector

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

procedure(read_dat), public :: read_data

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

procedure(write_dat), public :: write_data

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

procedure(t_step), public :: time_step

Calculates the appropriate time step for integration.

procedure(assign_ice), private :: 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, public :: assignment(=) => assign
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.

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.


Subroutines

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

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