ice_sheet_mod Module

Provides a concrete implementation of the glacier type, using a vertically integrated model of an ice sheet.


Uses

  • module~~ice_sheet_mod~~UsesGraph module~ice_sheet_mod ice_sheet_mod logger_mod logger_mod module~ice_sheet_mod->logger_mod factual_mod factual_mod module~ice_sheet_mod->factual_mod hdf5 hdf5 module~ice_sheet_mod->hdf5 module~glacier_mod glacier_mod module~ice_sheet_mod->module~glacier_mod module~jacobian_block_mod jacobian_block_mod module~ice_sheet_mod->module~jacobian_block_mod iso_fortran_env iso_fortran_env module~ice_sheet_mod->iso_fortran_env module~viscosity_mod viscosity_mod module~ice_sheet_mod->module~viscosity_mod module~glacier_mod->logger_mod module~glacier_mod->factual_mod module~glacier_mod->hdf5 module~glacier_mod->iso_fortran_env penf penf module~glacier_mod->penf module~nitsol_mod nitsol_mod module~glacier_mod->module~nitsol_mod module~jacobian_block_mod->logger_mod module~jacobian_block_mod->factual_mod module~jacobian_block_mod->iso_fortran_env module~jacobian_block_mod->penf f95_lapack f95_lapack module~jacobian_block_mod->f95_lapack module~boundary_types_mod boundary_types_mod module~jacobian_block_mod->module~boundary_types_mod module~viscosity_mod->factual_mod module~viscosity_mod->iso_fortran_env module~nitsol_mod->iso_fortran_env

Contents


Interfaces

public interface ice_sheet

  • private function constructor(domain, resolution, thickness, velocity, viscosity_law, lambda, chi) result(this)

    Author
    Christopher MacMackin
    Date
    April 2016

    Creates a new ice_sheet object with initial conditions provided by the arguments. At present only a 1D model is supported. If information is provided for higher dimensions then it will be ignored.

    Arguments

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

    An array containing the upper and lower limits of the domain for the ice sheet. The first index represents the dimension for which the boundaries apply. If the second index is 1 then it corresponds to the lower bound. If the second index is 2 then it corresponds to the upper bound.

    integer, intent(in), dimension(:):: resolution

    The number of data points in each dimension

    procedure(thickness_func) :: thickness

    A function which calculates the initial value of the thickness of the ice sheet at a given location.

    procedure(velocity_func) :: velocity

    A function which calculates the initial value of the velocity (vector) of the ice at a given location in an ice sheet.

    class(abstract_viscosity), intent(in), optional :: viscosity_law

    An object which calculates the viscosity of the ice.

    real(kind=r8), intent(in), optional :: lambda

    The dimensionless ratio $\lambda \equiv \frac{\rho_0m_0x_0}{\rho_iH-0u_0}$

    real(kind=r8), intent(in), optional :: chi

    The dimensionless ratio $\chi \equiv \frac{\rho_igh_0x_x}{2\eta_0u_0}$

    Return Value type(ice_sheet)

    An ice sheet object with its domain and initial conditions set according to the arguments of the constructor function.


Derived Types

type, public, extends(glacier) :: ice_sheet

A concrete implementation of the glacier type, using a vertically integrated model of an ice sheet. This model is 1-dimensional only.

Components

TypeVisibility AttributesNameInitial
type(cheb1d_scalar_field), private :: thickness

Thickness of ice sheet, $h$

type(cheb1d_vector_field), private :: velocity

Flow velocity of ice sheet, $\vec{u}$

real(kind=r8), private :: lambda

The dimensionless ratio $\lambda \equiv \frac{\rho_0m_0x_0}{\rho_iH-0u_0}$

real(kind=r8), private :: chi

The dimensionless ratio $\chi \equiv \frac{\rho_igh_0x_x}{2\eta_0u_0}$

class(abstract_viscosity), private, allocatable:: viscosity_law

An object representing the model used for ice viscosity.

real(kind=r8), private :: time

The time at which the ice sheet is in this state

Constructor

private function constructor(domain, resolution, thickness, velocity, viscosity_law, lambda, chi)

Creates a new ice_sheet object with initial conditions provided by the arguments. At present only a 1D model is supported. If information is provided for higher dimensions then it will be ignored.

Type-Bound Procedures

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.

procedure, public :: ice_thickness => sheet_thickness
procedure, public :: ice_density => sheet_density
procedure, public :: ice_temperature => sheet_temperature
procedure, public :: residual => sheet_residual
procedure, public :: update => sheet_update
procedure, public :: precondition => sheet_precondition
procedure, public :: set_time => sheet_set_time
procedure, public :: data_size => sheet_data_size
procedure, public :: state_vector => sheet_state_vector
procedure, public :: solve_velocity => sheet_solve_velocity
procedure, public :: read_data => sheet_read_data
procedure, public :: write_data => sheet_write_data
procedure, public :: time_step => sheet_time_step
procedure, private :: assign => sheet_assign

Functions

private function constructor(domain, resolution, thickness, velocity, viscosity_law, lambda, chi) result(this)

Author
Christopher MacMackin
Date
April 2016

Creates a new ice_sheet object with initial conditions provided by the arguments. At present only a 1D model is supported. If information is provided for higher dimensions then it will be ignored.

Arguments

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

An array containing the upper and lower limits of the domain for the ice sheet. The first index represents the dimension for which the boundaries apply. If the second index is 1 then it corresponds to the lower bound. If the second index is 2 then it corresponds to the upper bound.

integer, intent(in), dimension(:):: resolution

The number of data points in each dimension

procedure(thickness_func) :: thickness

A function which calculates the initial value of the thickness of the ice sheet at a given location.

procedure(velocity_func) :: velocity

A function which calculates the initial value of the velocity (vector) of the ice at a given location in an ice sheet.

class(abstract_viscosity), intent(in), optional :: viscosity_law

An object which calculates the viscosity of the ice.

real(kind=r8), intent(in), optional :: lambda

The dimensionless ratio $\lambda \equiv \frac{\rho_0m_0x_0}{\rho_iH-0u_0}$

real(kind=r8), intent(in), optional :: chi

The dimensionless ratio $\chi \equiv \frac{\rho_igh_0x_x}{2\eta_0u_0}$

Return Value type(ice_sheet)

An ice sheet object with its domain and initial conditions set according to the arguments of the constructor function.

private pure function sheet_thickness(this) result(thickness)

Author
Christopher MacMackin
Date
April 2016

Returns the thickness of the ice sheet across its domain.

Arguments

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

Return Value class(scalar_field), pointer

The ice thickness.

private function sheet_velocity(this) result(velocity)

Author
Christopher MacMackin
Date
July 2016

Returns the velocity of the ice sheet across its domain.

Arguments

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

Return Value class(vector_field), allocatable

The ice velocity.

private pure function sheet_density(this) result(density)

Author
Christopher MacMackin
Date
April 2016

Returns the density of the ice in the sheet, which is assumed to be uniform across its domain.

Arguments

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

Return Value real(kind=r8)

The ice density.

private pure function sheet_temperature(this) result(temperature)

Author
Christopher MacMackin
Date
April 2016

Returns the density of the ice in the sheet, which is assumed to be uniform across its domain.

Arguments

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

Return Value real(kind=r8)

The ice density.

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

Author
Christopher MacMackin
Date
April 2016

Returns the residual when the current state of the glacier is run through the system of equations describing it. The residual takes the form of a 1D array, with each element respresenting the residual for one of the equations in the system.

Arguments

Type IntentOptional AttributesName
class(ice_sheet), 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.

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

Author
Chris MacMackin
Date
January 2016

Provides a preconditioner for the nonlinear solver trying to bring the residual to zero. The Jacobian is approximated as a block matrix, where each block is a tridiagonal matrix using a finite difference method for differentiation.

Arguments

Type IntentOptional AttributesName
class(ice_sheet), 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.

private pure function sheet_data_size(this)

Author
Christopher MacMackin
Date
August 2016

Returns the number of elements in the ice sheet's state vector. This is the size of the vector returned by residual and state_vector and taken as an argument by update.

Arguments

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

Return Value integer

The number of elements in the ice sheet's state vector.

private pure function sheet_state_vector(this) result(state_vector)

Author
Christopher MacMackin
Date
April 2016

Returns the state vector for the current state of the ice sheet. This takes the form of a 1D array.

Arguments

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

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

The state vector describing the ice sheet.

private function sheet_time_step(this) result(dt)

Author
Chris MacMackin
Date
December 2016

Calculates the time step for integrating the ice sheet, using the CFL condition.

Arguments

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

Return Value real(kind=r8)

The time-step to use


Subroutines

private subroutine sheet_update(this, state_vector)

Author
Christopher MacMackin
Date
April 2016

Updates the state of the ice sheet from its state vector. The state vector is a real array containing the value of each of the ice sheet's properties at each of the locations on the grid used in descretization.

Arguments

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

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

private subroutine sheet_set_time(this, time)

Author
Christopher MacMackin
Date
November 2016

Sets the time information held by the ice sheet object. This is the time at which the ice sheet is in its current state.

Arguments

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

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

private subroutine sheet_solve_velocity(this, basal_drag, success)

Author
Chris MacMackin
Date
May 2017

Computes the ice sheet velocity at the current time with the current ice thickness.

Arguments

Type IntentOptional AttributesName
class(ice_sheet), 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

private subroutine sheet_read_data(this, file_id, group_name, error)

Author
Chris MacMackin
Date
April 2017

Reads the state of the ice shelf object from the specified group in an HDF5 file. This sets the thickness, the velocity, and parameter values.

Arguments

Type IntentOptional AttributesName
class(ice_sheet), intent(inout) :: 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 ice shelf's data.

integer, intent(out) :: error

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

private subroutine sheet_write_data(this, file_id, group_name, error)

Author
Chris MacMackin
Date
November 2016

Writes the state of the ice sheet object to an HDF file in the specified group. This will consist of a thickness and a velocity dataset.

Arguments

Type IntentOptional AttributesName
class(ice_sheet), 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 ice sheet's data.

integer, intent(out) :: error

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

private subroutine sheet_assign(this, rhs)

Author
Chris MacMackin
Date
February 2017

Copies the data from one ice sheet 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_sheet), intent(out) :: this
class(glacier), intent(in) :: rhs

The ice sheet to be assigned to this one.