Flow123d  JS_before_hm-984-g3a19f2f
assembly_dg.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 assembly_dg.hh
15  * @brief
16  */
17 
18 #ifndef ASSEMBLY_DG_HH_
19 #define ASSEMBLY_DG_HH_
20 
22 #include "fem/mapping_p1.hh"
23 #include "fem/fe_p.hh"
24 #include "fem/fe_rt.hh"
25 #include "fem/fe_values.hh"
27 #include "coupling/balance.hh"
28 #include "fields/eval_subset.hh"
29 #include "fields/eval_points.hh"
31 
32 
33 
34 /// Allow set mask of active integrals.
36  none = 0,
37  bulk = 0x0001,
38  edge = 0x0002,
39  coupling = 0x0004,
40  boundary = 0x0008
41 };
42 
43 
44 /// Set of all used integral necessary in assemblation
46  std::array<std::shared_ptr<BulkIntegral>, 3> bulk_; ///< Bulk integrals of elements of dimensions 1, 2, 3
47  std::array<std::shared_ptr<EdgeIntegral>, 3> edge_; ///< Edge integrals between elements of dimensions 1, 2, 3
48  std::array<std::shared_ptr<CouplingIntegral>, 2> coupling_; ///< Coupling integrals between elements of dimensions 1-2, 2-3
49  std::array<std::shared_ptr<BoundaryIntegral>, 3> boundary_; ///< Boundary integrals betwwen elements of dimensions 1, 2, 3 and boundaries
50 };
51 
52 
53 
54 /**
55  * @brief Generic class of assemblation.
56  *
57  * Class
58  * - holds assemblation structures (EvalPoints, Integral objects, Integral data table).
59  * - associates assemblation objects specified by dimension
60  * - provides general assemble method
61  * - provides methods that allow construction of element patches
62  */
63 template < template<IntDim...> class DimAssembly>
65 {
66 private:
69 
71  unsigned int subset_index;
72  };
73 
77 
79  unsigned int subset_index;
80  };
81 
84 
86  unsigned int bulk_subset_index;
88  unsigned int side_subset_index;
89  };
90 
93 
95  unsigned int subset_index;
96  };
97 
98 public:
99 
100  /// Constructor
101  GenericAssembly( typename DimAssembly<1>::EqDataDG *eq_data, int active_integrals )
102  : multidim_assembly_(eq_data),
103  active_integrals_(active_integrals), integrals_size_({0, 0, 0, 0})
104  {
105  eval_points_ = std::make_shared<EvalPoints>();
106  // first step - create integrals, then - initialize cache
107  multidim_assembly_[1_d]->create_integrals(eval_points_, integrals_, active_integrals_);
108  multidim_assembly_[2_d]->create_integrals(eval_points_, integrals_, active_integrals_);
109  multidim_assembly_[3_d]->create_integrals(eval_points_, integrals_, active_integrals_);
110  element_cache_map_.init(eval_points_);
111  }
112 
114  return multidim_assembly_;
115  }
116 
117  inline std::shared_ptr<EvalPoints> eval_points() const {
118  return eval_points_;
119  }
120 
121  /**
122  * @brief General assemble methods.
123  *
124  * Loops through local cells and calls assemble methods of assembly
125  * object of each cells over space dimension.
126  */
127  void assemble(std::shared_ptr<DOFHandlerMultiDim> dh) {
128  unsigned int i;
129  multidim_assembly_[1_d]->begin();
130  for (auto cell : dh->local_range() )
131  {
132  this->add_integrals_of_computing_step(cell);
133 
134  if ( cell.is_own() && (active_integrals_ & ActiveIntegrals::bulk) ) {
135  START_TIMER("assemble_volume_integrals");
136  for (i=0; i<integrals_size_[0]; ++i) { // volume integral
137  switch (bulk_integral_data_[i].cell.dim()) {
138  case 1:
139  multidim_assembly_[1_d]->assemble_volume_integrals(bulk_integral_data_[i].cell);
140  break;
141  case 2:
142  multidim_assembly_[2_d]->assemble_volume_integrals(bulk_integral_data_[i].cell);
143  break;
144  case 3:
145  multidim_assembly_[3_d]->assemble_volume_integrals(bulk_integral_data_[i].cell);
146  break;
147  }
148  }
149  END_TIMER("assemble_volume_integrals");
150  }
151 
152  if ( cell.is_own() && (active_integrals_ & ActiveIntegrals::boundary) ) {
153  START_TIMER("assemble_fluxes_boundary");
154  for (i=0; i<integrals_size_[3]; ++i) { // boundary integral
155  switch (boundary_integral_data_[i].side.dim()) {
156  case 1:
157  multidim_assembly_[1_d]->assemble_fluxes_boundary(boundary_integral_data_[i].side);
158  break;
159  case 2:
160  multidim_assembly_[2_d]->assemble_fluxes_boundary(boundary_integral_data_[i].side);
161  break;
162  case 3:
163  multidim_assembly_[3_d]->assemble_fluxes_boundary(boundary_integral_data_[i].side);
164  break;
165  }
166  }
167  END_TIMER("assemble_fluxes_boundary");
168  }
169 
170  if (active_integrals_ & ActiveIntegrals::edge) {
171  START_TIMER("assemble_fluxes_elem_elem");
172  for (i=0; i<integrals_size_[1]; ++i) { // edge integral
173  switch (edge_integral_data_[i].edge_side_range.begin()->dim()) {
174  case 1:
175  multidim_assembly_[1_d]->assemble_fluxes_element_element(edge_integral_data_[i].edge_side_range);
176  break;
177  case 2:
178  multidim_assembly_[2_d]->assemble_fluxes_element_element(edge_integral_data_[i].edge_side_range);
179  break;
180  case 3:
181  multidim_assembly_[3_d]->assemble_fluxes_element_element(edge_integral_data_[i].edge_side_range);
182  break;
183  }
184  }
185  END_TIMER("assemble_fluxes_elem_elem");
186  }
187 
188  if (active_integrals_ & ActiveIntegrals::coupling) {
189  START_TIMER("assemble_fluxes_elem_side");
190  for (i=0; i<integrals_size_[2]; ++i) { // coupling integral
191  switch (coupling_integral_data_[i].side.dim()) {
192  case 2:
193  multidim_assembly_[2_d]->assemble_fluxes_element_side(coupling_integral_data_[i].cell, coupling_integral_data_[i].side);
194  break;
195  case 3:
196  multidim_assembly_[3_d]->assemble_fluxes_element_side(coupling_integral_data_[i].cell, coupling_integral_data_[i].side);
197  break;
198  }
199  }
200  END_TIMER("assemble_fluxes_elem_side");
201  }
202  }
203  multidim_assembly_[1_d]->end();
204  }
205 
206 private:
207  /// Mark eval points in table of Element cache map.
209  for (unsigned int i=0; i<integrals_size_[0]; ++i) {
210  // add data to cache if there is free space, else return
211  unsigned int data_size = eval_points_->subset_size( bulk_integral_data_[i].cell.dim(), bulk_integral_data_[i].subset_index );
212  element_cache_map_.mark_used_eval_points(bulk_integral_data_[i].cell, bulk_integral_data_[i].subset_index, data_size);
213  }
214 
215  for (unsigned int i=0; i<integrals_size_[1]; ++i) {
216  // add data to cache if there is free space, else return
217  for (DHCellSide edge_side : edge_integral_data_[i].edge_side_range) {
218  unsigned int side_dim = edge_side.dim();
219  unsigned int data_size = eval_points_->subset_size( side_dim, edge_integral_data_[i].subset_index ) / (side_dim+1);
220  unsigned int start_point = data_size * edge_side.side_idx();
221  element_cache_map_.mark_used_eval_points(edge_side.cell(), edge_integral_data_[i].subset_index, data_size, start_point);
222  }
223  }
224 
225  for (unsigned int i=0; i<integrals_size_[2]; ++i) {
226  // add data to cache if there is free space, else return
227  unsigned int bulk_data_size = eval_points_->subset_size( coupling_integral_data_[i].cell.dim(), coupling_integral_data_[i].bulk_subset_index );
228  element_cache_map_.mark_used_eval_points(coupling_integral_data_[i].cell, coupling_integral_data_[i].bulk_subset_index, bulk_data_size);
229 
230  unsigned int side_dim = coupling_integral_data_[i].side.dim();
231  unsigned int side_data_size = eval_points_->subset_size( side_dim, coupling_integral_data_[i].side_subset_index ) / (side_dim+1);
232  unsigned int start_point = side_data_size * coupling_integral_data_[i].side.side_idx();
233  element_cache_map_.mark_used_eval_points(coupling_integral_data_[i].side.cell(), coupling_integral_data_[i].side_subset_index, side_data_size, start_point);
234  }
235 
236  for (unsigned int i=0; i<integrals_size_[3]; ++i) {
237  // add data to cache if there is free space, else return
238  unsigned int side_dim = boundary_integral_data_[i].side.dim();
239  unsigned int data_size = eval_points_->subset_size( side_dim, boundary_integral_data_[i].subset_index ) / (side_dim+1);
240  unsigned int start_point = data_size * boundary_integral_data_[i].side.side_idx();
241  element_cache_map_.mark_used_eval_points(boundary_integral_data_[i].side.cell(), boundary_integral_data_[i].subset_index, data_size, start_point);
242  }
243  }
244 
245  /**
246  * Add data of integrals to appropriate structure and register elements to ElementCacheMap.
247  *
248  * Types of used integrals must be set in data member \p active_integrals_.
249  */
251  for (unsigned int i=0; i<4; i++) integrals_size_[i] = 0; // clean integral data from previous step
252  element_cache_map_.start_elements_update();
253 
254  // generic_assembly.check_integral_data();
255  if (active_integrals_ & ActiveIntegrals::bulk)
256  if (cell.is_own()) { // Not ghost
257  this->add_volume_integral(cell);
258  element_cache_map_.add(cell);
259  }
260 
261  for( DHCellSide cell_side : cell.side_range() ) {
262  if (active_integrals_ & ActiveIntegrals::boundary)
263  if (cell.is_own()) // Not ghost
264  if ( (cell_side.side().edge().n_sides() == 1) && (cell_side.side().is_boundary()) ) {
265  this->add_boundary_integral(cell_side);
266  element_cache_map_.add(cell_side);
267  continue;
268  }
269  if (active_integrals_ & ActiveIntegrals::edge)
270  if ( (cell_side.n_edge_sides() >= 2) && (cell_side.edge_sides().begin()->element().idx() == cell.elm_idx())) {
271  this->add_edge_integral(cell_side.edge_sides());
272  for( DHCellSide edge_side : cell_side.edge_sides() ) {
273  element_cache_map_.add(edge_side);
274  }
275  }
276  }
277 
278  if (active_integrals_ & ActiveIntegrals::coupling)
279  for( DHCellSide neighb_side : cell.neighb_sides() ) { // cell -> elm lower dim, neighb_side -> elm higher dim
280  if (cell.dim() != neighb_side.dim()-1) continue;
281  this->add_coupling_integral(cell, neighb_side);
282  element_cache_map_.add(cell);
283  element_cache_map_.add(neighb_side);
284  }
285 
286  element_cache_map_.prepare_elements_to_update();
287  this->insert_eval_points_from_integral_data();
288  element_cache_map_.create_elements_points_map();
289  // not used yet: TODO need fix in MultiField, HeatModel ...; need better access to EqData
290  //multidim_assembly_[1]->data_->cache_update(element_cache_map_);
291  element_cache_map_.finish_elements_update();
292  }
293 
294  /// Add data of volume integral to appropriate data structure.
296  bulk_integral_data_[ integrals_size_[0] ].cell = cell;
297  bulk_integral_data_[ integrals_size_[0] ].subset_index = integrals_.bulk_[cell.dim()-1]->get_subset_idx();
298  integrals_size_[0]++;
299  }
300 
301  /// Add data of edge integral to appropriate data structure.
303  edge_integral_data_[ integrals_size_[1] ].edge_side_range = edge_side_range;
304  edge_integral_data_[ integrals_size_[1] ].subset_index = integrals_.edge_[edge_side_range.begin()->dim()-1]->get_subset_idx();
305  integrals_size_[1]++;
306  }
307 
308  /// Add data of coupling integral to appropriate data structure.
309  void add_coupling_integral(const DHCellAccessor &cell, const DHCellSide &ngh_side) {
310  coupling_integral_data_[ integrals_size_[2] ].cell = cell;
311  coupling_integral_data_[ integrals_size_[2] ].side = ngh_side;
312  coupling_integral_data_[ integrals_size_[2] ].bulk_subset_index = integrals_.coupling_[cell.dim()-1]->get_subset_low_idx();
313  coupling_integral_data_[ integrals_size_[2] ].side_subset_index = integrals_.coupling_[cell.dim()-1]->get_subset_high_idx();
314  integrals_size_[2]++;
315  }
316 
317  /// Add data of boundary integral to appropriate data structure.
318  void add_boundary_integral(const DHCellSide &bdr_side) {
319  boundary_integral_data_[ integrals_size_[3] ].side = bdr_side;
320  boundary_integral_data_[ integrals_size_[3] ].subset_index = integrals_.boundary_[bdr_side.dim()-1]->get_subset_idx();
321  integrals_size_[3]++;
322  }
323 
324 
325  /// Assembly object
327 
328  /// Holds mask of active integrals.
330 
331  AssemblyIntegrals integrals_; ///< Holds integral objects.
332  std::shared_ptr<EvalPoints> eval_points_; ///< EvalPoints object shared by all integrals
333  ElementCacheMap element_cache_map_; ///< ElementCacheMap according to EvalPoints
334 
335  // Following variables hold data of all integrals depending of actual computed element.
336  std::array<BulkIntegralData, 1> bulk_integral_data_; ///< Holds data for computing bulk integrals.
337  std::array<EdgeIntegralData, 4> edge_integral_data_; ///< Holds data for computing edge integrals.
338  std::array<CouplingIntegralData, 6> coupling_integral_data_; ///< Holds data for computing couplings integrals.
339  std::array<BoundaryIntegralData, 4> boundary_integral_data_; ///< Holds data for computing boundary integrals.
340  std::array<unsigned int, 4> integrals_size_; ///< Holds used sizes of previous integral data types
341 };
342 
343 
344 /**
345  * Base class define empty methods, these methods can be overwite in descendants.
346  */
347 template <unsigned int dim>
348 class AssemblyBase
349 {
350 public:
351  /// Constructor
352  AssemblyBase(unsigned int quad_order) {
353  quad_ = new QGauss(dim, 2*quad_order);
354  quad_low_ = new QGauss(dim-1, 2*quad_order);
355  }
356 
357  // Destructor
358  virtual ~AssemblyBase() {
359  delete quad_;
360  delete quad_low_;
361  }
362 
363  /// Assembles the volume integrals on cell.
365 
366  /// Assembles the fluxes on the boundary.
367  virtual void assemble_fluxes_boundary(FMT_UNUSED DHCellSide cell_side) {}
368 
369  /// Assembles the fluxes between sides on the edge.
371 
372  /// Assembles the fluxes between elements of different dimensions.
373  virtual void assemble_fluxes_element_side(FMT_UNUSED DHCellAccessor cell_lower_dim, FMT_UNUSED DHCellSide neighb_side) {}
374 
375  /// Method prepares object before assemblation (e.g. balance, ...).
376  virtual void begin() {}
377 
378  /// Method finishes object after assemblation (e.g. balance, ...).
379  virtual void end() {}
380 
381  /// Create integrals according to dim of assembly object
382  void create_integrals(std::shared_ptr<EvalPoints> eval_points, AssemblyIntegrals &integrals, int active_integrals) {
383  if (active_integrals & ActiveIntegrals::bulk)
384  integrals.bulk_[dim-1] = eval_points->add_bulk<dim>(*quad_);
385  if (active_integrals & ActiveIntegrals::edge)
386  integrals.edge_[dim-1] = eval_points->add_edge<dim>(*quad_low_);
387  if ((dim>1) && (active_integrals & ActiveIntegrals::coupling))
388  integrals.coupling_[dim-2] = eval_points->add_coupling<dim>(*quad_low_);
389  if (active_integrals & ActiveIntegrals::boundary)
390  integrals.boundary_[dim-1] = eval_points->add_boundary<dim>(*quad_low_);
391  }
392 
393 protected:
394  Quadrature *quad_; ///< Quadrature used in assembling methods.
395  Quadrature *quad_low_; ///< Quadrature used in assembling methods (dim-1).
396 };
397 
398 
399 /**
400  * Auxiliary container class for Finite element and related objects of given dimension.
401  */
402 template <unsigned int dim, class Model>
403 class MassAssemblyDG : public AssemblyBase<dim>
404 {
405 public:
407 
408  /// Constructor.
409  MassAssemblyDG(EqDataDG *data)
410  : AssemblyBase<dim>(data->dg_order), model_(nullptr), data_(data) {}
411 
412  /// Destructor.
414 
415  /// Initialize auxiliary vectors and other data members
417  this->model_ = &model;
418 
419  fe_ = std::make_shared< FE_P_disc<dim> >(data_->dg_order);
421  ndofs_ = fe_->n_dofs();
422  qsize_ = this->quad_->size();
423  dof_indices_.resize(ndofs_);
424  local_matrix_.resize(4*ndofs_*ndofs_);
425  local_retardation_balance_vector_.resize(ndofs_);
426  local_mass_balance_vector_.resize(ndofs_);
427 
428  mm_coef_.resize(qsize_);
429  ret_coef_.resize(model_->n_substances());
430  for (unsigned int sbi=0; sbi<model_->n_substances(); sbi++)
431  {
432  ret_coef_[sbi].resize(qsize_);
433  }
434  }
435 
436 
437  /// Assemble integral over element
439  {
440  ASSERT_EQ_DBG(cell.dim(), dim).error("Dimension of element mismatch!");
441  ElementAccessor<3> elm = cell.elm();
442 
443  fe_values_.reinit(elm);
444  cell.get_dof_indices(dof_indices_);
445 
446  model_->compute_mass_matrix_coefficient(fe_values_.point_list(), elm, mm_coef_);
447  model_->compute_retardation_coefficient(fe_values_.point_list(), elm, ret_coef_);
448 
449  for (unsigned int sbi=0; sbi<model_->n_substances(); ++sbi)
450  {
451  // assemble the local mass matrix
452  for (unsigned int i=0; i<ndofs_; i++)
453  {
454  for (unsigned int j=0; j<ndofs_; j++)
455  {
456  local_matrix_[i*ndofs_+j] = 0;
457  for (unsigned int k=0; k<qsize_; k++)
458  local_matrix_[i*ndofs_+j] += (mm_coef_[k]+ret_coef_[sbi][k])*fe_values_.shape_value(j,k)*fe_values_.shape_value(i,k)*fe_values_.JxW(k);
459  }
460  }
461 
462  for (unsigned int i=0; i<ndofs_; i++)
463  {
464  local_mass_balance_vector_[i] = 0;
465  local_retardation_balance_vector_[i] = 0;
466  for (unsigned int k=0; k<qsize_; k++)
467  {
468  local_mass_balance_vector_[i] += mm_coef_[k]*fe_values_.shape_value(i,k)*fe_values_.JxW(k);
469  local_retardation_balance_vector_[i] -= ret_coef_[sbi][k]*fe_values_.shape_value(i,k)*fe_values_.JxW(k);
470  }
471  }
472 
473  model_->balance()->add_mass_values(model_->get_subst_idx()[sbi], cell, cell.get_loc_dof_indices(),
474  local_mass_balance_vector_, 0);
475 
476  data_->ls_dt[sbi]->mat_set_values(ndofs_, &(dof_indices_[0]), ndofs_, &(dof_indices_[0]), &(local_matrix_[0]));
477  VecSetValues(data_->ret_vec[sbi], ndofs_, &(dof_indices_[0]), &(local_retardation_balance_vector_[0]), ADD_VALUES);
478  }
479  }
480 
481  /// Implements @p AssemblyBase::begin.
482  void begin() override
483  {
484  model_->balance()->start_mass_assembly( model_->subst_idx() );
485  }
486 
487  /// Implements @p AssemblyBase::end.
488  void end() override
489  {
490  model_->balance()->finish_mass_assembly( model_->subst_idx() );
491  }
492 
493 
494  private:
495  /**
496  * @brief Calculates the velocity field on a given cell.
497  *
498  * @param cell The cell.
499  * @param velocity The computed velocity field (at quadrature points).
500  * @param point_list The quadrature points.
501  */
503  const Armor::array &point_list)
504  {
505  velocity.resize(point_list.size());
506  model_->data().flow_flux->value_list(point_list, cell, velocity);
507  }
508 
509  shared_ptr<FiniteElement<dim>> fe_; ///< Finite element for the solution of the advection-diffusion equation.
510 
511  /// Pointer to model (we must use common ancestor of concentration and heat model)
513 
514  /// Data object shared with TransportDG
515  EqDataDG *data_;
516 
517  unsigned int ndofs_; ///< Number of dofs
518  unsigned int qsize_; ///< Size of quadrature of actual dim
519  FEValues<3> fe_values_; ///< FEValues of object (of P disc finite element type)
520 
521  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
522  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
523  vector<PetscScalar> local_retardation_balance_vector_; ///< Auxiliary vector for assemble mass matrix.
525 
526  /// Mass matrix coefficients.
528  /// Retardation coefficient due to sorption.
530 
531  friend class TransportDG<Model>;
532  template < template<IntDim...> class DimAssembly>
533  friend class GenericAssembly;
534 
535 };
536 
537 
538 /**
539  * Auxiliary container class for Finite element and related objects of given dimension.
540  */
541 template <unsigned int dim, class Model>
543 {
544 public:
546 
547  /// Constructor.
548  StiffnessAssemblyDG(EqDataDG *data)
549  : AssemblyBase<dim>(data->dg_order), fe_rt_(nullptr), model_(nullptr), data_(data) {}
550 
551  /// Destructor.
553  if (fe_rt_==nullptr) return; // uninitialized object
554 
555  delete fe_rt_;
556  delete fe_rt_low_;
557  }
558 
559  /// Initialize auxiliary vectors and other data members
561  this->model_ = &model;
562 
563  fe_ = std::make_shared< FE_P_disc<dim> >(data_->dg_order);
564  fe_low_ = std::make_shared< FE_P_disc<dim-1> >(data_->dg_order);
565  fe_rt_ = new FE_RT0<dim>();
566  fe_rt_low_ = new FE_RT0<dim-1>();
567  fv_rt_.initialize(*this->quad_, *fe_rt_, update_values | update_gradients | update_quadrature_points);
568  fe_values_.initialize(*this->quad_, *fe_, update_values | update_gradients | update_JxW_values | update_quadrature_points);
569  if (dim>1) {
570  fv_rt_vb_.initialize(*this->quad_low_, *fe_rt_low_, update_values | update_quadrature_points);
571  fe_values_vb_.initialize(*this->quad_low_, *fe_low_, update_values | update_gradients | update_JxW_values | update_quadrature_points);
572  }
573  fe_values_side_.initialize(*this->quad_low_, *fe_, update_values | update_gradients | update_side_JxW_values | update_normal_vectors | update_quadrature_points);
574  fsv_rt_.initialize(*this->quad_low_, *fe_rt_, update_values | update_quadrature_points);
575  ndofs_ = fe_->n_dofs();
576  qsize_ = this->quad_->size();
577  qsize_lower_dim_ = this->quad_low_->size();
578  dof_indices_.resize(ndofs_);
579  side_dof_indices_vb_.resize(2*ndofs_);
580  local_matrix_.resize(4*ndofs_*ndofs_);
581  local_retardation_balance_vector_.resize(ndofs_);
582  local_mass_balance_vector_.resize(ndofs_);
583  velocity_.resize(qsize_);
584  side_velocity_vec_.resize(data_->ad_coef_edg.size());
585  sources_sigma_.resize(model_->n_substances(), std::vector<double>(qsize_));
586  sigma_.resize(qsize_lower_dim_);
587  csection_.resize(qsize_lower_dim_);
588  csection_higher_.resize(qsize_lower_dim_);
589  dg_penalty_.resize(data_->ad_coef_edg.size());
590 
591  mm_coef_.resize(qsize_);
592  ret_coef_.resize(model_->n_substances());
593  for (unsigned int sbi=0; sbi<model_->n_substances(); sbi++)
594  {
595  ret_coef_[sbi].resize(qsize_);
596  }
597 
598  fe_values_vec_.resize(data_->ad_coef_edg.size());
599  for (unsigned int sid=0; sid<data_->ad_coef_edg.size(); sid++)
600  {
601  side_dof_indices_.push_back( vector<LongIdx>(ndofs_) );
602  fe_values_vec_[sid].initialize(*this->quad_low_, *fe_,
604  }
605 
606  // index 0 = element with lower dimension,
607  // index 1 = side of element with higher dimension
608  fv_sb_.resize(2);
609  fv_sb_[0] = &fe_values_vb_;
610  fv_sb_[1] = &fe_values_side_;
611  }
612 
613 
614  /// Assembles the volume integrals into the stiffness matrix.
616  {
617  ASSERT_EQ_DBG(cell.dim(), dim).error("Dimension of element mismatch!");
618  if (!cell.is_own()) return;
619 
620  ElementAccessor<3> elm = cell.elm();
621 
622  fe_values_.reinit(elm);
623  fv_rt_.reinit(elm);
624  cell.get_dof_indices(dof_indices_);
625 
626  calculate_velocity(elm, velocity_, fv_rt_.point_list());
627  model_->compute_advection_diffusion_coefficients(fe_values_.point_list(), velocity_, elm, data_->ad_coef, data_->dif_coef);
628  model_->compute_sources_sigma(fe_values_.point_list(), elm, sources_sigma_);
629 
630  // assemble the local stiffness matrix
631  for (unsigned int sbi=0; sbi<model_->n_substances(); sbi++)
632  {
633  for (unsigned int i=0; i<ndofs_; i++)
634  for (unsigned int j=0; j<ndofs_; j++)
635  local_matrix_[i*ndofs_+j] = 0;
636 
637  for (unsigned int k=0; k<qsize_; k++)
638  {
639  for (unsigned int i=0; i<ndofs_; i++)
640  {
641  arma::vec3 Kt_grad_i = data_->dif_coef[sbi][k].t()*fe_values_.shape_grad(i,k);
642  double ad_dot_grad_i = arma::dot(data_->ad_coef[sbi][k], fe_values_.shape_grad(i,k));
643 
644  for (unsigned int j=0; j<ndofs_; j++)
645  local_matrix_[i*ndofs_+j] += (arma::dot(Kt_grad_i, fe_values_.shape_grad(j,k))
646  -fe_values_.shape_value(j,k)*ad_dot_grad_i
647  +sources_sigma_[sbi][k]*fe_values_.shape_value(j,k)*fe_values_.shape_value(i,k))*fe_values_.JxW(k);
648  }
649  }
650  data_->ls[sbi]->mat_set_values(ndofs_, &(dof_indices_[0]), ndofs_, &(dof_indices_[0]), &(local_matrix_[0]));
651  }
652  }
653 
654 
655  /// Assembles the fluxes on the boundary.
656  void assemble_fluxes_boundary(DHCellSide cell_side) override
657  {
658  ASSERT_EQ_DBG(cell_side.dim(), dim).error("Dimension of element mismatch!");
659  if (!cell_side.cell().is_own()) return;
660 
661  Side side = cell_side.side();
662  const DHCellAccessor &cell = cell_side.cell();
663 
664  ElementAccessor<3> elm_acc = cell.elm();
665  cell.get_dof_indices(dof_indices_);
666  fe_values_side_.reinit(side);
667  fsv_rt_.reinit(side);
668 
669  calculate_velocity(elm_acc, velocity_, fsv_rt_.point_list());
670  model_->compute_advection_diffusion_coefficients(fe_values_side_.point_list(), velocity_, elm_acc, data_->ad_coef, data_->dif_coef);
671  arma::uvec bc_type;
672  model_->get_bc_type(side.cond().element_accessor(), bc_type);
673  data_->cross_section.value_list(fe_values_side_.point_list(), elm_acc, csection_);
674 
675  for (unsigned int sbi=0; sbi<model_->n_substances(); sbi++)
676  {
677  std::fill(local_matrix_.begin(), local_matrix_.end(), 0);
678 
679  // On Neumann boundaries we have only term from integrating by parts the advective term,
680  // on Dirichlet boundaries we additionally apply the penalty which enforces the prescribed value.
681  double side_flux = 0;
682  for (unsigned int k=0; k<qsize_lower_dim_; k++)
683  side_flux += arma::dot(data_->ad_coef[sbi][k], fe_values_side_.normal_vector(k))*fe_values_side_.JxW(k);
684  double transport_flux = side_flux/side.measure();
685 
686  if (bc_type[sbi] == AdvectionDiffusionModel::abc_dirichlet)
687  {
688  // set up the parameters for DG method
689  double gamma_l;
690  data_->set_DG_parameters_boundary(side, qsize_lower_dim_, data_->dif_coef[sbi], transport_flux, fe_values_side_.normal_vector(0), data_->dg_penalty[sbi].value(elm_acc.centre(), elm_acc), gamma_l);
691  data_->gamma[sbi][side.cond_idx()] = gamma_l;
692  transport_flux += gamma_l;
693  }
694 
695  // fluxes and penalty
696  for (unsigned int k=0; k<qsize_lower_dim_; k++)
697  {
698  double flux_times_JxW;
699  if (bc_type[sbi] == AdvectionDiffusionModel::abc_total_flux)
700  {
701  //sigma_ corresponds to robin_sigma
702  model_->get_flux_bc_sigma(sbi, fe_values_side_.point_list(), side.cond().element_accessor(), sigma_);
703  flux_times_JxW = csection_[k]*sigma_[k]*fe_values_side_.JxW(k);
704  }
705  else if (bc_type[sbi] == AdvectionDiffusionModel::abc_diffusive_flux)
706  {
707  model_->get_flux_bc_sigma(sbi, fe_values_side_.point_list(), side.cond().element_accessor(), sigma_);
708  flux_times_JxW = (transport_flux + csection_[k]*sigma_[k])*fe_values_side_.JxW(k);
709  }
710  else if (bc_type[sbi] == AdvectionDiffusionModel::abc_inflow && side_flux < 0)
711  flux_times_JxW = 0;
712  else
713  flux_times_JxW = transport_flux*fe_values_side_.JxW(k);
714 
715  for (unsigned int i=0; i<ndofs_; i++)
716  {
717  for (unsigned int j=0; j<ndofs_; j++)
718  {
719  // flux due to advection and penalty
720  local_matrix_[i*ndofs_+j] += flux_times_JxW*fe_values_side_.shape_value(i,k)*fe_values_side_.shape_value(j,k);
721 
722  // flux due to diffusion (only on dirichlet and inflow boundary)
723  if (bc_type[sbi] == AdvectionDiffusionModel::abc_dirichlet)
724  local_matrix_[i*ndofs_+j] -= (arma::dot(data_->dif_coef[sbi][k]*fe_values_side_.shape_grad(j,k),fe_values_side_.normal_vector(k))*fe_values_side_.shape_value(i,k)
725  + arma::dot(data_->dif_coef[sbi][k]*fe_values_side_.shape_grad(i,k),fe_values_side_.normal_vector(k))*fe_values_side_.shape_value(j,k)*data_->dg_variant
726  )*fe_values_side_.JxW(k);
727  }
728  }
729  }
730 
731  data_->ls[sbi]->mat_set_values(ndofs_, &(dof_indices_[0]), ndofs_, &(dof_indices_[0]), &(local_matrix_[0]));
732  }
733  }
734 
735 
736  /// Assembles the fluxes between elements of the same dimension.
738  ASSERT_EQ_DBG(edge_side_range.begin()->element().dim(), dim).error("Dimension of element mismatch!");
739 
740  sid=0;
741  for( DHCellSide edge_side : edge_side_range )
742  {
743  auto dh_edge_cell = data_->dh_->cell_accessor_from_element( edge_side.elem_idx() );
744  ElementAccessor<3> edg_elm = dh_edge_cell.elm();
745  dh_edge_cell.get_dof_indices(side_dof_indices_[sid]);
746  fe_values_vec_[sid].reinit(edge_side.side());
747  fsv_rt_.reinit(edge_side.side());
748  calculate_velocity(edg_elm, side_velocity_vec_[sid], fsv_rt_.point_list());
749  model_->compute_advection_diffusion_coefficients(fe_values_vec_[sid].point_list(), side_velocity_vec_[sid], edg_elm, data_->ad_coef_edg[sid], data_->dif_coef_edg[sid]);
750  dg_penalty_[sid].resize(model_->n_substances());
751  for (unsigned int sbi=0; sbi<model_->n_substances(); sbi++)
752  dg_penalty_[sid][sbi] = data_->dg_penalty[sbi].value(edg_elm.centre(), edg_elm);
753  ++sid;
754  }
755  arma::vec3 normal_vector = fe_values_vec_[0].normal_vector(0);
756 
757  // fluxes and penalty
758  for (unsigned int sbi=0; sbi<model_->n_substances(); sbi++)
759  {
760  vector<double> fluxes(edge_side_range.begin()->n_edge_sides());
761  double pflux = 0, nflux = 0; // calculate the total in- and out-flux through the edge
762  sid=0;
763  for( DHCellSide edge_side : edge_side_range )
764  {
765  fluxes[sid] = 0;
766  for (unsigned int k=0; k<qsize_lower_dim_; k++)
767  fluxes[sid] += arma::dot(data_->ad_coef_edg[sid][sbi][k], fe_values_vec_[sid].normal_vector(k))*fe_values_vec_[sid].JxW(k);
768  fluxes[sid] /= edge_side.measure();
769  if (fluxes[sid] > 0)
770  pflux += fluxes[sid];
771  else
772  nflux += fluxes[sid];
773  ++sid;
774  }
775 
776  s1=0;
777  for( DHCellSide edge_side1 : edge_side_range )
778  {
779  s2=-1; // need increment at begin of loop (see conditionally 'continue' directions)
780  for( DHCellSide edge_side2 : edge_side_range )
781  {
782  s2++;
783  if (s2<=s1) continue;
784  ASSERT(edge_side1.is_valid()).error("Invalid side of edge.");
785 
786  arma::vec3 nv = fe_values_vec_[s1].normal_vector(0);
787 
788  // set up the parameters for DG method
789  // calculate the flux from edge_side1 to edge_side2
790  if (fluxes[s2] > 0 && fluxes[s1] < 0)
791  transport_flux = fluxes[s1]*fabs(fluxes[s2]/pflux);
792  else if (fluxes[s2] < 0 && fluxes[s1] > 0)
793  transport_flux = fluxes[s1]*fabs(fluxes[s2]/nflux);
794  else
795  transport_flux = 0;
796 
797  gamma_l = 0.5*fabs(transport_flux);
798 
799  delta[0] = 0;
800  delta[1] = 0;
801  for (unsigned int k=0; k<qsize_lower_dim_; k++)
802  {
803  delta[0] += dot(data_->dif_coef_edg[s1][sbi][k]*normal_vector,normal_vector);
804  delta[1] += dot(data_->dif_coef_edg[s2][sbi][k]*normal_vector,normal_vector);
805  }
806  delta[0] /= qsize_lower_dim_;
807  delta[1] /= qsize_lower_dim_;
808 
809  delta_sum = delta[0] + delta[1];
810 
811 // if (delta_sum > numeric_limits<double>::epsilon())
812  if (fabs(delta_sum) > 0)
813  {
814  omega[0] = delta[1]/delta_sum;
815  omega[1] = delta[0]/delta_sum;
816  double local_alpha = max(dg_penalty_[s1][sbi], dg_penalty_[s2][sbi]);
817  double h = edge_side1.diameter();
818  aniso1 = data_->elem_anisotropy(edge_side1.element());
819  aniso2 = data_->elem_anisotropy(edge_side2.element());
820  gamma_l += local_alpha/h*aniso1*aniso2*(delta[0]*delta[1]/delta_sum);
821  }
822  else
823  for (int i=0; i<2; i++) omega[i] = 0;
824  // end of set up the parameters for DG method
825 
826  int sd[2]; bool is_side_own[2];
827  sd[0] = s1; is_side_own[0] = edge_side1.cell().is_own();
828  sd[1] = s2; is_side_own[1] = edge_side2.cell().is_own();
829 
830 #define AVERAGE(i,k,side_id) (fe_values_vec_[sd[side_id]].shape_value(i,k)*0.5)
831 #define WAVERAGE(i,k,side_id) (arma::dot(data_->dif_coef_edg[sd[side_id]][sbi][k]*fe_values_vec_[sd[side_id]].shape_grad(i,k),nv)*omega[side_id])
832 #define JUMP(i,k,side_id) ((side_id==0?1:-1)*fe_values_vec_[sd[side_id]].shape_value(i,k))
833 
834  // For selected pair of elements:
835  for (int n=0; n<2; n++)
836  {
837  if (!is_side_own[n]) continue;
838 
839  for (int m=0; m<2; m++)
840  {
841  for (unsigned int i=0; i<fe_values_vec_[sd[n]].n_dofs(); i++)
842  for (unsigned int j=0; j<fe_values_vec_[sd[m]].n_dofs(); j++)
843  local_matrix_[i*fe_values_vec_[sd[m]].n_dofs()+j] = 0;
844 
845  for (unsigned int k=0; k<qsize_lower_dim_; k++)
846  {
847  double flux_times_JxW = transport_flux*fe_values_vec_[0].JxW(k);
848  double gamma_times_JxW = gamma_l*fe_values_vec_[0].JxW(k);
849 
850  for (unsigned int i=0; i<fe_values_vec_[sd[n]].n_dofs(); i++)
851  {
852  double flux_JxW_jump_i = flux_times_JxW*JUMP(i,k,n);
853  double gamma_JxW_jump_i = gamma_times_JxW*JUMP(i,k,n);
854  double JxW_jump_i = fe_values_vec_[0].JxW(k)*JUMP(i,k,n);
855  double JxW_var_wavg_i = fe_values_vec_[0].JxW(k)*WAVERAGE(i,k,n)*data_->dg_variant;
856 
857  for (unsigned int j=0; j<fe_values_vec_[sd[m]].n_dofs(); j++)
858  {
859  int index = i*fe_values_vec_[sd[m]].n_dofs()+j;
860 
861  // flux due to transport (applied on interior edges) (average times jump)
862  local_matrix_[index] += flux_JxW_jump_i*AVERAGE(j,k,m);
863 
864  // penalty enforcing continuity across edges (applied on interior and Dirichlet edges) (jump times jump)
865  local_matrix_[index] += gamma_JxW_jump_i*JUMP(j,k,m);
866 
867  // terms due to diffusion
868  local_matrix_[index] -= WAVERAGE(j,k,m)*JxW_jump_i;
869  local_matrix_[index] -= JUMP(j,k,m)*JxW_var_wavg_i;
870  }
871  }
872  }
873  data_->ls[sbi]->mat_set_values(fe_values_vec_[sd[n]].n_dofs(), &(side_dof_indices_[sd[n]][0]), fe_values_vec_[sd[m]].n_dofs(), &(side_dof_indices_[sd[m]][0]), &(local_matrix_[0]));
874  }
875  }
876 #undef AVERAGE
877 #undef WAVERAGE
878 #undef JUMP
879  }
880  s1++;
881  }
882  }
883  }
884 
885 
886  /// Assembles the fluxes between elements of different dimensions.
887  void assemble_fluxes_element_side(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) override {
888  if (dim == 1) return;
889  ASSERT_EQ_DBG(cell_lower_dim.dim(), dim-1).error("Dimension of element mismatch!");
890 
891  // Note: use data members csection_ and velocity_ for appropriate quantities of lower dim element
892 
893  ElementAccessor<3> elm_lower_dim = cell_lower_dim.elm();
894  n_indices = cell_lower_dim.get_dof_indices(dof_indices_);
895  for(unsigned int i=0; i<n_indices; ++i) {
896  side_dof_indices_vb_[i] = dof_indices_[i];
897  }
898  fe_values_vb_.reinit(elm_lower_dim);
899  n_dofs[0] = fv_sb_[0]->n_dofs();
900 
901  DHCellAccessor cell_higher_dim = data_->dh_->cell_accessor_from_element( neighb_side.element().idx() );
902  ElementAccessor<3> elm_higher_dim = cell_higher_dim.elm();
903  n_indices = cell_higher_dim.get_dof_indices(dof_indices_);
904  for(unsigned int i=0; i<n_indices; ++i) {
905  side_dof_indices_vb_[i+n_dofs[0]] = dof_indices_[i];
906  }
907  fe_values_side_.reinit(neighb_side.side());
908  n_dofs[1] = fv_sb_[1]->n_dofs();
909 
910  // Testing element if they belong to local partition.
911  bool own_element_id[2];
912  own_element_id[0] = cell_lower_dim.is_own();
913  own_element_id[1] = cell_higher_dim.is_own();
914 
915  fsv_rt_.reinit(neighb_side.side());
916  fv_rt_vb_.reinit(elm_lower_dim);
917  calculate_velocity(elm_higher_dim, velocity_higher_, fsv_rt_.point_list());
918  calculate_velocity(elm_lower_dim, velocity_, fv_rt_vb_.point_list());
919  model_->compute_advection_diffusion_coefficients(fe_values_vb_.point_list(), velocity_, elm_lower_dim, data_->ad_coef_edg[0], data_->dif_coef_edg[0]);
920  model_->compute_advection_diffusion_coefficients(fe_values_vb_.point_list(), velocity_higher_, elm_higher_dim, data_->ad_coef_edg[1], data_->dif_coef_edg[1]);
921  data_->cross_section.value_list(fe_values_vb_.point_list(), elm_lower_dim, csection_);
922  data_->cross_section.value_list(fe_values_vb_.point_list(), elm_higher_dim, csection_higher_);
923 
924  for (unsigned int sbi=0; sbi<model_->n_substances(); sbi++) // Optimize: SWAP LOOPS
925  {
926  for (unsigned int i=0; i<n_dofs[0]+n_dofs[1]; i++)
927  for (unsigned int j=0; j<n_dofs[0]+n_dofs[1]; j++)
928  local_matrix_[i*(n_dofs[0]+n_dofs[1])+j] = 0;
929 
930  // sigma_ corresponds to frac_sigma
931  data_->fracture_sigma[sbi].value_list(fe_values_vb_.point_list(), elm_lower_dim, sigma_);
932 
933  // set transmission conditions
934  for (unsigned int k=0; k<qsize_lower_dim_; k++)
935  {
936  // The communication flux has two parts:
937  // - "diffusive" term containing sigma
938  // - "advective" term representing usual upwind
939  //
940  // The calculation differs from the reference manual, since ad_coef and dif_coef have different meaning
941  // than b and A in the manual.
942  // In calculation of sigma there appears one more csection_lower in the denominator.
943  double sigma = sigma_[k]*arma::dot(data_->dif_coef_edg[0][sbi][k]*fe_values_side_.normal_vector(k),fe_values_side_.normal_vector(k))*
944  2*csection_higher_[k]*csection_higher_[k]/(csection_[k]*csection_[k]);
945 
946  double transport_flux = arma::dot(data_->ad_coef_edg[1][sbi][k], fe_values_side_.normal_vector(k));
947 
948  comm_flux[0][0] = (sigma-min(0.,transport_flux))*fv_sb_[0]->JxW(k);
949  comm_flux[0][1] = -(sigma-min(0.,transport_flux))*fv_sb_[0]->JxW(k);
950  comm_flux[1][0] = -(sigma+max(0.,transport_flux))*fv_sb_[0]->JxW(k);
951  comm_flux[1][1] = (sigma+max(0.,transport_flux))*fv_sb_[0]->JxW(k);
952 
953  for (int n=0; n<2; n++)
954  {
955  if (!own_element_id[n]) continue;
956 
957  for (unsigned int i=0; i<n_dofs[n]; i++)
958  for (int m=0; m<2; m++)
959  for (unsigned int j=0; j<n_dofs[m]; j++)
960  local_matrix_[(i+n*n_dofs[0])*(n_dofs[0]+n_dofs[1]) + m*n_dofs[0] + j] +=
961  comm_flux[m][n]*fv_sb_[m]->shape_value(j,k)*fv_sb_[n]->shape_value(i,k);
962  }
963  }
964  data_->ls[sbi]->mat_set_values(n_dofs[0]+n_dofs[1], &(side_dof_indices_vb_[0]), n_dofs[0]+n_dofs[1], &(side_dof_indices_vb_[0]), &(local_matrix_[0]));
965  }
966  }
967 
968 
969 private:
970  /**
971  * @brief Calculates the velocity field on a given cell.
972  *
973  * @param cell The cell.
974  * @param velocity The computed velocity field (at quadrature points).
975  * @param point_list The quadrature points.
976  */
978  const Armor::array &point_list)
979  {
980  velocity.resize(point_list.size());
981  model_->data().flow_flux.value_list(point_list, cell, velocity);
982  }
983 
984  shared_ptr<FiniteElement<dim>> fe_; ///< Finite element for the solution of the advection-diffusion equation.
985  shared_ptr<FiniteElement<dim-1>> fe_low_; ///< Finite element for the solution of the advection-diffusion equation (dim-1).
986  FiniteElement<dim> *fe_rt_; ///< Finite element for the water velocity field.
987  FiniteElement<dim-1> *fe_rt_low_; ///< Finite element for the water velocity field (dim-1).
988 
989  /// Pointer to model (we must use common ancestor of concentration and heat model)
991 
992  /// Data object shared with TransportDG
993  EqDataDG *data_;
994 
995  unsigned int ndofs_; ///< Number of dofs
996  unsigned int qsize_; ///< Size of quadrature of actual dim
997  unsigned int qsize_lower_dim_; ///< Size of quadrature of dim-1
998  FEValues<3> fv_rt_; ///< FEValues of object (of RT0 finite element type)
999  FEValues<3> fe_values_; ///< FEValues of object (of P disc finite element type)
1000  FEValues<3> fv_rt_vb_; ///< FEValues of dim-1 object (of RT0 finite element type)
1001  FEValues<3> fe_values_vb_; ///< FEValues of dim-1 object (of P disc finite element type)
1002  FEValues<3> fe_values_side_; ///< FEValues of object (of P disc finite element type)
1003  FEValues<3> fsv_rt_; ///< FEValues of object (of RT0 finite element type)
1004  vector<FEValues<3>> fe_values_vec_; ///< Vector of FEValues of object (of P disc finite element types)
1005  vector<FEValues<3>*> fv_sb_; ///< Auxiliary vector, holds FEValues objects for assemble element-side
1006 
1007  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
1008  vector< vector<LongIdx> > side_dof_indices_; ///< Vector of vectors of side DOF indices
1009  vector<LongIdx> side_dof_indices_vb_; ///< Vector of side DOF indices (assemble element-side fluxex)
1010  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
1011  vector<PetscScalar> local_retardation_balance_vector_; ///< Auxiliary vector for assemble mass matrix.
1013  vector<arma::vec3> velocity_; ///< Auxiliary results.
1014  vector<arma::vec3> velocity_higher_; ///< Velocity results of higher dim element (element-side computation).
1015  vector<vector<arma::vec3> > side_velocity_vec_; ///< Vector of velocities results.
1016  vector<vector<double> > sources_sigma_; ///< Auxiliary vectors for assemble volume integrals and set_sources method.
1017  vector<double> sigma_; ///< Auxiliary vector for assemble boundary fluxes (robin sigma), element-side fluxes (frac sigma) and set boundary conditions method
1018  vector<double> csection_; ///< Auxiliary vector for assemble boundary fluxes, element-side fluxes and set boundary conditions
1019  vector<double> csection_higher_; ///< Auxiliary vector for assemble element-side fluxes
1020  vector<vector<double> > dg_penalty_; ///< Auxiliary vectors for assemble element-element fluxes
1021 
1022  /// Mass matrix coefficients.
1024  /// Retardation coefficient due to sorption.
1026 
1027  /// @name Auxiliary variables used during element-element assembly
1028  // @{
1029 
1030  double gamma_l, omega[2], transport_flux, delta[2], delta_sum;
1031  double aniso1, aniso2;
1032  int sid, s1, s2;
1033 
1034  // @}
1035 
1036  /// @name Auxiliary variables used during element-side assembly
1037  // @{
1038 
1039  unsigned int n_dofs[2], n_indices;
1040  double comm_flux[2][2];
1041 
1042  // @}
1043 
1044  friend class TransportDG<Model>;
1045  template < template<IntDim...> class DimAssembly>
1046  friend class GenericAssembly;
1047 
1048 };
1049 
1050 
1051 /**
1052  * Auxiliary container class for Finite element and related objects of given dimension.
1053  */
1054 template <unsigned int dim, class Model>
1055 class SourcesAssemblyDG : public AssemblyBase<dim>
1056 {
1057 public:
1059 
1060  /// Constructor.
1061  SourcesAssemblyDG(EqDataDG *data)
1062  : AssemblyBase<dim>(data->dg_order), model_(nullptr), data_(data) {}
1063 
1064  /// Destructor.
1066 
1067  /// Initialize auxiliary vectors and other data members
1069  this->model_ = &model;
1070 
1071  fe_ = std::make_shared< FE_P_disc<dim> >(data_->dg_order);
1073  ndofs_ = fe_->n_dofs();
1074  qsize_ = this->quad_->size();
1075  dof_indices_.resize(ndofs_);
1076  loc_dof_indices_.resize(ndofs_);
1077  local_rhs_.resize(ndofs_);
1078  local_source_balance_vector_.resize(ndofs_);
1079  local_source_balance_rhs_.resize(ndofs_);
1080  sources_conc_.resize(model_->n_substances(), std::vector<double>(qsize_));
1081  sources_density_.resize(model_->n_substances(), std::vector<double>(qsize_));
1082  sources_sigma_.resize(model_->n_substances(), std::vector<double>(qsize_));
1083  }
1084 
1085 
1086  /// Assemble integral over element
1088  {
1089  ASSERT_EQ_DBG(cell.dim(), dim).error("Dimension of element mismatch!");
1090 
1091  ElementAccessor<3> elm = cell.elm();
1092 
1093  fe_values_.reinit(elm);
1094  cell.get_dof_indices(dof_indices_);
1095 
1096  model_->compute_source_coefficients(fe_values_.point_list(), elm, sources_conc_, sources_density_, sources_sigma_);
1097 
1098  // assemble the local stiffness matrix
1099  for (unsigned int sbi=0; sbi<model_->n_substances(); sbi++)
1100  {
1101  fill_n( &(local_rhs_[0]), ndofs_, 0 );
1102  local_source_balance_vector_.assign(ndofs_, 0);
1103  local_source_balance_rhs_.assign(ndofs_, 0);
1104 
1105  // compute sources
1106  for (unsigned int k=0; k<qsize_; k++)
1107  {
1108  source = (sources_density_[sbi][k] + sources_conc_[sbi][k]*sources_sigma_[sbi][k])*fe_values_.JxW(k);
1109 
1110  for (unsigned int i=0; i<ndofs_; i++)
1111  local_rhs_[i] += source*fe_values_.shape_value(i,k);
1112  }
1113  data_->ls[sbi]->rhs_set_values(ndofs_, &(dof_indices_[0]), &(local_rhs_[0]));
1114 
1115  for (unsigned int i=0; i<ndofs_; i++)
1116  {
1117  for (unsigned int k=0; k<qsize_; k++)
1118  local_source_balance_vector_[i] -= sources_sigma_[sbi][k]*fe_values_.shape_value(i,k)*fe_values_.JxW(k);
1119 
1120  local_source_balance_rhs_[i] += local_rhs_[i];
1121  }
1122  model_->balance()->add_source_values(model_->get_subst_idx()[sbi], elm.region().bulk_idx(),
1123  cell.get_loc_dof_indices(),
1124  local_source_balance_vector_, local_source_balance_rhs_);
1125  }
1126  }
1127 
1128  /// Implements @p AssemblyBase::begin.
1129  void begin() override
1130  {
1131  model_->balance()->start_source_assembly( model_->subst_idx() );
1132  }
1133 
1134  /// Implements @p AssemblyBase::end.
1135  void end() override
1136  {
1137  model_->balance()->finish_source_assembly( model_->subst_idx() );
1138  }
1139 
1140 
1141  private:
1142  /**
1143  * @brief Calculates the velocity field on a given cell.
1144  *
1145  * @param cell The cell.
1146  * @param velocity The computed velocity field (at quadrature points).
1147  * @param point_list The quadrature points.
1148  */
1150  const Armor::array &point_list)
1151  {
1152  velocity.resize(point_list.size());
1153  model_->data().flow_flux.value_list(point_list, cell, velocity);
1154  }
1155 
1156  shared_ptr<FiniteElement<dim>> fe_; ///< Finite element for the solution of the advection-diffusion equation.
1157 
1158  /// Pointer to model (we must use common ancestor of concentration and heat model)
1160 
1161  /// Data object shared with TransportDG
1162  EqDataDG *data_;
1163 
1164  unsigned int ndofs_; ///< Number of dofs
1165  unsigned int qsize_; ///< Size of quadrature of actual dim
1166  FEValues<3> fe_values_; ///< FEValues of object (of P disc finite element type)
1167 
1168  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
1169  vector<LongIdx> loc_dof_indices_; ///< Vector of local DOF indices
1170  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for set_sources method.
1171  vector<PetscScalar> local_source_balance_vector_; ///< Auxiliary vector for set_sources method.
1172  vector<PetscScalar> local_source_balance_rhs_; ///< Auxiliary vector for set_sources method.
1173  vector<vector<double> > sources_conc_; ///< Auxiliary vectors for set_sources method.
1174  vector<vector<double> > sources_density_; ///< Auxiliary vectors for set_sources method.
1175  vector<vector<double> > sources_sigma_; ///< Auxiliary vectors for assemble volume integrals and set_sources method.
1176 
1177  /// @name Auxiliary variables used during set sources
1178  // @{
1179 
1180  double source;
1181 
1182  // @}
1183 
1184  friend class TransportDG<Model>;
1185  template < template<IntDim...> class DimAssembly>
1186  friend class GenericAssembly;
1187 
1188 };
1189 
1190 
1191 /**
1192  * Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions..
1193  */
1194 template <unsigned int dim, class Model>
1196 {
1197 public:
1199 
1200  /// Constructor.
1201  BdrConditionAssemblyDG(EqDataDG *data)
1202  : AssemblyBase<dim>(data->dg_order), fe_rt_(nullptr), model_(nullptr), data_(data) {}
1203 
1204  /// Destructor.
1206  if (fe_rt_==nullptr) return; // uninitialized object
1207 
1208  delete fe_rt_;
1209  }
1210 
1211  /// Initialize auxiliary vectors and other data members
1213  this->model_ = &model;
1214 
1215  fe_ = std::make_shared< FE_P_disc<dim> >(data_->dg_order);
1216  fe_rt_ = new FE_RT0<dim>();
1217  fe_values_side_.initialize(*this->quad_low_, *fe_, update_values | update_gradients | update_side_JxW_values | update_normal_vectors | update_quadrature_points);
1218  fsv_rt_.initialize(*this->quad_low_, *fe_rt_, update_values | update_quadrature_points);
1219  ndofs_ = fe_->n_dofs();
1220  qsize_ = this->quad_->size();
1221  qsize_lower_dim_ = this->quad_low_->size();
1222  dof_indices_.resize(ndofs_);
1223  local_rhs_.resize(ndofs_);
1224  local_flux_balance_vector_.resize(ndofs_);
1225  velocity_.resize(qsize_);
1226  sigma_.resize(qsize_lower_dim_);
1227  csection_.resize(qsize_lower_dim_);
1228  bc_values_.resize(qsize_lower_dim_);
1229  bc_fluxes_.resize(qsize_lower_dim_);
1230  bc_ref_values_.resize(qsize_lower_dim_);
1231  }
1232 
1233 
1234  /// Assemble integral over element
1236  {
1237  ElementAccessor<3> elm = cell.elm();
1238  if (elm->boundary_idx_ == nullptr) return;
1239 
1240  for (DHCellSide dh_side : cell.side_range())
1241  {
1242  if (dh_side.n_edge_sides() > 1) continue;
1243  // skip edges lying not on the boundary
1244  if (! dh_side.side().is_boundary()) continue;
1245 
1246  const unsigned int cond_idx = dh_side.side().cond_idx();
1247 
1248  ElementAccessor<3> bc_elm = dh_side.cond().element_accessor();
1249 
1250  arma::uvec bc_type;
1251  model_->get_bc_type(bc_elm, bc_type);
1252 
1253  fe_values_side_.reinit(dh_side.side());
1254  fsv_rt_.reinit(dh_side.side());
1255  calculate_velocity(elm, velocity_, fsv_rt_.point_list());
1256 
1257  cell.get_dof_indices(dof_indices_);
1258 
1259  model_->compute_advection_diffusion_coefficients(fe_values_side_.point_list(), velocity_, elm, data_->ad_coef, data_->dif_coef);
1260  data_->cross_section.value_list(fe_values_side_.point_list(), elm, csection_);
1261 
1262  for (unsigned int sbi=0; sbi<model_->n_substances(); sbi++)
1263  {
1264  fill_n(&(local_rhs_[0]), ndofs_, 0);
1265  local_flux_balance_vector_.assign(ndofs_, 0);
1266  local_flux_balance_rhs_ = 0;
1267 
1268  // The b.c. data are fetched for all possible b.c. types since we allow
1269  // different bc_type for each substance.
1270  data_->bc_dirichlet_value[sbi].value_list(fe_values_side_.point_list(), bc_elm, bc_values_);
1271 
1272  double side_flux = 0;
1273  for (unsigned int k=0; k<qsize_lower_dim_; k++)
1274  side_flux += arma::dot(data_->ad_coef[sbi][k], fe_values_side_.normal_vector(k))*fe_values_side_.JxW(k);
1275  double transport_flux = side_flux/dh_side.measure();
1276 
1277  if (bc_type[sbi] == AdvectionDiffusionModel::abc_inflow && side_flux < 0)
1278  {
1279  for (unsigned int k=0; k<qsize_lower_dim_; k++)
1280  {
1281  double bc_term = -transport_flux*bc_values_[k]*fe_values_side_.JxW(k);
1282  for (unsigned int i=0; i<ndofs_; i++)
1283  local_rhs_[i] += bc_term*fe_values_side_.shape_value(i,k);
1284  }
1285  for (unsigned int i=0; i<ndofs_; i++)
1286  local_flux_balance_rhs_ -= local_rhs_[i];
1287  }
1288  else if (bc_type[sbi] == AdvectionDiffusionModel::abc_dirichlet)
1289  {
1290  for (unsigned int k=0; k<qsize_lower_dim_; k++)
1291  {
1292  double bc_term = data_->gamma[sbi][cond_idx]*bc_values_[k]*fe_values_side_.JxW(k);
1293  arma::vec3 bc_grad = -bc_values_[k]*fe_values_side_.JxW(k)*data_->dg_variant*(arma::trans(data_->dif_coef[sbi][k])*fe_values_side_.normal_vector(k));
1294  for (unsigned int i=0; i<ndofs_; i++)
1295  local_rhs_[i] += bc_term*fe_values_side_.shape_value(i,k)
1296  + arma::dot(bc_grad,fe_values_side_.shape_grad(i,k));
1297  }
1298  for (unsigned int k=0; k<qsize_lower_dim_; k++)
1299  {
1300  for (unsigned int i=0; i<ndofs_; i++)
1301  {
1302  local_flux_balance_vector_[i] += (arma::dot(data_->ad_coef[sbi][k], fe_values_side_.normal_vector(k))*fe_values_side_.shape_value(i,k)
1303  - arma::dot(data_->dif_coef[sbi][k]*fe_values_side_.shape_grad(i,k),fe_values_side_.normal_vector(k))
1304  + data_->gamma[sbi][cond_idx]*fe_values_side_.shape_value(i,k))*fe_values_side_.JxW(k);
1305  }
1306  }
1307  if (model_->time().tlevel() > 0)
1308  for (unsigned int i=0; i<ndofs_; i++)
1309  local_flux_balance_rhs_ -= local_rhs_[i];
1310  }
1311  else if (bc_type[sbi] == AdvectionDiffusionModel::abc_total_flux)
1312  {
1313  model_->get_flux_bc_data(sbi, fe_values_side_.point_list(), bc_elm, bc_fluxes_, sigma_, bc_ref_values_);
1314  for (unsigned int k=0; k<qsize_lower_dim_; k++)
1315  {
1316  double bc_term = csection_[k]*(sigma_[k]*bc_ref_values_[k]+bc_fluxes_[k])*fe_values_side_.JxW(k);
1317  for (unsigned int i=0; i<ndofs_; i++)
1318  local_rhs_[i] += bc_term*fe_values_side_.shape_value(i,k);
1319  }
1320 
1321  for (unsigned int i=0; i<ndofs_; i++)
1322  {
1323  for (unsigned int k=0; k<qsize_lower_dim_; k++)
1324  local_flux_balance_vector_[i] += csection_[k]*sigma_[k]*fe_values_side_.JxW(k)*fe_values_side_.shape_value(i,k);
1325  local_flux_balance_rhs_ -= local_rhs_[i];
1326  }
1327  }
1328  else if (bc_type[sbi] == AdvectionDiffusionModel::abc_diffusive_flux)
1329  {
1330  model_->get_flux_bc_data(sbi, fe_values_side_.point_list(), bc_elm, bc_fluxes_, sigma_, bc_ref_values_);
1331  for (unsigned int k=0; k<qsize_lower_dim_; k++)
1332  {
1333  double bc_term = csection_[k]*(sigma_[k]*bc_ref_values_[k]+bc_fluxes_[k])*fe_values_side_.JxW(k);
1334  for (unsigned int i=0; i<ndofs_; i++)
1335  local_rhs_[i] += bc_term*fe_values_side_.shape_value(i,k);
1336  }
1337 
1338  for (unsigned int i=0; i<ndofs_; i++)
1339  {
1340  for (unsigned int k=0; k<qsize_lower_dim_; k++)
1341  local_flux_balance_vector_[i] += csection_[k]*(arma::dot(data_->ad_coef[sbi][k], fe_values_side_.normal_vector(k)) + sigma_[k])*fe_values_side_.JxW(k)*fe_values_side_.shape_value(i,k);
1342  local_flux_balance_rhs_ -= local_rhs_[i];
1343  }
1344  }
1345  else if (bc_type[sbi] == AdvectionDiffusionModel::abc_inflow && side_flux >= 0)
1346  {
1347  for (unsigned int k=0; k<qsize_lower_dim_; k++)
1348  {
1349  for (unsigned int i=0; i<ndofs_; i++)
1350  local_flux_balance_vector_[i] += arma::dot(data_->ad_coef[sbi][k], fe_values_side_.normal_vector(k))*fe_values_side_.JxW(k)*fe_values_side_.shape_value(i,k);
1351  }
1352  }
1353  data_->ls[sbi]->rhs_set_values(ndofs_, &(dof_indices_[0]), &(local_rhs_[0]));
1354 
1355  model_->balance()->add_flux_values(model_->get_subst_idx()[sbi], dh_side,
1356  cell.get_loc_dof_indices(),
1357  local_flux_balance_vector_, local_flux_balance_rhs_);
1358  }
1359  }
1360  }
1361 
1362  /// Implements @p AssemblyBase::begin.
1363  void begin() override
1364  {
1365  model_->balance()->start_flux_assembly( model_->subst_idx() );
1366  }
1367 
1368  /// Implements @p AssemblyBase::end.
1369  void end() override
1370  {
1371  model_->balance()->finish_flux_assembly( model_->subst_idx() );
1372  }
1373 
1374 
1375  private:
1376  /**
1377  * @brief Calculates the velocity field on a given cell.
1378  *
1379  * @param cell The cell.
1380  * @param velocity The computed velocity field (at quadrature points).
1381  * @param point_list The quadrature points.
1382  */
1384  const Armor::array &point_list)
1385  {
1386  velocity.resize(point_list.size());
1387  model_->data().flow_flux.value_list(point_list, cell, velocity);
1388  }
1389 
1390  shared_ptr<FiniteElement<dim>> fe_; ///< Finite element for the solution of the advection-diffusion equation.
1391  FiniteElement<dim> *fe_rt_; ///< Finite element for the water velocity field.
1392 
1393  /// Pointer to model (we must use common ancestor of concentration and heat model)
1395 
1396  /// Data object shared with TransportDG
1397  EqDataDG *data_;
1398 
1399  unsigned int ndofs_; ///< Number of dofs
1400  unsigned int qsize_; ///< Size of quadrature of actual dim
1401  unsigned int qsize_lower_dim_; ///< Size of quadrature of dim-1
1402  FEValues<3> fe_values_side_; ///< FEValues of object (of P disc finite element type)
1403  FEValues<3> fsv_rt_; ///< FEValues of object (of RT0 finite element type)
1404 
1405  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
1406  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for set_sources method.
1407  vector<PetscScalar> local_flux_balance_vector_; ///< Auxiliary vector for set_boundary_conditions method.
1408  PetscScalar local_flux_balance_rhs_; ///< Auxiliary variable for set_boundary_conditions method.
1409  vector<arma::vec3> velocity_; ///< Auxiliary results.
1410  vector<double> sigma_; ///< Auxiliary vector for assemble boundary fluxes (robin sigma), element-side fluxes (frac sigma) and set boundary conditions method
1411  vector<double> csection_; ///< Auxiliary vector for assemble boundary fluxes, element-side fluxes and set boundary conditions
1412  vector<double> bc_values_; ///< Auxiliary vector for set boundary conditions method
1413  vector<double> bc_fluxes_; ///< Same as previous
1414  vector<double> bc_ref_values_; ///< Same as previous
1415 
1416  friend class TransportDG<Model>;
1417  template < template<IntDim...> class DimAssembly>
1418  friend class GenericAssembly;
1419 
1420 };
1421 
1422 
1423 /**
1424  * Auxiliary container class sets the initial condition.
1425  */
1426 template <unsigned int dim, class Model>
1428 {
1429 public:
1431 
1432  /// Constructor.
1433  InitConditionAssemblyDG(EqDataDG *data)
1434  : AssemblyBase<dim>(data->dg_order), model_(nullptr), data_(data) {}
1435 
1436  /// Destructor.
1438 
1439  /// Initialize auxiliary vectors and other data members
1441  this->model_ = &model;
1442 
1443  fe_ = std::make_shared< FE_P_disc<dim> >(data_->dg_order);
1445  ndofs_ = fe_->n_dofs();
1446  qsize_ = this->quad_->size();
1447  dof_indices_.resize(ndofs_);
1448  local_matrix_.resize(4*ndofs_*ndofs_);
1449  local_rhs_.resize(ndofs_);
1450  init_values_.resize(model_->n_substances(), std::vector<double>(qsize_));
1451  }
1452 
1453 
1454  /// Assemble integral over element
1456  {
1457  ASSERT_EQ_DBG(cell.dim(), dim).error("Dimension of element mismatch!");
1458 
1459  ElementAccessor<3> elem = cell.elm();
1460  cell.get_dof_indices(dof_indices_);
1461  fe_values_.reinit(elem);
1462  model_->compute_init_cond(fe_values_.point_list(), elem, init_values_);
1463 
1464  for (unsigned int sbi=0; sbi<model_->n_substances(); sbi++)
1465  {
1466  for (unsigned int i=0; i<ndofs_; i++)
1467  {
1468  local_rhs_[i] = 0;
1469  for (unsigned int j=0; j<ndofs_; j++)
1470  local_matrix_[i*ndofs_+j] = 0;
1471  }
1472 
1473  for (unsigned int k=0; k<qsize_; k++)
1474  {
1475  double rhs_term = init_values_[sbi][k]*fe_values_.JxW(k);
1476 
1477  for (unsigned int i=0; i<ndofs_; i++)
1478  {
1479  for (unsigned int j=0; j<ndofs_; j++)
1480  local_matrix_[i*ndofs_+j] += fe_values_.shape_value(i,k)*fe_values_.shape_value(j,k)*fe_values_.JxW(k);
1481 
1482  local_rhs_[i] += fe_values_.shape_value(i,k)*rhs_term;
1483  }
1484  }
1485  data_->ls[sbi]->set_values(ndofs_, &(dof_indices_[0]), ndofs_, &(dof_indices_[0]), &(local_matrix_[0]), &(local_rhs_[0]));
1486  }
1487  }
1488 
1489 
1490  private:
1491  /**
1492  * @brief Calculates the velocity field on a given cell.
1493  *
1494  * @param cell The cell.
1495  * @param velocity The computed velocity field (at quadrature points).
1496  * @param point_list The quadrature points.
1497  */
1499  const Armor::array &point_list)
1500  {
1501  velocity.resize(point_list.size());
1502  model_->data().flow_flux.value_list(point_list, cell, velocity);
1503  }
1504 
1505  shared_ptr<FiniteElement<dim>> fe_; ///< Finite element for the solution of the advection-diffusion equation.
1506 
1507  /// Pointer to model (we must use common ancestor of concentration and heat model)
1509 
1510  /// Data object shared with TransportDG
1511  EqDataDG *data_;
1512 
1513  unsigned int ndofs_; ///< Number of dofs
1514  unsigned int qsize_; ///< Size of quadrature of actual dim
1515  FEValues<3> fe_values_; ///< FEValues of object (of P disc finite element type)
1516 
1517  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
1518  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
1519  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for set_sources method.
1520  std::vector<std::vector<double> > init_values_; ///< Auxiliary vectors for prepare initial condition
1521 
1522  friend class TransportDG<Model>;
1523  template < template<IntDim...> class DimAssembly>
1524  friend class GenericAssembly;
1525 
1526 };
1527 
1528 
1529 
1530 #endif /* ASSEMBLY_DG_HH_ */
1531 
TransportDG< Model >::EqData EqDataDG
unsigned int ndofs_
Number of dofs.
Definition: assembly_dg.hh:517
vector< vector< double > > ret_coef_
Retardation coefficient due to sorption.
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
virtual void assemble_fluxes_element_side(FMT_UNUSED DHCellAccessor cell_lower_dim, FMT_UNUSED DHCellSide neighb_side)
Assembles the fluxes between elements of different dimensions.
Definition: assembly_dg.hh:373
vector< PetscScalar > local_mass_balance_vector_
Same as previous.
std::array< CouplingIntegralData, 6 > coupling_integral_data_
Holds data for computing couplings integrals.
Definition: assembly_dg.hh:338
Shape function values.
Definition: update_flags.hh:87
std::array< EdgeIntegralData, 4 > edge_integral_data_
Holds data for computing edge integrals.
Definition: assembly_dg.hh:337
vector< double > csection_
Auxiliary vector for assemble boundary fluxes, element-side fluxes and set boundary conditions...
shared_ptr< FiniteElement< dim > > fe_
Finite element for the solution of the advection-diffusion equation.
Definition: assembly_dg.hh:509
Range< DHCellSide > side_range() const
Returns range of cell sides.
TransportDG< Model > * model_
Pointer to model (we must use common ancestor of concentration and heat model)
unsigned int n_dofs() const
Returns the number of degrees of freedom needed by the finite element.
unsigned int n_indices
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
vector< vector< double > > ret_coef_
Retardation coefficient due to sorption.
Definition: assembly_dg.hh:529
FEValues< 3 > fe_values_
FEValues of object (of P disc finite element type)
Definition: assembly_dg.hh:999
unsigned int qsize_
Size of quadrature of actual dim.
Definition: assembly_dg.hh:518
unsigned int size() const
Definition: armor.hh:718
Transformed quadrature weight for cell sides.
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:332
unsigned int elm_idx() const
Return serial idx to element of loc_ele_idx_.
Quadrature * quad_low_
Quadrature used in assembling methods (dim-1).
Definition: assembly_dg.hh:395
vector< arma::vec3 > velocity_
Auxiliary results.
std::shared_ptr< EvalPoints > eval_points() const
Definition: assembly_dg.hh:117
unsigned int * boundary_idx_
Definition: elements.h:80
std::array< unsigned int, 4 > integrals_size_
Holds used sizes of previous integral data types.
Definition: assembly_dg.hh:340
vector< LongIdx > dof_indices_
Vector of global DOF indices.
vector< double > sigma_
Auxiliary vector for assemble boundary fluxes (robin sigma), element-side fluxes (frac sigma) and set...
unsigned int qsize_
Size of quadrature of actual dim.
virtual void assemble_fluxes_boundary(FMT_UNUSED DHCellSide cell_side)
Assembles the fluxes on the boundary.
Definition: assembly_dg.hh:367
unsigned int dim() const
Return dimension of element appropriate to the side.
GenericAssembly(typename DimAssembly< 1 >::EqDataDG *eq_data, int active_integrals)
Constructor.
Definition: assembly_dg.hh:101
Transport with dispersion implemented using discontinuous Galerkin method.
#define AVERAGE(i, k, side_id)
unsigned int get_dof_indices(std::vector< LongIdx > &indices) const
Fill vector of the global indices of dofs associated to the cell.
void initialize(TransportDG< Model > &model)
Initialize auxiliary vectors and other data members.
unsigned int ndofs_
Number of dofs.
FEValues< 3 > fe_values_
FEValues of object (of P disc finite element type)
Definition: assembly_dg.hh:519
std::vector< std::vector< double > > init_values_
Auxiliary vectors for prepare initial condition.
~InitConditionAssemblyDG()
Destructor.
int active_integrals_
Holds mask of active integrals.
Definition: assembly_dg.hh:329
shared_ptr< FiniteElement< dim > > fe_
Finite element for the solution of the advection-diffusion equation.
FEValues< 3 > fv_rt_vb_
FEValues of dim-1 object (of RT0 finite element type)
std::array< std::shared_ptr< CouplingIntegral >, 2 > coupling_
Coupling integrals between elements of dimensions 1-2, 2-3.
Definition: assembly_dg.hh:48
void assemble(std::shared_ptr< DOFHandlerMultiDim > dh)
General assemble methods.
Definition: assembly_dg.hh:127
ActiveIntegrals
Allow set mask of active integrals.
Definition: assembly_dg.hh:35
vector< double > csection_
Auxiliary vector for assemble boundary fluxes, element-side fluxes and set boundary conditions...
vector< vector< LongIdx > > side_dof_indices_
Vector of vectors of side DOF indices.
unsigned int dim() const
Return dimension of element appropriate to cell.
vector< LongIdx > dof_indices_
Vector of global DOF indices.
FEValues< 3 > fe_values_vb_
FEValues of dim-1 object (of P disc finite element type)
void assemble_fluxes_boundary(DHCellSide cell_side) override
Assembles the fluxes on the boundary.
Definition: assembly_dg.hh:656
unsigned int qsize_lower_dim_
Size of quadrature of dim-1.
#define WAVERAGE(i, k, side_id)
MixedPtr< DimAssembly, 1 > multidim_assembly() const
Definition: assembly_dg.hh:113
StiffnessAssemblyDG(EqDataDG *data)
Constructor.
Definition: assembly_dg.hh:548
Directing class of FieldValueCache.
Iter< Object > make_iter(Object obj)
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
Cell accessor allow iterate over DOF handler cells.
vector< arma::vec3 > velocity_
Auxiliary results.
Class FEValues calculates finite element data on the actual cells such as shape function values...
void add_edge_integral(RangeConvert< DHEdgeSide, DHCellSide > edge_side_range)
Add data of edge integral to appropriate data structure.
Definition: assembly_dg.hh:302
vector< LongIdx > dof_indices_
Vector of global DOF indices.
vector< PetscScalar > local_source_balance_rhs_
Auxiliary vector for set_sources method.
virtual void begin()
Method prepares object before assemblation (e.g. balance, ...).
Definition: assembly_dg.hh:376
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
SideIter side(const unsigned int loc_index)
~SourcesAssemblyDG()
Destructor.
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:347
shared_ptr< FiniteElement< dim > > fe_
Finite element for the solution of the advection-diffusion equation.
EqDataDG * data_
Data object shared with TransportDG.
Definition: assembly_dg.hh:993
void assemble_volume_integrals(DHCellAccessor cell) override
Assemble integral over element.
virtual void end()
Method finishes object after assemblation (e.g. balance, ...).
Definition: assembly_dg.hh:379
IterConvert< ObjectIn, ObjectOut > begin()
Iterator to begin item of range.
EqDataDG * data_
Data object shared with TransportDG.
bool is_own() const
Return true if accessor represents own element (false for ghost element)
std::array< BoundaryIntegralData, 4 > boundary_integral_data_
Holds data for computing boundary integrals.
Definition: assembly_dg.hh:339
FEValues< 3 > fsv_rt_
FEValues of object (of RT0 finite element type)
FEValues< 3 > fe_values_side_
FEValues of object (of P disc finite element type)
Side side() const
Return Side of given cell and side_idx.
TransportDG< Model > * model_
Pointer to model (we must use common ancestor of concentration and heat model)
void end() override
Implements AssemblyBase::end.
Definition: assembly_dg.hh:488
void begin() override
Implements AssemblyBase::begin.
void assemble_volume_integrals(DHCellAccessor cell) override
Assemble integral over element.
vector< PetscScalar > local_mass_balance_vector_
Same as previous.
Definition: assembly_dg.hh:524
void begin() override
Implements AssemblyBase::begin.
Definition: assembly_dg.hh:482
vector< double > bc_fluxes_
Same as previous.
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
Symmetric Gauss-Legendre quadrature formulae on simplices.
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
vector< double > sigma_
Auxiliary vector for assemble boundary fluxes (robin sigma), element-side fluxes (frac sigma) and set...
void assemble_fluxes_element_element(RangeConvert< DHEdgeSide, DHCellSide > edge_side_range) override
Assembles the fluxes between elements of the same dimension.
Definition: assembly_dg.hh:737
MassAssemblyDG(EqDataDG *data)
Constructor.
Definition: assembly_dg.hh:409
void initialize() override
std::array< std::shared_ptr< EdgeIntegral >, 3 > edge_
Edge integrals between elements of dimensions 1, 2, 3.
Definition: assembly_dg.hh:47
Discontinuous Lagrangean finite element on dim dimensional simplex.
Definition: fe_p.hh:108
Transformed quadrature points.
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
FEValues< 3 > fv_rt_
FEValues of object (of RT0 finite element type)
Definition: assembly_dg.hh:998
unsigned int ndofs_
Number of dofs.
void initialize(TransportDG< Model > &model)
Initialize auxiliary vectors and other data members.
Definition: assembly_dg.hh:416
FEValues< 3 > fe_values_
FEValues of object (of P disc finite element type)
void calculate_velocity(const ElementAccessor< 3 > &cell, vector< arma::vec3 > &velocity, const Armor::array &point_list)
Calculates the velocity field on a given cell.
Definition: assembly_dg.hh:977
FEValues< 3 > fe_values_side_
FEValues of object (of P disc finite element type)
#define FMT_UNUSED
Definition: posix.h:75
InitConditionAssemblyDG(EqDataDG *data)
Constructor.
void end() override
Implements AssemblyBase::end.
void assemble_fluxes_element_side(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) override
Assembles the fluxes between elements of different dimensions.
Definition: assembly_dg.hh:887
TransportDG< Model > * model_
Pointer to model (we must use common ancestor of concentration and heat model)
TransportDG< Model >::EqData EqDataDG
vector< vector< double > > sources_conc_
Auxiliary vectors for set_sources method.
void calculate_velocity(const ElementAccessor< 3 > &cell, vector< arma::vec3 > &velocity, const Armor::array &point_list)
Calculates the velocity field on a given cell.
TransportDG< Model > * model_
Pointer to model (we must use common ancestor of concentration and heat model)
Definition: assembly_dg.hh:990
RangeConvert< DHEdgeSide, DHCellSide > edge_sides() const
Returns range of all sides looped over common Edge.
#define JUMP(i, k, side_id)
vector< FEValues< 3 > * > fv_sb_
Auxiliary vector, holds FEValues objects for assemble element-side.
Raviart-Thomas element of order 0.
Definition: fe_rt.hh:60
Definitions of basic Lagrangean finite elements with polynomial shape functions.
FEValues< 3 > fe_values_
FEValues of object (of P disc finite element type)
void add_boundary_integral(const DHCellSide &bdr_side)
Add data of boundary integral to appropriate data structure.
Definition: assembly_dg.hh:318
void calculate_velocity(const ElementAccessor< 3 > &cell, vector< arma::vec3 > &velocity, const Armor::array &point_list)
Calculates the velocity field on a given cell.
Definition: assembly_dg.hh:502
RangeConvert< DHEdgeSide, DHCellSide > edge_side_range
Definition: assembly_dg.hh:78
ElementAccessor< 3 > element() const
RangeConvert< DHNeighbSide, DHCellSide > neighb_sides() const
Returns range of neighbour cell of lower dimension corresponding to cell of higher dimension...
vector< PetscScalar > local_rhs_
Auxiliary vector for set_sources method.
#define START_TIMER(tag)
Starts a timer with specified tag.
void calculate_velocity(const ElementAccessor< 3 > &cell, vector< arma::vec3 > &velocity, const Armor::array &point_list)
Calculates the velocity field on a given cell.
void initialize(TransportDG< Model > &model)
Initialize auxiliary vectors and other data members.
BdrConditionAssemblyDG(EqDataDG *data)
Constructor.
unsigned int qsize_
Size of quadrature of actual dim.
Definition: assembly_dg.hh:996
vector< double > bc_values_
Auxiliary vector for set boundary conditions method.
vector< LongIdx > loc_dof_indices_
Vector of local DOF indices.
vector< PetscScalar > local_rhs_
Auxiliary vector for set_sources method.
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
Definition: assembly_dg.hh:522
vector< double > bc_ref_values_
Same as previous.
TransportDG< Model >::EqData EqDataDG
Definition: assembly_dg.hh:406
vector< double > mm_coef_
Mass matrix coefficients.
FiniteElement< dim > * fe_rt_
Finite element for the water velocity field.
Definition: assembly_dg.hh:986
unsigned int ndofs_
Number of dofs.
Definition: assembly_dg.hh:995
ElementCacheMap element_cache_map_
ElementCacheMap according to EvalPoints.
Definition: assembly_dg.hh:333
Shape function gradients.
Definition: update_flags.hh:95
vector< LongIdx > dof_indices_
Vector of global DOF indices.
unsigned int IntDim
A dimension index type.
Definition: mixed.hh:19
vector< double > mm_coef_
Mass matrix coefficients.
Definition: assembly_dg.hh:527
double measure() const
Calculate metrics of the side.
Definition: accessors.cc:27
void assemble_volume_integrals(DHCellAccessor cell) override
Assemble integral over element.
Definition: assembly_dg.hh:438
unsigned int qsize_lower_dim_
Size of quadrature of dim-1.
Definition: assembly_dg.hh:997
Region region() const
Definition: accessors.hh:165
Normal vectors.
std::array< std::shared_ptr< BulkIntegral >, 3 > bulk_
Bulk integrals of elements of dimensions 1, 2, 3.
Definition: assembly_dg.hh:46
vector< LongIdx > dof_indices_
Vector of global DOF indices.
Definition: assembly_dg.hh:521
vector< arma::vec3 > velocity_higher_
Velocity results of higher dim element (element-side computation).
Discontinuous Galerkin method for equation of transport with dispersion.
void assemble_volume_integrals(DHCellAccessor cell) override
Assemble integral over element.
vector< vector< arma::vec3 > > side_velocity_vec_
Vector of velocities results.
PetscScalar local_flux_balance_rhs_
Auxiliary variable for set_boundary_conditions method.
virtual ~AssemblyBase()
Definition: assembly_dg.hh:358
vector< double > csection_higher_
Auxiliary vector for assemble element-side fluxes.
const DHCellAccessor & cell() const
Return DHCellAccessor appropriate to the side.
vector< vector< double > > sources_sigma_
Auxiliary vectors for assemble volume integrals and set_sources method.
TransportDG< Model >::EqData EqDataDG
void assemble_volume_integrals(DHCellAccessor cell) override
Assembles the volume integrals into the stiffness matrix.
Definition: assembly_dg.hh:615
MixedPtr< DimAssembly, 1 > multidim_assembly_
Assembly object.
Definition: assembly_dg.hh:326
vector< vector< double > > dg_penalty_
Auxiliary vectors for assemble element-element fluxes.
unsigned int cond_idx() const
Returns global index of the prescribed boundary condition.
void create_integrals(std::shared_ptr< EvalPoints > eval_points, AssemblyIntegrals &integrals, int active_integrals)
Create integrals according to dim of assembly object.
Definition: assembly_dg.hh:382
unsigned int qsize_
Size of quadrature of actual dim.
vector< LongIdx > side_dof_indices_vb_
Vector of side DOF indices (assemble element-side fluxex)
shared_ptr< FiniteElement< dim > > fe_
Finite element for the solution of the advection-diffusion equation.
Definition: assembly_dg.hh:984
~BdrConditionAssemblyDG()
Destructor.
EqDataDG * data_
Data object shared with TransportDG.
Definition: assembly_dg.hh:515
void add_coupling_integral(const DHCellAccessor &cell, const DHCellSide &ngh_side)
Add data of coupling integral to appropriate data structure.
Definition: assembly_dg.hh:309
void end() override
Implements AssemblyBase::end.
vector< FEValues< 3 > > fe_values_vec_
Vector of FEValues of object (of P disc finite element types)
unsigned int ndofs_
Number of dofs.
Generic class of assemblation.
Definition: assembly_dg.hh:64
std::array< std::shared_ptr< BoundaryIntegral >, 3 > boundary_
Boundary integrals betwwen elements of dimensions 1, 2, 3 and boundaries.
Definition: assembly_dg.hh:49
TransportDG< Model >::EqData EqDataDG
Definition: assembly_dg.hh:545
FiniteElement< dim > * fe_rt_
Finite element for the water velocity field.
AssemblyBase(unsigned int quad_order)
Constructor.
Definition: assembly_dg.hh:352
vector< vector< double > > sources_density_
Auxiliary vectors for set_sources method.
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
Definition: region.hh:91
Definitions of particular quadrature rules on simplices.
void add_integrals_of_computing_step(DHCellAccessor cell)
Definition: assembly_dg.hh:250
vector< PetscScalar > local_rhs_
Auxiliary vector for set_sources method.
void initialize(TransportDG< Model > &model)
Initialize auxiliary vectors and other data members.
shared_ptr< FiniteElement< dim > > fe_
Finite element for the solution of the advection-diffusion equation.
#define END_TIMER(tag)
Ends a timer with specified tag.
shared_ptr< FiniteElement< dim-1 > > fe_low_
Finite element for the solution of the advection-diffusion equation (dim-1).
Definition: assembly_dg.hh:985
std::shared_ptr< EvalPoints > eval_points_
EvalPoints object shared by all integrals.
Definition: assembly_dg.hh:332
unsigned int qsize_
Size of quadrature of actual dim.
std::array< BulkIntegralData, 1 > bulk_integral_data_
Holds data for computing bulk integrals.
Definition: assembly_dg.hh:336
EqDataDG * data_
Data object shared with TransportDG.
void add_volume_integral(const DHCellAccessor &cell)
Add data of volume integral to appropriate data structure.
Definition: assembly_dg.hh:295
virtual void assemble_volume_integrals(FMT_UNUSED DHCellAccessor cell)
Assembles the volume integrals on cell.
Definition: assembly_dg.hh:364
EqDataDG * data_
Data object shared with TransportDG.
virtual void assemble_fluxes_element_element(FMT_UNUSED RangeConvert< DHEdgeSide, DHCellSide > edge_side_range)
Assembles the fluxes between sides on the edge.
Definition: assembly_dg.hh:370
SourcesAssemblyDG(EqDataDG *data)
Constructor.
void insert_eval_points_from_integral_data()
Mark eval points in table of Element cache map.
Definition: assembly_dg.hh:208
vector< PetscScalar > local_flux_balance_vector_
Auxiliary vector for set_boundary_conditions method.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
unsigned int dim() const
Definition: accessors.hh:157
Boundary cond() const
TransportDG< Model > * model_
Pointer to model (we must use common ancestor of concentration and heat model)
Definition: assembly_dg.hh:512
vector< vector< double > > sources_sigma_
Auxiliary vectors for assemble volume integrals and set_sources method.
~StiffnessAssemblyDG()
Destructor.
Definition: assembly_dg.hh:552
vector< PetscScalar > local_source_balance_vector_
Auxiliary vector for set_sources method.
vector< PetscScalar > local_retardation_balance_vector_
Auxiliary vector for assemble mass matrix.
Definition: assembly_dg.hh:523
FiniteElement< dim-1 > * fe_rt_low_
Finite element for the water velocity field (dim-1).
Definition: assembly_dg.hh:987
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
Definition: accessors.hh:181
~MassAssemblyDG()
Destructor.
Definition: assembly_dg.hh:413
void initialize(TransportDG< Model > &model)
Initialize auxiliary vectors and other data members.
Definition: assembly_dg.hh:560
FEValues< 3 > fsv_rt_
FEValues of object (of RT0 finite element type)
Quadrature * quad_
Quadrature used in assembling methods.
Definition: assembly_dg.hh:394
vector< PetscScalar > local_retardation_balance_vector_
Auxiliary vector for assemble mass matrix.
Definitions of Raviart-Thomas finite elements.
Side accessor allows to iterate over sides of DOF handler cell.
void begin() override
Implements AssemblyBase::begin.
void calculate_velocity(const ElementAccessor< 3 > &cell, vector< arma::vec3 > &velocity, const Armor::array &point_list)
Calculates the velocity field on a given cell.
ElementAccessor< 3 > element_accessor()
Set of all used integral necessary in assemblation.
Definition: assembly_dg.hh:45
Transformed quadrature weights.
Class allows to iterate over sides of edge.
AssemblyIntegrals integrals_
Holds integral objects.
Definition: assembly_dg.hh:331