Flow123d  DF_patch_fe_data_tables-9d4017b
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 
24 #include "fem/fe_p.hh"
25 #include "fem/patch_fe_values.hh"
27 #include "coupling/balance.hh"
29 
30 
31 
32 /**
33  * Auxiliary container class for Finite element and related objects of given dimension.
34  */
35 template <unsigned int dim, class Model>
36 class MassAssemblyDG : public AssemblyBasePatch<dim>
37 {
38 public:
41 
42  static constexpr const char * name() { return "MassAssemblyDG"; }
43 
44  /// Constructor.
45  MassAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues<3> *fe_values)
46  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
47  JxW_( this->bulk_values().JxW() ),
48  conc_shape_( this->bulk_values().scalar_shape() ) {
50  this->used_fields_ += eq_fields_->mass_matrix_coef;
51  this->used_fields_ += eq_fields_->retardation_coef;
52  }
53 
54  /// Destructor.
56 
57  /// Initialize auxiliary vectors and other data members
58  void initialize(ElementCacheMap *element_cache_map) {
59  this->element_cache_map_ = element_cache_map;
60 
61  this->fe_values_->template initialize<dim>(*this->quad_);
62  ndofs_ = this->n_dofs();
63  dof_indices_.resize(ndofs_);
64  local_matrix_.resize(4*ndofs_*ndofs_);
65  local_retardation_balance_vector_.resize(ndofs_);
67 
68  }
69 
70 
71  /// Assemble integral over element
72  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
73  {
74  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
75 
77 
78  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); ++sbi)
79  {
80  // assemble the local mass matrix
81  for (unsigned int i=0; i<ndofs_; i++)
82  {
83  for (unsigned int j=0; j<ndofs_; j++)
84  {
85  local_matrix_[i*ndofs_+j] = 0;
86  for (auto p : this->bulk_points(element_patch_idx) )
87  {
88  local_matrix_[i*ndofs_+j] += (eq_fields_->mass_matrix_coef(p)+eq_fields_->retardation_coef[sbi](p)) *
89  conc_shape_(j,p)*conc_shape_(i,p)*JxW_(p);
90  }
91  }
92  }
93 
94  for (unsigned int i=0; i<ndofs_; i++)
95  {
98  for (auto p : this->bulk_points(element_patch_idx) )
99  {
100  local_mass_balance_vector_[i] += eq_fields_->mass_matrix_coef(p)*conc_shape_(i,p)*JxW_(p);
101  local_retardation_balance_vector_[i] -= eq_fields_->retardation_coef[sbi](p)*conc_shape_(i,p)*JxW_(p);
102  }
103  }
104 
105  eq_data_->balance_->add_mass_values(eq_data_->subst_idx()[sbi], cell, cell.get_loc_dof_indices(),
107 
109  VecSetValues(eq_data_->ret_vec[sbi], ndofs_, &(dof_indices_[0]), &(local_retardation_balance_vector_[0]), ADD_VALUES);
110  }
111  }
112 
113  /// Implements @p AssemblyBase::begin.
114  void begin() override
115  {
116  eq_data_->balance_->start_mass_assembly( eq_data_->subst_idx() );
117  }
118 
119  /// Implements @p AssemblyBase::end.
120  void end() override
121  {
122  eq_data_->balance_->finish_mass_assembly( eq_data_->subst_idx() );
123  }
124 
125  private:
126  /// Data objects shared with TransportDG
129 
130  /// Sub field set contains fields used in calculation.
132 
133  unsigned int ndofs_; ///< Number of dofs
134  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
135  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
136  vector<PetscScalar> local_retardation_balance_vector_; ///< Auxiliary vector for assemble mass matrix.
138 
141 
142  template < template<IntDim...> class DimAssembly>
143  friend class GenericAssembly;
144 
145 };
146 
147 
148 // return the ratio of longest and shortest edge
150 {
151  double h_max = 0, h_min = numeric_limits<double>::infinity();
152  for (unsigned int i=0; i<e->n_nodes(); i++)
153  for (unsigned int j=i+1; j<e->n_nodes(); j++)
154  {
155  double dist = arma::norm(*e.node(i) - *e.node(j));
156  h_max = max(h_max, dist);
157  h_min = min(h_min, dist);
158  }
159  return h_max/h_min;
160 }
161 
162 /**
163  * @brief Computes average normal diffusivity over a set of points
164  *
165  * @param diff_coef Diffusion tensor.
166  * @param pts Points.
167  * @param nv Normal vector.
168  * @return double
169  */
171 {
172  double delta = 0;
173  unsigned int n = 0;
174  for (auto p : pts )
175  {
176  delta += dot(diff_coef(p)*nv, nv);
177  n++;
178  }
179  return n == 0 ? 0 : (delta/n);
180 }
181 
182 
183 /**
184  * @brief Computes the penalty parameter of the DG method on a given boundary edge.
185  *
186  * Assumption is that the edge consists of only 1 side.
187  * @param side The boundary side.
188  * @param diff_delta Average normal dispersivity K*n.n computed by diffusion_delta()
189  * @param ad_vector Advection vector.
190  * @param alpha Penalty parameter that influences the continuity
191  * of the solution (large value=more continuity).
192  */
194  const double &diff_delta,
195  const double flux,
196  const double alpha)
197 {
198  return 0.5*fabs(flux) + alpha/side.diameter()*diff_delta*elem_anisotropy(side.element());
199 }
200 
201 
202 /**
203  * @brief Computes advective flux.
204  *
205  * @param advection_coef Advection vector.
206  * @param pts Quadrature points.
207  * @param JxW JxW accessor.
208  * @param normal normalW accessor.
209  * @return double
210  */
211 template <class PointType>
213 {
214  double side_flux = 0;
215  for (auto p : pts) {
216  side_flux += arma::dot(advection_coef(p), normal(p))*JxW(p);
217  }
218  return side_flux;
219 }
220 
221 
222 /**
223  * Auxiliary container class for Finite element and related objects of given dimension.
224  */
225 template <unsigned int dim, class Model>
227 {
228 public:
231 
232  static constexpr const char * name() { return "StiffnessAssemblyDG"; }
233 
234  /// Constructor.
235  StiffnessAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues<3> *fe_values)
236  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
237  JxW_( this->bulk_values().JxW() ),
238  JxW_side_( this->side_values().JxW() ),
239  normal_( this->side_values().normal_vector() ),
240  conc_shape_( this->bulk_values().scalar_shape() ),
241  conc_shape_side_( this->side_values().scalar_shape() ),
242  conc_grad_( this->bulk_values().grad_scalar_shape() ),
243  conc_grad_sidw_( this->side_values().grad_scalar_shape() ),
244  conc_join_shape_( Range< JoinShapeAccessor<Scalar> >( this->join_values().scalar_join_shape() ) ) {
246  this->used_fields_ += eq_fields_->advection_coef;
247  this->used_fields_ += eq_fields_->diffusion_coef;
248  this->used_fields_ += eq_fields_->cross_section;
250  this->used_fields_ += eq_fields_->sources_sigma_out;
252  this->used_fields_ += eq_fields_->bc_type;
253  this->used_fields_ += eq_fields_->bc_robin_sigma;
254  }
255 
256  /// Destructor.
258  for (auto a : averages) if (a != nullptr) delete[] a;
259  for (auto a : waverages) if (a != nullptr) delete[] a;
260  for (auto a : jumps) if (a != nullptr) delete[] a;
261  }
262 
263  /// Initialize auxiliary vectors and other data members
264  void initialize(ElementCacheMap *element_cache_map) {
265  this->element_cache_map_ = element_cache_map;
266 
267  this->fe_values_->template initialize<dim>(*this->quad_);
268  this->fe_values_->template initialize<dim>(*this->quad_low_);
269  ndofs_ = this->n_dofs();
270  qsize_lower_dim_ = this->quad_low_->size();
271  dof_indices_.resize(ndofs_);
272  side_dof_indices_vb_.resize(2*ndofs_);
273  local_matrix_.resize(4*ndofs_*ndofs_);
274 
275  for (unsigned int sid=0; sid<eq_data_->max_edg_sides; sid++)
276  {
278  }
279 
281  for (unsigned int s=0; s<eq_data_->max_edg_sides; s++)
282  averages[s] = new double[qsize_lower_dim_*ndofs_];
283  waverages.resize(2);
284  jumps.resize(2);
285  for (unsigned int s=0; s<2; s++)
286  {
287  waverages[s] = new double[qsize_lower_dim_*ndofs_];
288  jumps[s] = new double[qsize_lower_dim_*ndofs_];
289  }
290  }
291 
292 
293  /// Assembles the cell (volume) integral into the stiffness matrix.
294  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
295  {
296  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
297  if (!cell.is_own()) return;
298 
300 
301  // assemble the local stiffness matrix
302  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
303  {
304  for (unsigned int i=0; i<ndofs_; i++)
305  for (unsigned int j=0; j<ndofs_; j++)
306  local_matrix_[i*ndofs_+j] = 0;
307 
308  for (auto p : this->bulk_points(element_patch_idx) )
309  {
310  for (unsigned int i=0; i<ndofs_; i++)
311  {
312  arma::vec3 Kt_grad_i = eq_fields_->diffusion_coef[sbi](p).t()*conc_grad_(i,p);
313  double ad_dot_grad_i = arma::dot(eq_fields_->advection_coef[sbi](p), conc_grad_(i,p));
314 
315  for (unsigned int j=0; j<ndofs_; j++)
316  local_matrix_[i*ndofs_+j] += (arma::dot(Kt_grad_i, conc_grad_(j,p))
317  -conc_shape_(j,p)*ad_dot_grad_i
318  +eq_fields_->sources_sigma_out[sbi](p)*conc_shape_(j,p)*conc_shape_(i,p))*JxW_(p);
319  }
320  }
322  }
323  }
324 
325 
326  /// Assembles the fluxes on the boundary.
327  inline void boundary_side_integral(DHCellSide cell_side)
328  {
329  ASSERT_EQ(cell_side.dim(), dim).error("Dimension of element mismatch!");
330  if (!cell_side.cell().is_own()) return;
331 
332  Side side = cell_side.side();
333  const DHCellAccessor &cell = cell_side.cell();
334 
336  unsigned int k;
337  double gamma_l;
338 
339  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
340  {
341  std::fill(local_matrix_.begin(), local_matrix_.end(), 0);
342 
343  double side_flux = advective_flux(eq_fields_->advection_coef[sbi], this->boundary_points(cell_side), JxW_side_, normal_);
344  double transport_flux = side_flux/side.measure();
345 
346  // On Neumann boundaries we have only term from integrating by parts the advective term,
347  // on Dirichlet boundaries we additionally apply the penalty which enforces the prescribed value.
348  auto p_side = *( this->boundary_points(cell_side).begin() );
349  auto p_bdr = p_side.point_bdr( side.cond().element_accessor() );
350  unsigned int bc_type = eq_fields_->bc_type[sbi](p_bdr);
352  {
353  // set up the parameters for DG method
354  auto p = *( this->boundary_points(cell_side).begin() );
355  gamma_l = DG_penalty_boundary(side,
356  diffusion_delta(eq_fields_->diffusion_coef[sbi], this->boundary_points(cell_side), normal_(p)),
357  transport_flux,
358  eq_fields_->dg_penalty[sbi](p_side));
359  transport_flux += gamma_l;
360  }
361 
362  // fluxes and penalty
363  k=0;
364  for (auto p : this->boundary_points(cell_side) )
365  {
366  double flux_times_JxW;
368  {
369  //sigma_ corresponds to robin_sigma
370  auto p_bdr = p.point_bdr(side.cond().element_accessor());
371  flux_times_JxW = eq_fields_->cross_section(p)*eq_fields_->bc_robin_sigma[sbi](p_bdr)*JxW_side_(p);
372  }
374  {
375  auto p_bdr = p.point_bdr(side.cond().element_accessor());
376  flux_times_JxW = (transport_flux + eq_fields_->cross_section(p)*eq_fields_->bc_robin_sigma[sbi](p_bdr))*JxW_side_(p);
377  }
378  else if (bc_type == AdvectionDiffusionModel::abc_inflow && side_flux < 0)
379  flux_times_JxW = 0;
380  else
381  flux_times_JxW = transport_flux*JxW_side_(p);
382 
383  for (unsigned int i=0; i<ndofs_; i++)
384  {
385  for (unsigned int j=0; j<ndofs_; j++)
386  {
387  // flux due to advection and penalty
388  local_matrix_[i*ndofs_+j] += flux_times_JxW*conc_shape_side_(i,p)*conc_shape_side_(j,p);
389 
390  // flux due to diffusion (only on dirichlet and inflow boundary)
392  local_matrix_[i*ndofs_+j] -= (arma::dot(eq_fields_->diffusion_coef[sbi](p)*conc_grad_sidw_(j,p),normal_(p))*conc_shape_side_(i,p)
393  + arma::dot(eq_fields_->diffusion_coef[sbi](p)*conc_grad_sidw_(i,p),normal_(p))*conc_shape_side_(j,p)*eq_data_->dg_variant
394  )*JxW_side_(p);
395  }
396  }
397  k++;
398  }
399 
401  }
402  }
403 
404 
405  /// Assembles the fluxes between sides of elements of the same dimension.
407  ASSERT_EQ(edge_side_range.begin()->element().dim(), dim).error("Dimension of element mismatch!");
408 
409  unsigned int k;
410  double gamma_l, omega[2], transport_flux, delta[2], delta_sum;
411  double aniso1, aniso2;
412  double local_alpha=0.0;
413  int sid=0, s1, s2;
414  for( DHCellSide edge_side : edge_side_range )
415  {
416  auto dh_edge_cell = eq_data_->dh_->cell_accessor_from_element( edge_side.elem_idx() );
417  dh_edge_cell.get_dof_indices(side_dof_indices_[sid]);
418  ++sid;
419  }
420  auto zero_edge_side = *edge_side_range.begin();
421  auto p = *( this->edge_points(zero_edge_side).begin() );
422  arma::vec3 normal_vector = normal_(p);
423 
424  // fluxes and penalty
425  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
426  {
427  vector<double> fluxes(edge_side_range.begin()->n_edge_sides());
428  double pflux = 0, nflux = 0; // calculate the total in- and out-flux through the edge
429  sid=0;
430  for( DHCellSide edge_side : edge_side_range )
431  {
432  fluxes[sid] = advective_flux(eq_fields_->advection_coef[sbi], this->edge_points(edge_side), JxW_side_, normal_) / edge_side.measure();
433  if (fluxes[sid] > 0)
434  pflux += fluxes[sid];
435  else
436  nflux += fluxes[sid];
437  ++sid;
438  }
439 
440  // precompute averages of shape functions over pairs of sides
441  s1=0;
442  for (DHCellSide edge_side : edge_side_range)
443  {
444  k=0;
445  for (auto p : this->edge_points(edge_side) )
446  {
447  for (unsigned int i=0; i<ndofs_; i++)
448  averages[s1][k*ndofs_+i] = conc_shape_side_(i,p)*0.5;
449  k++;
450  }
451  s1++;
452  }
453 
454 
455 
456  s1=0;
457  for( DHCellSide edge_side1 : edge_side_range )
458  {
459  s2=-1; // need increment at begin of loop (see conditionally 'continue' directions)
460  for( DHCellSide edge_side2 : edge_side_range )
461  {
462  s2++;
463  if (s2<=s1) continue;
464  ASSERT(edge_side1.is_valid()).error("Invalid side of edge.");
465 
466  auto p = *( this->edge_points(edge_side1).begin() );
467  arma::vec3 nv = normal_(p);
468 
469  // set up the parameters for DG method
470  // calculate the flux from edge_side1 to edge_side2
471  if (fluxes[s2] > 0 && fluxes[s1] < 0)
472  transport_flux = fluxes[s1]*fabs(fluxes[s2]/pflux);
473  else if (fluxes[s2] < 0 && fluxes[s1] > 0)
474  transport_flux = fluxes[s1]*fabs(fluxes[s2]/nflux);
475  else
476  transport_flux = 0;
477 
478  gamma_l = 0.5*fabs(transport_flux);
479 
480  delta[0] = 0;
481  delta[1] = 0;
482  for (auto p1 : this->edge_points(edge_side1) )
483  {
484  auto p2 = p1.point_on(edge_side2);
485  delta[0] += dot(eq_fields_->diffusion_coef[sbi](p1)*normal_vector,normal_vector);
486  delta[1] += dot(eq_fields_->diffusion_coef[sbi](p2)*normal_vector,normal_vector);
487  local_alpha = max(eq_fields_->dg_penalty[sbi](p1), eq_fields_->dg_penalty[sbi](p2));
488  }
489  delta[0] /= qsize_lower_dim_;
490  delta[1] /= qsize_lower_dim_;
491 
492  delta_sum = delta[0] + delta[1];
493 
494 // if (delta_sum > numeric_limits<double>::epsilon())
495  if (fabs(delta_sum) > 0)
496  {
497  omega[0] = delta[1]/delta_sum;
498  omega[1] = delta[0]/delta_sum;
499  double h = edge_side1.diameter();
500  aniso1 = elem_anisotropy(edge_side1.element());
501  aniso2 = elem_anisotropy(edge_side2.element());
502  gamma_l += local_alpha/h*aniso1*aniso2*(delta[0]*delta[1]/delta_sum);
503  }
504  else
505  for (int i=0; i<2; i++) omega[i] = 0;
506  // end of set up the parameters for DG method
507 
508  int sd[2]; bool is_side_own[2];
509  sd[0] = s1; is_side_own[0] = edge_side1.cell().is_own();
510  sd[1] = s2; is_side_own[1] = edge_side2.cell().is_own();
511 
512  // precompute jumps and weighted averages of shape functions over the pair of sides (s1,s2)
513  k=0;
514  for (auto p1 : this->edge_points(edge_side1) )
515  {
516  auto p2 = p1.point_on(edge_side2);
517  for (unsigned int i=0; i<ndofs_; i++)
518  {
519  jumps[0][k*ndofs_+i] = conc_shape_side_(i,p1);
520  jumps[1][k*ndofs_+i] = - conc_shape_side_(i,p2);
521  waverages[0][k*ndofs_+i] = arma::dot(eq_fields_->diffusion_coef[sbi](p1)*conc_grad_sidw_(i,p1),nv)*omega[0];
522  waverages[1][k*ndofs_+i] = arma::dot(eq_fields_->diffusion_coef[sbi](p2)*conc_grad_sidw_(i,p2),nv)*omega[1];
523  }
524  k++;
525  }
526 
527  // For selected pair of elements:
528  for (int n=0; n<2; n++)
529  {
530  if (!is_side_own[n]) continue;
531 
532  for (int m=0; m<2; m++)
533  {
534  for (unsigned int i=0; i<ndofs_; i++)
535  for (unsigned int j=0; j<ndofs_; j++)
536  local_matrix_[i*ndofs_+j] = 0;
537 
538  k=0;
539  for (auto p1 : this->edge_points(zero_edge_side) )
540  //for (k=0; k<this->quad_low_->size(); ++k)
541  {
542  for (unsigned int i=0; i<ndofs_; i++)
543  {
544  for (unsigned int j=0; j<ndofs_; j++)
545  {
546  int index = i*ndofs_+j;
547 
548  local_matrix_[index] += (
549  // flux due to transport (applied on interior edges) (average times jump)
550  transport_flux*jumps[n][k*ndofs_+i]*averages[sd[m]][k*ndofs_+j]
551 
552  // penalty enforcing continuity across edges (applied on interior and Dirichlet edges) (jump times jump)
553  + gamma_l*jumps[n][k*ndofs_+i]*jumps[m][k*ndofs_+j]
554 
555  // terms due to diffusion
556  - jumps[n][k*ndofs_+i]*waverages[m][k*ndofs_+j]
557  - eq_data_->dg_variant*waverages[n][k*ndofs_+i]*jumps[m][k*ndofs_+j]
559  }
560  }
561  k++;
562  }
563  eq_data_->ls[sbi]->mat_set_values(ndofs_, &(side_dof_indices_[sd[n]][0]), ndofs_, &(side_dof_indices_[sd[m]][0]), &(local_matrix_[0]));
564  }
565  }
566  }
567  s1++;
568  }
569  }
570  }
571 
572 
573  /// Assembles the fluxes between elements of different dimensions.
574  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
575  if (dim == 1) return;
576  ASSERT_EQ(cell_lower_dim.dim(), dim-1).error("Dimension of element mismatch!");
577 
578  // Note: use data members csection_ and velocity_ for appropriate quantities of lower dim element
579 
580  unsigned int n_dofs[2];
581  unsigned int n_indices = cell_lower_dim.get_dof_indices(dof_indices_);
582  for(unsigned int i=0; i<n_indices; ++i) {
584  }
585  n_dofs[0] = conc_join_shape_.begin()->n_dofs_low();
586 
587  DHCellAccessor cell_higher_dim = eq_data_->dh_->cell_accessor_from_element( neighb_side.element().idx() );
588  n_indices = cell_higher_dim.get_dof_indices(dof_indices_);
589  for(unsigned int i=0; i<n_indices; ++i) {
591  }
592  n_dofs[1] = conc_join_shape_.begin()->n_dofs_high();
593 
594  // Testing element if they belong to local partition.
595  bool own_element_id[2];
596  own_element_id[0] = cell_lower_dim.is_own();
597  own_element_id[1] = cell_higher_dim.is_own();
598 
599  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++) // Optimize: SWAP LOOPS
600  {
601  for (unsigned int i=0; i<n_dofs[0]+n_dofs[1]; i++)
602  for (unsigned int j=0; j<n_dofs[0]+n_dofs[1]; j++)
603  local_matrix_[i*(n_dofs[0]+n_dofs[1])+j] = 0;
604 
605  // set transmission conditions
606  for (auto p_high : this->coupling_points(neighb_side) )
607  {
608  auto p_low = p_high.lower_dim(cell_lower_dim);
609  // The communication flux has two parts:
610  // - "diffusive" term containing sigma
611  // - "advective" term representing usual upwind
612  //
613  // The calculation differs from the reference manual, since ad_coef and dif_coef have different meaning
614  // than b and A in the manual.
615  // In calculation of sigma there appears one more csection_lower in the denominator.
616  double sigma = eq_fields_->fracture_sigma[sbi](p_low)*arma::dot(eq_fields_->diffusion_coef[sbi](p_low)*normal_(p_high),normal_(p_high))*
617  2*eq_fields_->cross_section(p_high)*eq_fields_->cross_section(p_high)/(eq_fields_->cross_section(p_low)*eq_fields_->cross_section(p_low));
618 
619  double transport_flux = arma::dot(eq_fields_->advection_coef[sbi](p_high), normal_(p_high));
620 
621  for( auto conc_shape_i : conc_join_shape_) {
622  uint is_high_i = conc_shape_i.is_high_dim();
623  if (!own_element_id[is_high_i]) continue;
624  uint i_mat_idx = conc_shape_i.join_idx(); // i + is_high * n_dofs_low
625  double diff_shape_i = conc_shape_i(p_high) - conc_shape_i(p_low);
626  for( auto conc_shape_j : conc_join_shape_) {
627  uint j_mat_idx = conc_shape_j.join_idx();
628  local_matrix_[i_mat_idx * (n_dofs[0]+n_dofs[1]) + j_mat_idx] += (
629  sigma * diff_shape_i * (conc_shape_j(p_high) - conc_shape_j(p_low))
630  + diff_shape_i * ( max(0.,transport_flux) * conc_shape_j(p_high) + min(0.,transport_flux) * conc_shape_j(p_low))
632  }
633  }
634  }
636  }
637  }
638 
639 
640 private:
641  /// Data objects shared with TransportDG
644 
645  /// Sub field set contains fields used in calculation.
647 
648  unsigned int ndofs_; ///< Number of dofs
649  unsigned int qsize_lower_dim_; ///< Size of quadrature of dim-1
650 
651  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
652  vector< vector<LongIdx> > side_dof_indices_; ///< Vector of vectors of side DOF indices
653  vector<LongIdx> side_dof_indices_vb_; ///< Vector of side DOF indices (assemble element-side fluxex)
654  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
655 
656  vector<double*> averages; ///< Auxiliary storage for averages of shape functions.
657  vector<double*> waverages; ///< Auxiliary storage for weighted averages of shape functions.
658  vector<double*> jumps; ///< Auxiliary storage for jumps of shape functions.
659 
668 
669  template < template<IntDim...> class DimAssembly>
670  friend class GenericAssembly;
671 
672 };
673 
674 
675 /**
676  * Auxiliary container class for Finite element and related objects of given dimension.
677  */
678 template <unsigned int dim, class Model>
680 {
681 public:
684 
685  static constexpr const char * name() { return "SourcesAssemblyDG"; }
686 
687  /// Constructor.
688  SourcesAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues<3> *fe_values)
689  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
690  JxW_( this->bulk_values().JxW() ),
691  conc_shape_( this->bulk_values().scalar_shape() ) {
693  this->used_fields_ += eq_fields_->sources_density_out;
694  this->used_fields_ += eq_fields_->sources_conc_out;
695  this->used_fields_ += eq_fields_->sources_sigma_out;
696  }
697 
698  /// Destructor.
700 
701  /// Initialize auxiliary vectors and other data members
702  void initialize(ElementCacheMap *element_cache_map) {
703  this->element_cache_map_ = element_cache_map;
704 
705  this->fe_values_->template initialize<dim>(*this->quad_);
706  ndofs_ = this->n_dofs();
707  dof_indices_.resize(ndofs_);
708  local_rhs_.resize(ndofs_);
711  }
712 
713 
714  /// Assemble integral over element
715  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
716  {
717  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
718 
719  ElementAccessor<3> elm = cell.elm();
720  double source;
721 
723 
724  // assemble the local stiffness matrix
725  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
726  {
727  fill_n( &(local_rhs_[0]), ndofs_, 0 );
730 
731  for (auto p : this->bulk_points(element_patch_idx) )
732  {
733  source = (eq_fields_->sources_density_out[sbi](p) + eq_fields_->sources_conc_out[sbi](p)*eq_fields_->sources_sigma_out[sbi](p))*JxW_(p);
734 
735  for (unsigned int i=0; i<ndofs_; i++)
736  local_rhs_[i] += source*conc_shape_(i,p);
737  }
738  eq_data_->ls[sbi]->rhs_set_values(ndofs_, &(dof_indices_[0]), &(local_rhs_[0]));
739 
740  for (unsigned int i=0; i<ndofs_; i++)
741  {
742  for (auto p : this->bulk_points(element_patch_idx) )
743  {
744  local_source_balance_vector_[i] -= eq_fields_->sources_sigma_out[sbi](p)*conc_shape_(i,p)*JxW_(p);
745  }
746 
748  }
749  eq_data_->balance_->add_source_values(eq_data_->subst_idx()[sbi], elm.region().bulk_idx(),
750  cell.get_loc_dof_indices(),
752  }
753  }
754 
755  /// Implements @p AssemblyBase::begin.
756  void begin() override
757  {
758  eq_data_->balance_->start_source_assembly( eq_data_->subst_idx() );
759  }
760 
761  /// Implements @p AssemblyBase::end.
762  void end() override
763  {
764  eq_data_->balance_->finish_source_assembly( eq_data_->subst_idx() );
765  }
766 
767 
768  private:
769  /// Data objects shared with TransportDG
772 
773  /// Sub field set contains fields used in calculation.
775 
776  unsigned int ndofs_; ///< Number of dofs
777 
778  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
779  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for set_sources method.
780  vector<PetscScalar> local_source_balance_vector_; ///< Auxiliary vector for set_sources method.
781  vector<PetscScalar> local_source_balance_rhs_; ///< Auxiliary vector for set_sources method.
782 
785 
786  template < template<IntDim...> class DimAssembly>
787  friend class GenericAssembly;
788 
789 };
790 
791 
792 /**
793  * Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions..
794  */
795 template <unsigned int dim, class Model>
797 {
798 public:
801 
802  static constexpr const char * name() { return "BdrConditionAssemblyDG"; }
803 
804  /// Constructor.
805  BdrConditionAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues<3> *fe_values)
806  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
807  JxW_( this->side_values().JxW() ),
808  normal_( this->side_values().normal_vector() ),
809  conc_shape_( this->side_values().scalar_shape() ),
810  conc_grad_( this->side_values().grad_scalar_shape() ) {
812  this->used_fields_ += eq_fields_->advection_coef;
813  this->used_fields_ += eq_fields_->diffusion_coef;
814  this->used_fields_ += eq_fields_->cross_section;
815  this->used_fields_ += eq_fields_->bc_type;
816  this->used_fields_ += eq_fields_->bc_dirichlet_value;
817  this->used_fields_ += eq_fields_->bc_robin_sigma;
818  this->used_fields_ += eq_fields_->bc_flux;
819  }
820 
821  /// Destructor.
823 
824  /// Initialize auxiliary vectors and other data members
825  void initialize(ElementCacheMap *element_cache_map) {
826  this->element_cache_map_ = element_cache_map;
827 
828  this->fe_values_->template initialize<dim>(*this->quad_low_);
829  ndofs_ = this->n_dofs();
830  dof_indices_.resize(ndofs_);
831  local_rhs_.resize(ndofs_);
833  }
834 
835 
836  /// Implements @p AssemblyBase::boundary_side_integral.
837  inline void boundary_side_integral(DHCellSide cell_side)
838  {
839 
840  ElementAccessor<3> bc_elm = cell_side.cond().element_accessor();
841 
842  const DHCellAccessor &cell = cell_side.cell();
844 
845  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
846  {
847  fill_n(&(local_rhs_[0]), ndofs_, 0);
850 
851  double side_flux = advective_flux(eq_fields_->advection_coef[sbi], this->boundary_points(cell_side), JxW_, normal_);
852  double transport_flux = side_flux/cell_side.measure();
853 
854  auto p_side = *( this->boundary_points(cell_side)).begin();
855  auto p_bdr = p_side.point_bdr(cell_side.cond().element_accessor() );
856  unsigned int bc_type = eq_fields_->bc_type[sbi](p_bdr);
857  if (bc_type == AdvectionDiffusionModel::abc_inflow && side_flux < 0)
858  {
859  for (auto p : this->boundary_points(cell_side) )
860  {
861  auto p_bdr = p.point_bdr(bc_elm);
862  double bc_term = -transport_flux*eq_fields_->bc_dirichlet_value[sbi](p_bdr)*JxW_(p);
863  for (unsigned int i=0; i<ndofs_; i++)
864  local_rhs_[i] += bc_term*conc_shape_(i,p);
865  }
866  for (unsigned int i=0; i<ndofs_; i++)
868  }
869  else if (bc_type == AdvectionDiffusionModel::abc_dirichlet)
870  {
871  double side_flux = advective_flux(eq_fields_->advection_coef[sbi], this->boundary_points(cell_side), JxW_, normal_);
872  double transport_flux = side_flux/cell_side.measure();
873 
874  auto p = *( this->boundary_points(cell_side).begin() );
875  double gamma_l = DG_penalty_boundary(cell_side.side(),
876  diffusion_delta(eq_fields_->diffusion_coef[sbi], this->boundary_points(cell_side), normal_(p)),
877  transport_flux,
878  eq_fields_->dg_penalty[sbi](p_bdr));
879 
880  for (auto p : this->boundary_points(cell_side) )
881  {
882  auto p_bdr = p.point_bdr(bc_elm);
883  double bc_term = gamma_l*eq_fields_->bc_dirichlet_value[sbi](p_bdr)*JxW_(p);
884  arma::vec3 bc_grad = -eq_fields_->bc_dirichlet_value[sbi](p_bdr)*JxW_(p)*eq_data_->dg_variant*(arma::trans(eq_fields_->diffusion_coef[sbi](p))*normal_(p));
885  for (unsigned int i=0; i<ndofs_; i++)
886  local_rhs_[i] += bc_term*conc_shape_(i,p)
887  + arma::dot(bc_grad,conc_grad_(i,p));
888  }
889  for (auto p : this->boundary_points(cell_side) )
890  {
891  for (unsigned int i=0; i<ndofs_; i++)
892  {
893  local_flux_balance_vector_[i] += (arma::dot(eq_fields_->advection_coef[sbi](p), normal_(p))*conc_shape_(i,p)
894  - arma::dot(eq_fields_->diffusion_coef[sbi](p)*conc_grad_(i,p),normal_(p))
895  + gamma_l*conc_shape_(i,p))*JxW_(p);
896  }
897  }
898  if (eq_data_->time_->step().index() > 0)
899  for (unsigned int i=0; i<ndofs_; i++)
901  }
902  else if (bc_type == AdvectionDiffusionModel::abc_total_flux)
903  {
904  for (auto p : this->boundary_points(cell_side) )
905  {
906  auto p_bdr = p.point_bdr(bc_elm);
907  double bc_term = eq_fields_->cross_section(p) * (eq_fields_->bc_robin_sigma[sbi](p_bdr)*eq_fields_->bc_dirichlet_value[sbi](p_bdr) +
908  eq_fields_->bc_flux[sbi](p_bdr)) * JxW_(p);
909  for (unsigned int i=0; i<ndofs_; i++)
910  local_rhs_[i] += bc_term*conc_shape_(i,p);
911  }
912 
913  for (unsigned int i=0; i<ndofs_; i++)
914  {
915  for (auto p : this->boundary_points(cell_side) ) {
916  auto p_bdr = p.point_bdr(bc_elm);
917  local_flux_balance_vector_[i] += eq_fields_->cross_section(p) * eq_fields_->bc_robin_sigma[sbi](p_bdr) *
918  JxW_(p) * conc_shape_(i,p);
919  }
921  }
922  }
924  {
925  for (auto p : this->boundary_points(cell_side) )
926  {
927  auto p_bdr = p.point_bdr(bc_elm);
928  double bc_term = eq_fields_->cross_section(p) * (eq_fields_->bc_robin_sigma[sbi](p_bdr)*eq_fields_->bc_dirichlet_value[sbi](p_bdr) +
929  eq_fields_->bc_flux[sbi](p_bdr)) * JxW_(p);
930  for (unsigned int i=0; i<ndofs_; i++)
931  local_rhs_[i] += bc_term*conc_shape_(i,p);
932  }
933 
934  for (unsigned int i=0; i<ndofs_; i++)
935  {
936  for (auto p : this->boundary_points(cell_side) ) {
937  auto p_bdr = p.point_bdr(bc_elm);
938  local_flux_balance_vector_[i] += eq_fields_->cross_section(p)*(arma::dot(eq_fields_->advection_coef[sbi](p), normal_(p)) +
939  eq_fields_->bc_robin_sigma[sbi](p_bdr))*JxW_(p)*conc_shape_(i,p);
940  }
942  }
943  }
944  else if (bc_type == AdvectionDiffusionModel::abc_inflow && side_flux >= 0)
945  {
946  for (auto p : this->boundary_points(cell_side) )
947  {
948  for (unsigned int i=0; i<ndofs_; i++)
949  local_flux_balance_vector_[i] += arma::dot(eq_fields_->advection_coef[sbi](p), normal_(p))*JxW_(p)*conc_shape_(i,p);
950  }
951  }
952  eq_data_->ls[sbi]->rhs_set_values(ndofs_, &(dof_indices_[0]), &(local_rhs_[0]));
953 
954  eq_data_->balance_->add_flux_values(eq_data_->subst_idx()[sbi], cell_side,
955  cell.get_loc_dof_indices(),
957  }
958  }
959 
960  /// Implements @p AssemblyBase::begin.
961  void begin() override
962  {
963  eq_data_->balance_->start_flux_assembly( eq_data_->subst_idx() );
964  }
965 
966  /// Implements @p AssemblyBase::end.
967  void end() override
968  {
969  eq_data_->balance_->finish_flux_assembly( eq_data_->subst_idx() );
970  }
971 
972 
973  private:
974  /// Data objects shared with TransportDG
977 
978  unsigned int ndofs_; ///< Number of dofs
979 
980  /// Sub field set contains fields used in calculation.
982 
983  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
984  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for set_sources method.
985  vector<PetscScalar> local_flux_balance_vector_; ///< Auxiliary vector for set_boundary_conditions method.
986  PetscScalar local_flux_balance_rhs_; ///< Auxiliary variable for set_boundary_conditions method.
987 
992 
993  template < template<IntDim...> class DimAssembly>
994  friend class GenericAssembly;
995 
996 };
997 
998 
999 /**
1000  * Auxiliary container class sets the initial condition.
1001  */
1002 template <unsigned int dim, class Model>
1004 {
1005 public:
1008 
1009  static constexpr const char * name() { return "InitProjectionAssemblyDG"; }
1010 
1011  /// Constructor.
1012  InitProjectionAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues<3> *fe_values)
1013  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
1014  JxW_( this->bulk_values().JxW() ),
1015  init_shape_( this->bulk_values().scalar_shape() ) {
1017  this->used_fields_ += eq_fields_->init_condition;
1018  }
1019 
1020  /// Destructor.
1022 
1023  /// Initialize auxiliary vectors and other data members
1024  void initialize(ElementCacheMap *element_cache_map) {
1025  this->element_cache_map_ = element_cache_map;
1026 
1027  this->fe_values_->template initialize<dim>(*this->quad_);
1028  ndofs_ = this->n_dofs();
1029  dof_indices_.resize(ndofs_);
1030  local_matrix_.resize(4*ndofs_*ndofs_);
1031  local_rhs_.resize(ndofs_);
1032  }
1033 
1034 
1035  /// Assemble integral over element
1036  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
1037  {
1038  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
1039 
1041 
1042  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
1043  {
1044  for (unsigned int i=0; i<ndofs_; i++)
1045  {
1046  local_rhs_[i] = 0;
1047  for (unsigned int j=0; j<ndofs_; j++)
1048  local_matrix_[i*ndofs_+j] = 0;
1049  }
1050 
1051  for (auto p : this->bulk_points(element_patch_idx) )
1052  {
1053  double rhs_term = eq_fields_->init_condition[sbi](p)*JxW_(p);
1054 
1055  for (unsigned int i=0; i<ndofs_; i++)
1056  {
1057  for (unsigned int j=0; j<ndofs_; j++)
1058  local_matrix_[i*ndofs_+j] += init_shape_(i,p)*init_shape_(j,p)*JxW_(p);
1059 
1060  local_rhs_[i] += init_shape_(i,p)*rhs_term;
1061  }
1062  }
1063  eq_data_->ls[sbi]->set_values(ndofs_, &(dof_indices_[0]), ndofs_, &(dof_indices_[0]), &(local_matrix_[0]), &(local_rhs_[0]));
1064  }
1065  }
1066 
1067 
1068  private:
1069  /// Data objects shared with TransportDG
1072 
1073  /// Sub field set contains fields used in calculation.
1075 
1076  unsigned int ndofs_; ///< Number of dofs
1077  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
1078  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
1079  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for set_sources method.
1080 
1083 
1084  template < template<IntDim...> class DimAssembly>
1085  friend class GenericAssembly;
1086 
1087 };
1088 
1089 
1090 /**
1091  * Auxiliary container class sets the initial condition.
1092  */
1093 template <unsigned int dim, class Model>
1095 {
1096 public:
1099 
1100  static constexpr const char * name() { return "InitConditionAssemblyDG"; }
1101 
1102  /// Constructor.
1104  : AssemblyBase<dim>(), eq_fields_(eq_fields), eq_data_(eq_data) {
1106  this->used_fields_ += eq_fields_->init_condition;
1107 
1108  this->quad_ = new Quadrature(dim, RefElement<dim>::n_nodes);
1109  for(unsigned int i = 0; i<RefElement<dim>::n_nodes; i++)
1110  {
1111  this->quad_->weight(i) = 1.0;
1112  this->quad_->set(i) = RefElement<dim>::node_coords(i);
1113  }
1114  }
1115 
1116  /// Destructor.
1118 
1119  /// Initialize auxiliary vectors and other data members
1120  void initialize(ElementCacheMap *element_cache_map) {
1121  this->element_cache_map_ = element_cache_map;
1122  }
1123 
1124 
1125  /// Assemble integral over element
1126  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
1127  {
1128  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
1129 
1130  unsigned int k;
1131 
1132  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
1133  {
1134  k=0;
1135  for (auto p : this->bulk_points(element_patch_idx) )
1136  {
1137  double val = eq_fields_->init_condition[sbi](p);
1138 
1139  eq_data_->output_vec[sbi].set(cell.get_loc_dof_indices()[k], val);
1140  k++;
1141  }
1142  }
1143  }
1144 
1145 
1146  private:
1147  /// Data objects shared with TransportDG
1150 
1151  /// Sub field set contains fields used in calculation.
1153 
1154  template < template<IntDim...> class DimAssembly>
1155  friend class GenericAssembly;
1156 
1157 };
1158 
1159 #endif /* ASSEMBLY_DG_HH_ */
1160 
double elem_anisotropy(ElementAccessor< 3 > e)
Definition: assembly_dg.hh:149
double DG_penalty_boundary(Side side, const double &diff_delta, const double flux, const double alpha)
Computes the penalty parameter of the DG method on a given boundary edge.
Definition: assembly_dg.hh:193
double diffusion_delta(Field< 3, FieldValue< 3 >::TensorFixed > &diff_coef, Range< BoundaryPoint > pts, const arma::vec3 &nv)
Computes average normal diffusivity over a set of points.
Definition: assembly_dg.hh:170
double advective_flux(Field< 3, FieldValue< 3 >::VectorFixed > &advection_coef, Range< PointType > pts, ElQ< Scalar > &JxW, ElQ< Vector > normal)
Computes advective flux.
Definition: assembly_dg.hh:212
#define ASSERT(expr)
Definition: asserts.hh:351
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:333
JoinValues< dim > join_values()
Return JoinValues object.
PatchFEValues< 3 > * fe_values_
Common FEValues object over all dimensions.
BulkValues< dim > bulk_values()
Return BulkValues object.
unsigned int n_dofs()
Return BulkValues object.
SideValues< dim > side_values()
Return SideValues object.
Range< BulkPoint > bulk_points(unsigned int element_patch_idx) const
Return BulkPoint range of appropriate dimension.
Range< EdgePoint > edge_points(const DHCellSide &cell_side) const
Return EdgePoint range of appropriate dimension.
Quadrature * quad_
Quadrature used in assembling methods.
Quadrature * quad_low_
Quadrature used in assembling methods (dim-1).
Range< CouplingPoint > coupling_points(const DHCellSide &cell_side) const
Return CouplingPoint range of appropriate dimension.
int active_integrals_
Holds mask of active integrals.
Range< BoundaryPoint > boundary_points(const DHCellSide &cell_side) const
Return BoundaryPoint range of appropriate dimension.
ElementCacheMap * element_cache_map_
ElementCacheMap shared with GenericAssembly object.
static constexpr const char * name()
Definition: assembly_dg.hh:802
ElQ< Vector > normal_
Definition: assembly_dg.hh:989
FieldSet used_fields_
Sub field set contains fields used in calculation.
Definition: assembly_dg.hh:981
vector< PetscScalar > local_flux_balance_vector_
Auxiliary vector for set_boundary_conditions method.
Definition: assembly_dg.hh:985
void boundary_side_integral(DHCellSide cell_side)
Implements AssemblyBase::boundary_side_integral.
Definition: assembly_dg.hh:837
BdrConditionAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
Definition: assembly_dg.hh:805
TransportDG< Model >::EqData EqData
Definition: assembly_dg.hh:800
TransportDG< Model >::EqFields EqFields
Definition: assembly_dg.hh:799
void begin() override
Implements AssemblyBase::begin.
Definition: assembly_dg.hh:961
~BdrConditionAssemblyDG()
Destructor.
Definition: assembly_dg.hh:822
PetscScalar local_flux_balance_rhs_
Auxiliary variable for set_boundary_conditions method.
Definition: assembly_dg.hh:986
FeQ< Scalar > conc_shape_
Definition: assembly_dg.hh:990
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
Definition: assembly_dg.hh:825
vector< LongIdx > dof_indices_
Vector of global DOF indices.
Definition: assembly_dg.hh:983
EqFields * eq_fields_
Data objects shared with TransportDG.
Definition: assembly_dg.hh:975
void end() override
Implements AssemblyBase::end.
Definition: assembly_dg.hh:967
vector< PetscScalar > local_rhs_
Auxiliary vector for set_sources method.
Definition: assembly_dg.hh:984
unsigned int ndofs_
Number of dofs.
Definition: assembly_dg.hh:978
FeQ< Vector > conc_grad_
Definition: assembly_dg.hh:991
ElementAccessor< 3 > element_accessor()
Cell accessor allow iterate over DOF handler cells.
bool is_own() const
Return true if accessor represents own element (false for ghost element)
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
unsigned int get_dof_indices(std::vector< LongIdx > &indices) const
Fill vector of the global indices of dofs associated to the cell.
unsigned int dim() const
Return dimension of element appropriate to cell.
ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
Side accessor allows to iterate over sides of DOF handler cell.
Side side() const
Return Side of given cell and side_idx.
Boundary cond() const
const DHCellAccessor & cell() const
Return DHCellAccessor appropriate to the side.
unsigned int dim() const
Return dimension of element appropriate to the side.
double measure() const
ElementAccessor< 3 > element() const
NodeAccessor< 3 > node(unsigned int ni) const
Definition: accessors.hh:230
Region region() const
Definition: accessors.hh:198
unsigned int idx() const
We need this method after replacing Region by RegionIdx, and movinf RegionDB instance into particular...
Definition: accessors.hh:215
Directing class of FieldValueCache.
unsigned int n_nodes() const
Definition: elements.h:125
Container for various descendants of FieldCommonBase.
Definition: field_set.hh:159
Class template representing a field with values dependent on: point, element, and region.
Definition: field.hh:92
Generic class of assemblation.
FieldSet used_fields_
Sub field set contains fields used in calculation.
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
static constexpr const char * name()
~InitConditionAssemblyDG()
Destructor.
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
InitConditionAssemblyDG(EqFields *eq_fields, EqData *eq_data)
Constructor.
EqFields * eq_fields_
Data objects shared with TransportDG.
TransportDG< Model >::EqData EqData
TransportDG< Model >::EqFields EqFields
EqFields * eq_fields_
Data objects shared with TransportDG.
vector< PetscScalar > local_rhs_
Auxiliary vector for set_sources method.
vector< LongIdx > dof_indices_
Vector of global DOF indices.
TransportDG< Model >::EqData EqData
FieldSet used_fields_
Sub field set contains fields used in calculation.
unsigned int ndofs_
Number of dofs.
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
InitProjectionAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
static constexpr const char * name()
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
~InitProjectionAssemblyDG()
Destructor.
TransportDG< Model >::EqFields EqFields
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
virtual void rhs_set_values(int nrow, int *rows, double *vals)=0
virtual void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)=0
void set_values(int nrow, int *rows, int ncol, int *cols, PetscScalar *mat_vals, PetscScalar *rhs_vals)
Set values in the system matrix and values in the right-hand side vector on corresponding rows.
Definition: linsys.hh:386
static constexpr double almost_zero
vector< LongIdx > dof_indices_
Vector of global DOF indices.
Definition: assembly_dg.hh:134
vector< PetscScalar > local_mass_balance_vector_
Same as previous.
Definition: assembly_dg.hh:137
TransportDG< Model >::EqFields EqFields
Definition: assembly_dg.hh:39
MassAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
Definition: assembly_dg.hh:45
FieldSet used_fields_
Sub field set contains fields used in calculation.
Definition: assembly_dg.hh:131
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
Definition: assembly_dg.hh:135
void end() override
Implements AssemblyBase::end.
Definition: assembly_dg.hh:120
~MassAssemblyDG()
Destructor.
Definition: assembly_dg.hh:55
ElQ< Scalar > JxW_
Definition: assembly_dg.hh:139
EqData * eq_data_
Definition: assembly_dg.hh:128
void begin() override
Implements AssemblyBase::begin.
Definition: assembly_dg.hh:114
static constexpr const char * name()
Definition: assembly_dg.hh:42
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
Definition: assembly_dg.hh:58
FeQ< Scalar > conc_shape_
Definition: assembly_dg.hh:140
unsigned int ndofs_
Number of dofs.
Definition: assembly_dg.hh:133
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
Definition: assembly_dg.hh:72
TransportDG< Model >::EqData EqData
Definition: assembly_dg.hh:40
EqFields * eq_fields_
Data objects shared with TransportDG.
Definition: assembly_dg.hh:127
vector< PetscScalar > local_retardation_balance_vector_
Auxiliary vector for assemble mass matrix.
Definition: assembly_dg.hh:136
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:86
Armor::Array< double >::ArrayMatSet set(uint i)
Definition: quadrature.hh:93
double weight(unsigned int i) const
Returns the ith weight.
Definition: quadrature.hh:103
IterConvert< ObjectIn, ObjectOut > begin()
Iterator to begin item of range.
Range helper class.
Iter< Object > begin()
Iterator to begin item of range.
static LocalPoint node_coords(unsigned int nid)
Definition: ref_element.hh:631
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
Definition: region.hh:90
double measure() const
Calculate metrics of the side.
Definition: accessors.cc:27
Boundary cond() const
ElementAccessor< 3 > element() const
Returns iterator to the element of the side.
double diameter() const
Calculate the side diameter.
Definition: accessors.cc:132
vector< PetscScalar > local_rhs_
Auxiliary vector for set_sources method.
Definition: assembly_dg.hh:779
TransportDG< Model >::EqData EqData
Definition: assembly_dg.hh:683
TransportDG< Model >::EqFields EqFields
Definition: assembly_dg.hh:682
ElQ< Scalar > JxW_
Definition: assembly_dg.hh:783
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
Definition: assembly_dg.hh:715
~SourcesAssemblyDG()
Destructor.
Definition: assembly_dg.hh:699
SourcesAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
Definition: assembly_dg.hh:688
vector< PetscScalar > local_source_balance_vector_
Auxiliary vector for set_sources method.
Definition: assembly_dg.hh:780
void end() override
Implements AssemblyBase::end.
Definition: assembly_dg.hh:762
static constexpr const char * name()
Definition: assembly_dg.hh:685
EqFields * eq_fields_
Data objects shared with TransportDG.
Definition: assembly_dg.hh:770
vector< PetscScalar > local_source_balance_rhs_
Auxiliary vector for set_sources method.
Definition: assembly_dg.hh:781
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
Definition: assembly_dg.hh:702
FieldSet used_fields_
Sub field set contains fields used in calculation.
Definition: assembly_dg.hh:774
unsigned int ndofs_
Number of dofs.
Definition: assembly_dg.hh:776
vector< LongIdx > dof_indices_
Vector of global DOF indices.
Definition: assembly_dg.hh:778
void begin() override
Implements AssemblyBase::begin.
Definition: assembly_dg.hh:756
FeQ< Scalar > conc_shape_
Definition: assembly_dg.hh:784
TransportDG< Model >::EqData EqData
Definition: assembly_dg.hh:230
FeQ< Scalar > conc_shape_
Definition: assembly_dg.hh:663
void edge_integral(RangeConvert< DHEdgeSide, DHCellSide > edge_side_range)
Assembles the fluxes between sides of elements of the same dimension.
Definition: assembly_dg.hh:406
FeQ< Vector > conc_grad_sidw_
Definition: assembly_dg.hh:666
FieldSet used_fields_
Sub field set contains fields used in calculation.
Definition: assembly_dg.hh:646
unsigned int ndofs_
Number of dofs.
Definition: assembly_dg.hh:648
vector< vector< LongIdx > > side_dof_indices_
Vector of vectors of side DOF indices.
Definition: assembly_dg.hh:652
Range< JoinShapeAccessor< Scalar > > conc_join_shape_
Definition: assembly_dg.hh:667
StiffnessAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
Definition: assembly_dg.hh:235
void boundary_side_integral(DHCellSide cell_side)
Assembles the fluxes on the boundary.
Definition: assembly_dg.hh:327
vector< double * > jumps
Auxiliary storage for jumps of shape functions.
Definition: assembly_dg.hh:658
vector< double * > waverages
Auxiliary storage for weighted averages of shape functions.
Definition: assembly_dg.hh:657
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
Assembles the fluxes between elements of different dimensions.
Definition: assembly_dg.hh:574
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
Definition: assembly_dg.hh:264
vector< LongIdx > side_dof_indices_vb_
Vector of side DOF indices (assemble element-side fluxex)
Definition: assembly_dg.hh:653
FeQ< Vector > conc_grad_
Definition: assembly_dg.hh:665
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assembles the cell (volume) integral into the stiffness matrix.
Definition: assembly_dg.hh:294
TransportDG< Model >::EqFields EqFields
Definition: assembly_dg.hh:229
vector< double * > averages
Auxiliary storage for averages of shape functions.
Definition: assembly_dg.hh:656
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
Definition: assembly_dg.hh:654
~StiffnessAssemblyDG()
Destructor.
Definition: assembly_dg.hh:257
ElQ< Scalar > JxW_
Definition: assembly_dg.hh:660
static constexpr const char * name()
Definition: assembly_dg.hh:232
ElQ< Vector > normal_
Definition: assembly_dg.hh:662
vector< LongIdx > dof_indices_
Vector of global DOF indices.
Definition: assembly_dg.hh:651
ElQ< Scalar > JxW_side_
Definition: assembly_dg.hh:661
FeQ< Scalar > conc_shape_side_
Definition: assembly_dg.hh:664
EqFields * eq_fields_
Data objects shared with TransportDG.
Definition: assembly_dg.hh:642
unsigned int qsize_lower_dim_
Size of quadrature of dim-1.
Definition: assembly_dg.hh:649
const TimeStep & step(int index=-1) const
unsigned int index() const
std::shared_ptr< DOFHandlerMultiDim > dh_
Object for distribution of dofs.
std::vector< VectorMPI > output_vec
Vector of solution data.
LinSys ** ls_dt
Linear algebra system for the time derivative (actually it is used only for handling the matrix struc...
int dg_variant
DG variant ((non-)symmetric/incomplete.
TimeGovernor * time_
std::shared_ptr< Balance > balance_
LinSys ** ls
Linear algebra system for the transport equation.
std::vector< Vec > ret_vec
Auxiliary vectors for calculation of sources in balance due to retardation (e.g. sorption).
unsigned int max_edg_sides
Maximal number of edge sides (evaluate from dim 1,2,3)
MultiField< 3, FieldValue< 3 >::Scalar > dg_penalty
Penalty enforcing inter-element continuity of solution (for each substance).
MultiField< 3, FieldValue< 3 >::Scalar > fracture_sigma
Transition parameter for diffusive transfer on fractures (for each substance).
Definitions of basic Lagrangean finite elements with polynomial shape functions.
@ coupling
@ boundary
@ bulk
@ edge
unsigned int uint
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,...
double Scalar
Definitions of particular quadrature rules on simplices.
Discontinuous Galerkin method for equation of transport with dispersion.