Flow123d  master-921db77
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 #include "fem/patch_internals.hh"
29 
30 
31 
32 /**
33  * Common interface class for all Assembly classes.
34  */
36 {
37 public:
39  {}
40 
42  : patch_internals_(fe)
43  {}
44 
46 
47  virtual void assemble(std::shared_ptr<DOFHandlerMultiDim> dh) = 0;
48 
49  /// Getter to EvalPoints object
50  inline std::shared_ptr<EvalPoints> eval_points() const {
52  }
53 
54 protected:
55  PatchInternals patch_internals_; ///< Holds shared internals data
56 };
57 
58 
59 /**
60  * @brief Generic class of assemblation.
61  *
62  * Class
63  * - holds assemblation structures (EvalPoints, Integral objects, Integral data table).
64  * - associates assemblation objects specified by dimension
65  * - provides general assemble method
66  * - provides methods that allow construction of element patches
67  */
68 template < template<IntDim...> class DimAssembly>
70 {
71 public:
72  /**
73  * Constructor.
74  *
75  * Used in equations that don't need PatchFeValues objects in evaluation.
76  * @param eq_fields Descendant of FieldSet declared in equation
77  * @param eq_data Object defined in equation containing shared data of eqation and assembly class.
78  */
79  GenericAssembly( typename DimAssembly<1>::EqData *eq_data)
81  multidim_assembly_(eq_data, &this->patch_internals_)
82  {
83  initialize();
84  }
85 
86  /**
87  * Constructor.
88  *
89  * Used in equations working with 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  * @param dh DOF handler object
93  */
94  GenericAssembly( typename DimAssembly<1>::EqData *eq_data, DOFHandlerMultiDim* dh)
95  : GenericAssemblyBase(dh->ds()->fe()),
96  multidim_assembly_(eq_data, &this->patch_internals_)
97  {
98  initialize();
99  }
100 
101  /// Getter to set of assembly objects
103  return multidim_assembly_;
104  }
105 
106  /// Allows to rewrite number of minimal edge sides.
107  void set_min_edge_sides(unsigned int val) {
108  multidim_assembly_[1_d]->set_min_edge_sides(val);
109  multidim_assembly_[2_d]->set_min_edge_sides(val);
110  multidim_assembly_[3_d]->set_min_edge_sides(val);
111  }
112 
113  /**
114  * @brief General assemble method.
115  *
116  * Loops through local cells and calls assemble methods of assembly
117  * object of each cells over space dimension.
118  *
119  * TODO:
120  * - make estimate of the cache fill for combination of (integral_type x element dimension)
121  * - add next cell to patch if current_patch_size + next_element_size <= fixed_cache_size
122  * - avoid reverting the integral data lists.
123  */
124  void assemble(std::shared_ptr<DOFHandlerMultiDim> dh) override {
125  START_TIMER( DimAssembly<1>::name() );
126  multidim_assembly_[1_d]->eq_fields_->cache_reallocate(patch_internals_, multidim_assembly_[1_d]->used_fields_);
127  multidim_assembly_[1_d]->begin();
128 
129  bool add_into_patch = false; // control variable
130  for(auto cell_it = dh->local_range().begin(); cell_it != dh->local_range().end(); )
131  {
132  unsigned int cell_dim = cell_it->dim();
133  if (!add_into_patch) {
135  add_into_patch = true;
136  }
137 
138  START_TIMER("add_integrals_to_patch");
139  bool is_patch_full = false;
140  switch ( cell_dim ) {
141  case 1:
142  is_patch_full = multidim_assembly_[1_d]->add_integrals_of_computing_step(*cell_it);
143  break;
144  case 2:
145  is_patch_full = multidim_assembly_[2_d]->add_integrals_of_computing_step(*cell_it);
146  break;
147  case 3:
148  is_patch_full = multidim_assembly_[3_d]->add_integrals_of_computing_step(*cell_it);
149  break;
150  default:
151  ASSERT(false).error("Should not happen!");
152  }
153  END_TIMER("add_integrals_to_patch");
154 
155  if (is_patch_full) {
157  this->assemble_integrals();
158  add_into_patch = false;
159  } else {
163  this->assemble_integrals();
164  add_into_patch = false;
165  }
166  ++cell_it;
167  }
168  }
169  if (add_into_patch) {
170  this->assemble_integrals();
171  }
172 
173  multidim_assembly_[1_d]->end();
174  END_TIMER( DimAssembly<1>::name() );
175  }
176 
177  /// Return ElementCacheMap
178  inline const ElementCacheMap &cache_map() const {
180  }
181 
182 private:
183  /// Common part of GenericAssemblz constructors.
184  void initialize() {
186  multidim_assembly_[1_d]->initialize();
187  multidim_assembly_[2_d]->initialize();
188  multidim_assembly_[3_d]->initialize();
190  }
191 
192  /// Call assemblations when patch is filled
194  START_TIMER("create_patch");
196  END_TIMER("create_patch");
197 
198  START_TIMER("cache_update");
199  multidim_assembly_[1_d]->eq_fields_->cache_update(patch_internals_.element_cache_map_); // TODO replace with sub FieldSet
200  END_TIMER("cache_update");
201 
202  START_TIMER("patch_reinit");
203  patch_reinit(); // reinit PatchFeValues
204  END_TIMER("patch_reinit");
205 
207 
208  {
209  START_TIMER("assemble_volume_integrals");
210  multidim_assembly_[1_d]->assemble_cell_integrals();
211  multidim_assembly_[2_d]->assemble_cell_integrals();
212  multidim_assembly_[3_d]->assemble_cell_integrals();
213  END_TIMER("assemble_volume_integrals");
214  }
215 
216  {
217  START_TIMER("assemble_fluxes_boundary");
218  multidim_assembly_[1_d]->assemble_boundary_side_integrals();
219  multidim_assembly_[2_d]->assemble_boundary_side_integrals();
220  multidim_assembly_[3_d]->assemble_boundary_side_integrals();
221  END_TIMER("assemble_fluxes_boundary");
222  }
223 
224  {
225  START_TIMER("assemble_fluxes_elem_elem");
226  multidim_assembly_[1_d]->assemble_edge_integrals();
227  multidim_assembly_[2_d]->assemble_edge_integrals();
228  multidim_assembly_[3_d]->assemble_edge_integrals();
229  END_TIMER("assemble_fluxes_elem_elem");
230  }
231 
232  {
233  START_TIMER("assemble_fluxes_elem_side");
234  multidim_assembly_[1_d]->assemble_neighbour_integrals();
235  multidim_assembly_[2_d]->assemble_neighbour_integrals();
236  END_TIMER("assemble_fluxes_elem_side");
237  }
238  // clean integral data
239  multidim_assembly_[1_d]->clean_integral_data();
240  multidim_assembly_[2_d]->clean_integral_data();
241  multidim_assembly_[3_d]->clean_integral_data();
244  }
245 
246  /// Reinit PatchFeValues object during construction of patch
247  void patch_reinit() {
252 
254  }
255 
256 
258 };
259 
260 
261 #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.
virtual ~GenericAssemblyBase()
GenericAssemblyBase(MixedPtr< FiniteElement > fe)
PatchInternals patch_internals_
Holds shared internals data.
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.
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.
std::shared_ptr< EvalPoints > eval_points_
EvalPoints object shared by all integrals.
ElementCacheMap element_cache_map_
ElementCacheMap according to EvalPoints.
PatchFEValues< 3 > fe_values_
Common FEValues object over all dimensions.
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.