Next: , Previous: , Up: Reference manual   [Contents][Index]


7.18 Package FL.DISCRETIZATION

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.

Variable: *QUADRATURE-ORDER*

Quadrature order to be used. NIL determines the order automatically.

Variable: *SUGGESTED-DISCRETIZATION-ORDER*

The suggested order of discretization. In non-nil, this value should be taken into account by methods to select-discretization.

Class: <ANSATZ-SPACE-AUTOMORPHISM>

A automorphism of an ansatz space.

Superclasses: <ANSATZ-SPACE-MORPHISM> <ANSATZ-SPACE-OBJECT>

Class: <ANSATZ-SPACE-MORPHISM>

A mapping between two ansatz-spaces.

Class: <ANSATZ-SPACE-VECTOR>

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>

Class: <ANSATZ-SPACE>

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:

Class: <DISCRETIZATION>

Discretization base class.

Superclasses: PROPERTY-MIXIN

Class: <FE-DISCRETIZATION>

FE discretization base class.

Superclasses: <DISCRETIZATION>

Class: <SCALAR-FE-DISCRETIZATION>

Class for scalar fe discretizations.

Superclasses: <STANDARD-FE-DISCRETIZATION>

Direct slots:

Class: <SCALAR-FE>

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:

Class: <STANDARD-FE-DISCRETIZATION>

Finite element discretization where the finite element depends only on the cell (usually via its reference cell).

Superclasses: <FE-DISCRETIZATION>

Direct slots:

Class: <VECTOR-FE-DISCRETIZATION>

Vector FE discretization class.

Superclasses: <STANDARD-FE-DISCRETIZATION>

Direct slots:

Class: <VECTOR-FE>

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:

Function: ASSEMBLE-INTERIOR ANSATZ-SPACE WHERE &KEY LEVEL MATRIX MASS-MATRIX RHS PARALLEL-CLUSTERING &ALLOW-OTHER-KEYS

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.

Function: CELL-INTEGRATE CELL OBJECT &KEY INITIAL-VALUE COMBINER KEY COEFF-FUNC NO-WEIGHTS-P QUADRATURE-ORDER &ALLOW-OTHER-KEYS

Integrates object on cell.

Function: CELL-KEY CELL MESH

If cell is identified, its identification is the key.

Function: CELL-LAGRANGE-FE CELL ORDER TYPE &OPTIONAL DISC

Returns a Lagrange fe depending on reference cell, an order (which can be number or vector), and a type symbol.

Function: CHOOSE-START-VECTOR ANSATZ-SPACE

Choose a reasonable start vector for some strategy.

Function: COMPONENT FE COMP

Returns the comp-th component of the (vector) finite element fe.

Function: CONSTRAINED-INTERPOLATION-MATRIX ANSATZ-SPACE &KEY LEVEL WHERE IMAT (TYPE LOCAL)

The multigrid algorithm needs an interpolation which satisfies the constraints like essential or periodic boundary conditions.

Function: CONSTRUCT-COEFF-INPUT CELL GLOBAL DPHI VALUES GRADIENTS FE-PARAMETERS PROBLEM

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.

Function: DISCRETIZATION-ORDER DISC

Returns the order of the discretization.

Function: DISCRETIZE DISCRETIZATION PROBLEM BLACKBOARD

Calculates matrix and right-hand side for the discretization and the problem. Further parameters are provided inside the blackboard.

Function: DISCRETIZE-GLOBALLY PROBLEM H-MESH FE-CLASS

Discretize problem on the hierarchical mesh h-mesh using finite elments given by fe-class.

Function: DISCRETIZE-LOCALLY PROBLEM COEFFS FE QRULE FE-GEOMETRY &KEY MATRIX MASS-MATRIX RHS LOCAL-U LOCAL-V FE-PARAMETERS RESIDUAL-P &ALLOW-OTHER-KEYS

Computes a local stiffness matrix and right-hand side. The algorithm will usually work as follows:

  1. Get coefficient functions for the patch of the cell.
  2. Compute geometry information for all ips (values and gradients of the shape functions).
  3. Loop over integration points ip:
    1. If necessary, compute input for coefficient functions. This input may contain values of finite element function in the property list fe-parameters.
    2. Evaluate coefficient functions at ips.
    3. Add the contributions for matrix and right-hand side to local-mat and local-rhs.

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.

Class: DOF

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:

Function: ELIMINATE-CONSTRAINTS MAT RHS CONSTRAINTS-P CONSTRAINTS-Q CONSTRAINTS-R &KEY ASSEMBLE-LOCALLY INCLUDE-CONSTRAINTS

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.

Function: ESSENTIAL-BOUNDARY-CONSTRAINTS PROBLEM ANSATZ-SPACE &KEY LEVEL WHERE INTERFACE &ALLOW-OTHER-KEYS

Computation of essential constraints. Should probably be incorporated into the ansatz-space definition.

Function: EXTENDED-EXTRACT MAT KEYS &KEY ROW? COL?

Extract a sub-matrix from a sparse matrix for the given keys.

Function: FE-CELL-GEOMETRY CELL SAMPLE-POINTS &KEY WEIGHTS METRIC VOLUME &AUX (N (LENGTH SAMPLE-POINTS))

Collects cell geometry information at sample-points inside a property list.

Function: FE-DISCRETIZE BLACKBOARD

Finite element discretization for an ansatz space provided on the blackboard.

Function: FE-EXTRACTION-INFORMATION FE INDICES

Computes information for extracting components out of a vector finite element.

Function: FE-EXTREME-VALUES ASV &KEY CELLS SKELETON COMPONENT

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.

Function: FE-GRADIENT ASV GLOBAL-POS

Evaluates the gradient of the FE ansatz-space-vector asv at global-pos.

Function: FE-INTEGRATE ASV &KEY CELLS SKELETON INITIAL-VALUE COMBINER KEY COEFF-FUNC NO-WEIGHTS-P

Integrates a finite element function over the domain. key is a transformer function, as always (e.g. #’abs if you want the L1-norm).

Function: FE-LOCAL-GRADIENT ASV CELL LOCAL-POS

Evaluates the gradient of a FE ansatz-space-vector on cell at local-pos.

Function: FE-LOCAL-VALUE ASV CELL LOCAL-POS

Evaluates a FE ansatz-space-vector in cell at local-pos.

Function: FE-VALUE ASV GLOBAL-POS

Evaluates a FE ansatz-space-vector at global-pos.

Function: GAUSS-RULE FACTOR-DIMS S

Returns an s-point Gauss integration rule.

Function: GET-FE FE-DISC CELL

Returns the finite element for the given discretization and reference cell.

Function: GET-LOCAL-FROM-GLOBAL-MAT GLOBAL-MAT CELL &OPTIONAL DOMAIN-CELL

Maps the region in the global stiffness matrix determined by cell to a local matrix array.

Function: GET-LOCAL-FROM-GLOBAL-VEC CELL GLOBAL-VEC

Maps the region in global-vec determined by cell to a local vector.

Function: HIERARCHICAL-MESH ASO

The hierarchical mesh for the given ansatz-space or ansatz-space object.

Function: INCREMENT-GLOBAL-BY-LOCAL-MAT GLOBAL-MAT LOCAL-MAT CELL &OPTIONAL DOMAIN-CELL

Increments the region in global-mat determined by cell to the values of local-mat.

Function: INCREMENT-GLOBAL-BY-LOCAL-VEC CELL GLOBAL-VEC LOCAL-VEC

Increments the region in global-vec determined by cell to the values of the local vector array.

Function: INTEGRATION-RULE CELL ORDER

Returns an integration rule for cell of the given order.

Function: INTERPOLATE-ON-REFCELL FE FUNCTION

Interpolates function on the reference cell of the finite element fe. Returns a standard-matrix corresponding to the block in the sparse vector.

Function: INTERPOLATION-MATRIX ANSATZ-SPACE &KEY LEVEL REGION IMAT (TYPE *INTERPOLATION-TYPE*)

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*.

Function: IP-GRADIENTS FE OBJ

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.

Function: IP-VALUES FE OBJ

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.

Function: LAGRANGE-ANSATZ-SPACE PROBLEM MESH &KEY (ORDER 1) (TYPE UNIFORM)

A constructor for a problem-dependent Lagrange fe. Here, the number of components may vary with the respective patch.

Function: LAGRANGE-FE ORDER &KEY (NR-COMPS 1) (TYPE UNIFORM)

Constructor for Lagrange fe. nr-comps=nil builds a scalar fe-discretization, otherwise a vector-fe-discretization is built.

Function: LAGRANGE-MAPPING ORDER &OPTIONAL (TYPE UNIFORM)

Returns a function which maps a cell by a polynomial which is obtained by interpolating the boundary map via Lagrange interpolation.

Function: LOCAL-IMATRIX RULE FE &OPTIONAL (TYPE LOCAL)

Memoized call of compute-local-imatrix.

Function: LOCAL-PMATRIX RULE FE

Memoized call of compute-local-pmatrix.

Function: LOCAL-TRANSFER-MATRIX FE-FROM FE-TO

Computes a local transfer matrix between different FE spaces.

Function: MAKE-FE-ANSATZ-SPACE FE-CLASS PROBLEM MESH

Constructor of <ansatz-space>.

Function: MAKE-LOCAL-MAT AS1 CELL1 &OPTIONAL AS2 CELL2

Generates a local matrix discretization for the given ansatz-space(s).

Function: MAKE-LOCAL-VEC ANSATZ-SPACE CELL &KEY MULTIPLICITY

Generates a local vector for local discretization. We allow overriding the default multiplicity for the ansatz space for special use cases.

Function: MULTIPLE-EVALUATE-LOCAL-FE LOCAL-VEC SHAPE-VALUES

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.

Function: NR-OF-DOFS FE

Number of degrees of freedom for the (vector) finite element.

Function: NR-OF-INNER-DOFS FE

Number of inner degrees of freedom for the (vector) finite element.

Function: PRODUCT-RULE &REST QUADRATURE-RULES

Computes a product rule for several lower-dimensional quadrature rules.

Function: PROJECTION-MATRIX ANSATZ-SPACE &KEY LEVEL REGION PMAT

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.

Function: QUADRATURE-RULE FE

Computes the quadrature rule to be used for the finite element fe.

Function: RANDOM-ANSATZ-SPACE-VECTOR ANSATZ-SPACE

Returns a ansatz space vector for ansatz-space filled with random entries. Essential constraints are satisfied.

Function: REPRESENTATIVE OBJ

Returns a representative for this object.

Function: SELECT-DISCRETIZATION PROBLEM BLACKBOARD

Select a discretization for the given problem depending on the parameters on the blackboard.

Function: SET-GLOBAL-TO-LOCAL-MAT GLOBAL-MAT LOCAL-MAT CELL &OPTIONAL DOMAIN-CELL

Sets the region in global-mat determined by cell to the values of the local matrix array.

Function: SET-GLOBAL-TO-LOCAL-VEC CELL GLOBAL-VEC LOCAL-VEC

Sets the region in global-vec determined by cell to the values of the local vector array.

Function: SPECIAL-ANSATZ-SPACE-VECTOR ANSATZ-SPACE &OPTIONAL (TYPE RANDOM) (VALUE 1.0) (CONSTRAINTS-P T)

Returns a ansatz space vector for ansatz-space filled with constant or random entries. Essential constraints are satisfied.

Function: SUBCELL-OFFSETS FE

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.

Function: TRANSFER-MATRIX DOMAIN-AS IMAGE-AS &KEY NO-SLAVES

Builds a transfer matrix from domain-as to image-as.

Function: WEIGHTS-FOR-IPS BETA IPS

Determines weights for the integration points ips such that they integrate int_{-1}^{1} (1+y)^beta ... dy with optimal order.


Next: , Previous: , Up: Reference manual   [Contents][Index]