dallaston2015_melt_mod Module

Provides an implementation of abstract_melt_relationship which mimics the simple model used by Dallaston, Hewitt, and Wells (2015) for an ice shelf melting into a vertically integrated plume.


Uses

  • module~~dallaston2015_melt_mod~~UsesGraph module~dallaston2015_melt_mod dallaston2015_melt_mod factual_mod factual_mod module~dallaston2015_melt_mod->factual_mod iso_fortran_env iso_fortran_env module~dallaston2015_melt_mod->iso_fortran_env module~melt_relationship_mod melt_relationship_mod module~dallaston2015_melt_mod->module~melt_relationship_mod module~melt_relationship_mod->factual_mod module~melt_relationship_mod->iso_fortran_env

Used by

  • module~~dallaston2015_melt_mod~~UsedByGraph module~dallaston2015_melt_mod dallaston2015_melt_mod module~asymmetric_plume_mod asymmetric_plume_mod module~asymmetric_plume_mod->module~dallaston2015_melt_mod module~static_plume_mod static_plume_mod module~static_plume_mod->module~dallaston2015_melt_mod module~plume_mod plume_mod module~plume_mod->module~dallaston2015_melt_mod

Contents


Interfaces

public interface dallaston2015_melt

  • private pure function constructor(beta, melt_conversion, salinity_denom) result(this)

    Arguments

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

    The inverse stefan number,

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

    The factor to convert between the scale for melt used by Dallaston et al. (2015) and that used in ISOFT, where is the melt scale used by Dalalston et al.

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

    The factor which, when used to divide the melt_conversion term, produces a conversion factor for the melt-terms in the salinity equation of Dallaston et al. (2015). It has the form where is the subglacial discharge across the grounding line.

    Return Value type(dallaston2015_melt)

    The newly created object representing the melt relationship.


Derived Types

A parameterisation of melting into a plume which comes from heavily simplifying the 3 equation model. It is taken from Dallaston, Hewitt, and Wells (2015). The melt rate, as well as effect on termperature and salinity, are calculated by calling solve_for_melt and then accessed using melt_rate, heat_equation_terms, salt_equation_terms.

Components

TypeVisibility AttributesNameInitial
class(scalar_field), public, allocatable:: melt_values

Stores the resulting melt rate

real(kind=r8), public :: coef =1449.29936

The coefficient by which the melt rate is multiplied in order to determine the contribution to the heat equation.

real(kind=r8), public :: melt_conversion =6.9e-4_r8

The factor to convert between the scale for melt used by Dallaston et al. (2015) and that used in ISOFT, where is the melt scale used by Dalalston et al.

real(kind=r8), public :: salinity_denom =1e100_r8

Constructor

private pure function constructor(beta, melt_conversion, salinity_denom)

Type-Bound Procedures

procedure, public :: solve_for_melt => dallaston2015_solve
procedure, public :: heat_equation_terms => dallaston2015_heat

Returns the terms this melt formulation contributes to the heat equation, after they have been solved for using solve_for_melt.

procedure, public :: salt_equation_terms => dallaston2015_salt

Returns the terms this melt formulation contributes to the salt equation, after they have been solved for using solve_for_melt.

procedure, public :: melt_rate => dallaston2015_melt_rate

Returns the melt rate calculated using this formulation, after it has been solved for using solve_for_melt.

procedure, public :: has_heat_terms => dallaston2015_has_heat

Whether this formulation of melting contributes any terms to a plume's heat equation.

procedure, public :: has_salt_terms => dallaston2015_has_salt

Whether this formulation of melting contributes any terms to a plume's salinity equation.


Functions

private pure function constructor(beta, melt_conversion, salinity_denom) result(this)

Arguments

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

The inverse stefan number,

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

The factor to convert between the scale for melt used by Dallaston et al. (2015) and that used in ISOFT, where is the melt scale used by Dalalston et al.

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

The factor which, when used to divide the melt_conversion term, produces a conversion factor for the melt-terms in the salinity equation of Dallaston et al. (2015). It has the form where is the subglacial discharge across the grounding line.

Return Value type(dallaston2015_melt)

The newly created object representing the melt relationship.

private function dallaston2015_heat(this) result(heat)

Arguments

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

Return Value class(scalar_field), pointer

The value of the contribution made by melting/thermal transfer to the heat equation for a plume

private function dallaston2015_salt(this) result(salt)

Arguments

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

Return Value class(scalar_field), pointer

The value of the contribution made by melting/thermal transfer to the salt equation for a plume

private function dallaston2015_melt_rate(this) result(melt)

Arguments

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

Return Value class(scalar_field), pointer

The melt rate from the ice into the plume water.

private pure function dallaston2015_has_heat(this) result(has_heat)

Arguments

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

Return Value logical

Whether this formulation of melting contributes terms to the heat equation of the plume.

private pure function dallaston2015_has_salt(this) result(has_salt)

Arguments

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

Return Value logical

Whether this formulation of melting contributes terms to the salinity equation of the plume.


Subroutines

private subroutine dallaston2015_solve(this, velocity, pressure, temperature, salinity, plume_thickness, time)

Arguments

Type IntentOptional AttributesName
class(dallaston2015_melt), intent(inout) :: this
class(vector_field), intent(in) :: velocity

The velocity field of the plume into which fluid is melting.

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

The water pressure at the interface where the melting occurs.

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

The temperature of the plume into which fluid is melting.

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

The salinity of the plume into which fluid is melting.

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

The thickness of the plume into which fluid is melting.

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

The time at which the melting is being solved for. If not present then assumed to be same as previous value passed.