Next: Package FL.ELLSYS-FE, Previous: Package FL.NAVIER-STOKES, Up: Reference manual [Contents][Index]
The FL.DISCRETIZATION
package defines
<discretization>
as an abstract class and
<fe-discretization>
as a concrete derived class.
The key for local assembly is given by the generic function
get-fe
, which yields a suitable finite element for a given cell.
The value of get-fe
is a class <fe>
for scalar problems
or <vector-fe>
for vector-valued problems which contains
information on base functions and node functionals. Another generic
function quadrature-rule
computes memoized quadrature rules for
those finite elements.
Obviously, also non-standard finite element discretizations like hp-methods
would fit into this scheme. The key for local assembly is given by the
generic function get-fe
, which yields a suitable finite element
for a given cell. The value of get-fe
is a class <fe>
for scalar problems or <vector-fe>
for vector-valued problems which
contains information on base functions and node functionals. Another
generic function quadrature-rule
computes memoized quadrature
rules for those finite elements.
The file lagrange.lisp provides Lagrange finite elements of arbitrary order. The evaluation points for the node functionals may be chosen either uniformly distributed or at the Gauss-Lobatto points. The latter choice of points yields better stability properties but is restricted to cube meshes. Also functions for constructing cell mappings by pointwise evaluation of the domain boundary are provided here, which may be used to construct isoparametric domain approximations.
In the file fedisc.lisp, the function fe-discretize
is
defined. This function performs the standard steps for finite element
discretization: interior assembly, assembly and elimination of hanging-node
and essential-boundary constraints. It works on a blackboard as explained
in Section Blackboards and can reuse already available matrix-vector
structure. There is a somewhat less flexible interface provided by the
funtion discretize-globally
which calls
fe-discretize
.
In the files ellsys-fe.lisp, elasticity-fe.lisp and
navier-stokes.lisp one can find methods for local assembly for the
different problems. They are defined in own packages which use both the
package FL.DISCRETIZATION
and the package for the particular
problem.
Quadrature order to be used. NIL determines the order automatically.
The suggested order of discretization. In non-nil, this value should be taken into account by methods to select-discretization.
A automorphism of an ansatz space.
Superclasses: <ANSATZ-SPACE-MORPHISM> <ANSATZ-SPACE-OBJECT>
A mapping between two ansatz-spaces.
A sparse vector which is interpreted as the ansatz-space for a specific fe-class on a given mesh.
Superclasses: <ANSATZ-SPACE-OBJECT> <HT-SPARSE-VECTOR>
A finite element ansatz space is determined by finite element discretization, mesh and problem. The constraints are stored in the slot properties.
Superclasses: PROPERTY-MIXIN
Direct slots:
Discretization base class.
Superclasses: PROPERTY-MIXIN
FE discretization base class.
Superclasses: <DISCRETIZATION>
Class for scalar fe discretizations.
Superclasses: <STANDARD-FE-DISCRETIZATION>
Direct slots:
A finite element <fe> is given for each reference cell, e.g. <2-simplex>. dofs are the degrees of freedom associated with the cell, basis is the dual basis to dofs in some polynomial space. subcell-ndofs is the number of ndofs on each subcell. subcell-indices is a list of indices for all subcells with dofs. Usually, the <scalar-fe> will occur as values of a procedure or as values in a hash-table with the reference cells as keys.
Superclasses: <FE>
Direct slots:
Finite element discretization where the finite element depends only on the cell (usually via its reference cell).
Superclasses: <FE-DISCRETIZATION>
Direct slots:
Vector FE discretization class.
Superclasses: <STANDARD-FE-DISCRETIZATION>
Direct slots:
Finite element for vector functions. Components is a vector of scalar finite elements. Subcell-offsets is an array consisting of arrays which yield such an offset for every subcell.
Superclasses: <FE>
Direct slots:
Assemble the interior, i.e. ignore constraints arising from boundaries
and hanging nodes. Discretization is done using the ansatz space
ansatz-space on level level. The level argument will usually
be NIL
when performing a global assembly, and be equal to some
number when assembling coarse level matrices for multigrid. The argument
where is a flag indicating where assembly is to be done. It should
be one of the keywords :surface
, :refined
, :all
. The
arguments matrix, rhs should contain vectors/matrices where the
local assembly is accumulated. Boundary conditions and constraints are not
taken into account within this routine.
In general, this function does most of the assembly work. Other steps like handling constraints are intricate, but usually of lower computational complexity.
Integrates object on cell.
If cell is identified, its identification is the key.
Returns a Lagrange fe depending on reference cell, an order (which can be number or vector), and a type symbol.
Choose a reasonable start vector for some strategy.
Returns the comp-th component of the (vector) finite element fe.
The multigrid algorithm needs an interpolation which satisfies the constraints like essential or periodic boundary conditions.
Constructs a coefficient input list from FE data cell is the cell, global is the global coordinate of the integration point, values and gradients the values and gradients of the shape functions at the ip, and fe-parameters are the corresponding data of fe-functions to be evalutated.
Returns the order of the discretization.
Calculates matrix and right-hand side for the discretization and the problem. Further parameters are provided inside the blackboard.
Discretize problem on the hierarchical mesh h-mesh using finite elments given by fe-class.
Computes a local stiffness matrix and right-hand side. The algorithm will usually work as follows:
If local-u and local-v are set, then local-v*local-mat*local-u and local-v*local-rhs is computed. This feature may be used later on for implementing matrixless computations.
Degree of freedom in a finite element. It is defined as a functional defined by integration over a sub-cell or by evaluation at a local coordinate of a sub-cell of a reference cell.
Direct slots:
Constraints are given by an equation: P x = Q x + r
Here x in V, P is an orthogonal projection on a subspace V_P of V, Q maps some other space which may have nonempty intersection with V_P to V_P. With S we denote the mapping Id-P. This function returns the matrix for a Galerkin method on the constrained space. It is used for treating hanging nodes and essential boundary conditions. When assemble-locally is t the sparse structure of mat is used instead of the sparse structure of the hanging node interface. When include-constraints is non-nil, the constraints are included in matrix and rhs.
Computation of essential constraints. Should probably be incorporated into the ansatz-space definition.
Extract a sub-matrix from a sparse matrix for the given keys.
Collects cell geometry information at sample-points inside a property list.
Finite element discretization for an ansatz space provided on the blackboard.
Computes information for extracting components out of a vector finite element.
Computes the extreme values of a finite element function over the domain or some region. The result is a pair, the car being the minimum values and the cdr the maximum values. Each part is a matrix of the format ncomps x multiplicity.
Evaluates the gradient of the FE ansatz-space-vector asv at global-pos.
Integrates a finite element function over the domain. key is a transformer function, as always (e.g. #’abs if you want the L1-norm).
Evaluates the gradient of a FE ansatz-space-vector on cell at local-pos.
Evaluates a FE ansatz-space-vector in cell at local-pos.
Evaluates a FE ansatz-space-vector at global-pos.
Returns an s-point Gauss integration rule.
Returns the finite element for the given discretization and reference cell.
Maps the region in the global stiffness matrix determined by cell to a local matrix array.
Maps the region in global-vec determined by cell to a local vector.
The hierarchical mesh for the given ansatz-space or ansatz-space object.
Increments the region in global-mat determined by cell to the values of local-mat.
Increments the region in global-vec determined by cell to the values of the local vector array.
Returns an integration rule for cell of the given order.
Interpolates function on the reference cell of the finite element fe. Returns a standard-matrix corresponding to the block in the sparse vector.
On each cell of the skeleton region or on all cells of level level of the mesh of ansatz-space, a local interpolation matrix is collected into an interpolation matrix. type is the interpolation type having a default value *interpolation-type*.
Returns a vector of local gradient matrices for obj which may be a vector of integration points or a quadrature rule. Note that this function is memoized using an :around method.
Returns a vector of ip values for obj which may be a vector of integration points or a quadrature rule. Note that this function is memoized using an :around method.
A constructor for a problem-dependent Lagrange fe. Here, the number of components may vary with the respective patch.
Constructor for Lagrange fe. nr-comps=nil builds a scalar fe-discretization, otherwise a vector-fe-discretization is built.
Returns a function which maps a cell by a polynomial which is obtained by interpolating the boundary map via Lagrange interpolation.
Memoized call of compute-local-imatrix.
Memoized call of compute-local-pmatrix.
Computes a local transfer matrix between different FE spaces.
Constructor of <ansatz-space>
.
Generates a local matrix discretization for the given ansatz-space(s).
Generates a local vector for local discretization. We allow overriding the default multiplicity for the ansatz space for special use cases.
Evaluates the vector given in local-vec at multiple points. Here
local-vec should be a data vector obtained with
get-local-from-global-vec
and ip-values should be a vector
obtained from ip-values
.
Number of degrees of freedom for the (vector) finite element.
Number of inner degrees of freedom for the (vector) finite element.
Computes a product rule for several lower-dimensional quadrature rules.
The algorithm works as follows: On each cell of the provided cell list or the whole refinement a local projection matrix computed on the reference finite element is copied into the global projection matrix.
Computes the quadrature rule to be used for the finite element fe.
Returns a ansatz space vector for ansatz-space filled with random entries. Essential constraints are satisfied.
Returns a representative for this object.
Select a discretization for the given problem depending on the parameters on the blackboard.
Sets the region in global-mat determined by cell to the values of the local matrix array.
Sets the region in global-vec determined by cell to the values of the local vector array.
Returns a ansatz space vector for ansatz-space filled with constant or random entries. Essential constraints are satisfied.
Reader for subcell-offsets. This is an array of length the number of components. Each component is an array giving the offset of this component in a sparse vector value block corresponding to the subcell.
Builds a transfer matrix from domain-as to image-as.
Determines weights for the integration points ips such that they integrate int_{-1}^{1} (1+y)^beta ... dy with optimal order.
Next: Package FL.ELLSYS-FE, Previous: Package FL.NAVIER-STOKES, Up: Reference manual [Contents][Index]