shelf_integrate Subroutine

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

Integrates the glacier's state forward to time. This is done using an explicit method for the thickness and a Newton's solver for velocity.

Arguments

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


Contents

Source Code


Source Code

  subroutine shelf_integrate(this, old_states, basal_melt, basal_drag, &
                             water_density, time, success)
    !* Author: Chris MacMackin
    !  Date: May 2017
    !
    ! Integrates the glacier's state forward to `time`. This is done
    ! using an explicit method for the thickness and a Newton's solver
    ! for velocity.
    !
    class(ice_shelf), intent(inout)          :: this
    class(glacier), dimension(:), intent(in) :: 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(r8), intent(in)                     :: water_density
      !! The density of the water below the glacier.
    real(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
    
    associate(h => this%thickness, uvec => this%velocity, m => basal_melt, &
              lambda => this%lambda, chi => this%chi, zeta => this%zeta,   &
              t_old => this%time)
      this%thickness = this%thickness - (time - t_old)*(lambda*m + .div.(h*uvec))
    end associate
  end subroutine shelf_integrate