Top level methods and types

CryoGrid.AlgebraicType
Algebraic{name,S,T,units,domain} <: Var{name,S,T,units,domain}

Defines an algebraic (implicit) state variable.

CryoGrid.BCKindType
BCKind

Trait that specifies the kind of boundary condition. This can be used to write generic implementations of interact! that are (relatively) agnostic to specific implementations of BoundaryProcess. A good example of this can be found in the boundaryflux method interface.

CryoGrid.BCKindMethod
BCKind(::Type{T})

Can be overriden by BoundaryProcess types to indicate the type of boundary condition, e.g:

BCKind(::Type{BP}) = Dirichlet()

where BP is a BoundaryProcess that provides the boundary conditions.

CryoGrid.BottomType
Bottom{TProc} <: Layer

Generic "bottom" layer that marks the lower boundary of the subsurface grid.

CryoGrid.BoundaryProcessType
BoundaryProcess{T<:SubSurfaceProcess}

Abstract base type for boundary processes, i.e. processes that operate at the boundaries of the subsurface. A BoundaryProcess represents the boundary conditions of one or more SubSurfaceProcesses but may include its own diagnostic (or even prognostic) variables, if necessary.

CryoGrid.CGEulerType
CGEuler <: SciMLBase.AbstractODEAlgorithm

Simple, lightweight implementation of the forward Euler integration algorithm. Does not include support for fancier features such as interpolation, adaptive timestepping, or event handling. In order to get these features, you must install the OrdinaryDiffEq package.

CryoGrid.ConstantBCType
ConstantBC{P,S,T} <: BoundaryProcess{P}

Constant boundary condition (of any type/unit) specified by value.

CryoGrid.Coupled2Type
Coupled2{P1,P2} = CoupledProcesses{Tuple{T1,T2}} where {T1,T2}

Type alias for coupled processes, i.e. CoupledProcesses{Tuple{P1,P2}}. Coupled provides a simple mechanism for defining new behaviors on multi-processes systems.

CryoGrid.CoupledProcessesType
CoupledProcesses{TProcs} <: Process

Represents an explicitly or implicitly coupled system of processes. TProcs is always a Tuple of other processes.

CryoGrid.CryoGridProblemType
CryoGridProblem(
    tile::Tile,
    u0::ComponentVector,
    tspan::NTuple{2,Float64},
    p=nothing;
    saveat=3600.0,
    savevars=(),
    save_everystep=false,
    save_start=true,
    save_end=true,
    step_limiter=timestep,
    safety_factor=1,
    max_step=true,
    callback=nothing,
    isoutofdomain=Tiles.domain(tile),
    specialization=SciMLBase.AutoSpecialize,
    function_kwargs=(),
    prob_kwargs...
)

Constructor for CryoGridProblem that automatically generates all necessary callbacks.

CryoGrid.CryoGridProblemType
CryoGridProblem{iip,Tu,Tt,Tp,TT,Tsv,Tsf,Tcb,Tdf,Tkw} <: SciMLBase.AbstractODEProblem{Tu,Tt,iip}

Represents a CryoGrid discretized PDE forward model configuration using the SciMLBase/DiffEqBase problem interface.

CryoGrid.CryoGridProblemMethod
CryoGridProblem(tile::Tile, u0::ComponentVector, tspan::NTuple{2,DateTime}, args...;kwargs...)
CryoGrid.DiagnosticType
Diagnostic{name,S,T,units,domain} <: Var{name,S,T,units,domain}

Defines a diagnostic variable which is allocated and cached per timestep but not integrated over time.

CryoGrid.DiagnosticVolumeType
DiagnosticVolume

Volume trait instance for layers with time varying volume where the volume should be treated as a diagnostic state variable.

CryoGrid.FixedVolumeType
FixedVolume

Volume trait instance for fixed volume layers. Default for all Layer types.

CryoGrid.PeriodicBCType
PeriodicBC{P,S,T1,T2,T3,T4} <: BoundaryProcess{P}

Periodic boundary condition (of any type/unit) specified by period, amplitude, and phaseshift.

CryoGrid.PrognosticType
Prognostic{name,S,T,units,domain} <: Var{name,S,T,units,domain}

Defines a prognostic (time-integrated) state variable.

CryoGrid.PrognosticVolumeType
PrognosticVolume

Volume trait instance for layers with time varying volume where the volume should be treated as a prognostic state variable.

CryoGrid.SubSurfaceType
SubSurface <: Layer

Abstract base type for layers in the stratigraphy, e.g. soil, snow, pond, etc.

CryoGrid.SubSurfaceProcessType
SubSurfaceProcess <: Process

Abstract base type for subsurface processes, i.e. processes that operate at or below the surface, such as heat conduction, water infiltration, etc.

CryoGrid.TopType
Top{TProc} <: Layer

Generic "top" layer that marks the upper boundary of the subsurface grid.

CryoGrid.VarType
Var{name,S<:VarDim,T,units,domain}

Base type for symbolic state variables in the model.

CryoGrid.VarInitializerType
VarInitializer{varname}

Base type for state variable initializers. Initializers of this type are dedicated to a specific state variable with name varname and are run before all other initializers.

CryoGrid.VolumeMethod
Volume(::Type{<:Layer})

Trait for layer types that determines whether its spatial volume is temporally invariant, FixedVolume, or varying with time, DiagnosticVolume or PrognosticVolume.

CryoGrid.CoupledMethod
Coupled(ps::SubSurfaceProcess...)

Constructs a composite/coupled process from one or more subsurface processes. Alias for CoupledProcesses(ps...).

CryoGrid.CoupledMethod
Coupled(types::Type{<:Process}...)

Convenince method which constructs a CoupledProcesses type corresponding to each type in types, e.g:

Coupled(SnowMassBalance, HeatBalance) = CoupledProcesses{Tuple{T1,T2}} where {T1<:SnowMassBalance, T2<:HeatBalance}

also equivalent to Coupled2{<:SnowMassBalance,<:HeatBalance}.

CryoGrid.boundaryfluxMethod
boundaryflux(bc::BoundaryProcess, b::Union{Top,Bottom}, p::SubSurfaceProcess, sub::SubSurface, sbc, ssub)
boundaryflux(s::BCKind, bc::BoundaryProcess, b::Union{Top,Bottom}, p::SubSurfaceProcess, sub::SubSurface, sbc, ssub)

Computes the flux dH/dt at the boundary layer. Calls boundaryflux(BCKind(B),...) to allow for generic implementations by boundary condition type. Note that this method uses a different argument order convention than interact!. This is intended to faciliate stratigraphy independent implementations of certain boundary conditions (e.g. a simple Dirichlet boundary could be applied in the same manner to both the upper and lower boundary).

CryoGrid.boundaryvalueMethod
boundaryvalue(bc::BoundaryProcess, state)

Computes the value of the boundary condition specified by bc for the given layer/process combinations. Note that this method uses a different argument order convention than interact!. This is intended to faciliate stratigraphy independent implementations of certain boundary conditions (e.g. a simple Dirichlet boundary could be applied in the same manner to both the upper and lower boundary).

CryoGrid.caninteractMethod
caninteract(layer1::Layer, layer2::Layer, state1, state2)
caninteract(l1::Layer, ::Process, l2::Layer, ::Process, state1, state2)

Returns true if and only if the given layer/process types are able to interact based on the current state. Defaults to checking whether both layers are currently active. This behavior should be overridden by subtypes where necessary.

CryoGrid.computediagnostic!Method
computediagnostic!(l::Layer, state)
computediagnostic!(l::Layer, p::Process, state)

Updates all diagnostic/non-flux state variables for the given Layer based on the current prognostic state.

CryoGrid.computefluxes!Method
computefluxes!(l::Layer, p::Process, state)

Calculates all internal fluxes for a given layer. Note that an instance of computefluxes! must be provided for all non-boundary (subsurface) processes/layers.

CryoGrid.criterion!Method
criterion!(out::AbstractArray, ev::GridContinuousEvent, ::Layer, ::Process, state)

Event criterion for on-grid (i.e. multi-valued) continuous events. The condition for each grid cell should be stored in out.

CryoGrid.criterionMethod
criterion(::Event, ::Layer, ::Process, state)

Event criterion/condition. Should return a Bool for discrete events. For continuous events, this should be a real-valued function where the event is fired at the zeros/roots.

CryoGrid.diagnosticstep!Method
diagnosticstep!(layer::Layer, state)

Optionally performs discontinuous/discrete-time updates to the layer state. Should return true if the prognostic state was modified and false otherwise. Defaults to returning false.

CryoGrid.initialcondition!Method
initialcondition!(::Layer, state)
initialcondition!(::Layer, ::Process, state)
initialcondition!(::VarInitializer, ::Layer, state)

Defines the initial condition for a given Layer and possibly an initializer. initialcondition! should compute initial values into all relevant state variables in state.

CryoGrid.initialcondition!Method
initialcondition!(layer1::Layer, layer2::Layer, state1, state2)
initialcondition!(::Layer, ::Process, ::Layer, ::Process, state1, state2)

Defines the initial condition for two processes on adjacent layers. initialcondition! should write initial values into all relevant state variables in state.

CryoGrid.initializerMethod
initializer(varname::Symbol, args...) = initializer(Val{varname}(), args...)
initializer(::Val{varname}, x::Number) => FunctionInitializer w/ constant
initializer(::Val{varname}, f::Function) => FunctionInitializer
initializer(::Val{varname}, profile::Profile, interp=Linear(), extrap=Flat()) => InterpInitializer

Convenience constructor for VarInitializer that selects the appropriate initializer type based on the arguments.

CryoGrid.initializersMethod
initializers(::Layer)
initializers(::Layer, ::Process)

Optional method that can be used to provide default initializers for state variables that will be run before user provided ones.

CryoGrid.interact!Method
interact!(::Layer, ::Process, ::Layer, ::Process, state1, state2)

Defines a boundary interaction between two processes on adjacent layers. For any interaction, the order of the arguments follows decreasing depth, i.e. the first layer/process is always on top of the second layer/process. This ordering matters and separate dispatches must be provided for interactions in reverse order.

CryoGrid.interactmaybe!Method
interactmaybe!(layer1::Layer, layer2::Layer, state1, state2)
interactmaybe!(layer1::Layer, p1::Process, layer2::Layer, p2::Process, state1, state2)

Conditionally invokes interact! if and only if caninteract is true.

CryoGrid.isactiveMethod
isactive(::Layer, state)

Returns a boolean whether or not this layer is currently active in the stratigraphy and should interact with other layers. Note that computediagnostic! and computefluxes! are always invoked regardless of the current state of isactive. The default implementation of isactive always returns true.

CryoGrid.processesMethod
processes(l::Layer)

Fetches the process(es) attached to this layer, if any. Returned value must be of type Process. If the layer has more than one process, they should be combined together with Coupled(procs...).

CryoGrid.resetfluxes!Method
resetfluxes!(layer::Layer, state)
resetfluxes!(layer::Layer, ::Process, state)

Resets all flux variables for the given layer/process to zero.

CryoGrid.timestepMethod
timestep(::Layer, ::Process, state)

Retrieves the recommended timestep for the given Process defined on the given Layer. The default implementation returns Inf which indicates no timestep restriction. The actual chosen timestep will depend on the integrator being used and other user configuration options.

CryoGrid.timestepMethod
timestep(l::Layer, ps::CoupledProcesses{P}, state) where {P}

Default implementation of timestep for coupled process types. Calls each process in sequence.

CryoGrid.trigger!Method
trigger!(::Event, ::Layer, ::Process, state)
trigger!(ev::ContinuousEvent, ::ContinuousTrigger, ::Layer, ::Process, state)
trigger!(ev::GridContinuousEvent, ::ContinuousTrigger, ::Layer, ::Process, state)

Event action executed when criterion is met.

CryoGrid.variablesMethod
variables(layer::Layer, process::Process)
variables(::Layer)
variables(::Any)

Defines variables for a given Layer, Process, or arbitrary user-defined type. Implementations should return a Tuple of Vars.

CryoGrid.volumetricfractionsMethod
volumetricfractions(::SubSurface, state)
volumetricfractions(::SubSurface, state, i)

Get the volumetric fractions of each constituent in the volume (at grid cell i, if specificed). All implementations of volumetricfractions are expected to obey a semi-consistent order in the returned Tuple of fractions; the first three consituents should always be θw,θi,θa, i.e. water, ice, and air, followed by any number of additional constituents which may be defined by the specific layer. There is no feasible way to verify that client code actually obeys this ordering, so be sure to double check your implementation, otherwise this can cause very subtle bugs!