Flow123d
master-3768d5dec
|
Go to the documentation of this file.
18 #ifndef GENERIC_ASSEMBLY_HH_
19 #define GENERIC_ASSEMBLY_HH_
42 std::array<std::shared_ptr<BulkIntegral>, 3>
bulk_;
43 std::array<std::shared_ptr<EdgeIntegral>, 3>
edge_;
44 std::array<std::shared_ptr<CouplingIntegral>, 2>
coupling_;
45 std::array<std::shared_ptr<BoundaryIntegral>, 3>
boundary_;
146 virtual void assemble(std::shared_ptr<DOFHandlerMultiDim> dh) = 0;
159 template <
template<
IntDim...>
class DimAssembly>
164 GenericAssembly(
typename DimAssembly<1>::EqFields *eq_fields,
typename DimAssembly<1>::EqData *eq_data)
204 void assemble(std::shared_ptr<DOFHandlerMultiDim> dh)
override {
209 bool add_into_patch =
false;
210 for(
auto cell_it = dh->local_range().begin(); cell_it != dh->local_range().end(); )
213 if (!add_into_patch) {
215 add_into_patch =
true;
229 add_into_patch =
false;
238 add_into_patch =
false;
243 if (add_into_patch) {
319 if ( (cell_side.side().edge().n_sides() == 1) && (cell_side.side().is_boundary()) ) {
324 if ( (cell_side.n_edge_sides() >=
min_edge_sides_) && (cell_side.edge_sides().begin()->element().idx() == cell.
elm_idx())) {
332 if (cell.
dim() != neighb_side.dim()-1)
continue;
359 unsigned int reg_idx = edge_side.element().region_idx().idx();
377 auto p_low = p.lower_dim(cell);
unsigned int side_subset_index
Index (order) of higher dim subset in EvalPoints object.
RevertableList< BulkIntegralData > bulk_integral_data_
Holds data for computing bulk integrals.
void add_boundary_integral(const DHCellSide &bdr_side)
Add data of boundary integral to appropriate data structure.
static unsigned int get()
Return number of stored elements.
BoundaryIntegralData(unsigned int bdr_idx, DHCellSide dhside, unsigned int side_idx)
Constructor with data mebers initialization.
std::array< std::shared_ptr< CouplingIntegral >, 2 > coupling_
Coupling integrals between elements of dimensions 1-2, 2-3.
void add_integrals_of_computing_step(DHCellAccessor cell)
std::shared_ptr< EvalPoints > eval_points_
EvalPoints object shared by all integrals.
Set of all used integral necessary in assemblation.
void create_patch()
Create patch of cached elements before reading data to cache.
std::size_t make_permanent()
Finalize temporary part of data.
BulkIntegralData(DHCellAccessor dhcell, unsigned int subset_idx)
Constructor with data mebers initialization.
unsigned int bulk_subset_index
Index (order) of lower dim subset in EvalPoints object.
void set_min_edge_sides(unsigned int val)
std::size_t emplace_back(Args &&... args)
RangeConvert< DHEdgeSide, DHCellSide > edge_sides() const
Returns range of all sides looped over common Edge.
unsigned int side_subset_index
Index (order) of subset on side of bulk element in EvalPoints object.
CouplingIntegralData()
Default constructor.
RevertableList< EdgeIntegralData > edge_integral_data_
Holds data for computing edge integrals.
Directing class of FieldValueCache.
bool is_own() const
Return true if accessor represents own element (false for ghost element)
std::array< std::shared_ptr< EdgeIntegral >, 3 > edge_
Edge integrals between elements of dimensions 1, 2, 3.
Base point accessor class.
void reallocate_cache()
Calls cache_reallocate method on.
unsigned int dim() const
Return dimension of element appropriate to cell.
void reset()
Clear the list.
unsigned int idx() const
Returns a global index of the region.
unsigned int dim() const
Return dimension of element appropriate to the side.
GenericAssembly(typename DimAssembly< 1 >::EqFields *eq_fields, typename DimAssembly< 1 >::EqData *eq_data)
Constructor.
RevertableList< EvalPointData > eval_point_data_
virtual void assemble(std::shared_ptr< DOFHandlerMultiDim > dh)=0
void assemble(std::shared_ptr< DOFHandlerMultiDim > dh) override
General assemble methods.
void add_edge_integral(const DHCellSide &cell_side)
Add data of edge integral to appropriate data structure.
void init(std::shared_ptr< EvalPoints > eval_points)
Init cache.
std::size_t temporary_size() const
Return temporary size of list (full size of stored data).
unsigned int subset_index
Index (order) of subset in EvalPoints object.
RevertableList< BoundaryIntegralData > boundary_integral_data_
Holds data for computing boundary integrals.
Side accessor allows to iterate over sides of DOF handler cell.
DHCellSide side
Specified cell side (bulk element)
Range< DHCellSide > side_range() const
Returns range of cell sides.
BoundaryIntegralData(const BoundaryIntegralData &other)
Copy constructor.
EdgeIntegralData()
Default constructor.
int active_integrals_
Holds mask of active integrals.
CouplingIntegralData(DHCellAccessor dhcell, unsigned int bulk_idx, DHCellSide dhside, unsigned int side_idx)
Constructor with data mebers initialization.
const DHCellAccessor & cell() const
Return DHCellAccessor appropriate to the side.
RangeConvert< DHNeighbSide, DHCellSide > neighb_sides() const
Returns range of neighbour cell of lower dimension corresponding to cell of higher dimension.
std::array< std::shared_ptr< BulkIntegral >, 3 > bulk_
Bulk integrals of elements of dimensions 1, 2, 3.
Definitions of particular quadrature rules on simplices.
DHCellSide side
Specified cell side (higher dim element)
void clear_element_eval_points_map()
Reset all items of elements_eval_points_map.
unsigned int local_idx() const
Return local index to element (index of DOF handler).
MixedPtr< DimAssembly, 1 > multidim_assembly() const
Getter to set of assembly objects.
AssemblyIntegrals integrals_
Holds integral objects.
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
void add_coupling_integral(const DHCellAccessor &cell, const DHCellSide &ngh_side, bool add_low)
Add data of coupling integral to appropriate data structure.
void add_volume_integral(const DHCellAccessor &cell)
Add data of volume integral to appropriate data structure.
ActiveIntegrals
Allow set mask of active integrals.
RevertableList< CouplingIntegralData > coupling_integral_data_
Holds data for computing couplings integrals.
unsigned int eval_point_idx() const
Return index in EvalPoints object.
ElementAccessor< 3 > element_accessor()
unsigned int elm_idx() const
Return serial idx to element of loc_ele_idx_.
BulkIntegralData()
Default constructor.
void finish_elements_update()
Finish update after reading data to cache.
ElementAccessor< 3 > element() const
CouplingIntegralData(const CouplingIntegralData &other)
Copy constructor.
unsigned int elem_idx() const
const ElementCacheMap & cache_map() const
Return ElementCacheMap.
void assemble_integrals()
Call assemblations when patch is filled.
unsigned int subset_index
Index (order) of subset in EvalPoints object.
Cell accessor allow iterate over DOF handler cells.
BoundaryIntegralData()
Default constructor.
RangeConvert< DHEdgeSide, DHCellSide > edge_side_range
Specified cell side (element)
std::size_t revert_temporary()
Erase temporary part of data.
Class allows to iterate over sides of edge.
unsigned int min_edge_sides_
Iter< Object > make_iter(Object obj)
std::shared_ptr< EvalPoints > eval_points() const
Geter to EvalPoints object.
EdgeIntegralData(const EdgeIntegralData &other)
Copy constructor.
DHCellAccessor cell
Specified cell (element)
ElementCacheMap element_cache_map_
ElementCacheMap according to EvalPoints.
void start_elements_update()
Start update of cache.
BulkIntegralData(const BulkIntegralData &other)
Copy constructor.
unsigned int bdr_subset_index
Index (order) of subset on boundary element in EvalPoints object.
virtual ~GenericAssemblyBase()
#define START_TIMER(tag)
Starts a timer with specified tag.
Generic class of assemblation.
EdgeIntegralData(RangeConvert< DHEdgeSide, DHCellSide > range, unsigned int subset_idx)
Constructor with data mebers initialization.
#define END_TIMER(tag)
Ends a timer with specified tag.
RegionIdx region_idx() const
MixedPtr< DimAssembly, 1 > multidim_assembly_
Assembly object.
std::array< std::shared_ptr< BoundaryIntegral >, 3 > boundary_
Boundary integrals betwwen elements of dimensions 1, 2, 3 and boundaries.
unsigned int IntDim
A dimension index type.