Flow123d  JS_constraints-e651b99
generic_assembly.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file generic_assembly.hh
15  * @brief
16  */
17 
18 #ifndef GENERIC_ASSEMBLY_HH_
19 #define GENERIC_ASSEMBLY_HH_
20 
22 #include "fem/integral_acc.hh"
23 #include "fem/eval_points.hh"
24 #include "fem/element_cache_map.hh"
25 #include "fem/patch_fe_values.hh"
26 #include "tools/revertable_list.hh"
27 #include "system/sys_profiler.hh"
28 
29 
30 
31 /// Holds common data shared between GenericAssemblz and Assembly<dim> classes.
33 public:
35  : eval_points_(std::make_shared<EvalPoints>()) {}
36 
38  : eval_points_(std::make_shared<EvalPoints>()), fe_values_(fe) {}
39 
40  std::shared_ptr<EvalPoints> eval_points_; ///< EvalPoints object shared by all integrals
41  ElementCacheMap element_cache_map_; ///< ElementCacheMap according to EvalPoints
42  PatchFEValues<3> fe_values_; ///< Common FEValues object over all dimensions
43 };
44 
45 
46 /**
47  * Common interface class for all Assembly classes.
48  */
50 {
51 public:
53  {}
54 
56  : asm_internals_(fe)
57  {}
58 
60 
61  virtual void assemble(std::shared_ptr<DOFHandlerMultiDim> dh) = 0;
62 
63  /// Getter to EvalPoints object
64  inline std::shared_ptr<EvalPoints> eval_points() const {
66  }
67 
68 protected:
69  AssemblyInternals asm_internals_; ///< Holds shared internals data
70 };
71 
72 
73 /**
74  * @brief Generic class of assemblation.
75  *
76  * Class
77  * - holds assemblation structures (EvalPoints, Integral objects, Integral data table).
78  * - associates assemblation objects specified by dimension
79  * - provides general assemble method
80  * - provides methods that allow construction of element patches
81  */
82 template < template<IntDim...> class DimAssembly>
84 {
85 public:
86  /**
87  * Constructor.
88  *
89  * Used in equations that don't need PatchFeValues objects in evaluation.
90  * @param eq_fields Descendant of FieldSet declared in equation
91  * @param eq_data Object defined in equation containing shared data of eqation and assembly class.
92  */
93  GenericAssembly( typename DimAssembly<1>::EqData *eq_data)
95  multidim_assembly_(eq_data, &this->asm_internals_)
96  {
97  initialize();
98  }
99 
100  /**
101  * Constructor.
102  *
103  * Used in equations working with PatchFeValues objects in evaluation.
104  * @param eq_fields Descendant of FieldSet declared in equation
105  * @param eq_data Object defined in equation containing shared data of eqation and assembly class.
106  * @param dh DOF handler object
107  */
108  GenericAssembly( typename DimAssembly<1>::EqData *eq_data, DOFHandlerMultiDim* dh)
109  : GenericAssemblyBase(dh->ds()->fe()),
110  multidim_assembly_(eq_data, &this->asm_internals_)
111  {
112  initialize();
113  }
114 
115  /// Getter to set of assembly objects
117  return multidim_assembly_;
118  }
119 
120  /// Allows to rewrite number of minimal edge sides.
121  void set_min_edge_sides(unsigned int val) {
122  multidim_assembly_[1_d]->set_min_edge_sides(val);
123  multidim_assembly_[2_d]->set_min_edge_sides(val);
124  multidim_assembly_[3_d]->set_min_edge_sides(val);
125  }
126 
127  /**
128  * @brief General assemble method.
129  *
130  * Loops through local cells and calls assemble methods of assembly
131  * object of each cells over space dimension.
132  *
133  * TODO:
134  * - make estimate of the cache fill for combination of (integral_type x element dimension)
135  * - add next cell to patch if current_patch_size + next_element_size <= fixed_cache_size
136  * - avoid reverting the integral data lists.
137  */
138  void assemble(std::shared_ptr<DOFHandlerMultiDim> dh) override {
139  START_TIMER( DimAssembly<1>::name() );
140  this->reallocate_cache();
141  multidim_assembly_[1_d]->begin();
142 
143  bool add_into_patch = false; // control variable
144  for(auto cell_it = dh->local_range().begin(); cell_it != dh->local_range().end(); )
145  {
146  unsigned int cell_dim = cell_it->dim();
147  if (!add_into_patch) {
149  add_into_patch = true;
150  }
151 
152  START_TIMER("add_integrals_to_patch");
153  bool is_patch_full = false;
154  switch ( cell_dim ) {
155  case 1:
156  is_patch_full = multidim_assembly_[1_d]->add_integrals_of_computing_step(*cell_it);
157  break;
158  case 2:
159  is_patch_full = multidim_assembly_[2_d]->add_integrals_of_computing_step(*cell_it);
160  break;
161  case 3:
162  is_patch_full = multidim_assembly_[3_d]->add_integrals_of_computing_step(*cell_it);
163  break;
164  default:
165  ASSERT(false).error("Should not happen!");
166  }
167  END_TIMER("add_integrals_to_patch");
168 
169  if (is_patch_full) {
171  this->assemble_integrals();
172  add_into_patch = false;
173  } else {
177  this->assemble_integrals();
178  add_into_patch = false;
179  }
180  ++cell_it;
181  }
182  }
183  if (add_into_patch) {
184  this->assemble_integrals();
185  }
186 
187  multidim_assembly_[1_d]->end();
188  END_TIMER( DimAssembly<1>::name() );
189  }
190 
191  /// Return ElementCacheMap
192  inline const ElementCacheMap &cache_map() const {
194  }
195 
196 private:
197  /// Common part of GenericAssemblz constructors.
198  void initialize() {
200  multidim_assembly_[1_d]->initialize();
201  multidim_assembly_[2_d]->initialize();
202  multidim_assembly_[3_d]->initialize();
204  }
205 
206  /// Call assemblations when patch is filled
208  START_TIMER("create_patch");
210  END_TIMER("create_patch");
211 
212  START_TIMER("cache_update");
213  multidim_assembly_[1_d]->eq_fields_->cache_update(asm_internals_.element_cache_map_); // TODO replace with sub FieldSet
214  END_TIMER("cache_update");
215 
216  START_TIMER("patch_reinit");
217  patch_reinit(); // reinit PatchFeValues
218  END_TIMER("patch_reinit");
219 
221 
222  {
223  START_TIMER("assemble_volume_integrals");
224  multidim_assembly_[1_d]->assemble_cell_integrals();
225  multidim_assembly_[2_d]->assemble_cell_integrals();
226  multidim_assembly_[3_d]->assemble_cell_integrals();
227  END_TIMER("assemble_volume_integrals");
228  }
229 
230  {
231  START_TIMER("assemble_fluxes_boundary");
232  multidim_assembly_[1_d]->assemble_boundary_side_integrals();
233  multidim_assembly_[2_d]->assemble_boundary_side_integrals();
234  multidim_assembly_[3_d]->assemble_boundary_side_integrals();
235  END_TIMER("assemble_fluxes_boundary");
236  }
237 
238  {
239  START_TIMER("assemble_fluxes_elem_elem");
240  multidim_assembly_[1_d]->assemble_edge_integrals();
241  multidim_assembly_[2_d]->assemble_edge_integrals();
242  multidim_assembly_[3_d]->assemble_edge_integrals();
243  END_TIMER("assemble_fluxes_elem_elem");
244  }
245 
246  {
247  START_TIMER("assemble_fluxes_elem_side");
248  multidim_assembly_[1_d]->assemble_neighbour_integrals();
249  multidim_assembly_[2_d]->assemble_neighbour_integrals();
250  END_TIMER("assemble_fluxes_elem_side");
251  }
252  // clean integral data
253  multidim_assembly_[1_d]->clean_integral_data();
254  multidim_assembly_[2_d]->clean_integral_data();
255  multidim_assembly_[3_d]->clean_integral_data();
258  }
259 
260  /// Reinit PatchFeValues object during construction of patch
261  void patch_reinit() {
266 
268  }
269 
270  /// Calls cache_reallocate method on
271  inline void reallocate_cache() {
272  multidim_assembly_[1_d]->eq_fields_->cache_reallocate(asm_internals_.element_cache_map_, multidim_assembly_[1_d]->used_fields_);
273  // DebugOut() << "Order of evaluated fields (" << DimAssembly<1>::name() << "):" << multidim_assembly_[1_d]->eq_fields_->print_dependency();
274  }
275 
276 
278 };
279 
280 
281 #endif /* GENERIC_ASSEMBLY_HH_ */
#define ASSERT(expr)
Definition: asserts.hh:351
static unsigned int get()
Return number of stored elements.
Provides the numbering of the finite element degrees of freedom on the computational mesh.
Definition: dofhandler.hh:151
Directing class of FieldValueCache.
unsigned int get_simd_rounded_size()
Returns number of eval. points with addition of max simd duplicates due to regions.
void start_elements_update()
Start update of cache.
void finish_elements_update()
Finish update after reading data to cache.
void clear_element_eval_points_map()
Reset all items of elements_eval_points_map.
void create_patch()
Create patch of cached elements before reading data to cache.
RevertableList< EvalPointData > eval_point_data_
void init(std::shared_ptr< EvalPoints > eval_points)
Init cache.
void make_paermanent_eval_points()
Mark eval_point_data_ as permanent.
Class holds local coordinations of evaluating points (bulk and sides) specified by element dimension.
Definition: eval_points.hh:54
AssemblyInternals asm_internals_
Holds shared internals data.
virtual ~GenericAssemblyBase()
GenericAssemblyBase(MixedPtr< FiniteElement > fe)
std::shared_ptr< EvalPoints > eval_points() const
Getter to EvalPoints object.
virtual void assemble(std::shared_ptr< DOFHandlerMultiDim > dh)=0
Generic class of assemblation.
GenericAssembly(typename DimAssembly< 1 >::EqData *eq_data)
MixedPtr< DimAssembly, 1 > multidim_assembly_
Assembly object.
void assemble_integrals()
Call assemblations when patch is filled.
void assemble(std::shared_ptr< DOFHandlerMultiDim > dh) override
General assemble method.
void patch_reinit()
Reinit PatchFeValues object during construction of patch.
const ElementCacheMap & cache_map() const
Return ElementCacheMap.
GenericAssembly(typename DimAssembly< 1 >::EqData *eq_data, DOFHandlerMultiDim *dh)
void initialize()
Common part of GenericAssemblz constructors.
void set_min_edge_sides(unsigned int val)
Allows to rewrite number of minimal edge sides.
void reallocate_cache()
Calls cache_reallocate method on.
MixedPtr< DimAssembly, 1 > multidim_assembly() const
Getter to set of assembly objects.
void init_finalize()
Finalize initialization, creates child (patch) arena and passes it to PatchPointValue objects.
void add_patch_points(const DimIntegrals< dim > &integrals, ElementCacheMap *element_cache_map, std::shared_ptr< EvalPoints > eval_points)
Add elements, sides and quadrature points registered on patch.
void make_permanent_ppv_data()
Marks data of last successfully added element to patch as permanent.
void clean_elements_map()
Clear elements_map, set values to (-1)
void reinit_patch()
Reinit data.
void reset()
Reset PatchpointValues structures.
unsigned int IntDim
A dimension index type.
Definition: mixed.hh:19
Class FEValues calculates finite element data on the actual cells such as shape function values,...
Definitions of particular quadrature rules on simplices.
Holds common data shared between GenericAssemblz and Assembly<dim> classes.
PatchFEValues< 3 > fe_values_
Common FEValues object over all dimensions.
ElementCacheMap element_cache_map_
ElementCacheMap according to EvalPoints.
AssemblyInternals(MixedPtr< FiniteElement > fe)
std::shared_ptr< EvalPoints > eval_points_
EvalPoints object shared by all integrals.
std::size_t revert_temporary()
Erase temporary part of data.
#define END_TIMER(tag)
Ends a timer with specified tag.
#define START_TIMER(tag)
Starts a timer with specified tag.