Documentation ¶
Overview ¶
Package deflect implements components to solve the partial differential equations of Newton mechanics in their weak form and approximated with the finite element method, in particular using classical truss and beam elements for structural analysis. The package is concerned with element interpolation formulations (beams, trusses, frames, springs, ...), the mesh, boundary conditions, or assembling and solving the system of equations.
Index ¶
- func NewInclinedSupport(from, to Index, angle float64) (Transformer, NodalValue)
- type CrossSection
- type Dof
- type Element
- func NewFrame2d(id string, n0, n1 *Node, material *Material, hinges map[Index]struct{}) (Element, error)
- func NewTruss2d(id string, n0, n1 *Node, material *Material, hinges map[Index]struct{}) (Element, error)
- func NewTruss3d(id string, n0, n1 *Node, material *Material, hinges map[Index]struct{}) (Element, error)
- type EqLayout
- type EquationSolver
- type Fct
- type Index
- type Interpolation
- type LinearElastic
- type Material
- type NeumannElementBC
- type NodalValue
- type Node
- type PolyPiece
- type PolySequence
- type Problem
- type ProblemResult
- type ProblemSolver
- type Transformer
Constants ¶
This section is empty.
Variables ¶
This section is empty.
Functions ¶
func NewInclinedSupport ¶
func NewInclinedSupport(from, to Index, angle float64) (Transformer, NodalValue)
NewInclinedSupport instantiates a transformer and a zero Dirichlet boundary condition. Together, these two objects represent an inclined support. They are tightly coupled, as the prescribed index determines how the Transformer acts on the assembled tangent and vectors, and hence returned by the same constructor(-like) function.
Types ¶
type CrossSection ¶
type CrossSection interface { // Area returns the cross section's area in m^2 Area() float64 Iyy() float64 }
CrossSection provides access to all relevant geometric properties of a beam, truss, or frame cross section.
func NewConstantsCrossSections ¶
func NewConstantsCrossSections(param map[string]float64) (CrossSection, error)
NewConstantsCrossSections instantiates a cross section with all parameters specified as constants. The following keys are expected: - A - Iyy Returns an error if anything is not positive, or if the required parameters can't be found in param.
func NewRectangularCrossSection ¶
func NewRectangularCrossSection(b, h float64) (CrossSection, error)
NewRectangularCrossSection instantiates a new rectangular cross section. Returns an error if the dimensions are not positive.
type Dof ¶
type Dof uint8
Dof describes a degree of freedom symbolically and with an intuitive ordering. The predefined constants shall be used.
type Element ¶
type Element interface { // Assemble adds entries to the given tangent k and residual r, using the current primary nodal // values in d. Only non-linear elements need to read from primary. The matrices are global // entities, and the element uses indices to map symbolic indices to plain matrix indices, which // can be used to access the global matrices. Assemble(indices EqLayout, k *mat.SymDense, r, d *mat.VecDense) // AddLoad stores the given Neumann boundary condition to be used by the [Element.Assemble] // implementation. The same load can be added multiple times. The boolean return value indicates // if the load could be applied. If the inability to apply a BC is an error is up to the caller. AddLoad(bc NeumannElementBC) bool // RemoveLoad removes the given load. If there are multiple instances of the same load instance, // they are all removed. Removing a load that has previously not been added has no effect. RemoveLoad(bc NeumannElementBC) // Interpolate returns a sequence of piecewise polynomials to describe the primary solution or // stress/internal force in local element coordinates. Post-processing should be performed using // [PolySequence.TrimTrailingZeros] and [PolySequence.CompactIdentical] in that order, using a // context-dependent zero/equality tolerance (that's why it is not performed automatically). If // the element doesn't relate to the given dof, nil is returned. Interpolate(indices EqLayout, which Fct, d *mat.VecDense) PolySequence // Indices adds the indices this element uses (pairing of degree of freedom and nodal id) to set. Indices(set map[Index]struct{}) // NumNodes returns the number of nodes this element is connected to. NumNodes() uint // ID returns the potentially user-facing id (hence it's a string). ID() string }
Element defines the common API of any finite element formulation implemented in this package.
func NewFrame2d ¶
func NewFrame2d( id string, n0, n1 *Node, material *Material, hinges map[Index]struct{}, ) (Element, error)
NewFrame2d returns a 2d beam element implementation.
type EqLayout ¶
type EqLayout struct {
// contains filtered or unexported fields
}
EqLayout manages the mapping between symbolic and integral matrix indices. This is used for elements and boundary conditions to know where in the global matrices/vectors they should add their coefficients during assembly. It collects lookup errors behind the scene. These errors must be manually checked at a reasonable point in time - the individual lookups will _not_ fail, but return zero indices instead. This means an out-of-bounds panic is avoided, but any matrices populated with the erroneous index mapping can't be used. In the context of tangent/residual assembly, the solver process should probably be aborted.
func NewEqLayout ¶
NewEqLayout creates a new index layout for the given boundary value problem.
type EquationSolver ¶
EquationSolver implements an algorithm to solve a linear system of equations with a symmetric positive definite coefficient matrix. Examples: Cholesky or LU/QR decomposition, or an iterative method.
func NewCholeskySolver ¶
func NewCholeskySolver() EquationSolver
NewCholeskySolver creates a Cholesky solver for symmetric positive definite coefficient matrices.
type Fct ¶
type Fct uint8
Fct denotes a function to be interpolated on the element level. These can be primary solution variables (nodal values) as well as stresses/internal forces.
type Index ¶
Index pairs a nodal id with a degree of freedom and is the complete symbolic representation of a quantity within an array of a boundary value problem. Dealing with raw integral matrix indices is error-prone, so Index is used as long as possible instead.
type Interpolation ¶
type Interpolation struct { Element string Quantity Fct Piecewise PolySequence }
Interpolation associates an element ID with a Fct descriptor for the quantity, and a variable number of piecewise polynomials. When there is only one polynomial, it spans the entire element domain (usually from zero to the length of the element). When there is more than one polynomial, there are sub-ranges/intervals. There will not be gaps in between these sub-ranges, and there will only be a single PolyPiece instance per sub-range.
type LinearElastic ¶
type LinearElastic struct { // Young's modulus in N/m^2. YoungsModulus float64 // Poisson's ratio, dimensionless. Only one, since we assume an isotropic material. PoissonsRatio float64 // Density given in kg/m^3. We assume the Density is constant across every element. Density float64 }
LinearElastic "implements" a linear-elastic material law by providing access to required material parameters. This does not encapsulate or abstract away the actual material law - an element formulation is expected to inline the material law when computing residual and tangent. This might change in the future when we add non-linear material laws.
type Material ¶
type Material struct { CrossSection LinearElastic }
Material defines the API to retrieve material parameters to be used in element formulations. We currently restrain ourselves to linear-elastic materials, and the constant material parameters are directly accessed by each element to compute its residual and tangent. The Material abstraction is hence intentionally thin for now, and it is fine that element formulations and material law are tightly coupled.
In the future, non-linear material laws can be an important feature. Then, the Material API must be reworked (e.g. return a function to compute residual/tangent given current strain), and we would probably try to decouple the element's interpolation functionality from the material law. Then, an element formulation can work with different injected material instances.
type NeumannElementBC ¶
type NeumannElementBC any
NeumannElementBC is an opaque handle to be downcast by element implementations. It is always instantiated with a pointer.
func NewElementConcentratedLoad ¶
func NewElementConcentratedLoad(kind Dof, pos, value float64) (NeumannElementBC, error)
NewElementConcentratedLoad instantiates an element Neumann boundary condition that applies a concentrated quantity at the given position, e.g. a force or torque.
func NewElementConstantLoad ¶
func NewElementConstantLoad(kind Dof, value float64) NeumannElementBC
NewElementConstantLoad instantiates an element Neumann boundary condition that applies a distributed, constant load.
func NewElementLinearLoad ¶
func NewElementLinearLoad(kind Dof, first, last float64) NeumannElementBC
NewElementLinearLoad instantiates an element Neumann boundary condition that applies a distributed, linear load.
type NodalValue ¶
NodalValue associates a scalar solution value with an Index.
type Node ¶
Node describes a mesh vertex by 3-dimensional coordinates and an identifier. Coordinates are always in meters.
type PolyPiece ¶
type PolyPiece struct {
// Start and end of the domain. Precisely, the interval is [X0, XE], so that both values are
// included in the domain.
X0, XE float64
// The index of each coefficient is the associated exponent. Coefficients can be zero, but never
// the last one except when there is only a single entry. For example, 3·x² - 2.5 would be
// represented as []float{-2.5, 0.0, 3.0}.
Coeff []float64
}
PolyPiece describes a one-dimensional piecewise polynomial with non-negative exponents.
type PolySequence ¶
type PolySequence []PolyPiece
PolySequence is a piecewise polynomial.
func (PolySequence) CompactIdentical ¶
func (ps PolySequence) CompactIdentical(tol float64) PolySequence
CompactIdentical squashes consecutive identical polynomials into a single one, i.e. it mutates the object in place if there are consecutive identical PolyPiece instances. Since the size of the underlying slice may change, callers must use the return value the same way they would when appending to a slice.
func (PolySequence) IntervalIndex ¶
func (ps PolySequence) IntervalIndex(x0, xE float64) int
IntervalIndex behaves like slices.Index. A match is determined by comparing x0 and xE with the domain of the polynomials in ps, where both x0 and xE must match.
func (PolySequence) TrimTrailingZeros ¶
func (ps PolySequence) TrimTrailingZeros(zeroTol float64)
TrimTrailingZeros strips trailing approximate zeros in place, except the very first one, which is truncated if approximately zero.
type Problem ¶
type Problem struct { Nodes []Node Elements []Element // Implementation note: nodal BCs are saved separately from basic.Node, since there is no coupling // between nodes and nodal boundary conditions (this is different from element loading, which is // tightly coupled). Neumann []NodalValue // Only nodal Neumann BCs, no element loading Dirichlet []NodalValue EqTransforms []Transformer }
Problem stores data to solve a user specified boundary value problem.
func ProblemFromJSON ¶
ProblemFromJSON parses the given JSON data and constructs a boundary value problem from it.
type ProblemResult ¶
type ProblemResult interface { Primary(i Index) (NodalValue, error) Reaction(i Index) (NodalValue, error) PrimaryAll() []NodalValue ReactionAll() []NodalValue Interpolate(elmtID string, quantity Fct, zeroTol float64) (Interpolation, error) InterpolateAll(zeroTol float64) []Interpolation Dimension() (total, net int) }
ProblemResult is the API to retrieve BVP results as needed. An implementation can choose to evaluate individual results lazily or in batches, whatever makes most sense for the representation of the data.
type ProblemSolver ¶
type ProblemSolver interface { Solve( p *Problem, idx EqLayout, strategy EquationSolver, ) (ProblemResult, error) }
ProblemSolver solves the given boundary value problem and returns the complete set of primary nodal values and reactions for nodes with Dirichlet BC.
func NewLinearProblemSolver ¶
func NewLinearProblemSolver() ProblemSolver
NewLinearProblemSolver creates a linear solver for boundary value problems.
type Transformer ¶
type Transformer interface { Pre(indices EqLayout, k *mat.SymDense, r, d *mat.VecDense) Post(indices EqLayout, r, d *mat.VecDense) }
Transformer mutate the system of equation before and after it is solved. For example, an inclined support can be split into a Transformer and a Dirichlet BC, such that they act independently to link degrees of freedom through a trigonometric relation. Transformer instances don't prescribe values in r or d (this is done by NodalBC instances).
Source Files ¶
- algorithm.go
- beam2d.go
- cholesky.go
- cross_sections.go
- deflect.go
- doc.go
- dof_string.go
- element_loads.go
- eq_layout.go
- fct_string.go
- frame2d.go
- geometry.go
- inclined_support.go
- linear_bvp_solver.go
- linear_elastic.go
- one_dim_element.go
- polynomial.go
- problem_from_json.go
- solver_result.go
- static_condensation.go
- truss2d.go
- truss3d.go