18 #ifndef ASSEMBLY_BASE_HH_
19 #define ASSEMBLY_BASE_HH_
33 template <
unsigned int dim>
114 ASSERT( cell_side.
dim() > 0 ).error(
"Invalid cell dimension, must be 1, 2 or 3!\n");
120 ASSERT( cell_side.
dim() > 1 ).error(
"Invalid cell dimension, must be 2 or 3!\n");
126 ASSERT( cell_side.
dim() > 0 ).error(
"Invalid cell dimension, must be 1, 2 or 3!\n");
132 for (
unsigned int i=0; i<bulk_integral_data.
permanent_size(); ++i) {
133 if (bulk_integral_data[i].cell.dim() != dim)
continue;
146 for (
unsigned int i=0; i<boundary_integral_data.
permanent_size(); ++i) {
147 if (boundary_integral_data[i].side.dim() != dim)
continue;
154 for (
unsigned int i=0; i<edge_integral_data.
permanent_size(); ++i) {
155 auto range = edge_integral_data[i].edge_side_range;
156 if (range.begin()->dim() != dim)
continue;
163 for (
unsigned int i=0; i<coupling_integral_data.
permanent_size(); ++i) {
164 if (coupling_integral_data[i].side.dim() != dim)
continue;
165 this->
dimjoin_intergral(coupling_integral_data[i].cell, coupling_integral_data[i].side);
185 std::shared_ptr<BulkIntegral>
bulk_;
186 std::shared_ptr<EdgeIntegral>
edge_;
214 template <
unsigned int dim>
231 for (
unsigned int i=0; i<bulk_integral_data.
permanent_size(); ++i) {
232 if (bulk_integral_data[i].cell.dim() != dim)
continue;
236 for (
auto p : this->
bulk_points(element_patch_idx) ) {
237 unsigned int value_cache_idx = p.elm_cache_map()->element_eval_point(p.elem_patch_idx(), p.eval_point_idx());
245 for (
unsigned int i=0; i<boundary_integral_data.
permanent_size(); ++i) {
246 if (boundary_integral_data[i].side.dim() != dim)
continue;
248 for (
auto p : this->
boundary_points(boundary_integral_data[i].side) ) {
249 unsigned int value_cache_idx = p.elm_cache_map()->element_eval_point(p.elem_patch_idx(), p.eval_point_idx());
257 for (
unsigned int i=0; i<edge_integral_data.
permanent_size(); ++i) {
258 auto range = edge_integral_data[i].edge_side_range;
259 if (range.begin()->dim() != dim)
continue;
264 unsigned int value_cache_idx = p.elm_cache_map()->element_eval_point(p.elem_patch_idx(), p.eval_point_idx());
273 uint side_pos, element_patch_idx, elm_pos=0;
274 uint last_element_idx = -1;
276 for (
unsigned int i=0; i<coupling_integral_data.
permanent_size(); ++i) {
277 if (coupling_integral_data[i].side.dim() != dim)
continue;
279 if (coupling_integral_data[i].cell.elm_idx() != last_element_idx) {
285 for (
auto p_high : this->
coupling_points(coupling_integral_data[i].side) )
287 unsigned int value_cache_idx = p_high.elm_cache_map()->element_eval_point(p_high.elem_patch_idx(), p_high.eval_point_idx());
289 if (coupling_integral_data[i].cell.elm_idx() != last_element_idx) {
290 auto p_low = p_high.lower_dim(coupling_integral_data[i].cell);
291 value_cache_idx = p_low.elm_cache_map()->element_eval_point(p_low.elem_patch_idx(), p_low.eval_point_idx());
295 last_element_idx = coupling_integral_data[i].cell.elm_idx();
306 return fe_values_->template bulk_values<dim>();
311 return fe_values_->template side_values<dim>();
316 return fe_values_->template join_values<dim>();
#define ASSERT_PERMANENT_PTR(ptr)
Definition of assert macro checking non-null pointer (PTR)
AssemblyBasePatch(PatchFEValues< 3 > *fe_values)
void add_patch_coupling_integrals(const RevertableList< CouplingIntegralData > &coupling_integral_data) override
Register bulk and side points of coupling integral.
GenericAssemblyBase::CouplingIntegralData CouplingIntegralData
void add_patch_bulk_points(const RevertableList< BulkIntegralData > &bulk_integral_data) override
Register cell points of volume integral.
void add_patch_edge_points(const RevertableList< EdgeIntegralData > &edge_integral_data) override
Register side points of edge integral.
JoinValues< dim > join_values()
Return JoinValues object.
PatchFEValues< 3 > * fe_values_
Common FEValues object over all dimensions.
void add_patch_bdr_side_points(const RevertableList< BoundaryIntegralData > &boundary_integral_data) override
Register side points of boundary side integral.
BulkValues< dim > bulk_values()
Return BulkValues object.
GenericAssemblyBase::BoundaryIntegralData BoundaryIntegralData
unsigned int n_dofs()
Return BulkValues object.
GenericAssemblyBase::EdgeIntegralData EdgeIntegralData
GenericAssemblyBase::BulkIntegralData BulkIntegralData
SideValues< dim > side_values()
Return SideValues object.
virtual void assemble_cell_integrals(const RevertableList< BulkIntegralData > &bulk_integral_data)
Assembles the cell integrals for the given dimension.
Range< BulkPoint > bulk_points(unsigned int element_patch_idx) const
Return BulkPoint range of appropriate dimension.
virtual void boundary_side_integral(FMT_UNUSED DHCellSide cell_side)
Assembles the fluxes on the boundary.
virtual void end()
Method finishes object after assemblation (e.g. balance, ...).
virtual void add_patch_edge_points(FMT_UNUSED const RevertableList< EdgeIntegralData > &edge_integral_data)
Register side points of edge integral.
virtual ~AssemblyBase()
Destructor.
virtual void cell_integral(FMT_UNUSED DHCellAccessor cell, FMT_UNUSED unsigned int element_patch_idx)
Assembles the volume integrals on cell.
Range< EdgePoint > edge_points(const DHCellSide &cell_side) const
Return EdgePoint range of appropriate dimension.
GenericAssemblyBase::CouplingIntegralData CouplingIntegralData
Quadrature * quad_
Quadrature used in assembling methods.
void create_integrals(std::shared_ptr< EvalPoints > eval_points, AssemblyIntegrals &integrals)
Create integrals according to dim of assembly object.
void assemble_edge_integrals(const RevertableList< EdgeIntegralData > &edge_integral_data)
Assembles the edge integrals for the given dimension.
virtual void add_patch_coupling_integrals(FMT_UNUSED const RevertableList< CouplingIntegralData > &coupling_integral_data)
Register bulk and side points of coupling integral.
virtual void patch_reinit(FMT_UNUSED std::array< PatchElementsList, 4 > &patch_elements)
DimIntegrals integrals_
Set of used integrals.
virtual void begin()
Method prepares object before assemblation (e.g. balance, ...).
std::string print_update_flags(UpdateFlags u) const
Print update flags to string format.
Quadrature * quad_low_
Quadrature used in assembling methods (dim-1).
virtual void add_patch_bdr_side_points(FMT_UNUSED const RevertableList< BoundaryIntegralData > &boundary_integral_data)
Register side points of boundary side integral.
GenericAssemblyBase::BulkIntegralData BulkIntegralData
void assemble_neighbour_integrals(const RevertableList< CouplingIntegralData > &coupling_integral_data)
Assembles the neighbours integrals for the given dimension.
Range< CouplingPoint > coupling_points(const DHCellSide &cell_side) const
Return CouplingPoint range of appropriate dimension.
virtual void dimjoin_intergral(FMT_UNUSED DHCellAccessor cell_lower_dim, FMT_UNUSED DHCellSide neighb_side)
Assembles the fluxes between elements of different dimensions.
virtual void edge_integral(FMT_UNUSED RangeConvert< DHEdgeSide, DHCellSide > edge_side_range)
Assembles the fluxes between sides on the edge.
int active_integrals_
Holds mask of active integrals.
Range< BoundaryPoint > boundary_points(const DHCellSide &cell_side) const
Return BoundaryPoint range of appropriate dimension.
AssemblyBase(unsigned int quad_order)
Constructor.
void assemble_boundary_side_integrals(const RevertableList< BoundaryIntegralData > &boundary_integral_data)
Assembles the boundary side integrals for the given dimension.
GenericAssemblyBase::EdgeIntegralData EdgeIntegralData
virtual void add_patch_bulk_points(FMT_UNUSED const RevertableList< BulkIntegralData > &bulk_integral_data)
Register cell points of volume integral.
int n_active_integrals() const
Getter of active_integrals.
ElementCacheMap * element_cache_map_
ElementCacheMap shared with GenericAssembly object.
GenericAssemblyBase::BoundaryIntegralData BoundaryIntegralData
Cell accessor allow iterate over DOF handler cells.
Side accessor allows to iterate over sides of DOF handler cell.
unsigned int dim() const
Return dimension of element appropriate to the side.
Directing class of FieldValueCache.
unsigned int position_in_cache(unsigned mesh_elm_idx, bool bdr=false) const
Return position of element stored in ElementCacheMap.
uint register_element(DHCellAccessor cell, uint element_patch_idx)
Register element to patch_point_vals_ table by dimension of element.
uint register_bulk_point(DHCellAccessor cell, uint elem_table_row, uint value_patch_idx, uint i_point_on_elem)
Register bulk point to patch_point_vals_ table by dimension of element.
uint register_side(DHCellSide cell_side)
Register side to patch_point_vals_ table by dimension of side.
Quadrature * get_quadrature(uint dim, bool is_bulk) const
Return bulk or side quadrature of given dimension.
uint register_side_point(DHCellSide cell_side, uint elem_table_row, uint value_patch_idx)
Register side point to patch_point_vals_ table by dimension of side.
Symmetric Gauss-Legendre quadrature formulae on simplices.
Base class for quadrature rules on simplices in arbitrary dimensions.
Definitions of particular quadrature rules on simplices.
Set of integral of given dimension necessary in assemblation.
std::shared_ptr< EdgeIntegral > edge_
Edge integrals between elements of same dimensions.
std::shared_ptr< CouplingIntegral > coupling_
Coupling integrals between elements of dimensions dim and dim-1.
std::shared_ptr< BulkIntegral > bulk_
Bulk integrals of elements.
std::shared_ptr< BoundaryIntegral > boundary_
Boundary integrals betwwen side and boundary element of dim-1.
Set of all used integral necessary in assemblation.
std::array< std::shared_ptr< BulkIntegral >, 3 > bulk_
Bulk integrals of elements of dimensions 1, 2, 3.
std::array< std::shared_ptr< CouplingIntegral >, 2 > coupling_
Coupling integrals between elements of dimensions 1-2, 2-3.
std::array< std::shared_ptr< BoundaryIntegral >, 3 > boundary_
Boundary integrals betwwen elements of dimensions 1, 2, 3 and boundaries.
std::array< std::shared_ptr< EdgeIntegral >, 3 > edge_
Edge integrals between elements of dimensions 1, 2, 3.
std::size_t permanent_size() const
Return permanent size of list.
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.