All subclasses of Cut will implement intersection routines and quadrature point generation within the cut in the intersection of two elements. Although, this class is designed to support MortarAssembler and ParMortarAssembler, it can be used for any problem requiring to perform Petrov-Galerkin formulations on non-matching elements.
More...
#include <cut.hpp>
|
virtual | ~Cut ()=default |
|
virtual void | SetIntegrationOrder (const int order)=0 |
| This method creates/updates the quadrature on the reference element based on the integration order.
|
|
virtual bool | BuildQuadrature (const FiniteElementSpace &from_space, const int from_elem_idx, const FiniteElementSpace &to_space, const int to_elem_idx, IntegrationRule &from_quadrature, IntegrationRule &to_quadrature)=0 |
| This computes the intersection of the finite element geometries and generates quadratures formulas for the two reference configurations.
|
|
virtual void | Describe () const |
| Method for printing information to the command line.
|
|
All subclasses of Cut will implement intersection routines and quadrature point generation within the cut in the intersection of two elements. Although, this class is designed to support MortarAssembler and ParMortarAssembler, it can be used for any problem requiring to perform Petrov-Galerkin formulations on non-matching elements.
Definition at line 28 of file cut.hpp.
◆ ~Cut()
virtual mfem::Cut::~Cut |
( |
| ) |
|
|
virtualdefault |
◆ BuildQuadrature()
This computes the intersection of the finite element geometries and generates quadratures formulas for the two reference configurations.
- Parameters
-
| from_space | the space from which we want to transfer a function with MortarAssembler and ParMortarAssembler |
| from_elem_idx | the index of the element of the from_space |
| to_space | the space to which we want to transfer a function with MortarAssembler and ParMortarAssembler |
| to_elem_idx | the index of the element of the to_space |
[out] | from_quadrature | the quadrature rule in the reference coordinate system of the from element |
[out] | to_quadrature | the quadrature rule in the reference coordinate system of the to element |
- Returns
- true if the two element intersected and if the output must be used or ignored.
◆ Describe()
virtual void mfem::Cut::Describe |
( |
| ) |
const |
|
inlinevirtual |
Method for printing information to the command line.
Definition at line 66 of file cut.hpp.
◆ SetIntegrationOrder()
virtual void mfem::Cut::SetIntegrationOrder |
( |
const int | order | ) |
|
|
pure virtual |
This method creates/updates the quadrature on the reference element based on the integration order.
- Parameters
-
order | the order of the quadrature rule |
◆ SetQuadratureRule()
The documentation for this class was generated from the following file: