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