Flow123d  DF_patch_fe_data_tables-22c8b57
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 
62  this->fe_values_->template initialize<dim>(*this->quad_, u);
63  if (dim==1) // print to log only one time
64  DebugOut() << "List of MassAssembly FEValues updates flags: " << this->print_update_flags(u);
65  ndofs_ = this->fe_values_->n_dofs(dim);
66  dof_indices_.resize(ndofs_);
67  local_matrix_.resize(4*ndofs_*ndofs_);
70 
71  }
72 
73 
74  /// Assemble integral over element
75  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
76  {
77  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
78 
80 
81  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); ++sbi)
82  {
83  // assemble the local mass matrix
84  for (unsigned int i=0; i<ndofs_; i++)
85  {
86  for (unsigned int j=0; j<ndofs_; j++)
87  {
88  local_matrix_[i*ndofs_+j] = 0;
89  for (auto p : this->bulk_points(element_patch_idx) )
90  {
91  local_matrix_[i*ndofs_+j] += (eq_fields_->mass_matrix_coef(p)+eq_fields_->retardation_coef[sbi](p)) *
92  conc_shape_(j,p)*conc_shape_(i,p)*JxW_(p);
93  }
94  }
95  }
96 
97  for (unsigned int i=0; i<ndofs_; i++)
98  {
101  for (auto p : this->bulk_points(element_patch_idx) )
102  {
103  local_mass_balance_vector_[i] += eq_fields_->mass_matrix_coef(p)*conc_shape_(i,p)*JxW_(p);
104  local_retardation_balance_vector_[i] -= eq_fields_->retardation_coef[sbi](p)*conc_shape_(i,p)*JxW_(p);
105  }
106  }
107 
108  eq_data_->balance_->add_mass_values(eq_data_->subst_idx()[sbi], cell, cell.get_loc_dof_indices(),
110 
112  VecSetValues(eq_data_->ret_vec[sbi], ndofs_, &(dof_indices_[0]), &(local_retardation_balance_vector_[0]), ADD_VALUES);
113  }
114  }
115 
116  /// Implements @p AssemblyBase::begin.
117  void begin() override
118  {
119  eq_data_->balance_->start_mass_assembly( eq_data_->subst_idx() );
120  }
121 
122  /// Implements @p AssemblyBase::end.
123  void end() override
124  {
125  eq_data_->balance_->finish_mass_assembly( eq_data_->subst_idx() );
126  }
127 
128  private:
129  /// Data objects shared with TransportDG
132 
133  /// Sub field set contains fields used in calculation.
135 
136  unsigned int ndofs_; ///< Number of dofs
137  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
138  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
139  vector<PetscScalar> local_retardation_balance_vector_; ///< Auxiliary vector for assemble mass matrix.
141 
144 
145  template < template<IntDim...> class DimAssembly>
146  friend class GenericAssembly;
147 
148 };
149 
150 
151 // return the ratio of longest and shortest edge
153 {
154  double h_max = 0, h_min = numeric_limits<double>::infinity();
155  for (unsigned int i=0; i<e->n_nodes(); i++)
156  for (unsigned int j=i+1; j<e->n_nodes(); j++)
157  {
158  double dist = arma::norm(*e.node(i) - *e.node(j));
159  h_max = max(h_max, dist);
160  h_min = min(h_min, dist);
161  }
162  return h_max/h_min;
163 }
164 
165 /**
166  * @brief Computes average normal diffusivity over a set of points
167  *
168  * @param diff_coef Diffusion tensor.
169  * @param pts Points.
170  * @param nv Normal vector.
171  * @return double
172  */
174 {
175  double delta = 0;
176  unsigned int n = 0;
177  for (auto p : pts )
178  {
179  delta += dot(diff_coef(p)*nv, nv);
180  n++;
181  }
182  return n == 0 ? 0 : (delta/n);
183 }
184 
185 
186 /**
187  * @brief Computes the penalty parameter of the DG method on a given boundary edge.
188  *
189  * Assumption is that the edge consists of only 1 side.
190  * @param side The boundary side.
191  * @param diff_delta Average normal dispersivity K*n.n computed by diffusion_delta()
192  * @param ad_vector Advection vector.
193  * @param alpha Penalty parameter that influences the continuity
194  * of the solution (large value=more continuity).
195  */
197  const double &diff_delta,
198  const double flux,
199  const double alpha)
200 {
201  return 0.5*fabs(flux) + alpha/side.diameter()*diff_delta*elem_anisotropy(side.element());
202 }
203 
204 
205 /**
206  * @brief Computes advective flux.
207  *
208  * @param advection_coef Advection vector.
209  * @param pts Quadrature points.
210  * @param JxW JxW accessor.
211  * @param normal normalW accessor.
212  * @return double
213  */
214 template <class PointType>
216 {
217  double side_flux = 0;
218  for (auto p : pts) {
219  side_flux += arma::dot(advection_coef(p), normal(p))*JxW(p);
220  }
221  return side_flux;
222 }
223 
224 
225 /**
226  * Auxiliary container class for Finite element and related objects of given dimension.
227  */
228 template <unsigned int dim, class Model>
230 {
231 public:
234 
235  static constexpr const char * name() { return "StiffnessAssemblyDG"; }
236 
237  /// Constructor.
238  StiffnessAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues<3> *fe_values)
239  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
240  JxW_( this->bulk_values().JxW() ),
241  JxW_side_( this->side_values().JxW() ),
242  normal_( this->side_values().normal_vector() ),
243  conc_shape_( this->bulk_values().scalar_shape() ),
244  conc_shape_side_( this->side_values().scalar_shape() ),
245  conc_grad_( this->bulk_values().grad_scalar_shape() ),
246  conc_grad_sidw_( this->side_values().grad_scalar_shape() ),
247  conc_join_shape_( Range< JoinShapeAccessor<Scalar> >( this->join_values().scalar_join_shape() ) ) {
249  this->used_fields_ += eq_fields_->advection_coef;
250  this->used_fields_ += eq_fields_->diffusion_coef;
251  this->used_fields_ += eq_fields_->cross_section;
253  this->used_fields_ += eq_fields_->sources_sigma_out;
255  this->used_fields_ += eq_fields_->bc_type;
256  this->used_fields_ += eq_fields_->bc_robin_sigma;
257  }
258 
259  /// Destructor.
261  for (auto a : averages) if (a != nullptr) delete[] a;
262  for (auto a : waverages) if (a != nullptr) delete[] a;
263  for (auto a : jumps) if (a != nullptr) delete[] a;
264  }
265 
266  /// Initialize auxiliary vectors and other data members
267  void initialize(ElementCacheMap *element_cache_map) {
268  this->element_cache_map_ = element_cache_map;
269 
272  this->fe_values_->template initialize<dim>(*this->quad_, u);
273  this->fe_values_->template initialize<dim>(*this->quad_low_, u_side);
274  if (dim==1) { // print to log only one time
275  DebugOut() << "List of StiffnessAssemblyDG FEValues (cell) updates flags: " << this->print_update_flags(u);
276  DebugOut() << "List of StiffnessAssemblyDG FEValues (side) updates flags: " << this->print_update_flags(u_side);
277  }
278  ndofs_ = this->fe_values_->n_dofs(dim);
279  qsize_lower_dim_ = this->quad_low_->size();
280  dof_indices_.resize(ndofs_);
281  side_dof_indices_vb_.resize(2*ndofs_);
282  local_matrix_.resize(4*ndofs_*ndofs_);
283 
284  for (unsigned int sid=0; sid<eq_data_->max_edg_sides; sid++)
285  {
287  }
288 
290  for (unsigned int s=0; s<eq_data_->max_edg_sides; s++)
291  averages[s] = new double[qsize_lower_dim_*ndofs_];
292  waverages.resize(2);
293  jumps.resize(2);
294  for (unsigned int s=0; s<2; s++)
295  {
296  waverages[s] = new double[qsize_lower_dim_*ndofs_];
297  jumps[s] = new double[qsize_lower_dim_*ndofs_];
298  }
299  }
300 
301 
302  /// Assembles the cell (volume) integral into the stiffness matrix.
303  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
304  {
305  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
306  if (!cell.is_own()) return;
307 
309 
310  // assemble the local stiffness matrix
311  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
312  {
313  for (unsigned int i=0; i<ndofs_; i++)
314  for (unsigned int j=0; j<ndofs_; j++)
315  local_matrix_[i*ndofs_+j] = 0;
316 
317  for (auto p : this->bulk_points(element_patch_idx) )
318  {
319  for (unsigned int i=0; i<ndofs_; i++)
320  {
321  arma::vec3 Kt_grad_i = eq_fields_->diffusion_coef[sbi](p).t()*conc_grad_(i,p);
322  double ad_dot_grad_i = arma::dot(eq_fields_->advection_coef[sbi](p), conc_grad_(i,p));
323 
324  for (unsigned int j=0; j<ndofs_; j++)
325  local_matrix_[i*ndofs_+j] += (arma::dot(Kt_grad_i, conc_grad_(j,p))
326  -conc_shape_(j,p)*ad_dot_grad_i
327  +eq_fields_->sources_sigma_out[sbi](p)*conc_shape_(j,p)*conc_shape_(i,p))*JxW_(p);
328  }
329  }
331  }
332  }
333 
334 
335  /// Assembles the fluxes on the boundary.
336  inline void boundary_side_integral(DHCellSide cell_side)
337  {
338  ASSERT_EQ(cell_side.dim(), dim).error("Dimension of element mismatch!");
339  if (!cell_side.cell().is_own()) return;
340 
341  Side side = cell_side.side();
342  const DHCellAccessor &cell = cell_side.cell();
343 
345  unsigned int k;
346  double gamma_l;
347 
348  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
349  {
350  std::fill(local_matrix_.begin(), local_matrix_.end(), 0);
351 
352  double side_flux = advective_flux(eq_fields_->advection_coef[sbi], this->boundary_points(cell_side), JxW_side_, normal_);
353  double transport_flux = side_flux/side.measure();
354 
355  // On Neumann boundaries we have only term from integrating by parts the advective term,
356  // on Dirichlet boundaries we additionally apply the penalty which enforces the prescribed value.
357  auto p_side = *( this->boundary_points(cell_side).begin() );
358  auto p_bdr = p_side.point_bdr( side.cond().element_accessor() );
359  unsigned int bc_type = eq_fields_->bc_type[sbi](p_bdr);
361  {
362  // set up the parameters for DG method
363  auto p = *( this->boundary_points(cell_side).begin() );
364  gamma_l = DG_penalty_boundary(side,
365  diffusion_delta(eq_fields_->diffusion_coef[sbi], this->boundary_points(cell_side), normal_(p)),
366  transport_flux,
367  eq_fields_->dg_penalty[sbi](p_side));
368  transport_flux += gamma_l;
369  }
370 
371  // fluxes and penalty
372  k=0;
373  for (auto p : this->boundary_points(cell_side) )
374  {
375  double flux_times_JxW;
377  {
378  //sigma_ corresponds to robin_sigma
379  auto p_bdr = p.point_bdr(side.cond().element_accessor());
380  flux_times_JxW = eq_fields_->cross_section(p)*eq_fields_->bc_robin_sigma[sbi](p_bdr)*JxW_side_(p);
381  }
383  {
384  auto p_bdr = p.point_bdr(side.cond().element_accessor());
385  flux_times_JxW = (transport_flux + eq_fields_->cross_section(p)*eq_fields_->bc_robin_sigma[sbi](p_bdr))*JxW_side_(p);
386  }
387  else if (bc_type == AdvectionDiffusionModel::abc_inflow && side_flux < 0)
388  flux_times_JxW = 0;
389  else
390  flux_times_JxW = transport_flux*JxW_side_(p);
391 
392  for (unsigned int i=0; i<ndofs_; i++)
393  {
394  for (unsigned int j=0; j<ndofs_; j++)
395  {
396  // flux due to advection and penalty
397  local_matrix_[i*ndofs_+j] += flux_times_JxW*conc_shape_side_(i,p)*conc_shape_side_(j,p);
398 
399  // flux due to diffusion (only on dirichlet and inflow boundary)
401  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)
402  + arma::dot(eq_fields_->diffusion_coef[sbi](p)*conc_grad_sidw_(i,p),normal_(p))*conc_shape_side_(j,p)*eq_data_->dg_variant
403  )*JxW_side_(p);
404  }
405  }
406  k++;
407  }
408 
410  }
411  }
412 
413 
414  /// Assembles the fluxes between sides of elements of the same dimension.
416  ASSERT_EQ(edge_side_range.begin()->element().dim(), dim).error("Dimension of element mismatch!");
417 
418  unsigned int k;
419  double gamma_l, omega[2], transport_flux, delta[2], delta_sum;
420  double aniso1, aniso2;
421  double local_alpha=0.0;
422  int sid=0, s1, s2;
423  for( DHCellSide edge_side : edge_side_range )
424  {
425  auto dh_edge_cell = eq_data_->dh_->cell_accessor_from_element( edge_side.elem_idx() );
426  dh_edge_cell.get_dof_indices(side_dof_indices_[sid]);
427  ++sid;
428  }
429  auto zero_edge_side = *edge_side_range.begin();
430  auto p = *( this->edge_points(zero_edge_side).begin() );
431  arma::vec3 normal_vector = normal_(p);
432 
433  // fluxes and penalty
434  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
435  {
436  vector<double> fluxes(edge_side_range.begin()->n_edge_sides());
437  double pflux = 0, nflux = 0; // calculate the total in- and out-flux through the edge
438  sid=0;
439  for( DHCellSide edge_side : edge_side_range )
440  {
441  fluxes[sid] = advective_flux(eq_fields_->advection_coef[sbi], this->edge_points(edge_side), JxW_side_, normal_) / edge_side.measure();
442  if (fluxes[sid] > 0)
443  pflux += fluxes[sid];
444  else
445  nflux += fluxes[sid];
446  ++sid;
447  }
448 
449  // precompute averages of shape functions over pairs of sides
450  s1=0;
451  for (DHCellSide edge_side : edge_side_range)
452  {
453  k=0;
454  for (auto p : this->edge_points(edge_side) )
455  {
456  for (unsigned int i=0; i<ndofs_; i++)
457  averages[s1][k*ndofs_+i] = conc_shape_side_(i,p)*0.5;
458  k++;
459  }
460  s1++;
461  }
462 
463 
464 
465  s1=0;
466  for( DHCellSide edge_side1 : edge_side_range )
467  {
468  s2=-1; // need increment at begin of loop (see conditionally 'continue' directions)
469  for( DHCellSide edge_side2 : edge_side_range )
470  {
471  s2++;
472  if (s2<=s1) continue;
473  ASSERT(edge_side1.is_valid()).error("Invalid side of edge.");
474 
475  auto p = *( this->edge_points(edge_side1).begin() );
476  arma::vec3 nv = normal_(p);
477 
478  // set up the parameters for DG method
479  // calculate the flux from edge_side1 to edge_side2
480  if (fluxes[s2] > 0 && fluxes[s1] < 0)
481  transport_flux = fluxes[s1]*fabs(fluxes[s2]/pflux);
482  else if (fluxes[s2] < 0 && fluxes[s1] > 0)
483  transport_flux = fluxes[s1]*fabs(fluxes[s2]/nflux);
484  else
485  transport_flux = 0;
486 
487  gamma_l = 0.5*fabs(transport_flux);
488 
489  delta[0] = 0;
490  delta[1] = 0;
491  for (auto p1 : this->edge_points(edge_side1) )
492  {
493  auto p2 = p1.point_on(edge_side2);
494  delta[0] += dot(eq_fields_->diffusion_coef[sbi](p1)*normal_vector,normal_vector);
495  delta[1] += dot(eq_fields_->diffusion_coef[sbi](p2)*normal_vector,normal_vector);
496  local_alpha = max(eq_fields_->dg_penalty[sbi](p1), eq_fields_->dg_penalty[sbi](p2));
497  }
498  delta[0] /= qsize_lower_dim_;
499  delta[1] /= qsize_lower_dim_;
500 
501  delta_sum = delta[0] + delta[1];
502 
503 // if (delta_sum > numeric_limits<double>::epsilon())
504  if (fabs(delta_sum) > 0)
505  {
506  omega[0] = delta[1]/delta_sum;
507  omega[1] = delta[0]/delta_sum;
508  double h = edge_side1.diameter();
509  aniso1 = elem_anisotropy(edge_side1.element());
510  aniso2 = elem_anisotropy(edge_side2.element());
511  gamma_l += local_alpha/h*aniso1*aniso2*(delta[0]*delta[1]/delta_sum);
512  }
513  else
514  for (int i=0; i<2; i++) omega[i] = 0;
515  // end of set up the parameters for DG method
516 
517  int sd[2]; bool is_side_own[2];
518  sd[0] = s1; is_side_own[0] = edge_side1.cell().is_own();
519  sd[1] = s2; is_side_own[1] = edge_side2.cell().is_own();
520 
521  // precompute jumps and weighted averages of shape functions over the pair of sides (s1,s2)
522  k=0;
523  for (auto p1 : this->edge_points(edge_side1) )
524  {
525  auto p2 = p1.point_on(edge_side2);
526  for (unsigned int i=0; i<ndofs_; i++)
527  {
528  jumps[0][k*ndofs_+i] = conc_shape_side_(i,p1);
529  jumps[1][k*ndofs_+i] = - conc_shape_side_(i,p2);
530  waverages[0][k*ndofs_+i] = arma::dot(eq_fields_->diffusion_coef[sbi](p1)*conc_grad_sidw_(i,p1),nv)*omega[0];
531  waverages[1][k*ndofs_+i] = arma::dot(eq_fields_->diffusion_coef[sbi](p2)*conc_grad_sidw_(i,p2),nv)*omega[1];
532  }
533  k++;
534  }
535 
536  // For selected pair of elements:
537  for (int n=0; n<2; n++)
538  {
539  if (!is_side_own[n]) continue;
540 
541  for (int m=0; m<2; m++)
542  {
543  for (unsigned int i=0; i<this->fe_values_->n_dofs(dim); i++)
544  for (unsigned int j=0; j<this->fe_values_->n_dofs(dim); j++)
545  local_matrix_[i*this->fe_values_->n_dofs(dim)+j] = 0;
546 
547  k=0;
548  for (auto p1 : this->edge_points(zero_edge_side) )
549  //for (k=0; k<this->quad_low_->size(); ++k)
550  {
551  for (unsigned int i=0; i<this->fe_values_->n_dofs(dim); i++)
552  {
553  for (unsigned int j=0; j<this->fe_values_->n_dofs(dim); j++)
554  {
555  int index = i*this->fe_values_->n_dofs(dim)+j;
556 
557  local_matrix_[index] += (
558  // flux due to transport (applied on interior edges) (average times jump)
559  transport_flux*jumps[n][k*ndofs_+i]*averages[sd[m]][k*ndofs_+j]
560 
561  // penalty enforcing continuity across edges (applied on interior and Dirichlet edges) (jump times jump)
562  + gamma_l*jumps[n][k*ndofs_+i]*jumps[m][k*ndofs_+j]
563 
564  // terms due to diffusion
565  - jumps[n][k*ndofs_+i]*waverages[m][k*ndofs_+j]
566  - eq_data_->dg_variant*waverages[n][k*ndofs_+i]*jumps[m][k*ndofs_+j]
568  }
569  }
570  k++;
571  }
572  eq_data_->ls[sbi]->mat_set_values(ndofs_, &(side_dof_indices_[sd[n]][0]), ndofs_, &(side_dof_indices_[sd[m]][0]), &(local_matrix_[0]));
573  }
574  }
575  }
576  s1++;
577  }
578  }
579  }
580 
581 
582  /// Assembles the fluxes between elements of different dimensions.
583  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
584  if (dim == 1) return;
585  ASSERT_EQ(cell_lower_dim.dim(), dim-1).error("Dimension of element mismatch!");
586 
587  // Note: use data members csection_ and velocity_ for appropriate quantities of lower dim element
588 
589  unsigned int n_dofs[2];
590  unsigned int n_indices = cell_lower_dim.get_dof_indices(dof_indices_);
591  for(unsigned int i=0; i<n_indices; ++i) {
593  }
594  n_dofs[0] = conc_join_shape_.begin()->n_dofs_low();
595 
596  DHCellAccessor cell_higher_dim = eq_data_->dh_->cell_accessor_from_element( neighb_side.element().idx() );
597  n_indices = cell_higher_dim.get_dof_indices(dof_indices_);
598  for(unsigned int i=0; i<n_indices; ++i) {
599  side_dof_indices_vb_[i+n_dofs[0]] = dof_indices_[i];
600  }
601  n_dofs[1] = conc_join_shape_.begin()->n_dofs_high();
602 
603  // Testing element if they belong to local partition.
604  bool own_element_id[2];
605  own_element_id[0] = cell_lower_dim.is_own();
606  own_element_id[1] = cell_higher_dim.is_own();
607 
608  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++) // Optimize: SWAP LOOPS
609  {
610  for (unsigned int i=0; i<n_dofs[0]+n_dofs[1]; i++)
611  for (unsigned int j=0; j<n_dofs[0]+n_dofs[1]; j++)
612  local_matrix_[i*(n_dofs[0]+n_dofs[1])+j] = 0;
613 
614  // set transmission conditions
615  for (auto p_high : this->coupling_points(neighb_side) )
616  {
617  auto p_low = p_high.lower_dim(cell_lower_dim);
618  auto p_low_conv = p_high.lower_dim_converted(cell_lower_dim);
619  // The communication flux has two parts:
620  // - "diffusive" term containing sigma
621  // - "advective" term representing usual upwind
622  //
623  // The calculation differs from the reference manual, since ad_coef and dif_coef have different meaning
624  // than b and A in the manual.
625  // In calculation of sigma there appears one more csection_lower in the denominator.
626  double sigma = eq_fields_->fracture_sigma[sbi](p_low)*arma::dot(eq_fields_->diffusion_coef[sbi](p_low)*normal_(p_high),normal_(p_high))*
627  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));
628 
629  double transport_flux = arma::dot(eq_fields_->advection_coef[sbi](p_high), normal_(p_high));
630 
631  for( auto conc_shape_i : conc_join_shape_) {
632  uint is_high_i = conc_shape_i.is_high_dim();
633  if (!own_element_id[is_high_i]) continue;
634  uint i_mat_idx = conc_shape_i.join_idx(); // i + is_high * n_dofs_low
635  double diff_shape_i = conc_shape_i(p_high) - conc_shape_i(p_low_conv);
636  for( auto conc_shape_j : conc_join_shape_) {
637  uint j_mat_idx = conc_shape_j.join_idx();
638  local_matrix_[i_mat_idx * (n_dofs[0]+n_dofs[1]) + j_mat_idx] += (
639  sigma * diff_shape_i * (conc_shape_j(p_high) - conc_shape_j(p_low_conv))
640  + diff_shape_i * ( max(0.,transport_flux) * conc_shape_j(p_high) + min(0.,transport_flux) * conc_shape_j(p_low_conv))
642  }
643  }
644  }
645  eq_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]));
646  }
647  }
648 
649 
650 private:
651  /// Data objects shared with TransportDG
654 
655  /// Sub field set contains fields used in calculation.
657 
658  unsigned int ndofs_; ///< Number of dofs
659  unsigned int qsize_lower_dim_; ///< Size of quadrature of dim-1
660 
661  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
662  vector< vector<LongIdx> > side_dof_indices_; ///< Vector of vectors of side DOF indices
663  vector<LongIdx> side_dof_indices_vb_; ///< Vector of side DOF indices (assemble element-side fluxex)
664  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
665 
666  vector<double*> averages; ///< Auxiliary storage for averages of shape functions.
667  vector<double*> waverages; ///< Auxiliary storage for weighted averages of shape functions.
668  vector<double*> jumps; ///< Auxiliary storage for jumps of shape functions.
669 
678 
679  template < template<IntDim...> class DimAssembly>
680  friend class GenericAssembly;
681 
682 };
683 
684 
685 /**
686  * Auxiliary container class for Finite element and related objects of given dimension.
687  */
688 template <unsigned int dim, class Model>
690 {
691 public:
694 
695  static constexpr const char * name() { return "SourcesAssemblyDG"; }
696 
697  /// Constructor.
698  SourcesAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues<3> *fe_values)
699  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
700  JxW_( this->bulk_values().JxW() ),
701  conc_shape_( this->bulk_values().scalar_shape() ) {
703  this->used_fields_ += eq_fields_->sources_density_out;
704  this->used_fields_ += eq_fields_->sources_conc_out;
705  this->used_fields_ += eq_fields_->sources_sigma_out;
706  }
707 
708  /// Destructor.
710 
711  /// Initialize auxiliary vectors and other data members
712  void initialize(ElementCacheMap *element_cache_map) {
713  this->element_cache_map_ = element_cache_map;
714 
716  this->fe_values_->template initialize<dim>(*this->quad_, u);
717  if (dim==1) // print to log only one time
718  DebugOut() << "List of SourcesAssemblyDG FEValues updates flags: " << this->print_update_flags(u);
719  ndofs_ = this->fe_values_->n_dofs(dim);
720  dof_indices_.resize(ndofs_);
721  local_rhs_.resize(ndofs_);
724  }
725 
726 
727  /// Assemble integral over element
728  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
729  {
730  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
731 
732  ElementAccessor<3> elm = cell.elm();
733  double source;
734 
736 
737  // assemble the local stiffness matrix
738  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
739  {
740  fill_n( &(local_rhs_[0]), ndofs_, 0 );
743 
744  for (auto p : this->bulk_points(element_patch_idx) )
745  {
746  source = (eq_fields_->sources_density_out[sbi](p) + eq_fields_->sources_conc_out[sbi](p)*eq_fields_->sources_sigma_out[sbi](p))*JxW_(p);
747 
748  for (unsigned int i=0; i<ndofs_; i++)
749  local_rhs_[i] += source*conc_shape_(i,p);
750  }
751  eq_data_->ls[sbi]->rhs_set_values(ndofs_, &(dof_indices_[0]), &(local_rhs_[0]));
752 
753  for (unsigned int i=0; i<ndofs_; i++)
754  {
755  for (auto p : this->bulk_points(element_patch_idx) )
756  {
757  local_source_balance_vector_[i] -= eq_fields_->sources_sigma_out[sbi](p)*conc_shape_(i,p)*JxW_(p);
758  }
759 
761  }
762  eq_data_->balance_->add_source_values(eq_data_->subst_idx()[sbi], elm.region().bulk_idx(),
763  cell.get_loc_dof_indices(),
765  }
766  }
767 
768  /// Implements @p AssemblyBase::begin.
769  void begin() override
770  {
771  eq_data_->balance_->start_source_assembly( eq_data_->subst_idx() );
772  }
773 
774  /// Implements @p AssemblyBase::end.
775  void end() override
776  {
777  eq_data_->balance_->finish_source_assembly( eq_data_->subst_idx() );
778  }
779 
780 
781  private:
782  /// Data objects shared with TransportDG
785 
786  /// Sub field set contains fields used in calculation.
788 
789  unsigned int ndofs_; ///< Number of dofs
790 
791  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
792  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for set_sources method.
793  vector<PetscScalar> local_source_balance_vector_; ///< Auxiliary vector for set_sources method.
794  vector<PetscScalar> local_source_balance_rhs_; ///< Auxiliary vector for set_sources method.
795 
798 
799  template < template<IntDim...> class DimAssembly>
800  friend class GenericAssembly;
801 
802 };
803 
804 
805 /**
806  * Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions..
807  */
808 template <unsigned int dim, class Model>
810 {
811 public:
814 
815  static constexpr const char * name() { return "BdrConditionAssemblyDG"; }
816 
817  /// Constructor.
818  BdrConditionAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues<3> *fe_values)
819  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
820  JxW_( this->side_values().JxW() ),
821  normal_( this->side_values().normal_vector() ),
822  conc_shape_( this->side_values().scalar_shape() ),
823  conc_grad_( this->side_values().grad_scalar_shape() ) {
825  this->used_fields_ += eq_fields_->advection_coef;
826  this->used_fields_ += eq_fields_->diffusion_coef;
827  this->used_fields_ += eq_fields_->cross_section;
828  this->used_fields_ += eq_fields_->bc_type;
829  this->used_fields_ += eq_fields_->bc_dirichlet_value;
830  this->used_fields_ += eq_fields_->bc_robin_sigma;
831  this->used_fields_ += eq_fields_->bc_flux;
832  }
833 
834  /// Destructor.
836 
837  /// Initialize auxiliary vectors and other data members
838  void initialize(ElementCacheMap *element_cache_map) {
839  this->element_cache_map_ = element_cache_map;
840 
841  fe_ = std::make_shared< FE_P_disc<dim> >(eq_data_->dg_order);
843  this->fe_values_->template initialize<dim>(*this->quad_low_, u);
844  if (dim==1) // print to log only one time
845  DebugOut() << "List of BdrConditionAssemblyDG FEValues updates flags: " << this->print_update_flags(u);
846  ndofs_ = fe_->n_dofs();
847  dof_indices_.resize(ndofs_);
848  local_rhs_.resize(ndofs_);
850  }
851 
852 
853  /// Implements @p AssemblyBase::boundary_side_integral.
854  inline void boundary_side_integral(DHCellSide cell_side)
855  {
856 
857  ElementAccessor<3> bc_elm = cell_side.cond().element_accessor();
858 
859  const DHCellAccessor &cell = cell_side.cell();
861 
862  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
863  {
864  fill_n(&(local_rhs_[0]), ndofs_, 0);
867 
868  double side_flux = advective_flux(eq_fields_->advection_coef[sbi], this->boundary_points(cell_side), JxW_, normal_);
869  double transport_flux = side_flux/cell_side.measure();
870 
871  auto p_side = *( this->boundary_points(cell_side)).begin();
872  auto p_bdr = p_side.point_bdr(cell_side.cond().element_accessor() );
873  unsigned int bc_type = eq_fields_->bc_type[sbi](p_bdr);
874  if (bc_type == AdvectionDiffusionModel::abc_inflow && side_flux < 0)
875  {
876  for (auto p : this->boundary_points(cell_side) )
877  {
878  auto p_bdr = p.point_bdr(bc_elm);
879  double bc_term = -transport_flux*eq_fields_->bc_dirichlet_value[sbi](p_bdr)*JxW_(p);
880  for (unsigned int i=0; i<ndofs_; i++)
881  local_rhs_[i] += bc_term*conc_shape_(i,p);
882  }
883  for (unsigned int i=0; i<ndofs_; i++)
885  }
886  else if (bc_type == AdvectionDiffusionModel::abc_dirichlet)
887  {
888  double side_flux = advective_flux(eq_fields_->advection_coef[sbi], this->boundary_points(cell_side), JxW_, normal_);
889  double transport_flux = side_flux/cell_side.measure();
890 
891  auto p = *( this->boundary_points(cell_side).begin() );
892  double gamma_l = DG_penalty_boundary(cell_side.side(),
893  diffusion_delta(eq_fields_->diffusion_coef[sbi], this->boundary_points(cell_side), normal_(p)),
894  transport_flux,
895  eq_fields_->dg_penalty[sbi](p_bdr));
896 
897  for (auto p : this->boundary_points(cell_side) )
898  {
899  auto p_bdr = p.point_bdr(bc_elm);
900  double bc_term = gamma_l*eq_fields_->bc_dirichlet_value[sbi](p_bdr)*JxW_(p);
901  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));
902  for (unsigned int i=0; i<ndofs_; i++)
903  local_rhs_[i] += bc_term*conc_shape_(i,p)
904  + arma::dot(bc_grad,conc_grad_(i,p));
905  }
906  for (auto p : this->boundary_points(cell_side) )
907  {
908  for (unsigned int i=0; i<ndofs_; i++)
909  {
910  local_flux_balance_vector_[i] += (arma::dot(eq_fields_->advection_coef[sbi](p), normal_(p))*conc_shape_(i,p)
911  - arma::dot(eq_fields_->diffusion_coef[sbi](p)*conc_grad_(i,p),normal_(p))
912  + gamma_l*conc_shape_(i,p))*JxW_(p);
913  }
914  }
915  if (eq_data_->time_->step().index() > 0)
916  for (unsigned int i=0; i<ndofs_; i++)
918  }
919  else if (bc_type == AdvectionDiffusionModel::abc_total_flux)
920  {
921  for (auto p : this->boundary_points(cell_side) )
922  {
923  auto p_bdr = p.point_bdr(bc_elm);
924  double bc_term = eq_fields_->cross_section(p) * (eq_fields_->bc_robin_sigma[sbi](p_bdr)*eq_fields_->bc_dirichlet_value[sbi](p_bdr) +
925  eq_fields_->bc_flux[sbi](p_bdr)) * JxW_(p);
926  for (unsigned int i=0; i<ndofs_; i++)
927  local_rhs_[i] += bc_term*conc_shape_(i,p);
928  }
929 
930  for (unsigned int i=0; i<ndofs_; i++)
931  {
932  for (auto p : this->boundary_points(cell_side) ) {
933  auto p_bdr = p.point_bdr(bc_elm);
934  local_flux_balance_vector_[i] += eq_fields_->cross_section(p) * eq_fields_->bc_robin_sigma[sbi](p_bdr) *
935  JxW_(p) * conc_shape_(i,p);
936  }
938  }
939  }
941  {
942  for (auto p : this->boundary_points(cell_side) )
943  {
944  auto p_bdr = p.point_bdr(bc_elm);
945  double bc_term = eq_fields_->cross_section(p) * (eq_fields_->bc_robin_sigma[sbi](p_bdr)*eq_fields_->bc_dirichlet_value[sbi](p_bdr) +
946  eq_fields_->bc_flux[sbi](p_bdr)) * JxW_(p);
947  for (unsigned int i=0; i<ndofs_; i++)
948  local_rhs_[i] += bc_term*conc_shape_(i,p);
949  }
950 
951  for (unsigned int i=0; i<ndofs_; i++)
952  {
953  for (auto p : this->boundary_points(cell_side) ) {
954  auto p_bdr = p.point_bdr(bc_elm);
955  local_flux_balance_vector_[i] += eq_fields_->cross_section(p)*(arma::dot(eq_fields_->advection_coef[sbi](p), normal_(p)) +
956  eq_fields_->bc_robin_sigma[sbi](p_bdr))*JxW_(p)*conc_shape_(i,p);
957  }
959  }
960  }
961  else if (bc_type == AdvectionDiffusionModel::abc_inflow && side_flux >= 0)
962  {
963  for (auto p : this->boundary_points(cell_side) )
964  {
965  for (unsigned int i=0; i<ndofs_; i++)
966  local_flux_balance_vector_[i] += arma::dot(eq_fields_->advection_coef[sbi](p), normal_(p))*JxW_(p)*conc_shape_(i,p);
967  }
968  }
969  eq_data_->ls[sbi]->rhs_set_values(ndofs_, &(dof_indices_[0]), &(local_rhs_[0]));
970 
971  eq_data_->balance_->add_flux_values(eq_data_->subst_idx()[sbi], cell_side,
972  cell.get_loc_dof_indices(),
974  }
975  }
976 
977  /// Implements @p AssemblyBase::begin.
978  void begin() override
979  {
980  eq_data_->balance_->start_flux_assembly( eq_data_->subst_idx() );
981  }
982 
983  /// Implements @p AssemblyBase::end.
984  void end() override
985  {
986  eq_data_->balance_->finish_flux_assembly( eq_data_->subst_idx() );
987  }
988 
989 
990  private:
991  /// Data objects shared with TransportDG
994 
995  shared_ptr<FiniteElement<dim>> fe_; ///< Finite element for the solution of the advection-diffusion equation.
996  unsigned int ndofs_; ///< Number of dofs
997 
998  /// Sub field set contains fields used in calculation.
1000 
1001  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
1002  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for set_sources method.
1003  vector<PetscScalar> local_flux_balance_vector_; ///< Auxiliary vector for set_boundary_conditions method.
1004  PetscScalar local_flux_balance_rhs_; ///< Auxiliary variable for set_boundary_conditions method.
1005 
1010 
1011  template < template<IntDim...> class DimAssembly>
1012  friend class GenericAssembly;
1013 
1014 };
1015 
1016 
1017 /**
1018  * Auxiliary container class sets the initial condition.
1019  */
1020 template <unsigned int dim, class Model>
1022 {
1023 public:
1026 
1027  static constexpr const char * name() { return "InitProjectionAssemblyDG"; }
1028 
1029  /// Constructor.
1030  InitProjectionAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues<3> *fe_values)
1031  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
1032  JxW_( this->bulk_values().JxW() ),
1033  init_shape_( this->bulk_values().scalar_shape() ) {
1035  this->used_fields_ += eq_fields_->init_condition;
1036  }
1037 
1038  /// Destructor.
1040 
1041  /// Initialize auxiliary vectors and other data members
1042  void initialize(ElementCacheMap *element_cache_map) {
1043  this->element_cache_map_ = element_cache_map;
1044 
1046  this->fe_values_->template initialize<dim>(*this->quad_, u);
1047  // if (dim==1) // print to log only one time
1048  // DebugOut() << "List of InitProjectionAssemblyDG FEValues updates flags: " << this->print_update_flags(u);
1049  ndofs_ = this->fe_values_->n_dofs(dim);
1050  dof_indices_.resize(ndofs_);
1051  local_matrix_.resize(4*ndofs_*ndofs_);
1052  local_rhs_.resize(ndofs_);
1053  }
1054 
1055 
1056  /// Assemble integral over element
1057  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
1058  {
1059  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
1060 
1062 
1063  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
1064  {
1065  for (unsigned int i=0; i<ndofs_; i++)
1066  {
1067  local_rhs_[i] = 0;
1068  for (unsigned int j=0; j<ndofs_; j++)
1069  local_matrix_[i*ndofs_+j] = 0;
1070  }
1071 
1072  for (auto p : this->bulk_points(element_patch_idx) )
1073  {
1074  double rhs_term = eq_fields_->init_condition[sbi](p)*JxW_(p);
1075 
1076  for (unsigned int i=0; i<ndofs_; i++)
1077  {
1078  for (unsigned int j=0; j<ndofs_; j++)
1079  local_matrix_[i*ndofs_+j] += init_shape_(i,p)*init_shape_(j,p)*JxW_(p);
1080 
1081  local_rhs_[i] += init_shape_(i,p)*rhs_term;
1082  }
1083  }
1084  eq_data_->ls[sbi]->set_values(ndofs_, &(dof_indices_[0]), ndofs_, &(dof_indices_[0]), &(local_matrix_[0]), &(local_rhs_[0]));
1085  }
1086  }
1087 
1088 
1089  private:
1090  /// Data objects shared with TransportDG
1093 
1094  /// Sub field set contains fields used in calculation.
1096 
1097  unsigned int ndofs_; ///< Number of dofs
1098  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
1099  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
1100  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for set_sources method.
1101 
1104 
1105  template < template<IntDim...> class DimAssembly>
1106  friend class GenericAssembly;
1107 
1108 };
1109 
1110 
1111 /**
1112  * Auxiliary container class sets the initial condition.
1113  */
1114 template <unsigned int dim, class Model>
1116 {
1117 public:
1120 
1121  static constexpr const char * name() { return "InitConditionAssemblyDG"; }
1122 
1123  /// Constructor.
1125  : AssemblyBase<dim>(), eq_fields_(eq_fields), eq_data_(eq_data) {
1127  this->used_fields_ += eq_fields_->init_condition;
1128 
1129  this->quad_ = new Quadrature(dim, RefElement<dim>::n_nodes);
1130  for(unsigned int i = 0; i<RefElement<dim>::n_nodes; i++)
1131  {
1132  this->quad_->weight(i) = 1.0;
1133  this->quad_->set(i) = RefElement<dim>::node_coords(i);
1134  }
1135  }
1136 
1137  /// Destructor.
1139 
1140  /// Initialize auxiliary vectors and other data members
1141  void initialize(ElementCacheMap *element_cache_map) {
1142  this->element_cache_map_ = element_cache_map;
1143  }
1144 
1145 
1146  /// Assemble integral over element
1147  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
1148  {
1149  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
1150 
1151  unsigned int k;
1152 
1153  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
1154  {
1155  k=0;
1156  for (auto p : this->bulk_points(element_patch_idx) )
1157  {
1158  double val = eq_fields_->init_condition[sbi](p);
1159 
1160  eq_data_->output_vec[sbi].set(cell.get_loc_dof_indices()[k], val);
1161  k++;
1162  }
1163  }
1164  }
1165 
1166 
1167  private:
1168  /// Data objects shared with TransportDG
1171 
1172  /// Sub field set contains fields used in calculation.
1174 
1175  template < template<IntDim...> class DimAssembly>
1176  friend class GenericAssembly;
1177 
1178 };
1179 
1180 #endif /* ASSEMBLY_DG_HH_ */
1181 
double elem_anisotropy(ElementAccessor< 3 > e)
Definition: assembly_dg.hh:152
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:196
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:173
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:215
#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.
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.
std::string print_update_flags(UpdateFlags u) const
Print update flags to string format.
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:815
FieldSet used_fields_
Sub field set contains fields used in calculation.
Definition: assembly_dg.hh:999
vector< PetscScalar > local_flux_balance_vector_
Auxiliary vector for set_boundary_conditions method.
void boundary_side_integral(DHCellSide cell_side)
Implements AssemblyBase::boundary_side_integral.
Definition: assembly_dg.hh:854
BdrConditionAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
Definition: assembly_dg.hh:818
TransportDG< Model >::EqData EqData
Definition: assembly_dg.hh:813
TransportDG< Model >::EqFields EqFields
Definition: assembly_dg.hh:812
void begin() override
Implements AssemblyBase::begin.
Definition: assembly_dg.hh:978
~BdrConditionAssemblyDG()
Destructor.
Definition: assembly_dg.hh:835
PetscScalar local_flux_balance_rhs_
Auxiliary variable for set_boundary_conditions method.
FeQ< Scalar > conc_shape_
shared_ptr< FiniteElement< dim > > fe_
Finite element for the solution of the advection-diffusion equation.
Definition: assembly_dg.hh:995
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
Definition: assembly_dg.hh:838
vector< LongIdx > dof_indices_
Vector of global DOF indices.
EqFields * eq_fields_
Data objects shared with TransportDG.
Definition: assembly_dg.hh:992
void end() override
Implements AssemblyBase::end.
Definition: assembly_dg.hh:984
vector< PetscScalar > local_rhs_
Auxiliary vector for set_sources method.
unsigned int ndofs_
Number of dofs.
Definition: assembly_dg.hh:996
FeQ< Vector > conc_grad_
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:137
vector< PetscScalar > local_mass_balance_vector_
Same as previous.
Definition: assembly_dg.hh:140
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:134
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
Definition: assembly_dg.hh:138
void end() override
Implements AssemblyBase::end.
Definition: assembly_dg.hh:123
~MassAssemblyDG()
Destructor.
Definition: assembly_dg.hh:55
ElQ< Scalar > JxW_
Definition: assembly_dg.hh:142
EqData * eq_data_
Definition: assembly_dg.hh:131
void begin() override
Implements AssemblyBase::begin.
Definition: assembly_dg.hh:117
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:143
unsigned int ndofs_
Number of dofs.
Definition: assembly_dg.hh:136
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
Definition: assembly_dg.hh:75
TransportDG< Model >::EqData EqData
Definition: assembly_dg.hh:40
EqFields * eq_fields_
Data objects shared with TransportDG.
Definition: assembly_dg.hh:130
vector< PetscScalar > local_retardation_balance_vector_
Auxiliary vector for assemble mass matrix.
Definition: assembly_dg.hh:139
unsigned int n_dofs(unsigned int dim, FMT_UNUSED uint component_idx=0) const
Returns the number of shape functions.
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:792
TransportDG< Model >::EqData EqData
Definition: assembly_dg.hh:693
TransportDG< Model >::EqFields EqFields
Definition: assembly_dg.hh:692
ElQ< Scalar > JxW_
Definition: assembly_dg.hh:796
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
Definition: assembly_dg.hh:728
~SourcesAssemblyDG()
Destructor.
Definition: assembly_dg.hh:709
SourcesAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
Definition: assembly_dg.hh:698
vector< PetscScalar > local_source_balance_vector_
Auxiliary vector for set_sources method.
Definition: assembly_dg.hh:793
void end() override
Implements AssemblyBase::end.
Definition: assembly_dg.hh:775
static constexpr const char * name()
Definition: assembly_dg.hh:695
EqFields * eq_fields_
Data objects shared with TransportDG.
Definition: assembly_dg.hh:783
vector< PetscScalar > local_source_balance_rhs_
Auxiliary vector for set_sources method.
Definition: assembly_dg.hh:794
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
Definition: assembly_dg.hh:712
FieldSet used_fields_
Sub field set contains fields used in calculation.
Definition: assembly_dg.hh:787
unsigned int ndofs_
Number of dofs.
Definition: assembly_dg.hh:789
vector< LongIdx > dof_indices_
Vector of global DOF indices.
Definition: assembly_dg.hh:791
void begin() override
Implements AssemblyBase::begin.
Definition: assembly_dg.hh:769
FeQ< Scalar > conc_shape_
Definition: assembly_dg.hh:797
TransportDG< Model >::EqData EqData
Definition: assembly_dg.hh:233
FeQ< Scalar > conc_shape_
Definition: assembly_dg.hh:673
void edge_integral(RangeConvert< DHEdgeSide, DHCellSide > edge_side_range)
Assembles the fluxes between sides of elements of the same dimension.
Definition: assembly_dg.hh:415
FeQ< Vector > conc_grad_sidw_
Definition: assembly_dg.hh:676
FieldSet used_fields_
Sub field set contains fields used in calculation.
Definition: assembly_dg.hh:656
unsigned int ndofs_
Number of dofs.
Definition: assembly_dg.hh:658
vector< vector< LongIdx > > side_dof_indices_
Vector of vectors of side DOF indices.
Definition: assembly_dg.hh:662
Range< JoinShapeAccessor< Scalar > > conc_join_shape_
Definition: assembly_dg.hh:677
StiffnessAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
Definition: assembly_dg.hh:238
void boundary_side_integral(DHCellSide cell_side)
Assembles the fluxes on the boundary.
Definition: assembly_dg.hh:336
vector< double * > jumps
Auxiliary storage for jumps of shape functions.
Definition: assembly_dg.hh:668
vector< double * > waverages
Auxiliary storage for weighted averages of shape functions.
Definition: assembly_dg.hh:667
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
Assembles the fluxes between elements of different dimensions.
Definition: assembly_dg.hh:583
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
Definition: assembly_dg.hh:267
vector< LongIdx > side_dof_indices_vb_
Vector of side DOF indices (assemble element-side fluxex)
Definition: assembly_dg.hh:663
FeQ< Vector > conc_grad_
Definition: assembly_dg.hh:675
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assembles the cell (volume) integral into the stiffness matrix.
Definition: assembly_dg.hh:303
TransportDG< Model >::EqFields EqFields
Definition: assembly_dg.hh:232
vector< double * > averages
Auxiliary storage for averages of shape functions.
Definition: assembly_dg.hh:666
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
Definition: assembly_dg.hh:664
~StiffnessAssemblyDG()
Destructor.
Definition: assembly_dg.hh:260
ElQ< Scalar > JxW_
Definition: assembly_dg.hh:670
static constexpr const char * name()
Definition: assembly_dg.hh:235
ElQ< Vector > normal_
Definition: assembly_dg.hh:672
vector< LongIdx > dof_indices_
Vector of global DOF indices.
Definition: assembly_dg.hh:661
ElQ< Scalar > JxW_side_
Definition: assembly_dg.hh:671
FeQ< Scalar > conc_shape_side_
Definition: assembly_dg.hh:674
EqFields * eq_fields_
Data objects shared with TransportDG.
Definition: assembly_dg.hh:652
unsigned int qsize_lower_dim_
Size of quadrature of dim-1.
Definition: assembly_dg.hh:659
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_
unsigned int dg_order
Polynomial order of finite elements.
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
#define DebugOut()
Macro defining 'debug' record of log.
Definition: logger.hh:284
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.
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
Definition: update_flags.hh:68
@ update_values
Shape function values.
Definition: update_flags.hh:87
@ update_normal_vectors
Normal vectors.
@ update_JxW_values
Transformed quadrature weights.
@ update_side_JxW_values
Transformed quadrature weight for cell sides.
@ update_gradients
Shape function gradients.
Definition: update_flags.hh:95
@ update_quadrature_points
Transformed quadrature points.