Reference
Constants
Muscade.noFB
— ConstantnoFB
A constant, used by elements' residual
or lagrangian
as their 3rd output if they do provide any feedback to the solver (for example, on the reduction of the barrier parameter in interior point method).
Example: return L,noFB
See also: noFB
Muscade.ϵ
— Constantϵ (\epsilon)
an alias for Base.eps(𝕣)
.
Muscade.∞
— Constant∞ (\infty)
an alias for Base.inf
.
Types
Muscade.AbstractElement
— TypeAbstractElement
An abstract data type. An element type MyElement
must be declared as a subtype of AbstractElement
.
MyELement
must provide a constructor with interface
`eleobj = MyElement(nod::Vector{Node}; kwargs...)`
See also: coord
, Node
, lagrangian
Muscade.DirectXUA
— TypeDirectXUA{OX,OU,IA}
A non-linear direct solver for optimisation FEM.
An analysis is carried out by a call with the following syntax:
initialstate = initialize!(model)
stateXUA = solve(DirectXUA{OX,OU,IA};initialstate,time=0:1.:5)
The solver does not yet support interior point methods.
Parameters
OX
0 for static analysis 1 for first order problems in time (viscosity, friction, measurement of velocity) 2 for second order problems in time (inertia, measurement of acceleration)OU
0 for white noise prior to the unknown load process 2 otherwiseIA
0 for XU problems (variables of class A will be unchanged) 1 for XUA problems
Named arguments
dbg=(;)
a named tuple to trace the call tree (for debugging).verbose=true
set to false to suppress printed output (for testing).silenterror=false
set to true to suppress print out of error (for testing) .initialstate
aState
.time
anAbstractRange
of times at which to compute the steps. Example: 0:0.1:5.maxiter=50
maximum number of Newton-Raphson iterations.maxΔλ=1e-5
convergence criteria: a norm of the scaledΛ
increment.maxΔx=1e-5
convergence criteria: a norm of the scaledX
increment.maxΔu=1e-5
convergence criteria: a norm of the scaledU
increment.maxΔa=1e-5
convergence criteria: a norm of the scaledA
increment.saveiter=false
set to true so that the outputstate
is a vector (over the iterations) of vectors (over the steps) ofState
s of the model (for debugging non-convergence).
Setting the following flags to true
will improve the sparsity of the system. But setting a flag to true
when the condition isn't met causes the Hessian to be wrong, which is detrimental for convergence.
Xwhite=false
true
if response measurement error is a white noise process.XUindep=false
true
if response measurement error is independant ofU
UAindep=false
true
ifU
is independant ofA
XAindep=false
true
if response measurement error is independant ofA
Output
A vector of length equal to that of time
containing the state of the optimized model at each of these steps.
Muscade.DofConstraint
— TypeDofConstraint{λclass,Nx,Nu,Na,xinod,xfield,uinod,ufield,ainod,
afield,λinod,λfield,Tg,Tmode} <: AbstractElement
An element to apply physical/optimisation equality/inequality constraints on dofs.
The constraints are holonomic, i.e. they apply to the values, not the time derivatives, of the involved dofs. This element is very general but not very user-friendly to construct, factory functions are provided for better useability. The sign convention is that the gap g≥0
and the Lagrange multiplier λ≥0
.
This element can generate three classes of constraints, depending on the input argument λclass
.
λclass=:X
Physical constraint. In mechanics, the Lagrange multiplier dof is a generalized force, dual of the gap. The gap function must be of the formgap(x,t,gargs...)
.λclass=:U
Time varying optimisation constraint. For example: findA
-parameters so that at all times, the response does not exceed a given criteria. The gap function must be of the formgap(x,u,a,t,gargs...)
.λclass=:A
Time invariant optimisation constraint. For example: findA
-parameters such thatA[1]+A[2]=gargs.somevalue
. The gap function must be of the formgap(a,gargs...)
.
Named arguments to the constructor
xinod::NTuple{Nx,𝕫}=()
For each X-dof to be constrained, its element-node number.xfield::NTuple{Nx,Symbol}=()
For each X-dof to be constrained, its field.uinod::NTuple{Nu,𝕫}=()
For each U-dof to be constrained, its element-node number.ufield::NTuple{Nu,Symbol}=()
For each U-dof to be constrained, its field.ainod::NTuple{Na,𝕫}=()
For each A-dof to be constrained, its element-node number.afield::NTuple{Na,Symbol}=()
For each A-dof to be constrained, its field.λinod::𝕫
The element-node number of the Lagrange multiplier.λclass::Symbol
The class (:X
,:U
or:A
) of the Lagrange multiplier. See the explanation above for classes of constraintsλfield::Symbol
The field of the Lagrange multiplier.gap::Function
The gap function.gargs::NTuple
Additional inputs to the gap function.mode::Function
wheremode(t::ℝ) -> Symbol
, with value:equal
,:positive
or:off
at any time. An:off
constraint will set the Lagrange multiplier to zero.
Example
using Muscade
model = Model(:TestModel)
n1 = addnode!(model,𝕣[0])
e1 = addelement!(model,DofConstraint,[n1],xinod=(1,),xfield=(:t1,),
λinod=1, λclass=:X, λfield=:λ1,gap=(x,t)->x[1]+.1,
mode=positive)
e2 = addelement!(model,QuickFix ,[n1],inod=(1,),field=(:t1,),
res=(x,u,a,t)->0.4x.+.08+.5x.^2)
initialstate = initialize!(model)
setdof!(initialstate,1.;field=:λ1)
state = solve(SweepX{0};initialstate,time=[0.],verbose=false)
X = state[1].X[1]
See also: Hold
, ElementConstraint
, off
, equal
, positive
Muscade.DofCost
— TypeDofCost{Class,Nx,Nu,Na,xinod,xfield,uinod,ufield,ainod,
afield,Tcost,Tcostargs} <: AbstractElement
An element to apply costs on combinations of dofs.
Named arguments to the constructor
xinod::NTuple{Nx,𝕫}=()
For each X-dof to entercost
, its element-node number.xfield::NTuple{Nx,Symbol}=()
For each X-dof to entercost
, its field.uinod::NTuple{Nu,𝕫}=()
For each U-dof to entercost
, its element-node number.ufield::NTuple{Nu,Symbol}=()
For each U-dof to entercost
, its field.ainod::NTuple{Na,𝕫}=()
For each A-dof to entercost
, its element-node number.afield::NTuple{Na,Symbol}=()
For each A-dof to entercost
, its field.class:Symbol
:A
for cost on A-dofs only,:I
("instant") otherwise.cost::Function
ifclass==:I
,cost(X,U,A,t,costargs...)→ℝ
ifclass==:A
,cost(A,costargs...)→ℝ
X
andU
are tuples (derivates of dofs...), and∂0(X)
,∂1(X)
,∂2(X)
must be used bycost
to access the value and derivatives ofX
(resp.U
)[costargs::NTuple=() or NamedTuple] of additional arguments passed to
cost``
Requestable internal variables
cost
, the value of the cost.
Example
ele1 = addelement!(model,DofCost,[nod1],xinod=(1,),field=(:tx1,),
class=:I,cost=(X,U,A,t;X0)->(X[1]-X0)^2,costargs=(;X0=0.27)
See also: SingleDofCost
, ElementCost
, addelement!
Muscade.DofLoad
— TypeDofLoad{Tvalue,Field} <: AbstractElement
An element to apply a loading term to a single X-dof.
Named arguments to the constructor
field::Symbol
.value::Function
, wherevalue(t::ℝ) → ℝ
.
Requestable internal variables
F
, the value of the load.
Examples
using Muscade
model = Model(:TestModel)
node = addnode!(model,𝕣[0,0])
e = addelement!(model,DofLoad,[node];field=:tx,value=t->3t-1)
Muscade.ElementConstraint
— TypeElementConstraint{Teleobj,λinod,λfield,Nu,Treq,Tg,Tgargs,Tmode} <: AbstractElement
An element to apply optimisation equality/inequality constraints on the element-results of another "target" element. The target element must not be added separately to the model. Instead, the ElementType
, and the named arguments to the target element are provided as input to the ElementConstraint
constructor.
This element generates a time varying optimisation constraint. For example: find A
-parameters so that at all times, the element-result von-Mises stress does not exceed a given value.
The Lagrangian multiplier introduced by this optimisation constraint is of class :U
Named arguments to the constructor
λinod::𝕫
The element-node number of the Lagrange multiplier.λfield::Symbol
The field of the Lagrange multiplier.req
A request for element-results to be extracted from the target element, see@request
. The request is formulated as if adressed directly to the target element.gₛ::𝕣=1.
A scale for the gap.λₛ::𝕣=1.
A scale for the Lagrange multiplier.gap
A gap functiongap(eleres,X,U,A,t,gargs...)→ℝ
.eleres
is the output of the above-mentionned request to the target element.X
andU
are tuples (derivates of dofs...), and∂0(X)
,∂1(X)
,∂2(X)
must be used bycost
to access the value and derivatives ofX
(resp.U
).X
,U
andA
are the degrees of freedom of the elementElementType
.gargs::NTuple
Additional inputs to the gap function.mode::Function
wheremode(t::ℝ) -> Symbol
, with value:equal
,:positive
or:off
at any time. An:off
constraint will set the Lagrange multiplier to zero.ElementType
The named of the constructor for the relevant elementelementkwargs
A named tuple containing the named arguments of theElementType
constructor.
Requestable internal variables
Not to be confused with the req
provided as input to addelement!
when adding an ElementConstraint
, one can, after the analysis, request results from ElementConstraint
λ
The constraints Lagrange multipliergap
The constraints gap functioneleres(...)
where...
is the list of requestables wanted from the target element. The "prefix"eleres
is there to prevent possible confusion with variables requestable fromElementConstraint
. For example@request gap
would extract the value of theElementConstraint
's functiongap
, while@request eleres(gap)
refers to the value of a variable calledgap
in the target element.
Example
@once gap(eleres,X,U,A,t) = eleres.Fh^2
ele1 = addelement!(model,ElementCoonstraint,[nod1];req=@request(Fh),
gap,λinod=1,λfield=:λ,mode=equal,
ElementType=AnchorLine,
elementkwargs=(Δxₘtop=[5.,0,0], xₘbot=[250.,0],L=290., buoyancy=-5e3))
See also: Hold
, DofConstraint
, off
, equal
, positive
, @request
, @once
Muscade.ElementCost
— TypeElementCost{Teleobj,Treq,Tcost,Tcostargs} <: AbstractElement
An element to apply costs on another "target" element's dofs and element-results. The target element must not be added separately to the model. Instead, the ElementType
, and the named arguments to the target element are provided as input to the ElementCost
constructor.
Named arguments to the constructor
req
A request for element-results to be extracted from the target element, see@request
. The request is formulated as if adressed directly to the target element.cost
A cost functioncost(eleres,X,U,A,t,costargs...)→ℝ
.eleres
is the output of the above-mentionned request to the target element.X
andU
are tuples (derivates of dofs...), and∂0(X)
,∂1(X)
,∂2(X)
must be used bycost
to access the value and derivatives ofX
(resp.U
).X
,U
andA
are the degrees of freedom of the elementElementType
.[costargs::NTuple=() or NamedTuple]
Additional arguments passed tocost
.ElementType
The named of the constructor for the relevant element.elementkwargs
A named tuple containing the named arguments of theElementType
constructor.
Requestable internal variables
Not to be confused with the req
provided as input to addelement!
when adding an ElementCost
, one can, after the analysis, request results from ElementConstraint
cost
The value ofcost(eleres,X,U,A,t,costargs...)
.eleres(...)
where...
is the list of requestables wanted from the target element. The "prefix"eleres
is there to prevent possible confusion with variables requestable fromElementConstraint
. For example@request cost
would extract the value of theElementConstraint
's functioncost
, while@request eleres(cost)
refers to the value of a variable calledcost
in the target element.
Example
@once cost(eleres,X,U,A,t;Fh0) = (eleres.Fh-Fh0)^2
ele1 = addelement!(model,ElementCost,[nod1];req=@request(Fh),
cost=cost,
costargs = (;Fh0=0.27)
ElementType=AnchorLine,
elementkwargs=(Λₘtop=[5.,0,0], xₘbot=[250.,0], L=290., buoyancy=-5e3))
See also: SingleDofCost
, DofCost
, @request
, @once
Muscade.Hold
— TypeHold <: AbstractElement
An element to set a single X-dof to zero.
Named arguments to the constructor
field::Symbol
. The field of the X-dof to constraint.λfield::Symbol=Symbol(:λ,field)
. The field of the Lagrange multiplier.
Example
using Muscade
model = Model(:TestModel)
node = addnode!(model,𝕣[0,0])
e = addelement!(model,Hold,[node];field=:tx)
See also: DofConstraint
, DofLoad
, DofCost
Muscade.Model
— Typemodel = Model([ID=:my_model])
Construct a blank model
, which will be mutated to create a finite element [optimization] problem.
See also: addnode!
, addelement!
, describe
, solve
Muscade.Monitor
— TypeMonitor <: AbstractElement
Debbuging tool: An element for for monitoring inputs to and outputs from another element, during an analysis.
Instead of adding the element to be monitored directly into the model, add this element with the element to be monitored as argument.
Inputs and outs get @show'n to screen
Named arguments to the constructor
ElementType
The the type of element to be monitored-trigger
A function that takesdbg
as an input and returns a boolean (true
) to printout.elementkwargs
aNamedTuple
containing the named arguments of theElementType
constructor.
Muscade.Node
— TypeNode
The eltype of vectors handed by Muscade as first argument to element constructors.
Example: function SingleDofCost(nod::Vector{Node};class::Symbol, ... )
See also: coord
Muscade.QuickFix
— TypeQuickFix <: AbstractElement
An element for creating simple elements with "one line" of code. Elements thus created have several limitations:
- physical elements with only X-dofs.
- only
R
can be espied.
The element is intended for testing. Muscade-based applications should not include this in their API.
Named arguments to the constructor
inod::NTuple{Nx,𝕫}
. The element-node numbers of the X-dofs.field::NTuple{Nx,Symbol}
. The fields of the X-dofs.res::Function
, whereres(X::ℝ1,X′::ℝ1,X″::ℝ1,t::ℝ) → ℝ1
, the residual.
Examples
A one-dimensional linear elastic spring with stiffness 2.
using Muscade
model = Model(:TestModel)
node1 = addnode!(model,𝕣[0])
node2 = addnode!(model,𝕣[1])
e = addelement!(model,QuickFix,[node1,node2];inod=(1,2),field=(:tx1,:tx1),
res=(X,X′,X″,t)->Svector{2}(2*(X[1]-X[2]),2*(X[2]-X[1])) )
Muscade.SingleDofCost
— TypeSingleDofCost{Derivative,Class,Field,Tcost} <: AbstractElement
An element with a single node, for adding a cost to a given dof.
Named arguments to the constructor
class::Symbol
, either:X
,:U
or:A
.field::Symbol
.cost::Function
, wherecost(x::ℝ,t::ℝ[,costargs...]) → ℝ
ifclass
is:X
or:U
, andcost(x::ℝ, [,costargs...]) → ℝ
ifclass
is:A
.
[costargs::NTuple=() or NamedTuple] of additional arguments passed to
cost``derivative::Int=0
0, 1 or 2 - which time derivative of the dof enters the cost.
Requestable internal variables
cost
, the value of the cost.
Example
using Muscade
model = Model(:TestModel)
node = addnode!(model,𝕣[0,0])
e = addelement!(model,SingleDofCost,[node];class=:X,field=:tx,
costargs=(3.,),cost=(x,t,three)->(x/three)^2)
See also: DofCost
, ElementCost
Muscade.SingleUdof
— TypeSingleUdof{XField,Ufield,Tcost} <: AbstractElement
An element that creates a Udof, and associates a cost to its value. The value of the Udof is applied as a load to a Xdof on the same node.
Named arguments to the constructor
Xfield::Symbol
.Ufield::Symbol
.cost::Function
, wherecost(u::ℝ,t::ℝ[,costargs...]) → ℝ
[costargs::NTuple=() or NamedTuple] of additional arguments passed to
cost``
Requestable internal variables
cost
, the value of the cost.
Example
using Muscade
model = Model(:TestModel)
node = addnode!(model,𝕣[0,0])
e = addelement!(model,SingleUdof,[node];Xfield=:tx,Ufield=:utx,
costargs=(3.,),cost=(x,t,three)->(x/three)^2)
See also: DofCost
, ElementCost
Muscade.SweepX
— TypeSweepX{ORDER}
A non-linear, time domain solver, that solves the problem time-step by time-step. Only the X
-dofs of the model are solved for, while U
-dofs and A
-dofs are unchanged.
SweepX{0}
is Newton-Raphson, with feasibility line-search, to handle inequality constraints.SweepX{1}
is implicit Euler, with feasibility line-search.SweepX{2}
is Newmark-β, with Newton-Raphson iterations and feasibility line search
IMPORTANT NOTE: Muscade does not allow elements to have state variables, for example, plastic strain, or shear-free position for dry friction. Where the element implements such physics, this is implemented by introducing the state as a degree of freedom of the element, and solving for its evolution, even in a static problem, requires the use of ORDER≥1
An analysis is carried out by a call with the following syntax:
initialstate = initialize!(model)
setdof!(initialstate,1.;class=:U,field=:λcsr)
states = solve(SweepX{2};initialstate=initialstate,time=0:10)
Named arguments to solve
:
dbg=(;)
a named tuple to trace the call tree (for debugging)verbose=true
set to false to suppress printed output (for testing)silenterror=false
set to true to suppress print out of error (for testing)initialstate
aState
, obtain fromìnitialize!
orSweepX
.time
maximum number of Newton-Raphson iterationsβ=1/4
,γ=1/2
parameters to the Newmark-β algorithm. Dummy ifORDER<2
maxiter=50
maximum number of equilibrium iterations at each step.maxΔx=1e-5
convergence criteria: norm ofX
. DmaxLλ=∞
convergence criteria: norm of the residual.saveiter=false
set to true so that outputstates
contains the state at the iteration of the last step analysed. Useful to study a step that fails to converge.maxLineIter=50
Maximum number of iteration in the feasibility line search. set to 0 to skip the line search (not recommended for models with inequality constraints).sfac=0.5
Parameter in the line search for a feasible point. If a tentative result is not feasible, backtrack by a factorsfac
. If still not feasible, backtrack what is left by a factorsfac
, and so forth, up tomaxLineIter
times.γfac=0.5
Parameter for feasibility. For an inequality constraintg(X)
with reaction forceλ
, requireg(X)*λ==γ
, and multiplyγ *= γfac
at each iteration.
Output
A vector of length equal to that of the named input argument time
containing the states at the time steps.
See also: solve
, initialize!
, findlastassigned
, studysingular
Muscade.acceleration
— MethodP = constants(X,U,A,t)
N = length(X)
x = Muscade.motion{P,N}(X)
E = f(x)
̈ε = Muscade.acceleration{P,N}(E)
Extract the velocity or rate from a variable that is a function of the output of Muscade.motion
.
Muscade.default
— Typedefault{:fieldname}(namedtuple,defval)
attempt to get a field fieldname
from a NamedTuple
. If namedtuple
does not have such a field - or is not a NamedTuple
, return defval
.
Muscade.motion
— MethodP = constants(X,U,A,t)
x = Muscade.motion{P}(X)
Transform a NTuple
of SVector
s, for example the vector X
provided as an input to residual
or Lagrangian
into a SVector
of ∂ℝ
. This can be used by an element to compute time derivatives, for example Euler, Coriolis and centrifugal accelerations, or strain rates.
Some principles of safe automatic differentiation must be adhered to:
- the function that uses
Muscade.motion
must also 'unpack' : no variable that is touched by the output ofMuscade.motion
must be returned by the function without having been unpacked byMuscade.position
,Muscade.velocity
orMuscade.acceleration
. - The precendence
P
must be calculated usingconstants
with all variables that are input to the function and may be differentiated. - If other levels of automatic differentiation are introduced within the function, unpack in reverse order of packing.
See Muscade.position
, Muscade.velocity
, Muscade.acceleration
Muscade.position
— MethodP = constants(X,U,A,t)
N = length(X)
x = Muscade.motion{P,ND}(X)
E = f(x)
ε = Muscade.position{P,N}(E)
Extract the position from a variable that is a function of the output of Muscade.motion
.
Muscade.value
— Methody = value{P}(Y)
Extract the value of an automatic differentiation object, or SArray
of such objects.
Muscade.value_∂
— Methody,yₓ = value_∂{P,N}(Y)
y,y′ = value_∂{P }(Y)
Get value and derivative in one operation.
Muscade.variate
— MethodX = variate{P,N}(x)
where typeof(x)<:SVector{N}
, create a SVector
of automatic differentiation objects of precedence P
.
X = variate{P}(x)
where typeof(x)<:Real
, create an object of precedence P
.
Muscade.velocity
— MethodP = constants(X,U,A,t)
N = length(X)
x = Muscade.motion{P,N}(X)
E = f(x)
̇ε = Muscade.velocity{P,N}(E)
Extract the velocity or rate from a variable that is a function of the output of Muscade.motion
.
Muscade.δ
— MethodX = δ{P,N,R}()
create a SVector
of automatic differentiation objects of precedence P
and value zero
.
X = δ{P}()
Create automatic differentiation object of precedence P
and value zero
.
Muscade.ℕ
— Typeℕ (\bbN)
an alias for UInt64
. For use in dispatching. ℕ1
... ℕ4
are AbstractArrays
of dimensions 1 to 4.
Muscade.ℝ
— Typeℝ (\bbR)
an alias for abstract type Real
. For use in dispatching. ℝ1
... ℝ4
are AbstractArrays
of dimensions 1 to 4. ℝ11
is an AbstractVector
of AbstractVector
.
Muscade.ℤ
— Typeℤ (\bbZ)
an alias for abstract type Integer
. For use in dispatching. ℤ1
... ℤ4
are AbstractArrays
of dimensions 1 to 4. ℤ11
is an AbstractVector
of AbstractVector
.
Muscade.∂
— Methodyₓ = ∂{P,N}(Y)
Extract the gradient of an automatic differentiation object. If Y
is a SArray
, the index of the partial derivative is appended to the indices of Y
.
y′ = ∂{P}(Y)
Extract the derivative of an automatic differentiation object (or SArray
of such), where the variation was created by the syntax variate{P}
.
Muscade.𝔹
— Type𝔹 (\bbB)
an alias for Bool
. For use in dispatching. 𝔹1
... 𝔹4
are AbstractArrays
of dimensions 1 to 4.
Muscade.𝕓
— Type𝕓 (\bbb)
an alias for Bool
. For use in struct
definitions. 𝕓1
... 𝕓4
are Arrays
of dimensions 1 to 4.
Muscade.𝕟
— Type𝕟 (\bbn)
an alias for UInt64
. For use in struct
definitions. 𝕟1
... 𝕟4
are Arrays
of dimensions 1 to 4.
Muscade.𝕣
— Type𝕣 (\bbr)
an alias for Float64
. For use in struct
definitions. 𝕣1
... 𝕣4
are Arrays
of dimensions 1 to 4. 𝕣11
is a Vector
of Vector
.
Muscade.𝕫
— Type𝕫 (\bbz)
an alias for Int64
. For use in struct
definitions. 𝕫1
... 𝕫4
are Arrays
of dimensions 1 to 4. 𝕫11
is a Vector
of Vector
.
Functions
Muscade.:∘₁
— Methodc = a∘₁b
Compute the single-dot product of two arrays, so that cᵢⱼ=Σₖ aᵢₖ bₖⱼ
where i
and j
can be multiple indices.
Muscade.:∘₂
— Methodc = a∘₂b
Compute the double-dot product of two arrays, so that cᵢⱼ=Σₖₗ aᵢₖₗ bₖₗⱼ
where i
and j
can be multiple indices.
Muscade.:⊗
— Methodc = a⊗b
Compute the exterior product of two arrays, so that cᵢⱼ=aᵢ bⱼ
where i
and j
can be multiple indices.
Muscade.VALUE
— Method@show VALUE(Y)
Completely strip Y
of partial derivatives. Use only for debugging purpose.
Muscade.addelement!
— Methodeleid = addelement!(model,ElType,nodid;kwargs...)
Add one or several elements to model
, connecting them to the nodes specified by nodid
.
If nodid
is an AbstractVector
: add a single element to the model. eleid
is then a single element identifier.
If nodid
is an AbstractMatrix
: add multiple elements to the model. Each row of nodid
identifies the node of a single element. eleid
is then a vector of element identifiers.
For each element, addelement!
will call eleobj = ElType(nodes;kwargs...)
where nodes
is a vector of nodes of the element.
Muscade.addin!
— Methodaddin!(asm,global,block,ibr,ibc,factor=1.)
Add a sparse block
into a large out
sparse matrix, at block-row and -column ibr
and ibc
. Use prepare
to allocate memory for global
and build the assembler asm
.
Muscade.addin!
— Methodaddin!(asm,outvec,blockvec,ibr)
Add a full block
vector into a large outvec
full vector. at block-row ibr
. Use prepare
to create asm
.
See also: prepare
Muscade.addnode!
— Methodnodid = addnode!(model,coord)
If coord
is an AbstractVector
of Real
: add a single node to the model. Muscade does not prescribe what coordinate system to use. Muscade will handle to each element the coord
of the nodes of the element, and the element constructor must be able to make sense of it. nodid
is a node identifier, that is used as input to addelement!
.
If coord
is an AbstractMatrix
, its rows are treated as vectors of coordinates. nodid
is then a vector of node identifiers.
See also: addelement!
, coord
, describe
Muscade.constants
— MethodP = constants(a,b,c)
Generate a precedence P
that is higher than the precedence of the arguments.
Muscade.coord
— Methodc = coord(node)
Used by element constructors to obtain the coordinates of a vector of Nodes handed by Muscade to the constructor. c
is accessed as
c[inod][icoord]
where inod
is the element-node number and icoord
an index into a vector of coordinates.
Note that c[inod]
points at the same memory as nod[inod].coord
: do not mutate c[inod]
!
See also: addnode!
, addelement!
, describe
, solve
Muscade.describe
— Methoddescribe(model,spec)
Print out information about model
. spec
can be
- an
EleID
to describe an element, - a
DofID
to describe a dof. - a
NodID
to describe a node, :doftyp
to obtain a list of doftypes,:dof
to obtain a list of dofs or:eletyp
for a list of element types.
See also: addelement!
, addnode!
Muscade.describe
— Methoddescribe(state;class=:all)
Provide a description of the dofs stored in state
. class
can be either :all
, :Λ
, :ΛX
, :X
, :U
, :A
or :scale
See also: solve
Muscade.doflist
— MethodMuscade.doflist(::Type{E<:AbstractElement})
Elements must overload Muscade's doflist
function. The method must take the element type as only input, and return a NamedTuple
with fieldnames inod
,class
and field
. The tuple-fields are NTuple
s of the same length. For example
Muscade.doflist( ::Type{<:Turbine}) = (inod =(1 ,1 ,2 ,2 ),
class=(:X ,:X ,:A ,:A ),
field=(:tx1,:tx2,:Δseadrag,:Δskydrag))
In Λ
, X
, U
and A
handed by Muscade to residual
or lagrangian
, the dofs in the vectors will follow the order in the doflist. Element developers are free to number their dofs by node, by field, or in any other way.
See also: lagrangian
, residual
∂
Muscade.draw
— MethodMuscade.draw(ElementType,axe,o, Λ,X,U,A,t,SP,dbg;kwargs...)
Inputs are:
ElementType
, the method must dispatch on thisDataType
.axe
, aGLMakie
axeo
aAbstractVector
of element objects, of lengthnel
.Λ
a matrix of size(nXdof,nel)
X
aNTuple
(over the derivatives) of matrices of size(nXdof,nel)
U
aNTuple
(over the derivatives) of matrices of size(nUdof,nel)
A
a matrix of size(nAdof,nel)
t
timeSP
solver parametersdbg
debuging informationkwargs
aNamedTuple
containing the keyword arguments provided by the user. Seedefault
.
See also: lagrangian
, residual
, doflist
Muscade.draw
— Methoddraw(state[,els];kwargs...)
state
a single State
. els
specifies which elements to draw and can be either
- a vector of
EleID
s (obtained fromaddelement!
`), all corresponding to the same concrete element type - a concrete element type (see
eletyp
). - omitted: all the element of the model are drawn.
kwargs...
is any additional key words arguments that will be passed to the draw
method of each element, for example to specify colors, etc. See the elements' documentation.
See also: getdof
, @request
, @espy
, addelement!
, solve
Muscade.eletyp
— Methodet = eletyp(model)
Return a vector of the concrete types of elements in the model.
Muscade.equal
— Methodequal(t) → :equal
A function which for any value t
returns the symbol equal
. Usefull for specifying the keyword argument mode=equal
in adding an element of type `DofConstraint
to a Model
.
See also: DofConstraint
, ElementConstraint
, off
, positive
Muscade.findlastassigned
— Methodilast = findlastassigned(state)
Find the index ilast
of the element before the first non assigment element in a vector state
.
In multistep analyses, solve
returns a vector state
of length equal to the number of steps requested by the user. If the analysis is aborted, solve
still returns any available results at the begining of state
, and the vector state[1:ilast]
is fully assigned.
See also: solve
Muscade.getdof
— Methoddofres = getdof(state;[class=:X],field=:somefield,nodID=[nodids...],[orders=[0]])
Obtain the value of dofs of the same class and field, at various nodes and for various states.
If state
is a vector, the output dofres
has size (ndof,nder+1,nstate)
. If state
is a scalar, the output dofres
has size (ndof,nder+1)
.
Muscade.getndof
— Methodgetndof(model|Element)
getndof(model|Element,class)
getndof(model|Element,(class1,class2,[,...]))
where class
can be any of :X
, :U
, :A
: get the number of dofs of the specified dof-classes (default: all classes) for the variable model
or the type Element
.
See also: describe
Muscade.getresult
— Methodeleres = getresult(state,req,els)
Obtain an array of nested NamedTuples
and NTuples
of element results. req
is a request defined using @request
. state
a vector of State
s or a single State
. els
can be either
- a vector of
EleID
s (obtained fromaddelement!
) all corresponding to the same concrete element type - a concrete element type (see
eletyp
).
If state
is a vector, the output dofres
has size (nele,nstate)
. If state
is a scalar, the output dofres
has size (nele)
.
See also: getdof
, @request
, @espy
, addelement!
, solve
, eletyp
Muscade.getsomedofs
— Methodrotations = getsomedofs(X,[3,6])
Used by elements' residual
or lagrangian
to some degrees of freedom, and their time derivatives, from the variables X
and U
.
Muscade.gradient
— MethodL,Lλ,Lx,Lu,La = gradient(eleobj,Λ,X,U,A,t,SP,dbg)
Compute the Lagrangian, its gradients, and the memory of an element. For element debugging and testing.
See also: residual
,lagrangian
Muscade.initialize!
— Methodinitialstate = initialize!(model)
Return an initial State
for the model with all dofs set to zero
Modifying a model (invoquing addnode!
and addelement!
after initialize!
will result in an error)
Optional keyword arguments: nΛder=1
, nXder=1
, nUder=1
to specify the number of "derivatives" to store in the State
. note that "n⋅der==order+1
", that is, for a dynamic problem (with accelerations), nXder=3
so that order==2
. Setting n⋅der
is only required if setdof!
will be used. A dynamic solver handles a "static" initial state perfectly well. time=-∞
the time associated to the initial state. Note that for example setting time=0.
, and then calling a solver with a first time step also at 0.
causes an error.
See also: setdof!
, Model
, addnode!
, addelement!
, solve
Muscade.lagrangian
— Method@espy function Muscade.lagrangian(eleobj::MyElement,Λ,X,U,A,t,SP,dbg)
...
return L,FB
end
Inputs
eleobj
an element objectΛ
aSVector{nXdof,R} where{R<:Real}
, Lagrange multipliers (akaδX
virtual displacements).X
aNTuple
ofSVector{nXdof,R} where{R<:Real}
, containing the Xdofs and, depending on the solver, their time derivatives. Usex=∂0(X)
,v=∂1(X)
anda=∂2(X)
to safely obtain vectors of zeros where the solver leaves time derivatives undefined.U
aNTuple
ofSVector{nUdof,R} where{R<:Real}
, containing the Udofs and, depending on the solver, their time derivatives. Useu=∂0(U)
,̇u=∂1(U)
and̈u=∂2(U)
to safely obtain vectors of zeros where the solver leaves time derivatives undefined.A
aSVector{nAdof,R} where{R<:Real}
.t
a `Real
containing the time.SP
solver parameters (for example: the barrier parameterγ
for interior point methods).dbg
aNamedTuple
to be used only for debugging purposes.
Outputs
L
the lagrangianFB
feedback from the element to the solver (for example: canγ
be reduced?). ReturnnoFB
of the element has no feedback to provide.
Muscade.merge
— Methodreq = Muscade.merge(o.req)
Elements like ElementCost
and ElementConstraint
use requests to apply a cost or a constraint to requestables from another "target" element. These cost/constraint elements must be coded carefully so that getresult
can be used to extracted both requestable internal results from the cost/constraint and from the target element.
merge
(which is not exported by Muscade) is used to merge the requests for the request needed to enforce a cost or constraint, and the user's request for element to be obtained from the analysis. The call to merge
, to be inserted in the code of lagrange
for the cost/constraint element will be modified by @espy
.
See the code of ElementCost
's constructor and lagrange
method for an example.
See also: ElementCost
,@request
,getresult
Muscade.muscadeerror
— Methodmuscadeerror([[dbg,]msg])
Throw a MuscadeException
, where
dbg
is aNamedTuple
that contains "location information"
(for example: solver, step, iteration, element, quadrature point) that will be displayed with the error message.
msg
is aString
describing the problem.
Muscade.off
— Methodoff(t) → :off
A function which for any value t
returns the symbol off
. Usefull for specifying the keyword argument mode=off
in adding an element of type `DofConstraint
to a Model
.
See also: DofConstraint
, ElementConstraint
, equal
, positive
Muscade.positive
— Methodpositive(t) → :positive
A function which for any value t
returns the symbol positive
. Usefull for specifying the keyword argument mode=positive
in adding an element of type `DofConstraint
to a Model
.
See also: DofConstraint
, ElementConstraint
, off
, equal
Muscade.prepare
— Methodbigmat,bigmatasm,bigvecasm,bigvecdis = prepare(pattern)
Prepare for the assembly of sparse blocks into a large sparse matrix. bigmat
is allocated, with the correct sparsity structure, but its nzval
undef'ed. Where some blocks share the same sparsity structure, blocks
in pattern
can have ===
elements.
pattern
is a SparseMatrixCSC{<:SparseMatrixCSC}
, where empty blocks are structuraly zero
See also: addin!
Muscade.residual
— Method@espy function Muscade.residual(eleobj::MyElement,X,U,A,t,SP,dbg) ... return R,FB end
Inputs
eleobj
an element objectX
aNTuple
ofSVector{nXdof,R} where{R<:Real}
, containing the Xdofs and, depending on the solver, their time derivatives. Usex=∂0(X)
,v=∂1(X)
anda=∂2(X)
to safely obtain vectors of zeros where the solver leaves time derivatives undefined.U
aNTuple
ofSVector{nUdof,R} where{R<:Real}
, containing the Udofs and, depending on the solver, their time derivatives. Useu=∂0(U)
,̇u=∂1(U)
and̈u=∂2(U)
to safely obtain vectors of zeros where the solver leaves time derivatives undefined.A
aSVector{nAdof,R} where{R<:Real}
.t
a `Real
containing the time.SP
solver parameters (for example: the barrier parameterγ
for interior point methods).dbg
aNamedTuple
to be used only for debugging purposes.
Outputs
R
the residualFB
feedback from the element to the solver (for example: canγ
be reduced?). ReturnnoFB
of the element has no feedback to provide.
Muscade.setdof!
— Methodstate = setdof!(state,value ;[class=:X],field=:somefield, [order=0])
state = setdof!(state,value::Vector;[class=:X],field=:somefield,nodID=[nodids...],[order=0])
Set the value of dofs of the same class and field, at various nodes and for various states. There are two methods:
- A single
value
is applied to all relevant nodes in the model value
andnodID
are vectors of the same lengths, and each element invalue
is applied to the corresponding node.
setdof!
is peculiar in that it modifies its input state
variable, but must be used as a function. A call like stateout = setdof!(statein,value;class=:X,field=:somefield,order=1)
can turn out in two ways: If the state
already stores derivatives in X
to order 1, then statein
is mutated and statein===stateout
. Otherwise, statein
is unchanged, stateout
is a new object, sharing as much memory as possible with statein
. To avoid confusion, always use the syntax shown above.
Muscade.setscale!
— Methodsetscale!(model;scale=nothing,Λscale=nothing)
Provide scale value for each type (class and field) of dof in the model. This is usued to improve the conditioning of the incremental problems and for convergence criteria. scale
is a NamedTuple
with fieldnames within X
, U
and A
. Each field is itself a NamedTuple
with fieldnames being dof fields, and value being the expected order of magnitude.
For example scale = (X=(tx=10,rx=1),A=(drag=3.))
should be read as: X-dofs of field :tx
are expected to be of the order of magnitude of 10m in the solution, :rx
to be of the order of 1 radian, and A-dofs of field drag
of the order of 3. All other degrees of freedom are of the order of 1.
Determining scaling coefficients that improve the condition number of incremental matrices is a hard problem.
Λscale
is a scalar. The scale of a Λ-dof will be deemed to be the scale of the corresponding X-dof, times Λscale
.
See also: addnode!
, describe
, coord
, studyscale
Muscade.solve
— Methodsolve(Solver;dbg=(;),verbose=true,silenterror=false,kwargs...)
Execute an analysis using Solver
, and safeguard partial results in the case of error.
Named arguments
dbg=(;)
a named tuple to trace the call tree (for debugging)verbose=true
set to false to suppress printed output (for testing)silenterror=false
set to true to suppress print out of error (for testing)kwargs...
Further arguments passed on to the methodsolve
provided by the solver
This will call the method solve
provided by the solver with solve(Solver,pstate,verbose,(dbg...,solver=Symbol(Solver));kwargs...)
See also: SweepX
, DirectXUA
, initialize!
Muscade.studyscale
— Methodscale = studyscale(state;[verbose=false],[dbg=(;)])
Returns a named tuple of named tuples for scaling the model, accessed as scaled.myclass.myfield
, for example scale.X.tx1
.
Currently, the format of scale
is not identical to the input expected by setscale!
: work in progress
If verbose=true
, prints out a report of the analysis underlying the proposed scale
. The proposed scaling depends on the state
passed as input - as it is computed for a given incremental matrix.
See also: setscale!
Muscade.studysingular
— Methodmatrix = studysingular(state;SP,[iclasses=(Λ,:X,:U,:A)],[jclasses=iclasses],[verbose::𝕓=true],[dbg=(;)])
Generates an incremental matrix for state
(no time derivatives) corresponding to the classes required, and report on the null space of the matrix.
To do so, the incremental matrix is converted to full format, limiting the applicability to small models.
The function returns the incremental matrix.
Muscade.test_static_element
— Methodtest_static_element(eleobj;Λ,X,U,A,t=0,SP=nothing,verbose=true,dbg=(;))
Compute the Lagrangian, its gradients, and the memory of an element. For element debugging and testing.
See also: residual
,lagrangian
,gradient
Muscade.toggle
— Methodtoggle(condition,a,b)
Typestable equivalent of condition ? a : b
. Returns a value converted to promote_type(typeof(a),typeof(b))
Muscade.∂0
— Methodposition = ∂0(X)
Used by elements' residual
or lagrangian
to extract the zero-th order time derivative from the variables X
and U
.
See also: ∂1
,∂2
,getsomedofs
Muscade.∂1
— Methodvelocity = ∂1(X)
Used by elements' residual
or lagrangian
to extract the first order time derivative from the variables X
and U
. Where the solver does not provide this derivative (e.g. a static solver), the output is a vector of zeros.
See also: ∂0
,∂2
,getsomedofs
Muscade.∂2
— Methodposition = ∂2(X)
Used by elements' residual
or lagrangian
to extract the zero-th order time derivative from the variables X
and U
. Where the solver does not provide this derivative (e.g. a static solver), the output is a vector of zeros.
See also: ∂0
,∂1
,getsomedofs
Macros
Muscade.@espy
— Macro@espy function ... end
From an anotated function code, generate - "clean" code, in which the anotations have been deleted, and with the call syntax argout... = foo(argin...)
- "espying" code, with added input and ouput arguments argout...,res = foo(argin...,req)
where req
has been generated using @request
and res
is a nested structure of NamedTuple
s and NTuple
s containing the requested data.
The macro is not general: it is designed for residual
and lagrangian
, which for performance have to be programmed in "immutable" style: they must never mutate variables (this implies in particular, no adding into an array in a loop over Gauss points). So @espy
only supports the specific programming constructs needed in this style.
The following is an example of anotated code:
@espy function residual(x::Vector{R},y) where{R<:Real}
ngp=2
accum = ntuple(ngp) do igp
☼z = x[igp]+y[igp]
☼s,☼t = ☼material(z)
♢square = s^2
@named(s)
end
r = sum(i->accum[i].s,ngp)
return r,nothing,nothing
end
- The keyword
function
is preceded by the macro-call@espy
. - The name of requestable variables is preceded by
☼
(\sun
). Such anotation must always appear on the left of an assigment. - If the name of a variable is preceded by
♢
(\diamond
), then the variable is evaluated only if requested. Such a notation can only be used if there is only one variable left of the assignement. - The name of a function being called must be preceded by
☼
if the definition of the function is itself preceeded by the macro-call@espy
. for
-loops are not supported.do
-loops must be used: to be efficient,residual
andlagrangian
must not allocate and thus use immutables.- The keyword
return
must be explicitly used, and if must be followed the a comma separated list of output variables. Syntaxes likereturn if...
are not supported.
Muscade.@espydbg
— Macro@espydbg function ... end
Generate the same code as @espy
and print it (for debug purposes).
Muscade.@once
— Macro@once f(x)= x^2
do not parse the definition of function f
again if not modified. Using in a script, this prevents recompilations in Muscade
or applications based on it when they receive such functions as argument
Muscade.@request
— Macroreq = @request expr
Create a request of internal results wanted from a function. Considering the function presented as example for @espy
, examples of possible syntax include
req = @request gp(s,z,material(a,b))
req = @request gp(s)
req = @request gp(material(a))
The first expression can be read as follows: "In the function, there is a do
loop over variable igp
, and within this loop a call to a function material
. Results s
and z
are wanted from within the loop, and results a
and b
from within material
.
The corresponding datastructure containing the results for each element is a nesting of NTuples
and NamedTuples
, and can be accessed as out.gp[igp].material.a
and so forth.
Index
Muscade.noFB
Muscade.ϵ
Muscade.∞
Muscade.AbstractElement
Muscade.DirectXUA
Muscade.DofConstraint
Muscade.DofCost
Muscade.DofLoad
Muscade.ElementConstraint
Muscade.ElementCost
Muscade.Hold
Muscade.Model
Muscade.Monitor
Muscade.Node
Muscade.QuickFix
Muscade.SingleDofCost
Muscade.SingleUdof
Muscade.SweepX
Muscade.acceleration
Muscade.default
Muscade.motion
Muscade.position
Muscade.value
Muscade.value_∂
Muscade.variate
Muscade.velocity
Muscade.δ
Muscade.ℕ
Muscade.ℝ
Muscade.ℤ
Muscade.∂
Muscade.𝔹
Muscade.𝕓
Muscade.𝕟
Muscade.𝕣
Muscade.𝕫
Muscade.:∘₁
Muscade.:∘₂
Muscade.:⊗
Muscade.VALUE
Muscade.addelement!
Muscade.addin!
Muscade.addin!
Muscade.addnode!
Muscade.constants
Muscade.coord
Muscade.describe
Muscade.describe
Muscade.doflist
Muscade.draw
Muscade.draw
Muscade.eletyp
Muscade.equal
Muscade.findlastassigned
Muscade.getdof
Muscade.getndof
Muscade.getresult
Muscade.getsomedofs
Muscade.gradient
Muscade.initialize!
Muscade.lagrangian
Muscade.merge
Muscade.muscadeerror
Muscade.off
Muscade.positive
Muscade.prepare
Muscade.residual
Muscade.setdof!
Muscade.setscale!
Muscade.solve
Muscade.studyscale
Muscade.studysingular
Muscade.test_static_element
Muscade.toggle
Muscade.∂0
Muscade.∂1
Muscade.∂2
Muscade.@espy
Muscade.@espydbg
Muscade.@once
Muscade.@request