GridapROMs.ParamDataStructures
GridapROMs.ParamDataStructures.SparseSnapshots
— TypeGridapROMs.ParamDataStructures.TransientSparseSnapshots
— TypeGridapROMs.ParamDataStructures.AbstractParamArray
— Typeabstract type AbstractParamArray{T,N,A<:AbstractArray{T,N}} <: AbstractParamData{A,N} end
Type representing parametric abstract arrays of type A. Subtypes:
GridapROMs.ParamDataStructures.AbstractParamArray3D
— Typeconst AbstractParamArray3D{T} = AbstractParamArray{T,3,<:AbstractArray{T,3}}
GridapROMs.ParamDataStructures.AbstractParamData
— Typeabstract type AbstractParamData{T,N} <: AbstractArray{T,N} end
Type representing generic parametric quantities. Subtypes:
GridapROMs.ParamDataStructures.AbstractParamFunction
— Typeabstract type AbstractParamFunction{P<:Realization} <: Function end
Representation of parametric functions with domain a parametric space. Subtypes:
GridapROMs.ParamDataStructures.AbstractParamMatrix
— Typeconst AbstractParamMatrix{T} = AbstractParamArray{T,2,<:AbstractMatrix{T}}
GridapROMs.ParamDataStructures.AbstractParamVector
— Typeconst AbstractParamVector{T} = AbstractParamArray{T,1,<:AbstractVector{T}}
GridapROMs.ParamDataStructures.AbstractRealization
— Typeabstract type AbstractRealization end
Type representing parametric realizations, i.e. samples extracted from a given parameter space. Subtypes:
GridapROMs.ParamDataStructures.AbstractSnapshots
— Typeabstract type AbstractSnapshots{T,N} <: AbstractParamData{T,N} end
Type representing N-dimensional arrays of snapshots. Subtypes must contain the following information:
- data: a (parametric) array
- realization: a subtype of
AbstractRealization
, representing the points in the parameter space used to compute the arraydata
- dof map: a subtype of
AbstractDofMap
, representing a reindexing strategy for the arraydata
Subtypes:
GridapROMs.ParamDataStructures.ArrayOfArrays
— Typestruct ArrayOfArrays{T,N,A<:AbstractArray{T,N}} <: ParamArray{T,N}
data::Vector{A}
end
Parametric array with entries stored non-consecutively in memory. It is characterized by an inner size equal to size(data[1])
, and parametric length equal to length(data)
, where data
is a Vector{<:AbstractArray}
GridapROMs.ParamDataStructures.BlockConsecutiveParamMatrix
— Typeconst BlockConsecutiveParamMatrix{T,A<:Matrix{<:ConsecutiveParamMatrix{T}},B} = BlockParamMatrix{T,A,B}
GridapROMs.ParamDataStructures.BlockConsecutiveParamVector
— Typeconst BlockConsecutiveParamVector{T,A<:Vector{<:ConsecutiveParamVector{T}},B} = BlockParamVector{T,A,B}
GridapROMs.ParamDataStructures.BlockParamArray
— Typestruct BlockParamArray{T,N,A<:AbstractArray{<:AbstractParamArray{T,N},N},B<:NTuple{N,AbstractUnitRange{Int}}} <: ParamArray{T,N}
data::A
axes::B
end
Is to a ParamArray
as a BlockArray
is to a regular AbstractArray
GridapROMs.ParamDataStructures.BlockParamMatrix
— Typeconst BlockParamMatrix{T,A,B} = BlockParamArray{T,2,A,B}
GridapROMs.ParamDataStructures.BlockParamVector
— Typeconst BlockParamVector{T,A,B} = BlockParamArray{T,1,A,B}
GridapROMs.ParamDataStructures.BlockSnapshots
— Typestruct BlockSnapshots{S<:Snapshots,N} <: AbstractSnapshots{S,N}
array::Array{S,N}
touched::Array{Bool,N}
end
Block container for Snapshots of type S
in a MultiField
setting. This type is conceived similarly to ArrayBlock
in Gridap
GridapROMs.ParamDataStructures.ConsecutiveParamArray
— Typestruct ConsecutiveParamArray{T,N,M,A<:AbstractArray{T,M}} <: ParamArray{T,N}
data::A
end
Parametric array with entries stored consecutively in memory. It is characterized by an inner size equal to size(data)[1:N]
, and parametric length equal to size(data,N+1)
, where data
is an AbstractArray
of dimension M = N+1
GridapROMs.ParamDataStructures.ConsecutiveParamMatrix
— Typeconst ConsecutiveParamMatrix{T,A<:AbstractArray{T,3}} = ConsecutiveParamArray{T,2,3,A}
GridapROMs.ParamDataStructures.ConsecutiveParamSparseMatrixCSC
— Typestruct ConsecutiveParamSparseMatrixCSC{Tv,Ti<:Integer} <: ParamSparseMatrixCSC{Tv,Ti}
m::Int64
n::Int64
colptr::Vector{Ti}
rowval::Vector{Ti}
data::Matrix{Tv}
end
Represents a vector of sparse matrices in CSC format, with entries stored consecutively in memory. For sake of coherence, an instance of ConsecutiveParamSparseMatrixCSC
inherits from AbstractMatrix{<:SparseMatrixCSC{Tv,Ti}
rather than AbstractVector{<:SparseMatrixCSC{Tv,Ti}
, but should conceptually be thought as an AbstractVector{<:SparseMatrixCSC{Tv,Ti}
.
GridapROMs.ParamDataStructures.ConsecutiveParamSparseMatrixCSR
— Typestruct ConsecutiveParamSparseMatrixCSR{Bi,Tv,Ti<:Integer} <: ParamSparseMatrixCSR{Bi,Tv,Ti}
m::Int64
n::Int64
rowptr::Vector{Ti}
colval::Vector{Ti}
data::Matrix{Tv}
end
Represents a vector of sparse matrices in CSR format, with entries stored consecutively in memory. For sake of coherence, an instance of ConsecutiveParamSparseMatrixCSR
inherits from AbstractMatrix{<:SparseMatrixCSR{Bi,Tv,Ti}
rather than AbstractVector{<:SparseMatrixCSR{Bi,Tv,Ti}
, but should conceptually be thought as an AbstractVector{<:SparseMatrixCSR{Bi,Tv,Ti}
.
GridapROMs.ParamDataStructures.ConsecutiveParamVector
— Typeconst ConsecutiveParamVector{T,A<:AbstractArray{T,2}} = ConsecutiveParamArray{T,1,2,A}
GridapROMs.ParamDataStructures.GenericParamBlock
— Typestruct GenericParamBlock{A} <: ParamBlock{A}
data::Vector{A}
end
Most standard implementation of a ParamBlock
GridapROMs.ParamDataStructures.GenericParamMatrix
— Typestruct GenericParamMatrix{Tv,Ti} <: ParamArray{Tv,2}
data::Vector{Tv}
ptrs::Vector{Ti}
nrows::Vector{Ti}
end
Parametric matrix with entries stored non-consecutively in memory
GridapROMs.ParamDataStructures.GenericParamSparseMatrixCSC
— Typestruct GenericParamSparseMatrixCSC{Tv,Ti<:Integer} <: ParamSparseMatrixCSC{Tv,Ti}
m::Int64
n::Int64
colptr::Vector{Ti}
rowval::Vector{Ti}
data::Vector{Tv}
ptrs::Vector{Ti}
end
Represents a vector of sparse matrices in CSC format, with entries stored non-consecutively in memory. For sake of coherence, an instance of GenericParamSparseMatrixCSC
inherits from AbstractMatrix{<:SparseMatrixCSC{Tv,Ti}
rather than AbstractVector{<:SparseMatrixCSC{Tv,Ti}
, but should conceptually be thought as an AbstractVector{<:SparseMatrixCSC{Tv,Ti}
.
GridapROMs.ParamDataStructures.GenericParamSparseMatrixCSR
— Typestruct GenericParamSparseMatrixCSR{Bi,Tv,Ti<:Integer} <: ParamSparseMatrixCSR{Bi,Tv,Ti}
m::Int64
n::Int64
rowptr::Vector{Ti}
colval::Vector{Ti}
data::Vector{Tv}
ptrs::Vector{Ti}
end
Represents a vector of sparse matrices in CSR format, with entries stored non-consecutively in memory. For sake of coherence, an instance of GenericParamSparseMatrixCSR
inherits from AbstractMatrix{<:SparseMatrixCSR{Bi,Tv,Ti}
rather than AbstractVector{<:SparseMatrixCSR{Bi,Tv,Ti}
, but should conceptually be thought as an AbstractVector{<:SparseMatrixCSR{Bi,Tv,Ti}
.
GridapROMs.ParamDataStructures.GenericParamVector
— Typestruct GenericParamVector{Tv,Ti} <: ParamArray{Tv,1}
data::Vector{Tv}
ptrs::Vector{Ti}
end
Parametric vector with entries stored non-consecutively in memory
GridapROMs.ParamDataStructures.GenericSnapshots
— Typestruct GenericSnapshots{T,N,D,I,R,A} <: SteadySnapshots{T,N,D,I,R,A}
data::A
dof_map::I
realization::R
end
Most standard implementation of a SteadySnapshots
GridapROMs.ParamDataStructures.GenericTransientRealization
— Typestruct GenericTransientRealization{P,T,A} <: TransientRealization{P,T}
params::P
times::A
t0::T
end
Most standard implementation of a TransientRealization
.
GridapROMs.ParamDataStructures.HaltonSampling
— Typestruct HaltonSampling <: SamplingStyle end
Sampling according to a Halton sequence
Halton is a sequence, not a distribution, hence this sampling strategy repeats realizations since the draws are not randomized; to draw different parameters, one needs to provide a starting point in the sequence (start = 1 by default)
GridapROMs.ParamDataStructures.LatinHypercubeSampling
— Typestruct LatinHypercubeSampling <: SamplingStyle end
Sampling according to a Latin HyperCube distribution
GridapROMs.ParamDataStructures.ModeTransientSnapshots
— Typestruct ModeTransientSnapshots{M<:ModeAxes,T,I,R,A<:UnfoldingTransientSnapshots{T,I,R}} <: TransientSnapshots{T,2,1,I,R,A}
snaps::A
mode::M
end
Represents a TransientSnapshots
with a TrivialDofMap
indexing strategy as an AbstractMatrix
with a system of mode-unfolding representations. Only two mode-unfolding representations are considered:
Mode1Axes:
[u(x1,t1,μ1) ⋯ u(x1,t1,μP) u(x1,t2,μ1) ⋯ u(x1,t2,μP) u(x1,t3,μ1) ⋯ ⋯ u(x1,tT,μ1) ⋯ u(x1,tT,μP)] ⋮ ⋮ ⋮ ⋮ ⋮ ⋮ ⋮ u(xN,t1,μ1) ⋯ u(xN,t1,μP) u(xN,t2,μ1) ⋯ u(xN,t2,μP) u(xN,t3,μ1) ⋯ ⋯ u(xN,tT,μ1) ⋯ u(xN,tT,μP)]
Mode2Axes:
[u(x1,t1,μ1) ⋯ u(x1,t1,μP) u(x2,t1,μ1) ⋯ u(x2,t1,μP) u(x3,t1,μ1) ⋯ ⋯ u(xN,t1,μ1) ⋯ u(xN,t1,μP)] ⋮ ⋮ ⋮ ⋮ ⋮ ⋮ ⋮ u(x1,tT,μ1) ⋯ u(x1,tT,μP) u(x2,tT,μ1) ⋯ u(x2,tT,μP) u(x3,tT,μ1) ⋯ ⋯ u(xN,tT,μ1) ⋯ u(xN,tT,μP)]
GridapROMs.ParamDataStructures.NormalSampling
— Typestruct NormalSampling <: SamplingStyle end
Sampling according to a normal distribution
GridapROMs.ParamDataStructures.ParamArray
— Typeabstract type ParamArray{T,N} <: AbstractParamArray{T,N,Array{T,N}} end
Type representing parametric arrays of type A. Subtypes:
TrivialParamArray
ConsecutiveParamArray
GenericParamVector
GenericParamMatrix
ArrayOfArrays
BlockParamArray
Also acts as a constructor according to the following rules:
- ParamArray(A::AbstractArray{<:Number}) -> ParamNumber
- ParamArray(A::AbstractArray{<:Number},plength::Int) -> TrivialParamArray
- ParamArray(A::AbstractVector{<:AbstractArray}) -> ParamArray
- ParamArray(A::AbstractVector{<:AbstractSparseMatrix}) -> ParamSparseMatrix
- ParamArray(A::AbstractArray{<:ParamArray}) -> BlockParamArray
GridapROMs.ParamDataStructures.ParamBlock
— Typeabstract type ParamBlock{A} end
Type representing cell-wise quantities defined during the integration routine. They are primarily used when lazily evaluating parametric quantities on the mesh. The implementation of the lazy interface mimics that of ArrayBlock
in Gridap
. Subtypes: -GenericParamBlock
-TrivialParamBlock
GridapROMs.ParamDataStructures.ParamFunction
— Typestruct ParamFunction{F,P} <: AbstractParamFunction{P}
fun::F
params::P
end
Representation of parametric functions with domain a parametric space. Given a function f
: Ω₁ × ... × Ωₙ × D, where D is a ParamSpace
, the evaluation of f
in μ ∈ D
returns the restriction of f
to Ω₁ × ... × Ωₙ
GridapROMs.ParamDataStructures.ParamMatrix
— Typeconst ParamMatrix{T} = ParamArray{T,2}
GridapROMs.ParamDataStructures.ParamNumber
— Typestruct ParamNumber{T<:Number,A<:AbstractVector{T}} <: AbstractParamData{T,1}
data::A
end
Used as a wrapper for non-array structures, e.g. factorizations or numbers
GridapROMs.ParamDataStructures.ParamSpace
— Typestruct ParamSpace{P<:AbstractVector{<:AbstractVector},S<:SamplingStyle} <: AbstractSet{Realization}
param_domain::P
sampling_style::S
end
Fields:
param_domain
: domain of definition of the parameterssampling_style
: distribution onparam_domain
according to which we can sample the parameters (by default it is set toHaltonSampling
)
GridapROMs.ParamDataStructures.ParamSparseMatrix
— Typeabstract type ParamSparseMatrix{Tv,Ti,A<:AbstractSparseMatrix{Tv,Ti}
} <: AbstractParamArray{Tv,2,A} end
Type representing parametric abstract sparse matrices of type A. Subtypes:
GridapROMs.ParamDataStructures.ParamSparseMatrixCSC
— Typeabstract type ParamSparseMatrixCSC{Tv,Ti} <: ParamSparseMatrix{Tv,Ti,SparseMatrixCSC{Tv,Ti}} end
Type representing parametric sparse matrices in CSC format. Subtypes:
GridapROMs.ParamDataStructures.ParamSparseMatrixCSR
— Typeabstract type ParamSparseMatrixCSR{Bi,Tv,Ti} <: ParamSparseMatrix{Tv,Ti,SparseMatrixCSR{Bi,Tv,Ti}} end
Type representing parametric sparse matrices in CSR format. Subtypes:
GridapROMs.ParamDataStructures.ParamVector
— Typeconst ParamVector{T} = ParamArray{T,1}
GridapROMs.ParamDataStructures.ParamVectorWithEntryInserted
— Typestruct ParamVectorWithEntryInserted{T,A} <: ParamVector{T}
a::A
index::Int
value::Vector{T}
end
GridapROMs.ParamDataStructures.ParamVectorWithEntryRemoved
— Typestruct ParamVectorWithEntryRemoved{T,A} <: ParamVector{T}
a::A
index::Int
end
GridapROMs.ParamDataStructures.Realization
— Typestruct Realization{P<:AbstractVector} <: AbstractRealization
params::P
end
Represents standard parametric realizations, i.e. samples extracted from a given parameter space. The field params
is most commonly a vector of vectors. When the parameters are scalars, they still need to be defined as vectors of vectors of unit length. In other words, we treat the case in which params
is a vector of numbers as the case in which params
is a vector of one vector.
GridapROMs.ParamDataStructures.ReshapedSnapshots
— Typestruct ReshapedSnapshots{T,N,N′,D,I,R,A<:SteadySnapshots{T,N′,D,I,R},B} <: SteadySnapshots{T,N,D,I,R,A}
snaps::A
size::NTuple{N,Int}
mi::B
end
Represents a SteadySnapshots snaps
whose size is resized to size
. This struct is equivalent to ReshapedArray
, and is only used to make sure the result of this operation is still a subtype of SteadySnapshots
GridapROMs.ParamDataStructures.SamplingStyle
— Typeabstract type SamplingStyle end
Subtypes:
GridapROMs.ParamDataStructures.Snapshots
— Typeabstract type Snapshots{T,N,D,I<:AbstractDofMap{D},R<:AbstractRealization,A}
<: AbstractSnapshots{T,N} end
Type representing a collection of parametric abstract arrays of eltype T
, that are associated with a realization of type R
. The (spatial) entries of any instance of Snapshots
are indexed according to an index map of type I
:AbstractDofMap{D
}, where D
encodes the spatial dimension. Note that, as opposed to subtypes of AbstractParamArray
, which are arrays of arrays, subtypes of Snapshots
are arrays of numbers.
Subtypes:
GridapROMs.ParamDataStructures.SnapshotsAtIndices
— Typestruct SnapshotsAtIndices{T,N,D,I,R,A<:SteadySnapshots{T,N,D,I,R},B} <: SteadySnapshots{T,N,D,I,R,A}
snaps::A
prange::B
end
Represents a SteadySnapshots snaps
whose parametric range is restricted to the indices in prange
. This type essentially acts as a view for suptypes of SteadySnapshots, at every space location, on a selected number of parameter indices. An instance of SnapshotsAtIndices is created by calling the function select_snapshots
GridapROMs.ParamDataStructures.SteadySnapshots
— Typeabstract type SteadySnapshots{T,N,D,I,A} <: Snapshots{T,N,D,I,<:Realization,A} end
Spatial specialization of an Snapshots
. The dimension N
of a SteadySnapshots is equal to D
+ 1, where D
represents the number of spatial axes, to which a parametric dimension is added.
Subtypes:
GridapROMs.ParamDataStructures.TensorialUniformSampling
— Typestruct TensorialUniformSampling <: SamplingStyle end
Sampling according to a tensorial uniform distribution
GridapROMs.ParamDataStructures.TransientGenericSnapshots
— Typestruct TransientGenericSnapshots{T,N,D,I,R,A} <: TransientSnapshots{T,N,D,I,R,A}
data::A
dof_map::I
realization::R
end
Most standard implementation of a TransientSnapshots
GridapROMs.ParamDataStructures.TransientParamFunction
— Typestruct TransientParamFunction{F,P,T} <: AbstractParamFunction{P}
fun::F
params::P
times::T
end
Representation of parametric functions with domain a transient parametric space. Given a function f : Ω₁ × ... × Ωₙ × D × [t₁,t₂]
, where [t₁,t₂]
is a temporal domain and D
is a ParamSpace
, or equivalently f : Ω₁ × ... × Ωₙ × D × [t₁,t₂]
, where D
is a TransientParamSpace
, the evaluation of f
in (μ,t) ∈ D × [t₁,t₂]
returns the restriction of f
to Ω₁ × ... × Ωₙ
GridapROMs.ParamDataStructures.TransientParamSpace
— Typestruct TransientParamSpace{P<:ParamSpace,T} <: AbstractSet{TransientRealization}
parametric_space::P
temporal_domain::T
end
Fields:
parametric_space
: underlying parameter spacetemporal_domain
: underlying temporal space
It represents, in essence, the set of tuples (p,t) in parametric_space
× temporal_domain
GridapROMs.ParamDataStructures.TransientRealization
— Typeabstract type TransientRealization{P<:Realization,T<:Real} <: AbstractRealization end
Represents temporal parametric realizations, i.e. samples extracted from a given parameter space for every time step in a temporal range. The most obvious application of this type are transient PDEs, where an initial condition is given. Following this convention, the initial time instant is kept separate from the other time steps.
GridapROMs.ParamDataStructures.TransientRealizationAt
— Typestruct TransientRealizationAt{P,T} <: TransientRealization{P,T}
params::P
t::Base.RefValue{T}
end
Represents a GenericTransientRealization{P,T} at a certain time instant t
. To avoid making it a mutable struct, the time instant t
is stored as a Base.RefValue{T}.
GridapROMs.ParamDataStructures.TransientReshapedSnapshots
— Typestruct TransientReshapedSnapshots{T,N,N′,D,I,R,A<:TransientSnapshots{T,N′,D,I,R},B} <: TransientSnapshots{T,N,D,I,R,A}
snaps::A
size::NTuple{N,Int}
mi::B
end
Represents a TransientSnapshots
snaps
whose size is resized to size
. This struct is equivalent to ReshapedArray
, and is only used to make sure the result of this operation is still a subtype of TransientSnapshots
GridapROMs.ParamDataStructures.TransientSnapshots
— Typeabstract type TransientSnapshots{T,N,D,I,R<:TransientRealization,A} <: Snapshots{T,N,D,I,R,A} end
Transient specialization of a Snapshots
. The dimension N
of a SteadySnapshots is equal to D
+ 2, where D
represents the number of spatial axes, to which a temporal and a parametric dimension are added.
Subtypes:
GridapROMs.ParamDataStructures.TransientSnapshotsAtIndices
— Typestruct TransientSnapshotsAtIndices{T,N,D,I,R,A<:TransientSnapshots{T,N,D,I,R},B,C} <: TransientSnapshots{T,N,D,I,R,A}
snaps::A
trange::B
prange::C
end
Represents a TransientSnapshots
snaps
whose parametric and temporal ranges are restricted to the indices in prange
and trange
. This type essentially acts as a view for suptypes of TransientSnapshots
, at every space location, on a selected number of parameter/time indices. An instance of TransientSnapshotsAtIndices
is created by calling the function select_snapshots
GridapROMs.ParamDataStructures.TransientSnapshotsWithIC
— Typestruct TransientSnapshotsWithIC{T,N,D,I,R,A,B<:TransientSnapshots{T,N,D,I,R,A}} <: TransientSnapshots{T,N,D,I,R,A}
initial_data::A
snaps::B
end
Stores a TransientSnapshots
snaps
alongside a parametric initial condition initial_data
GridapROMs.ParamDataStructures.TrivialParamArray
— Typestruct TrivialParamArray{T<:Number,N,A<:AbstractArray{T,N}} <: ParamArray{T,N}
data::A
plength::Int
end
Wrapper for a non-parametric array data
that we wish assumed a parametric length
GridapROMs.ParamDataStructures.TrivialParamBlock
— Typestruct TrivialParamBlock{A} <: ParamBlock{A}
data::A
plength::Int
end
Wrapper for a non-paramentric quantity data
that we wish assumed a parametric length plength
GridapROMs.ParamDataStructures.UniformSampling
— Typestruct UniformSampling <: SamplingStyle end
Sampling according to a uniform distribution
GridapROMs.ParamDataStructures.change_mode
— Methodchange_mode(s::ModeTransientSnapshots) -> ModeTransientSnapshots
Returns the snapshots obtained by opposing the mode of s
. The result is a subtype of AbstractMatrix
with entries equal to those of s
, but with swapped spatial and temporal axes
GridapROMs.ParamDataStructures.eltype2
— Methodeltype2(a) -> Type
Returns the eltype of eltype(a)
, i.e. it extracts the eltype of a parametric entry of a
GridapROMs.ParamDataStructures.find_param_length
— Methodfind_param_length(a...) -> Int
Returns the parametric length of all parametric quantities. An error is thrown if there are no parametric quantities or if at least two quantities have different parametric length
GridapROMs.ParamDataStructures.get_all_data
— Methodget_all_data(A::ParamArray) -> AbstractArray{<:Any}
Returns all the entries stored in A
, assuming A
stores its entries consecutively
GridapROMs.ParamDataStructures.get_at_time
— Functionget_at_time(r::GenericTransientRealization,time) -> TransientRealizationAt
Returns a TransientRealizationAt
from a GenericTransientRealization
at a time instant specified by time
GridapROMs.ParamDataStructures.get_final_time
— Methodget_final_time(r::GenericTransientRealization) -> Real
GridapROMs.ParamDataStructures.get_indexed_data
— Methodget_indexed_data(s::Snapshots) -> AbstractArray
Returns the data in s
reindexed according to the indexing strategy provided in s
.
This function is not lazy, and should be used with parsimony
GridapROMs.ParamDataStructures.get_initial_time
— Methodget_initial_time(r::GenericTransientRealization) -> Real
GridapROMs.ParamDataStructures.get_param_data
— Methodget_param_data(a) -> Any
Returns the parametric data of a
, usually in the form of a AbstractVector
or a NTuple
GridapROMs.ParamDataStructures.get_param_entry
— Methodget_param_entry(A::AbstractParamArray{T},i...) where T -> Vector{eltype(T)}
Returns a vector of the entries of A
at index i
, for every parameter. The length of the output is equals to param_length(A)
GridapROMs.ParamDataStructures.get_params
— Methodget_params(r::AbstractRealization) -> Realization
GridapROMs.ParamDataStructures.get_realization
— Methodget_realization(s::AbstractSnapshots) -> AbstractRealization
Returns the realizations associated to the snapshots s
GridapROMs.ParamDataStructures.get_times
— Methodget_times(r::TransientRealization) -> Any
GridapROMs.ParamDataStructures.global_parameterize
— Methodglobal_parameterize(a,plength::Integer) -> AbstractParamArray
Returns a AbstractParamArray
with parametric length plength
from a
. This parameterization involves quantities defined at the global (or assembled) level. For local parameterizations, see the function local_parameterize
GridapROMs.ParamDataStructures.inneraxes
— Methodinneraxes(A::AbstractParamArray) -> Tuple{Vararg{Base.OneTo}}
Returns the axes of A
for a single parameter
GridapROMs.ParamDataStructures.innerlength
— Methodinnerlength(A::AbstractParamArray) -> Int
Returns the length of A
for a single parameter. Thus, the total entries of A
is equals to param_length(A)*innerlength(A)
GridapROMs.ParamDataStructures.lazy_parameterize
— Methodlazy_parameterize(a,plength::Integer) -> Any
Lazy version of parameterize
, does not allocate
GridapROMs.ParamDataStructures.local_parameterize
— Methodlocal_parameterize(a,plength::Integer) -> Any
Returns a quantity with parametric length plength
from a
. This parameterization involves quantities defined at the local (or cell) level. For global parameterizations, see the function global_parameterize
GridapROMs.ParamDataStructures.num_params
— Methodnum_params(r::AbstractRealization) -> Int
GridapROMs.ParamDataStructures.num_times
— Methodget_times(r::TransientRealization) -> Int
GridapROMs.ParamDataStructures.param_eachindex
— Methodparam_eachindex(a,i::Integer) -> Any
Returns the parametric range of a
1:param_length(a)
GridapROMs.ParamDataStructures.param_getindex
— Methodparam_getindex(a,i::Integer) -> Any
Returns the parametric entry of a
at the index i
∈ {1,...,param_length(a)
}
GridapROMs.ParamDataStructures.param_length
— Methodparam_length(a) -> Int
Returns the parametric length of a
GridapROMs.ParamDataStructures.param_setindex!
— Methodparam_setindex!(a,v,i::Integer) -> Any
Sets the parametric entry of a
to v
at the index i
∈ {1,...,param_length(a)
}
GridapROMs.ParamDataStructures.parameterize
— Methodparameterize(a,plength::Integer) -> Any
Returns a quantity with parametric length plength
from a
. When a
already possesses a parametric length, i.e. it is a parametrized quantity, it returns a
GridapROMs.ParamDataStructures.parameterize
— Methodparameterize(f::Function,r::Realization) -> ParamFunction
parameterize(f::Function,r::TransientRealization) -> TransientParamFunction
Method that parameterizes an input quantity by a realization r
GridapROMs.ParamDataStructures.realization
— Methodrealization(p::ParamSpace;nparams=1,sampling=get_sampling_style(p),kwargs...) -> Realization
realization(p::TransientParamSpace;nparams=1,sampling=get_sampling_style(p),kwargs...) -> TransientRealization
Extraction of a set of nparams
parameters from a given parametric space, by default according to the sampling strategy specified in p
.
GridapROMs.ParamDataStructures.select_snapshots
— Methodselect_snapshots(s::SteadySnapshots,prange) -> SnapshotsAtIndices
select_snapshots(s::TransientSnapshots,trange,prange) -> TransientSnapshotsAtIndices
Restricts the parametric range of s
to the indices prange
steady cases, to the indices trange
and prange
in transient cases, while leaving the spatial entries intact. The restriction operation is lazy.
GridapROMs.ParamDataStructures.shift!
— Methodshift!(r::TransientRealization,δ::Real) -> Nothing
In-place uniform shifting by a constant δ
of the temporal domain in the realization r
GridapROMs.ParamDataStructures.space_dofs
— Methodspace_dofs(s::Snapshots{T,N,D}) where {T,N,D} -> NTuple{D,Integer}
Returns the spatial size of the snapshots