Flow123d  DF_patch_fe_data_tables-522f93f
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->n_dofs();
66  dof_indices_.resize(ndofs_);
67  local_matrix_.resize(4*ndofs_*ndofs_);
68  local_retardation_balance_vector_.resize(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->n_dofs();
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<ndofs_; i++)
544  for (unsigned int j=0; j<ndofs_; j++)
545  local_matrix_[i*ndofs_+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<ndofs_; i++)
552  {
553  for (unsigned int j=0; j<ndofs_; j++)
554  {
555  int index = i*ndofs_+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) {
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  // The communication flux has two parts:
619  // - "diffusive" term containing sigma
620  // - "advective" term representing usual upwind
621  //
622  // The calculation differs from the reference manual, since ad_coef and dif_coef have different meaning
623  // than b and A in the manual.
624  // In calculation of sigma there appears one more csection_lower in the denominator.
625  double sigma = eq_fields_->fracture_sigma[sbi](p_low)*arma::dot(eq_fields_->diffusion_coef[sbi](p_low)*normal_(p_high),normal_(p_high))*
626  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));
627 
628  double transport_flux = arma::dot(eq_fields_->advection_coef[sbi](p_high), normal_(p_high));
629 
630  for( auto conc_shape_i : conc_join_shape_) {
631  uint is_high_i = conc_shape_i.is_high_dim();
632  if (!own_element_id[is_high_i]) continue;
633  uint i_mat_idx = conc_shape_i.join_idx(); // i + is_high * n_dofs_low
634  double diff_shape_i = conc_shape_i(p_high) - conc_shape_i(p_low);
635  for( auto conc_shape_j : conc_join_shape_) {
636  uint j_mat_idx = conc_shape_j.join_idx();
637  local_matrix_[i_mat_idx * (n_dofs[0]+n_dofs[1]) + j_mat_idx] += (
638  sigma * diff_shape_i * (conc_shape_j(p_high) - conc_shape_j(p_low))
639  + diff_shape_i * ( max(0.,transport_flux) * conc_shape_j(p_high) + min(0.,transport_flux) * conc_shape_j(p_low))
641  }
642  }
643  }
645  }
646  }
647 
648 
649 private:
650  /// Data objects shared with TransportDG
653 
654  /// Sub field set contains fields used in calculation.
656 
657  unsigned int ndofs_; ///< Number of dofs
658  unsigned int qsize_lower_dim_; ///< Size of quadrature of dim-1
659 
660  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
661  vector< vector<LongIdx> > side_dof_indices_; ///< Vector of vectors of side DOF indices
662  vector<LongIdx> side_dof_indices_vb_; ///< Vector of side DOF indices (assemble element-side fluxex)
663  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
664 
665  vector<double*> averages; ///< Auxiliary storage for averages of shape functions.
666  vector<double*> waverages; ///< Auxiliary storage for weighted averages of shape functions.
667  vector<double*> jumps; ///< Auxiliary storage for jumps of shape functions.
668 
677 
678  template < template<IntDim...> class DimAssembly>
679  friend class GenericAssembly;
680 
681 };
682 
683 
684 /**
685  * Auxiliary container class for Finite element and related objects of given dimension.
686  */
687 template <unsigned int dim, class Model>
689 {
690 public:
693 
694  static constexpr const char * name() { return "SourcesAssemblyDG"; }
695 
696  /// Constructor.
697  SourcesAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues<3> *fe_values)
698  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
699  JxW_( this->bulk_values().JxW() ),
700  conc_shape_( this->bulk_values().scalar_shape() ) {
702  this->used_fields_ += eq_fields_->sources_density_out;
703  this->used_fields_ += eq_fields_->sources_conc_out;
704  this->used_fields_ += eq_fields_->sources_sigma_out;
705  }
706 
707  /// Destructor.
709 
710  /// Initialize auxiliary vectors and other data members
711  void initialize(ElementCacheMap *element_cache_map) {
712  this->element_cache_map_ = element_cache_map;
713 
715  this->fe_values_->template initialize<dim>(*this->quad_, u);
716  if (dim==1) // print to log only one time
717  DebugOut() << "List of SourcesAssemblyDG FEValues updates flags: " << this->print_update_flags(u);
718  ndofs_ = this->n_dofs();
719  dof_indices_.resize(ndofs_);
720  local_rhs_.resize(ndofs_);
723  }
724 
725 
726  /// Assemble integral over element
727  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
728  {
729  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
730 
731  ElementAccessor<3> elm = cell.elm();
732  double source;
733 
735 
736  // assemble the local stiffness matrix
737  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
738  {
739  fill_n( &(local_rhs_[0]), ndofs_, 0 );
742 
743  for (auto p : this->bulk_points(element_patch_idx) )
744  {
745  source = (eq_fields_->sources_density_out[sbi](p) + eq_fields_->sources_conc_out[sbi](p)*eq_fields_->sources_sigma_out[sbi](p))*JxW_(p);
746 
747  for (unsigned int i=0; i<ndofs_; i++)
748  local_rhs_[i] += source*conc_shape_(i,p);
749  }
750  eq_data_->ls[sbi]->rhs_set_values(ndofs_, &(dof_indices_[0]), &(local_rhs_[0]));
751 
752  for (unsigned int i=0; i<ndofs_; i++)
753  {
754  for (auto p : this->bulk_points(element_patch_idx) )
755  {
756  local_source_balance_vector_[i] -= eq_fields_->sources_sigma_out[sbi](p)*conc_shape_(i,p)*JxW_(p);
757  }
758 
760  }
761  eq_data_->balance_->add_source_values(eq_data_->subst_idx()[sbi], elm.region().bulk_idx(),
762  cell.get_loc_dof_indices(),
764  }
765  }
766 
767  /// Implements @p AssemblyBase::begin.
768  void begin() override
769  {
770  eq_data_->balance_->start_source_assembly( eq_data_->subst_idx() );
771  }
772 
773  /// Implements @p AssemblyBase::end.
774  void end() override
775  {
776  eq_data_->balance_->finish_source_assembly( eq_data_->subst_idx() );
777  }
778 
779 
780  private:
781  /// Data objects shared with TransportDG
784 
785  /// Sub field set contains fields used in calculation.
787 
788  unsigned int ndofs_; ///< Number of dofs
789 
790  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
791  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for set_sources method.
792  vector<PetscScalar> local_source_balance_vector_; ///< Auxiliary vector for set_sources method.
793  vector<PetscScalar> local_source_balance_rhs_; ///< Auxiliary vector for set_sources method.
794 
797 
798  template < template<IntDim...> class DimAssembly>
799  friend class GenericAssembly;
800 
801 };
802 
803 
804 /**
805  * Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions..
806  */
807 template <unsigned int dim, class Model>
809 {
810 public:
813 
814  static constexpr const char * name() { return "BdrConditionAssemblyDG"; }
815 
816  /// Constructor.
817  BdrConditionAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues<3> *fe_values)
818  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
819  JxW_( this->side_values().JxW() ),
820  normal_( this->side_values().normal_vector() ),
821  conc_shape_( this->side_values().scalar_shape() ),
822  conc_grad_( this->side_values().grad_scalar_shape() ) {
824  this->used_fields_ += eq_fields_->advection_coef;
825  this->used_fields_ += eq_fields_->diffusion_coef;
826  this->used_fields_ += eq_fields_->cross_section;
827  this->used_fields_ += eq_fields_->bc_type;
828  this->used_fields_ += eq_fields_->bc_dirichlet_value;
829  this->used_fields_ += eq_fields_->bc_robin_sigma;
830  this->used_fields_ += eq_fields_->bc_flux;
831  }
832 
833  /// Destructor.
835 
836  /// Initialize auxiliary vectors and other data members
837  void initialize(ElementCacheMap *element_cache_map) {
838  this->element_cache_map_ = element_cache_map;
839 
841  this->fe_values_->template initialize<dim>(*this->quad_low_, u);
842  if (dim==1) // print to log only one time
843  DebugOut() << "List of BdrConditionAssemblyDG FEValues updates flags: " << this->print_update_flags(u);
844  ndofs_ = this->n_dofs();
845  dof_indices_.resize(ndofs_);
846  local_rhs_.resize(ndofs_);
848  }
849 
850 
851  /// Implements @p AssemblyBase::boundary_side_integral.
852  inline void boundary_side_integral(DHCellSide cell_side)
853  {
854 
855  ElementAccessor<3> bc_elm = cell_side.cond().element_accessor();
856 
857  const DHCellAccessor &cell = cell_side.cell();
859 
860  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
861  {
862  fill_n(&(local_rhs_[0]), ndofs_, 0);
865 
866  double side_flux = advective_flux(eq_fields_->advection_coef[sbi], this->boundary_points(cell_side), JxW_, normal_);
867  double transport_flux = side_flux/cell_side.measure();
868 
869  auto p_side = *( this->boundary_points(cell_side)).begin();
870  auto p_bdr = p_side.point_bdr(cell_side.cond().element_accessor() );
871  unsigned int bc_type = eq_fields_->bc_type[sbi](p_bdr);
872  if (bc_type == AdvectionDiffusionModel::abc_inflow && side_flux < 0)
873  {
874  for (auto p : this->boundary_points(cell_side) )
875  {
876  auto p_bdr = p.point_bdr(bc_elm);
877  double bc_term = -transport_flux*eq_fields_->bc_dirichlet_value[sbi](p_bdr)*JxW_(p);
878  for (unsigned int i=0; i<ndofs_; i++)
879  local_rhs_[i] += bc_term*conc_shape_(i,p);
880  }
881  for (unsigned int i=0; i<ndofs_; i++)
883  }
884  else if (bc_type == AdvectionDiffusionModel::abc_dirichlet)
885  {
886  double side_flux = advective_flux(eq_fields_->advection_coef[sbi], this->boundary_points(cell_side), JxW_, normal_);
887  double transport_flux = side_flux/cell_side.measure();
888 
889  auto p = *( this->boundary_points(cell_side).begin() );
890  double gamma_l = DG_penalty_boundary(cell_side.side(),
891  diffusion_delta(eq_fields_->diffusion_coef[sbi], this->boundary_points(cell_side), normal_(p)),
892  transport_flux,
893  eq_fields_->dg_penalty[sbi](p_bdr));
894 
895  for (auto p : this->boundary_points(cell_side) )
896  {
897  auto p_bdr = p.point_bdr(bc_elm);
898  double bc_term = gamma_l*eq_fields_->bc_dirichlet_value[sbi](p_bdr)*JxW_(p);
899  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));
900  for (unsigned int i=0; i<ndofs_; i++)
901  local_rhs_[i] += bc_term*conc_shape_(i,p)
902  + arma::dot(bc_grad,conc_grad_(i,p));
903  }
904  for (auto p : this->boundary_points(cell_side) )
905  {
906  for (unsigned int i=0; i<ndofs_; i++)
907  {
908  local_flux_balance_vector_[i] += (arma::dot(eq_fields_->advection_coef[sbi](p), normal_(p))*conc_shape_(i,p)
909  - arma::dot(eq_fields_->diffusion_coef[sbi](p)*conc_grad_(i,p),normal_(p))
910  + gamma_l*conc_shape_(i,p))*JxW_(p);
911  }
912  }
913  if (eq_data_->time_->step().index() > 0)
914  for (unsigned int i=0; i<ndofs_; i++)
916  }
917  else if (bc_type == AdvectionDiffusionModel::abc_total_flux)
918  {
919  for (auto p : this->boundary_points(cell_side) )
920  {
921  auto p_bdr = p.point_bdr(bc_elm);
922  double bc_term = eq_fields_->cross_section(p) * (eq_fields_->bc_robin_sigma[sbi](p_bdr)*eq_fields_->bc_dirichlet_value[sbi](p_bdr) +
923  eq_fields_->bc_flux[sbi](p_bdr)) * JxW_(p);
924  for (unsigned int i=0; i<ndofs_; i++)
925  local_rhs_[i] += bc_term*conc_shape_(i,p);
926  }
927 
928  for (unsigned int i=0; i<ndofs_; i++)
929  {
930  for (auto p : this->boundary_points(cell_side) ) {
931  auto p_bdr = p.point_bdr(bc_elm);
932  local_flux_balance_vector_[i] += eq_fields_->cross_section(p) * eq_fields_->bc_robin_sigma[sbi](p_bdr) *
933  JxW_(p) * conc_shape_(i,p);
934  }
936  }
937  }
939  {
940  for (auto p : this->boundary_points(cell_side) )
941  {
942  auto p_bdr = p.point_bdr(bc_elm);
943  double bc_term = eq_fields_->cross_section(p) * (eq_fields_->bc_robin_sigma[sbi](p_bdr)*eq_fields_->bc_dirichlet_value[sbi](p_bdr) +
944  eq_fields_->bc_flux[sbi](p_bdr)) * JxW_(p);
945  for (unsigned int i=0; i<ndofs_; i++)
946  local_rhs_[i] += bc_term*conc_shape_(i,p);
947  }
948 
949  for (unsigned int i=0; i<ndofs_; i++)
950  {
951  for (auto p : this->boundary_points(cell_side) ) {
952  auto p_bdr = p.point_bdr(bc_elm);
953  local_flux_balance_vector_[i] += eq_fields_->cross_section(p)*(arma::dot(eq_fields_->advection_coef[sbi](p), normal_(p)) +
954  eq_fields_->bc_robin_sigma[sbi](p_bdr))*JxW_(p)*conc_shape_(i,p);
955  }
957  }
958  }
959  else if (bc_type == AdvectionDiffusionModel::abc_inflow && side_flux >= 0)
960  {
961  for (auto p : this->boundary_points(cell_side) )
962  {
963  for (unsigned int i=0; i<ndofs_; i++)
964  local_flux_balance_vector_[i] += arma::dot(eq_fields_->advection_coef[sbi](p), normal_(p))*JxW_(p)*conc_shape_(i,p);
965  }
966  }
967  eq_data_->ls[sbi]->rhs_set_values(ndofs_, &(dof_indices_[0]), &(local_rhs_[0]));
968 
969  eq_data_->balance_->add_flux_values(eq_data_->subst_idx()[sbi], cell_side,
970  cell.get_loc_dof_indices(),
972  }
973  }
974 
975  /// Implements @p AssemblyBase::begin.
976  void begin() override
977  {
978  eq_data_->balance_->start_flux_assembly( eq_data_->subst_idx() );
979  }
980 
981  /// Implements @p AssemblyBase::end.
982  void end() override
983  {
984  eq_data_->balance_->finish_flux_assembly( eq_data_->subst_idx() );
985  }
986 
987 
988  private:
989  /// Data objects shared with TransportDG
992 
993  unsigned int ndofs_; ///< Number of dofs
994 
995  /// Sub field set contains fields used in calculation.
997 
998  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
999  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for set_sources method.
1000  vector<PetscScalar> local_flux_balance_vector_; ///< Auxiliary vector for set_boundary_conditions method.
1001  PetscScalar local_flux_balance_rhs_; ///< Auxiliary variable for set_boundary_conditions method.
1002 
1007 
1008  template < template<IntDim...> class DimAssembly>
1009  friend class GenericAssembly;
1010 
1011 };
1012 
1013 
1014 /**
1015  * Auxiliary container class sets the initial condition.
1016  */
1017 template <unsigned int dim, class Model>
1019 {
1020 public:
1023 
1024  static constexpr const char * name() { return "InitProjectionAssemblyDG"; }
1025 
1026  /// Constructor.
1027  InitProjectionAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues<3> *fe_values)
1028  : AssemblyBasePatch<dim>(fe_values), eq_fields_(eq_fields), eq_data_(eq_data),
1029  JxW_( this->bulk_values().JxW() ),
1030  init_shape_( this->bulk_values().scalar_shape() ) {
1032  this->used_fields_ += eq_fields_->init_condition;
1033  }
1034 
1035  /// Destructor.
1037 
1038  /// Initialize auxiliary vectors and other data members
1039  void initialize(ElementCacheMap *element_cache_map) {
1040  this->element_cache_map_ = element_cache_map;
1041 
1043  this->fe_values_->template initialize<dim>(*this->quad_, u);
1044  // if (dim==1) // print to log only one time
1045  // DebugOut() << "List of InitProjectionAssemblyDG FEValues updates flags: " << this->print_update_flags(u);
1046  ndofs_ = this->n_dofs();
1047  dof_indices_.resize(ndofs_);
1048  local_matrix_.resize(4*ndofs_*ndofs_);
1049  local_rhs_.resize(ndofs_);
1050  }
1051 
1052 
1053  /// Assemble integral over element
1054  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
1055  {
1056  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
1057 
1059 
1060  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
1061  {
1062  for (unsigned int i=0; i<ndofs_; i++)
1063  {
1064  local_rhs_[i] = 0;
1065  for (unsigned int j=0; j<ndofs_; j++)
1066  local_matrix_[i*ndofs_+j] = 0;
1067  }
1068 
1069  for (auto p : this->bulk_points(element_patch_idx) )
1070  {
1071  double rhs_term = eq_fields_->init_condition[sbi](p)*JxW_(p);
1072 
1073  for (unsigned int i=0; i<ndofs_; i++)
1074  {
1075  for (unsigned int j=0; j<ndofs_; j++)
1076  local_matrix_[i*ndofs_+j] += init_shape_(i,p)*init_shape_(j,p)*JxW_(p);
1077 
1078  local_rhs_[i] += init_shape_(i,p)*rhs_term;
1079  }
1080  }
1081  eq_data_->ls[sbi]->set_values(ndofs_, &(dof_indices_[0]), ndofs_, &(dof_indices_[0]), &(local_matrix_[0]), &(local_rhs_[0]));
1082  }
1083  }
1084 
1085 
1086  private:
1087  /// Data objects shared with TransportDG
1090 
1091  /// Sub field set contains fields used in calculation.
1093 
1094  unsigned int ndofs_; ///< Number of dofs
1095  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
1096  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
1097  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for set_sources method.
1098 
1101 
1102  template < template<IntDim...> class DimAssembly>
1103  friend class GenericAssembly;
1104 
1105 };
1106 
1107 
1108 /**
1109  * Auxiliary container class sets the initial condition.
1110  */
1111 template <unsigned int dim, class Model>
1113 {
1114 public:
1117 
1118  static constexpr const char * name() { return "InitConditionAssemblyDG"; }
1119 
1120  /// Constructor.
1122  : AssemblyBase<dim>(), eq_fields_(eq_fields), eq_data_(eq_data) {
1124  this->used_fields_ += eq_fields_->init_condition;
1125 
1126  this->quad_ = new Quadrature(dim, RefElement<dim>::n_nodes);
1127  for(unsigned int i = 0; i<RefElement<dim>::n_nodes; i++)
1128  {
1129  this->quad_->weight(i) = 1.0;
1130  this->quad_->set(i) = RefElement<dim>::node_coords(i);
1131  }
1132  }
1133 
1134  /// Destructor.
1136 
1137  /// Initialize auxiliary vectors and other data members
1138  void initialize(ElementCacheMap *element_cache_map) {
1139  this->element_cache_map_ = element_cache_map;
1140  }
1141 
1142 
1143  /// Assemble integral over element
1144  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
1145  {
1146  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
1147 
1148  unsigned int k;
1149 
1150  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
1151  {
1152  k=0;
1153  for (auto p : this->bulk_points(element_patch_idx) )
1154  {
1155  double val = eq_fields_->init_condition[sbi](p);
1156 
1157  eq_data_->output_vec[sbi].set(cell.get_loc_dof_indices()[k], val);
1158  k++;
1159  }
1160  }
1161  }
1162 
1163 
1164  private:
1165  /// Data objects shared with TransportDG
1168 
1169  /// Sub field set contains fields used in calculation.
1171 
1172  template < template<IntDim...> class DimAssembly>
1173  friend class GenericAssembly;
1174 
1175 };
1176 
1177 #endif /* ASSEMBLY_DG_HH_ */
1178 
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.
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.
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:814
FieldSet used_fields_
Sub field set contains fields used in calculation.
Definition: assembly_dg.hh:996
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:852
BdrConditionAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
Definition: assembly_dg.hh:817
TransportDG< Model >::EqData EqData
Definition: assembly_dg.hh:812
TransportDG< Model >::EqFields EqFields
Definition: assembly_dg.hh:811
void begin() override
Implements AssemblyBase::begin.
Definition: assembly_dg.hh:976
~BdrConditionAssemblyDG()
Destructor.
Definition: assembly_dg.hh:834
PetscScalar local_flux_balance_rhs_
Auxiliary variable for set_boundary_conditions method.
FeQ< Scalar > conc_shape_
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
Definition: assembly_dg.hh:837
vector< LongIdx > dof_indices_
Vector of global DOF indices.
Definition: assembly_dg.hh:998
EqFields * eq_fields_
Data objects shared with TransportDG.
Definition: assembly_dg.hh:990
void end() override
Implements AssemblyBase::end.
Definition: assembly_dg.hh:982
vector< PetscScalar > local_rhs_
Auxiliary vector for set_sources method.
Definition: assembly_dg.hh:999
unsigned int ndofs_
Number of dofs.
Definition: assembly_dg.hh:993
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
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:791
TransportDG< Model >::EqData EqData
Definition: assembly_dg.hh:692
TransportDG< Model >::EqFields EqFields
Definition: assembly_dg.hh:691
ElQ< Scalar > JxW_
Definition: assembly_dg.hh:795
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
Definition: assembly_dg.hh:727
~SourcesAssemblyDG()
Destructor.
Definition: assembly_dg.hh:708
SourcesAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
Definition: assembly_dg.hh:697
vector< PetscScalar > local_source_balance_vector_
Auxiliary vector for set_sources method.
Definition: assembly_dg.hh:792
void end() override
Implements AssemblyBase::end.
Definition: assembly_dg.hh:774
static constexpr const char * name()
Definition: assembly_dg.hh:694
EqFields * eq_fields_
Data objects shared with TransportDG.
Definition: assembly_dg.hh:782
vector< PetscScalar > local_source_balance_rhs_
Auxiliary vector for set_sources method.
Definition: assembly_dg.hh:793
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
Definition: assembly_dg.hh:711
FieldSet used_fields_
Sub field set contains fields used in calculation.
Definition: assembly_dg.hh:786
unsigned int ndofs_
Number of dofs.
Definition: assembly_dg.hh:788
vector< LongIdx > dof_indices_
Vector of global DOF indices.
Definition: assembly_dg.hh:790
void begin() override
Implements AssemblyBase::begin.
Definition: assembly_dg.hh:768
FeQ< Scalar > conc_shape_
Definition: assembly_dg.hh:796
TransportDG< Model >::EqData EqData
Definition: assembly_dg.hh:233
FeQ< Scalar > conc_shape_
Definition: assembly_dg.hh:672
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:675
FieldSet used_fields_
Sub field set contains fields used in calculation.
Definition: assembly_dg.hh:655
unsigned int ndofs_
Number of dofs.
Definition: assembly_dg.hh:657
vector< vector< LongIdx > > side_dof_indices_
Vector of vectors of side DOF indices.
Definition: assembly_dg.hh:661
Range< JoinShapeAccessor< Scalar > > conc_join_shape_
Definition: assembly_dg.hh:676
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:667
vector< double * > waverages
Auxiliary storage for weighted averages of shape functions.
Definition: assembly_dg.hh:666
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:662
FeQ< Vector > conc_grad_
Definition: assembly_dg.hh:674
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:665
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
Definition: assembly_dg.hh:663
~StiffnessAssemblyDG()
Destructor.
Definition: assembly_dg.hh:260
ElQ< Scalar > JxW_
Definition: assembly_dg.hh:669
static constexpr const char * name()
Definition: assembly_dg.hh:235
ElQ< Vector > normal_
Definition: assembly_dg.hh:671
vector< LongIdx > dof_indices_
Vector of global DOF indices.
Definition: assembly_dg.hh:660
ElQ< Scalar > JxW_side_
Definition: assembly_dg.hh:670
FeQ< Scalar > conc_shape_side_
Definition: assembly_dg.hh:673
EqFields * eq_fields_
Data objects shared with TransportDG.
Definition: assembly_dg.hh:651
unsigned int qsize_lower_dim_
Size of quadrature of dim-1.
Definition: assembly_dg.hh:658
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
#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.