Flow123d  DF_patch_fevalues-8016b85
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/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 AssemblyBase<dim>
37 {
38 public:
41 
42  static constexpr const char * name() { return "MassAssemblyDG"; }
43 
44  /// Constructor.
45  MassAssemblyDG(EqFields *eq_fields, EqData *eq_data)
46  : AssemblyBase<dim>(eq_data->dg_order), eq_fields_(eq_fields), eq_data_(eq_data),
49  this->used_fields_ += eq_fields_->mass_matrix_coef;
50  this->used_fields_ += eq_fields_->retardation_coef;
51  }
52 
53  /// Destructor.
55 
56  /// Initialize auxiliary vectors and other data members
57  void initialize(ElementCacheMap *element_cache_map) {
58  this->element_cache_map_ = element_cache_map;
59 
60  fe_ = std::make_shared< FE_P_disc<dim> >(eq_data_->dg_order);
62  fe_values_.initialize(*this->quad_, *fe_, 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_ = fe_->n_dofs();
66  dof_indices_.resize(ndofs_);
67  local_matrix_.resize(4*ndofs_*ndofs_);
70 
71  }
72 
73 
74  /// Reinit PatchFEValues objects (all computed elements in one step).
75  void patch_reinit(PatchElementsList patch_elements) override
76  {
77  fe_values_.reinit(patch_elements);
78  }
79 
80 
81  /// Assemble integral over element
82  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
83  {
84  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
85 
86  unsigned int k;
87 
88  fe_values_.get_cell(element_patch_idx);
90 
91  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); ++sbi)
92  {
93  // assemble the local mass matrix
94  for (unsigned int i=0; i<ndofs_; i++)
95  {
96  for (unsigned int j=0; j<ndofs_; j++)
97  {
98  local_matrix_[i*ndofs_+j] = 0;
99  k=0;
100  for (auto p : this->bulk_points(element_patch_idx) )
101  {
102  local_matrix_[i*ndofs_+j] += (eq_fields_->mass_matrix_coef(p)+eq_fields_->retardation_coef[sbi](p)) *
104  k++;
105  }
106  }
107  }
108 
109  for (unsigned int i=0; i<ndofs_; i++)
110  {
113  k=0;
114  for (auto p : this->bulk_points(element_patch_idx) )
115  {
116  local_mass_balance_vector_[i] += eq_fields_->mass_matrix_coef(p)*fe_values_.shape_value(i,k)*fe_values_.JxW(k);
117  local_retardation_balance_vector_[i] -= eq_fields_->retardation_coef[sbi](p)*fe_values_.shape_value(i,k)*fe_values_.JxW(k);
118  k++;
119  }
120  }
121 
122  eq_data_->balance_->add_mass_values(eq_data_->subst_idx()[sbi], cell, cell.get_loc_dof_indices(),
124 
126  VecSetValues(eq_data_->ret_vec[sbi], ndofs_, &(dof_indices_[0]), &(local_retardation_balance_vector_[0]), ADD_VALUES);
127  }
128  }
129 
130  /// Implements @p AssemblyBase::begin.
131  void begin() override
132  {
133  eq_data_->balance_->start_mass_assembly( eq_data_->subst_idx() );
134  }
135 
136  /// Implements @p AssemblyBase::end.
137  void end() override
138  {
139  eq_data_->balance_->finish_mass_assembly( eq_data_->subst_idx() );
140  }
141 
142  private:
143  shared_ptr<FiniteElement<dim>> fe_; ///< Finite element for the solution of the advection-diffusion equation.
144 
145  /// Data objects shared with TransportDG
148 
149  /// Sub field set contains fields used in calculation.
151 
152  unsigned int ndofs_; ///< Number of dofs
153  PatchFEValues<3> fe_values_; ///< FEValues of object (of P disc finite element type)
154 
155  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
156  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
157  vector<PetscScalar> local_retardation_balance_vector_; ///< Auxiliary vector for assemble mass matrix.
159 
160  template < template<IntDim...> class DimAssembly>
161  friend class GenericAssembly;
162 
163 };
164 
165 
166 // return the ratio of longest and shortest edge
168 {
169  double h_max = 0, h_min = numeric_limits<double>::infinity();
170  for (unsigned int i=0; i<e->n_nodes(); i++)
171  for (unsigned int j=i+1; j<e->n_nodes(); j++)
172  {
173  double dist = arma::norm(*e.node(i) - *e.node(j));
174  h_max = max(h_max, dist);
175  h_min = min(h_min, dist);
176  }
177  return h_max/h_min;
178 }
179 
180 /**
181  * @brief Computes average normal diffusivity over a set of points
182  *
183  * @param diff_coef Diffusion tensor.
184  * @param pts Points.
185  * @param nv Normal vector.
186  * @return double
187  */
189 {
190  double delta = 0;
191  unsigned int n = 0;
192  for (auto p : pts )
193  {
194  delta += dot(diff_coef(p)*nv, nv);
195  n++;
196  }
197  return n == 0 ? 0 : (delta/n);
198 }
199 
200 
201 /**
202  * @brief Computes the penalty parameter of the DG method on a given boundary edge.
203  *
204  * Assumption is that the edge consists of only 1 side.
205  * @param side The boundary side.
206  * @param diff_delta Average normal dispersivity K*n.n computed by diffusion_delta()
207  * @param ad_vector Advection vector.
208  * @param alpha Penalty parameter that influences the continuity
209  * of the solution (large value=more continuity).
210  */
212  const double &diff_delta,
213  const double flux,
214  const double alpha)
215 {
216  return 0.5*fabs(flux) + alpha/side.diameter()*diff_delta*elem_anisotropy(side.element());
217 }
218 
219 
220 /**
221  * @brief Computes advective flux.
222  *
223  * @param advection_coef Advection vector.
224  * @param pts Quadrature points.
225  * @param fv FE values object.
226  * @return double
227  */
228 template <class PointType>
230 {
231  double side_flux = 0;
232  unsigned int k=0;
233  for (auto p : pts) {
234  side_flux += arma::dot(advection_coef(p), fv.normal_vector(k))*fv.JxW(k);
235  k++;
236  }
237  return side_flux;
238 }
239 
240 
241 /**
242  * Auxiliary container class for Finite element and related objects of given dimension.
243  */
244 template <unsigned int dim, class Model>
246 {
247 public:
250 
251  static constexpr const char * name() { return "StiffnessAssemblyDG"; }
252 
253  /// Constructor.
254  StiffnessAssemblyDG(EqFields *eq_fields, EqData *eq_data)
255  : AssemblyBase<dim>(eq_data->dg_order), eq_fields_(eq_fields), eq_data_(eq_data) {
257  this->used_fields_ += eq_fields_->advection_coef;
258  this->used_fields_ += eq_fields_->diffusion_coef;
259  this->used_fields_ += eq_fields_->cross_section;
261  this->used_fields_ += eq_fields_->sources_sigma_out;
263  this->used_fields_ += eq_fields_->bc_type;
264  this->used_fields_ += eq_fields_->bc_robin_sigma;
265  }
266 
267  /// Destructor.
269  for (auto a : averages) if (a != nullptr) delete[] a;
270  for (auto a : waverages) if (a != nullptr) delete[] a;
271  for (auto a : jumps) if (a != nullptr) delete[] a;
272  }
273 
274  /// Initialize auxiliary vectors and other data members
275  void initialize(ElementCacheMap *element_cache_map) {
276  this->element_cache_map_ = element_cache_map;
277 
278  fe_ = std::make_shared< FE_P_disc<dim> >(eq_data_->dg_order);
279  fe_low_ = std::make_shared< FE_P_disc<dim-1> >(eq_data_->dg_order);
282  fe_values_.initialize(*this->quad_, *fe_, u);
283  if (dim>1) {
284  fe_values_vb_.initialize(*this->quad_low_, *fe_low_, u);
285  }
286  fe_values_side_.initialize(*this->quad_low_, *fe_, u_side);
287  if (dim==1) { // print to log only one time
288  DebugOut() << "List of StiffnessAssemblyDG FEValues (cell) updates flags: " << this->print_update_flags(u);
289  DebugOut() << "List of StiffnessAssemblyDG FEValues (side) updates flags: " << this->print_update_flags(u_side);
290  }
291  ndofs_ = fe_->n_dofs();
292  qsize_lower_dim_ = this->quad_low_->size();
293  dof_indices_.resize(ndofs_);
294  side_dof_indices_vb_.resize(2*ndofs_);
295  local_matrix_.resize(4*ndofs_*ndofs_);
296 
298  for (unsigned int sid=0; sid<eq_data_->max_edg_sides; sid++)
299  {
301  fe_values_vec_[sid].initialize(*this->quad_low_, *fe_,
303  }
304 
305  // index 0 = element with lower dimension,
306  // index 1 = side of element with higher dimension
307  fv_sb_.resize(2);
308  fv_sb_[0] = &fe_values_vb_;
309  fv_sb_[1] = &fe_values_side_;
310 
312  for (unsigned int s=0; s<eq_data_->max_edg_sides; s++)
313  averages[s] = new double[qsize_lower_dim_*fe_->n_dofs()];
314  waverages.resize(2);
315  jumps.resize(2);
316  for (unsigned int s=0; s<2; s++)
317  {
318  waverages[s] = new double[qsize_lower_dim_*fe_->n_dofs()];
319  jumps[s] = new double[qsize_lower_dim_*fe_->n_dofs()];
320  }
321  }
322 
323 
324  /// Assembles the cell (volume) integral into the stiffness matrix.
325  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
326  {
327  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
328  if (!cell.is_own()) return;
329 
330  ElementAccessor<3> elm = cell.elm();
331 
332  fe_values_.reinit(elm);
334  unsigned int k;
335 
336  // assemble the local stiffness matrix
337  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
338  {
339  for (unsigned int i=0; i<ndofs_; i++)
340  for (unsigned int j=0; j<ndofs_; j++)
341  local_matrix_[i*ndofs_+j] = 0;
342 
343  k=0;
344  for (auto p : this->bulk_points(element_patch_idx) )
345  {
346  for (unsigned int i=0; i<ndofs_; i++)
347  {
348  arma::vec3 Kt_grad_i = eq_fields_->diffusion_coef[sbi](p).t()*fe_values_.shape_grad(i,k);
349  double ad_dot_grad_i = arma::dot(eq_fields_->advection_coef[sbi](p), fe_values_.shape_grad(i,k));
350 
351  for (unsigned int j=0; j<ndofs_; j++)
352  local_matrix_[i*ndofs_+j] += (arma::dot(Kt_grad_i, fe_values_.shape_grad(j,k))
353  -fe_values_.shape_value(j,k)*ad_dot_grad_i
354  +eq_fields_->sources_sigma_out[sbi](p)*fe_values_.shape_value(j,k)*fe_values_.shape_value(i,k))*fe_values_.JxW(k);
355  }
356  k++;
357  }
359  }
360  }
361 
362 
363  /// Assembles the fluxes on the boundary.
364  inline void boundary_side_integral(DHCellSide cell_side)
365  {
366  ASSERT_EQ(cell_side.dim(), dim).error("Dimension of element mismatch!");
367  if (!cell_side.cell().is_own()) return;
368 
369  Side side = cell_side.side();
370  const DHCellAccessor &cell = cell_side.cell();
371 
373  fe_values_side_.reinit(side);
374  unsigned int k;
375  double gamma_l;
376 
377  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
378  {
379  std::fill(local_matrix_.begin(), local_matrix_.end(), 0);
380 
381  double side_flux = advective_flux(eq_fields_->advection_coef[sbi], this->boundary_points(cell_side), fe_values_side_);
382  double transport_flux = side_flux/side.measure();
383 
384  // On Neumann boundaries we have only term from integrating by parts the advective term,
385  // on Dirichlet boundaries we additionally apply the penalty which enforces the prescribed value.
386  auto p_side = *( this->boundary_points(cell_side).begin() );
387  auto p_bdr = p_side.point_bdr( side.cond().element_accessor() );
388  unsigned int bc_type = eq_fields_->bc_type[sbi](p_bdr);
390  {
391  // set up the parameters for DG method
392  gamma_l = DG_penalty_boundary(side,
393  diffusion_delta(eq_fields_->diffusion_coef[sbi], this->boundary_points(cell_side), fe_values_side_.normal_vector(0)),
394  transport_flux,
395  eq_fields_->dg_penalty[sbi](p_side));
396  transport_flux += gamma_l;
397  }
398 
399  // fluxes and penalty
400  k=0;
401  for (auto p : this->boundary_points(cell_side) )
402  {
403  double flux_times_JxW;
405  {
406  //sigma_ corresponds to robin_sigma
407  auto p_bdr = p.point_bdr(side.cond().element_accessor());
408  flux_times_JxW = eq_fields_->cross_section(p)*eq_fields_->bc_robin_sigma[sbi](p_bdr)*fe_values_side_.JxW(k);
409  }
411  {
412  auto p_bdr = p.point_bdr(side.cond().element_accessor());
413  flux_times_JxW = (transport_flux + eq_fields_->cross_section(p)*eq_fields_->bc_robin_sigma[sbi](p_bdr))*fe_values_side_.JxW(k);
414  }
415  else if (bc_type == AdvectionDiffusionModel::abc_inflow && side_flux < 0)
416  flux_times_JxW = 0;
417  else
418  flux_times_JxW = transport_flux*fe_values_side_.JxW(k);
419 
420  for (unsigned int i=0; i<ndofs_; i++)
421  {
422  for (unsigned int j=0; j<ndofs_; j++)
423  {
424  // flux due to advection and penalty
426 
427  // flux due to diffusion (only on dirichlet and inflow boundary)
429  local_matrix_[i*ndofs_+j] -= (arma::dot(eq_fields_->diffusion_coef[sbi](p)*fe_values_side_.shape_grad(j,k),fe_values_side_.normal_vector(k))*fe_values_side_.shape_value(i,k)
431  )*fe_values_side_.JxW(k);
432  }
433  }
434  k++;
435  }
436 
438  }
439  }
440 
441 
442  /// Assembles the fluxes between sides of elements of the same dimension.
444  ASSERT_EQ(edge_side_range.begin()->element().dim(), dim).error("Dimension of element mismatch!");
445 
446  unsigned int k;
447  double gamma_l, omega[2], transport_flux, delta[2], delta_sum;
448  double aniso1, aniso2;
449  double local_alpha=0.0;
450  int sid=0, s1, s2;
451  for( DHCellSide edge_side : edge_side_range )
452  {
453  auto dh_edge_cell = eq_data_->dh_->cell_accessor_from_element( edge_side.elem_idx() );
454  dh_edge_cell.get_dof_indices(side_dof_indices_[sid]);
455  fe_values_vec_[sid].reinit(edge_side.side());
456  ++sid;
457  }
458  arma::vec3 normal_vector = fe_values_vec_[0].normal_vector(0);
459 
460  // fluxes and penalty
461  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
462  {
463  vector<double> fluxes(edge_side_range.begin()->n_edge_sides());
464  double pflux = 0, nflux = 0; // calculate the total in- and out-flux through the edge
465  sid=0;
466  for( DHCellSide edge_side : edge_side_range )
467  {
468  fluxes[sid] = advective_flux(eq_fields_->advection_coef[sbi], this->edge_points(edge_side), fe_values_vec_[sid]) / edge_side.measure();
469  if (fluxes[sid] > 0)
470  pflux += fluxes[sid];
471  else
472  nflux += fluxes[sid];
473  ++sid;
474  }
475 
476  // precompute averages of shape functions over pairs of sides
477  s1=0;
478  for (DHCellSide edge_side : edge_side_range)
479  {
480  (void)edge_side;
481  for (unsigned int k=0; k<qsize_lower_dim_; k++)
482  {
483  for (unsigned int i=0; i<fe_->n_dofs(); i++)
484  averages[s1][k*fe_->n_dofs()+i] = fe_values_vec_[s1].shape_value(i,k)*0.5;
485  }
486  s1++;
487  }
488 
489 
490 
491  s1=0;
492  for( DHCellSide edge_side1 : edge_side_range )
493  {
494  s2=-1; // need increment at begin of loop (see conditionally 'continue' directions)
495  for( DHCellSide edge_side2 : edge_side_range )
496  {
497  s2++;
498  if (s2<=s1) continue;
499  ASSERT(edge_side1.is_valid()).error("Invalid side of edge.");
500 
501  arma::vec3 nv = fe_values_vec_[s1].normal_vector(0);
502 
503  // set up the parameters for DG method
504  // calculate the flux from edge_side1 to edge_side2
505  if (fluxes[s2] > 0 && fluxes[s1] < 0)
506  transport_flux = fluxes[s1]*fabs(fluxes[s2]/pflux);
507  else if (fluxes[s2] < 0 && fluxes[s1] > 0)
508  transport_flux = fluxes[s1]*fabs(fluxes[s2]/nflux);
509  else
510  transport_flux = 0;
511 
512  gamma_l = 0.5*fabs(transport_flux);
513 
514  delta[0] = 0;
515  delta[1] = 0;
516  for (auto p1 : this->edge_points(edge_side1) )
517  {
518  auto p2 = p1.point_on(edge_side2);
519  delta[0] += dot(eq_fields_->diffusion_coef[sbi](p1)*normal_vector,normal_vector);
520  delta[1] += dot(eq_fields_->diffusion_coef[sbi](p2)*normal_vector,normal_vector);
521  local_alpha = max(eq_fields_->dg_penalty[sbi](p1), eq_fields_->dg_penalty[sbi](p2));
522  }
523  delta[0] /= qsize_lower_dim_;
524  delta[1] /= qsize_lower_dim_;
525 
526  delta_sum = delta[0] + delta[1];
527 
528 // if (delta_sum > numeric_limits<double>::epsilon())
529  if (fabs(delta_sum) > 0)
530  {
531  omega[0] = delta[1]/delta_sum;
532  omega[1] = delta[0]/delta_sum;
533  double h = edge_side1.diameter();
534  aniso1 = elem_anisotropy(edge_side1.element());
535  aniso2 = elem_anisotropy(edge_side2.element());
536  gamma_l += local_alpha/h*aniso1*aniso2*(delta[0]*delta[1]/delta_sum);
537  }
538  else
539  for (int i=0; i<2; i++) omega[i] = 0;
540  // end of set up the parameters for DG method
541 
542  int sd[2]; bool is_side_own[2];
543  sd[0] = s1; is_side_own[0] = edge_side1.cell().is_own();
544  sd[1] = s2; is_side_own[1] = edge_side2.cell().is_own();
545 
546  // precompute jumps and weighted averages of shape functions over the pair of sides (s1,s2)
547  k=0;
548  for (auto p1 : this->edge_points(edge_side1) )
549  {
550  auto p2 = p1.point_on(edge_side2);
551  for (unsigned int i=0; i<fe_->n_dofs(); i++)
552  {
553  for (int n=0; n<2; n++)
554  {
555  jumps[n][k*fe_->n_dofs()+i] = (n==0)*fe_values_vec_[s1].shape_value(i,k) - (n==1)*fe_values_vec_[s2].shape_value(i,k);
556  waverages[n][k*fe_->n_dofs()+i] = arma::dot(eq_fields_->diffusion_coef[sbi]( (n==0 ? p1 : p2) )*fe_values_vec_[sd[n]].shape_grad(i,k),nv)*omega[n];
557  }
558  }
559  k++;
560  }
561 
562  // For selected pair of elements:
563  for (int n=0; n<2; n++)
564  {
565  if (!is_side_own[n]) continue;
566 
567  for (int m=0; m<2; m++)
568  {
569  for (unsigned int i=0; i<fe_values_vec_[sd[n]].n_dofs(); i++)
570  for (unsigned int j=0; j<fe_values_vec_[sd[m]].n_dofs(); j++)
571  local_matrix_[i*fe_values_vec_[sd[m]].n_dofs()+j] = 0;
572 
573  for (k=0; k<this->quad_low_->size(); ++k)
574  {
575  for (unsigned int i=0; i<fe_values_vec_[sd[n]].n_dofs(); i++)
576  {
577  for (unsigned int j=0; j<fe_values_vec_[sd[m]].n_dofs(); j++)
578  {
579  int index = i*fe_values_vec_[sd[m]].n_dofs()+j;
580 
581  local_matrix_[index] += (
582  // flux due to transport (applied on interior edges) (average times jump)
583  transport_flux*jumps[n][k*fe_->n_dofs()+i]*averages[sd[m]][k*fe_->n_dofs()+j]
584 
585  // penalty enforcing continuity across edges (applied on interior and Dirichlet edges) (jump times jump)
586  + gamma_l*jumps[n][k*fe_->n_dofs()+i]*jumps[m][k*fe_->n_dofs()+j]
587 
588  // terms due to diffusion
589  - jumps[n][k*fe_->n_dofs()+i]*waverages[m][k*fe_->n_dofs()+j]
590  - eq_data_->dg_variant*waverages[n][k*fe_->n_dofs()+i]*jumps[m][k*fe_->n_dofs()+j]
592  }
593  }
594  }
595  eq_data_->ls[sbi]->mat_set_values(fe_values_vec_[sd[n]].n_dofs(), &(side_dof_indices_[sd[n]][0]), fe_values_vec_[sd[m]].n_dofs(), &(side_dof_indices_[sd[m]][0]), &(local_matrix_[0]));
596  }
597  }
598  }
599  s1++;
600  }
601  }
602  }
603 
604 
605  /// Assembles the fluxes between elements of different dimensions.
606  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
607  if (dim == 1) return;
608  ASSERT_EQ(cell_lower_dim.dim(), dim-1).error("Dimension of element mismatch!");
609 
610  // Note: use data members csection_ and velocity_ for appropriate quantities of lower dim element
611 
612  double comm_flux[2][2];
613  unsigned int n_dofs[2];
614  ElementAccessor<3> elm_lower_dim = cell_lower_dim.elm();
615  unsigned int n_indices = cell_lower_dim.get_dof_indices(dof_indices_);
616  for(unsigned int i=0; i<n_indices; ++i) {
618  }
619  fe_values_vb_.reinit(elm_lower_dim);
620  n_dofs[0] = fv_sb_[0]->n_dofs();
621 
622  DHCellAccessor cell_higher_dim = eq_data_->dh_->cell_accessor_from_element( neighb_side.element().idx() );
623  n_indices = cell_higher_dim.get_dof_indices(dof_indices_);
624  for(unsigned int i=0; i<n_indices; ++i) {
625  side_dof_indices_vb_[i+n_dofs[0]] = dof_indices_[i];
626  }
627  fe_values_side_.reinit(neighb_side.side());
628  n_dofs[1] = fv_sb_[1]->n_dofs();
629 
630  // Testing element if they belong to local partition.
631  bool own_element_id[2];
632  own_element_id[0] = cell_lower_dim.is_own();
633  own_element_id[1] = cell_higher_dim.is_own();
634 
635  unsigned int k;
636  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++) // Optimize: SWAP LOOPS
637  {
638  for (unsigned int i=0; i<n_dofs[0]+n_dofs[1]; i++)
639  for (unsigned int j=0; j<n_dofs[0]+n_dofs[1]; j++)
640  local_matrix_[i*(n_dofs[0]+n_dofs[1])+j] = 0;
641 
642  // set transmission conditions
643  k=0;
644  for (auto p_high : this->coupling_points(neighb_side) )
645  {
646  auto p_low = p_high.lower_dim(cell_lower_dim);
647  // The communication flux has two parts:
648  // - "diffusive" term containing sigma
649  // - "advective" term representing usual upwind
650  //
651  // The calculation differs from the reference manual, since ad_coef and dif_coef have different meaning
652  // than b and A in the manual.
653  // In calculation of sigma there appears one more csection_lower in the denominator.
654  double sigma = eq_fields_->fracture_sigma[sbi](p_low)*arma::dot(eq_fields_->diffusion_coef[sbi](p_low)*fe_values_side_.normal_vector(k),fe_values_side_.normal_vector(k))*
655  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));
656 
657  double transport_flux = arma::dot(eq_fields_->advection_coef[sbi](p_high), fe_values_side_.normal_vector(k));
658 
659  comm_flux[0][0] = (sigma-min(0.,transport_flux))*fv_sb_[0]->JxW(k);
660  comm_flux[0][1] = -(sigma-min(0.,transport_flux))*fv_sb_[0]->JxW(k);
661  comm_flux[1][0] = -(sigma+max(0.,transport_flux))*fv_sb_[0]->JxW(k);
662  comm_flux[1][1] = (sigma+max(0.,transport_flux))*fv_sb_[0]->JxW(k);
663 
664  for (int n=0; n<2; n++)
665  {
666  if (!own_element_id[n]) continue;
667 
668  for (unsigned int i=0; i<n_dofs[n]; i++)
669  for (int m=0; m<2; m++)
670  for (unsigned int j=0; j<n_dofs[m]; j++)
671  local_matrix_[(i+n*n_dofs[0])*(n_dofs[0]+n_dofs[1]) + m*n_dofs[0] + j] +=
672  comm_flux[m][n]*fv_sb_[m]->shape_value(j,k)*fv_sb_[n]->shape_value(i,k) + LocalSystem::almost_zero;
673  }
674  k++;
675  }
676  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]));
677  }
678  }
679 
680 
681 private:
682  shared_ptr<FiniteElement<dim>> fe_; ///< Finite element for the solution of the advection-diffusion equation.
683  shared_ptr<FiniteElement<dim-1>> fe_low_; ///< Finite element for the solution of the advection-diffusion equation (dim-1).
684 
685  /// Data objects shared with TransportDG
688 
689  /// Sub field set contains fields used in calculation.
691 
692  unsigned int ndofs_; ///< Number of dofs
693  unsigned int qsize_lower_dim_; ///< Size of quadrature of dim-1
694  FEValues<3> fe_values_; ///< FEValues of object (of P disc finite element type)
695  FEValues<3> fe_values_vb_; ///< FEValues of dim-1 object (of P disc finite element type)
696  FEValues<3> fe_values_side_; ///< FEValues of object (of P disc finite element type)
697  vector<FEValues<3>> fe_values_vec_; ///< Vector of FEValues of object (of P disc finite element types)
698  vector<FEValues<3>*> fv_sb_; ///< Auxiliary vector, holds FEValues objects for assemble element-side
699 
700  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
701  vector< vector<LongIdx> > side_dof_indices_; ///< Vector of vectors of side DOF indices
702  vector<LongIdx> side_dof_indices_vb_; ///< Vector of side DOF indices (assemble element-side fluxex)
703  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
704 
705  vector<double*> averages; ///< Auxiliary storage for averages of shape functions.
706  vector<double*> waverages; ///< Auxiliary storage for weighted averages of shape functions.
707  vector<double*> jumps; ///< Auxiliary storage for jumps of shape functions.
708 
709  template < template<IntDim...> class DimAssembly>
710  friend class GenericAssembly;
711 
712 };
713 
714 
715 /**
716  * Auxiliary container class for Finite element and related objects of given dimension.
717  */
718 template <unsigned int dim, class Model>
719 class SourcesAssemblyDG : public AssemblyBase<dim>
720 {
721 public:
724 
725  static constexpr const char * name() { return "SourcesAssemblyDG"; }
726 
727  /// Constructor.
728  SourcesAssemblyDG(EqFields *eq_fields, EqData *eq_data)
729  : AssemblyBase<dim>(eq_data->dg_order), eq_fields_(eq_fields), eq_data_(eq_data) {
731  this->used_fields_ += eq_fields_->sources_density_out;
732  this->used_fields_ += eq_fields_->sources_conc_out;
733  this->used_fields_ += eq_fields_->sources_sigma_out;
734  }
735 
736  /// Destructor.
738 
739  /// Initialize auxiliary vectors and other data members
740  void initialize(ElementCacheMap *element_cache_map) {
741  this->element_cache_map_ = element_cache_map;
742 
743  fe_ = std::make_shared< FE_P_disc<dim> >(eq_data_->dg_order);
745  fe_values_.initialize(*this->quad_, *fe_, u);
746  if (dim==1) // print to log only one time
747  DebugOut() << "List of SourcesAssemblyDG FEValues updates flags: " << this->print_update_flags(u);
748  ndofs_ = fe_->n_dofs();
749  dof_indices_.resize(ndofs_);
750  local_rhs_.resize(ndofs_);
753  }
754 
755 
756  /// Assemble integral over element
757  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
758  {
759  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
760 
761  ElementAccessor<3> elm = cell.elm();
762  unsigned int k;
763  double source;
764 
765  fe_values_.reinit(elm);
767 
768  // assemble the local stiffness matrix
769  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
770  {
771  fill_n( &(local_rhs_[0]), ndofs_, 0 );
774 
775  k=0;
776  for (auto p : this->bulk_points(element_patch_idx) )
777  {
778  source = (eq_fields_->sources_density_out[sbi](p) + eq_fields_->sources_conc_out[sbi](p)*eq_fields_->sources_sigma_out[sbi](p))*fe_values_.JxW(k);
779 
780  for (unsigned int i=0; i<ndofs_; i++)
781  local_rhs_[i] += source*fe_values_.shape_value(i,k);
782  k++;
783  }
784  eq_data_->ls[sbi]->rhs_set_values(ndofs_, &(dof_indices_[0]), &(local_rhs_[0]));
785 
786  for (unsigned int i=0; i<ndofs_; i++)
787  {
788  k=0;
789  for (auto p : this->bulk_points(element_patch_idx) )
790  {
791  local_source_balance_vector_[i] -= eq_fields_->sources_sigma_out[sbi](p)*fe_values_.shape_value(i,k)*fe_values_.JxW(k);
792  k++;
793  }
794 
796  }
797  eq_data_->balance_->add_source_values(eq_data_->subst_idx()[sbi], elm.region().bulk_idx(),
798  cell.get_loc_dof_indices(),
800  }
801  }
802 
803  /// Implements @p AssemblyBase::begin.
804  void begin() override
805  {
806  eq_data_->balance_->start_source_assembly( eq_data_->subst_idx() );
807  }
808 
809  /// Implements @p AssemblyBase::end.
810  void end() override
811  {
812  eq_data_->balance_->finish_source_assembly( eq_data_->subst_idx() );
813  }
814 
815 
816  private:
817  shared_ptr<FiniteElement<dim>> fe_; ///< Finite element for the solution of the advection-diffusion equation.
818 
819  /// Data objects shared with TransportDG
822 
823  /// Sub field set contains fields used in calculation.
825 
826  unsigned int ndofs_; ///< Number of dofs
827  FEValues<3> fe_values_; ///< FEValues of object (of P disc finite element type)
828 
829  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
830  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for set_sources method.
831  vector<PetscScalar> local_source_balance_vector_; ///< Auxiliary vector for set_sources method.
832  vector<PetscScalar> local_source_balance_rhs_; ///< Auxiliary vector for set_sources method.
833 
834  template < template<IntDim...> class DimAssembly>
835  friend class GenericAssembly;
836 
837 };
838 
839 
840 /**
841  * Assembles the r.h.s. components corresponding to the Dirichlet boundary conditions..
842  */
843 template <unsigned int dim, class Model>
845 {
846 public:
849 
850  static constexpr const char * name() { return "BdrConditionAssemblyDG"; }
851 
852  /// Constructor.
853  BdrConditionAssemblyDG(EqFields *eq_fields, EqData *eq_data)
854  : AssemblyBase<dim>(eq_data->dg_order), eq_fields_(eq_fields), eq_data_(eq_data) {
856  this->used_fields_ += eq_fields_->advection_coef;
857  this->used_fields_ += eq_fields_->diffusion_coef;
858  this->used_fields_ += eq_fields_->cross_section;
859  this->used_fields_ += eq_fields_->bc_type;
860  this->used_fields_ += eq_fields_->bc_dirichlet_value;
861  this->used_fields_ += eq_fields_->bc_robin_sigma;
862  this->used_fields_ += eq_fields_->bc_flux;
863  }
864 
865  /// Destructor.
867 
868  /// Initialize auxiliary vectors and other data members
869  void initialize(ElementCacheMap *element_cache_map) {
870  this->element_cache_map_ = element_cache_map;
871 
872  fe_ = std::make_shared< FE_P_disc<dim> >(eq_data_->dg_order);
875  if (dim==1) // print to log only one time
876  DebugOut() << "List of BdrConditionAssemblyDG FEValues updates flags: " << this->print_update_flags(u);
877  ndofs_ = fe_->n_dofs();
878  dof_indices_.resize(ndofs_);
879  local_rhs_.resize(ndofs_);
881  }
882 
883 
884  /// Implements @p AssemblyBase::boundary_side_integral.
885  inline void boundary_side_integral(DHCellSide cell_side)
886  {
887  unsigned int k;
888 
889  ElementAccessor<3> bc_elm = cell_side.cond().element_accessor();
890 
891  fe_values_side_.reinit(cell_side.side());
892 
893  const DHCellAccessor &cell = cell_side.cell();
895 
896  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
897  {
898  fill_n(&(local_rhs_[0]), ndofs_, 0);
901 
902  double side_flux = advective_flux(eq_fields_->advection_coef[sbi], this->boundary_points(cell_side), fe_values_side_);
903  double transport_flux = side_flux/cell_side.measure();
904 
905  auto p_side = *( this->boundary_points(cell_side)).begin();
906  auto p_bdr = p_side.point_bdr(cell_side.cond().element_accessor() );
907  unsigned int bc_type = eq_fields_->bc_type[sbi](p_bdr);
908  if (bc_type == AdvectionDiffusionModel::abc_inflow && side_flux < 0)
909  {
910  k=0;
911  for (auto p : this->boundary_points(cell_side) )
912  {
913  auto p_bdr = p.point_bdr(bc_elm);
914  double bc_term = -transport_flux*eq_fields_->bc_dirichlet_value[sbi](p_bdr)*fe_values_side_.JxW(k);
915  for (unsigned int i=0; i<ndofs_; i++)
916  local_rhs_[i] += bc_term*fe_values_side_.shape_value(i,k);
917  k++;
918  }
919  for (unsigned int i=0; i<ndofs_; i++)
921  }
922  else if (bc_type == AdvectionDiffusionModel::abc_dirichlet)
923  {
924  double side_flux = advective_flux(eq_fields_->advection_coef[sbi], this->boundary_points(cell_side), fe_values_side_);
925  double transport_flux = side_flux/cell_side.measure();
926 
927  double gamma_l = DG_penalty_boundary(cell_side.side(),
928  diffusion_delta(eq_fields_->diffusion_coef[sbi], this->boundary_points(cell_side), fe_values_side_.normal_vector(0)),
929  transport_flux,
930  eq_fields_->dg_penalty[sbi](p_bdr));
931 
932  k=0;
933  for (auto p : this->boundary_points(cell_side) )
934  {
935  auto p_bdr = p.point_bdr(bc_elm);
936  double bc_term = gamma_l*eq_fields_->bc_dirichlet_value[sbi](p_bdr)*fe_values_side_.JxW(k);
937  arma::vec3 bc_grad = -eq_fields_->bc_dirichlet_value[sbi](p_bdr)*fe_values_side_.JxW(k)*eq_data_->dg_variant*(arma::trans(eq_fields_->diffusion_coef[sbi](p))*fe_values_side_.normal_vector(k));
938  for (unsigned int i=0; i<ndofs_; i++)
939  local_rhs_[i] += bc_term*fe_values_side_.shape_value(i,k)
940  + arma::dot(bc_grad,fe_values_side_.shape_grad(i,k));
941  k++;
942  }
943  k=0;
944  for (auto p : this->boundary_points(cell_side) )
945  {
946  for (unsigned int i=0; i<ndofs_; i++)
947  {
948  local_flux_balance_vector_[i] += (arma::dot(eq_fields_->advection_coef[sbi](p), fe_values_side_.normal_vector(k))*fe_values_side_.shape_value(i,k)
949  - arma::dot(eq_fields_->diffusion_coef[sbi](p)*fe_values_side_.shape_grad(i,k),fe_values_side_.normal_vector(k))
950  + gamma_l*fe_values_side_.shape_value(i,k))*fe_values_side_.JxW(k);
951  }
952  k++;
953  }
954  if (eq_data_->time_->step().index() > 0)
955  for (unsigned int i=0; i<ndofs_; i++)
957  }
958  else if (bc_type == AdvectionDiffusionModel::abc_total_flux)
959  {
960  k=0;
961  for (auto p : this->boundary_points(cell_side) )
962  {
963  auto p_bdr = p.point_bdr(bc_elm);
964  double bc_term = eq_fields_->cross_section(p) * (eq_fields_->bc_robin_sigma[sbi](p_bdr)*eq_fields_->bc_dirichlet_value[sbi](p_bdr) +
965  eq_fields_->bc_flux[sbi](p_bdr)) * fe_values_side_.JxW(k);
966  for (unsigned int i=0; i<ndofs_; i++)
967  local_rhs_[i] += bc_term*fe_values_side_.shape_value(i,k);
968  k++;
969  }
970 
971  for (unsigned int i=0; i<ndofs_; i++)
972  {
973  k=0;
974  for (auto p : this->boundary_points(cell_side) ) {
975  auto p_bdr = p.point_bdr(bc_elm);
976  local_flux_balance_vector_[i] += eq_fields_->cross_section(p) * eq_fields_->bc_robin_sigma[sbi](p_bdr) *
978  k++;
979  }
981  }
982  }
984  {
985  k=0;
986  for (auto p : this->boundary_points(cell_side) )
987  {
988  auto p_bdr = p.point_bdr(bc_elm);
989  double bc_term = eq_fields_->cross_section(p) * (eq_fields_->bc_robin_sigma[sbi](p_bdr)*eq_fields_->bc_dirichlet_value[sbi](p_bdr) +
990  eq_fields_->bc_flux[sbi](p_bdr)) * fe_values_side_.JxW(k);
991  for (unsigned int i=0; i<ndofs_; i++)
992  local_rhs_[i] += bc_term*fe_values_side_.shape_value(i,k);
993  k++;
994  }
995 
996  for (unsigned int i=0; i<ndofs_; i++)
997  {
998  k=0;
999  for (auto p : this->boundary_points(cell_side) ) {
1000  auto p_bdr = p.point_bdr(bc_elm);
1001  local_flux_balance_vector_[i] += eq_fields_->cross_section(p)*(arma::dot(eq_fields_->advection_coef[sbi](p), fe_values_side_.normal_vector(k)) +
1002  eq_fields_->bc_robin_sigma[sbi](p_bdr))*fe_values_side_.JxW(k)*fe_values_side_.shape_value(i,k);
1003  k++;
1004  }
1006  }
1007  }
1008  else if (bc_type == AdvectionDiffusionModel::abc_inflow && side_flux >= 0)
1009  {
1010  k=0;
1011  for (auto p : this->boundary_points(cell_side) )
1012  {
1013  for (unsigned int i=0; i<ndofs_; i++)
1015  k++;
1016  }
1017  }
1018  eq_data_->ls[sbi]->rhs_set_values(ndofs_, &(dof_indices_[0]), &(local_rhs_[0]));
1019 
1020  eq_data_->balance_->add_flux_values(eq_data_->subst_idx()[sbi], cell_side,
1021  cell.get_loc_dof_indices(),
1023  }
1024  }
1025 
1026  /// Implements @p AssemblyBase::begin.
1027  void begin() override
1028  {
1029  eq_data_->balance_->start_flux_assembly( eq_data_->subst_idx() );
1030  }
1031 
1032  /// Implements @p AssemblyBase::end.
1033  void end() override
1034  {
1035  eq_data_->balance_->finish_flux_assembly( eq_data_->subst_idx() );
1036  }
1037 
1038 
1039  private:
1040  shared_ptr<FiniteElement<dim>> fe_; ///< Finite element for the solution of the advection-diffusion equation.
1041 
1042  /// Data objects shared with TransportDG
1045 
1046  /// Sub field set contains fields used in calculation.
1048 
1049  unsigned int ndofs_; ///< Number of dofs
1050  FEValues<3> fe_values_side_; ///< FEValues of object (of P disc finite element type)
1051 
1052  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
1053  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for set_sources method.
1054  vector<PetscScalar> local_flux_balance_vector_; ///< Auxiliary vector for set_boundary_conditions method.
1055  PetscScalar local_flux_balance_rhs_; ///< Auxiliary variable for set_boundary_conditions method.
1056 
1057  template < template<IntDim...> class DimAssembly>
1058  friend class GenericAssembly;
1059 
1060 };
1061 
1062 
1063 /**
1064  * Auxiliary container class sets the initial condition.
1065  */
1066 template <unsigned int dim, class Model>
1068 {
1069 public:
1072 
1073  static constexpr const char * name() { return "InitProjectionAssemblyDG"; }
1074 
1075  /// Constructor.
1077  : AssemblyBase<dim>(eq_data->dg_order), eq_fields_(eq_fields), eq_data_(eq_data) {
1079  this->used_fields_ += eq_fields_->init_condition;
1080  }
1081 
1082  /// Destructor.
1084 
1085  /// Initialize auxiliary vectors and other data members
1086  void initialize(ElementCacheMap *element_cache_map) {
1087  this->element_cache_map_ = element_cache_map;
1088 
1089  fe_ = std::make_shared< FE_P_disc<dim> >(eq_data_->dg_order);
1091  fe_values_.initialize(*this->quad_, *fe_, u);
1092  // if (dim==1) // print to log only one time
1093  // DebugOut() << "List of InitProjectionAssemblyDG FEValues updates flags: " << this->print_update_flags(u);
1094  ndofs_ = fe_->n_dofs();
1095  dof_indices_.resize(ndofs_);
1096  local_matrix_.resize(4*ndofs_*ndofs_);
1097  local_rhs_.resize(ndofs_);
1098  }
1099 
1100 
1101  /// Assemble integral over element
1102  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
1103  {
1104  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
1105 
1106  unsigned int k;
1107  ElementAccessor<3> elem = cell.elm();
1109  fe_values_.reinit(elem);
1110 
1111  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
1112  {
1113  for (unsigned int i=0; i<ndofs_; i++)
1114  {
1115  local_rhs_[i] = 0;
1116  for (unsigned int j=0; j<ndofs_; j++)
1117  local_matrix_[i*ndofs_+j] = 0;
1118  }
1119 
1120  k=0;
1121  for (auto p : this->bulk_points(element_patch_idx) )
1122  {
1123  double rhs_term = eq_fields_->init_condition[sbi](p)*fe_values_.JxW(k);
1124 
1125  for (unsigned int i=0; i<ndofs_; i++)
1126  {
1127  for (unsigned int j=0; j<ndofs_; j++)
1129 
1130  local_rhs_[i] += fe_values_.shape_value(i,k)*rhs_term;
1131  }
1132  k++;
1133  }
1134  eq_data_->ls[sbi]->set_values(ndofs_, &(dof_indices_[0]), ndofs_, &(dof_indices_[0]), &(local_matrix_[0]), &(local_rhs_[0]));
1135  }
1136  }
1137 
1138 
1139  private:
1140  shared_ptr<FiniteElement<dim>> fe_; ///< Finite element for the solution of the advection-diffusion equation.
1141 
1142  /// Data objects shared with TransportDG
1145 
1146  /// Sub field set contains fields used in calculation.
1148 
1149  unsigned int ndofs_; ///< Number of dofs
1150  FEValues<3> fe_values_; ///< FEValues of object (of P disc finite element type)
1151 
1152  vector<LongIdx> dof_indices_; ///< Vector of global DOF indices
1153  vector<PetscScalar> local_matrix_; ///< Auxiliary vector for assemble methods
1154  vector<PetscScalar> local_rhs_; ///< Auxiliary vector for set_sources method.
1155 
1156  template < template<IntDim...> class DimAssembly>
1157  friend class GenericAssembly;
1158 
1159 };
1160 
1161 
1162 /**
1163  * Auxiliary container class sets the initial condition.
1164  */
1165 template <unsigned int dim, class Model>
1167 {
1168 public:
1171 
1172  static constexpr const char * name() { return "InitConditionAssemblyDG"; }
1173 
1174  /// Constructor.
1176  : AssemblyBase<dim>(), eq_fields_(eq_fields), eq_data_(eq_data) {
1178  this->used_fields_ += eq_fields_->init_condition;
1179 
1180  this->quad_ = new Quadrature(dim, RefElement<dim>::n_nodes);
1181  for(unsigned int i = 0; i<RefElement<dim>::n_nodes; i++)
1182  {
1183  this->quad_->weight(i) = 1.0;
1184  this->quad_->set(i) = RefElement<dim>::node_coords(i);
1185  }
1186  }
1187 
1188  /// Destructor.
1190 
1191  /// Initialize auxiliary vectors and other data members
1192  void initialize(ElementCacheMap *element_cache_map) {
1193  this->element_cache_map_ = element_cache_map;
1194  }
1195 
1196 
1197  /// Assemble integral over element
1198  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
1199  {
1200  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
1201 
1202  unsigned int k;
1203 
1204  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++)
1205  {
1206  k=0;
1207  for (auto p : this->bulk_points(element_patch_idx) )
1208  {
1209  double val = eq_fields_->init_condition[sbi](p);
1210 
1211  eq_data_->output_vec[sbi].set(cell.get_loc_dof_indices()[k], val);
1212  k++;
1213  }
1214  }
1215  }
1216 
1217 
1218  private:
1219  /// Data objects shared with TransportDG
1222 
1223  /// Sub field set contains fields used in calculation.
1225 
1226  template < template<IntDim...> class DimAssembly>
1227  friend class GenericAssembly;
1228 
1229 };
1230 
1231 #endif /* ASSEMBLY_DG_HH_ */
1232 
double elem_anisotropy(ElementAccessor< 3 > e)
Definition: assembly_dg.hh:167
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:211
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:188
double advective_flux(Field< 3, FieldValue< 3 >::VectorFixed > &advection_coef, Range< PointType > pts, FEValues< 3 > &fv)
Computes advective flux.
Definition: assembly_dg.hh:229
#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
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:850
FieldSet used_fields_
Sub field set contains fields used in calculation.
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:885
FEValues< 3 > fe_values_side_
FEValues of object (of P disc finite element type)
TransportDG< Model >::EqData EqData
Definition: assembly_dg.hh:848
TransportDG< Model >::EqFields EqFields
Definition: assembly_dg.hh:847
void begin() override
Implements AssemblyBase::begin.
~BdrConditionAssemblyDG()
Destructor.
Definition: assembly_dg.hh:866
PetscScalar local_flux_balance_rhs_
Auxiliary variable for set_boundary_conditions method.
shared_ptr< FiniteElement< dim > > fe_
Finite element for the solution of the advection-diffusion equation.
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
Definition: assembly_dg.hh:869
vector< LongIdx > dof_indices_
Vector of global DOF indices.
EqFields * eq_fields_
Data objects shared with TransportDG.
void end() override
Implements AssemblyBase::end.
vector< PetscScalar > local_rhs_
Auxiliary vector for set_sources method.
unsigned int ndofs_
Number of dofs.
BdrConditionAssemblyDG(EqFields *eq_fields, EqData *eq_data)
Constructor.
Definition: assembly_dg.hh:853
ElementAccessor< 3 > element_accessor()
Auxiliary data class holds number of elements in cache and allow to set this value explicitly (e....
Cell accessor allow iterate over DOF handler cells.
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
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.
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
void initialize(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags _flags)
Initialize structures and calculates cell-independent data.
Definition: fe_values.cc:132
double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
Definition: fe_values.hh:421
double shape_value(const unsigned int function_no, const unsigned int point_no) const
Return the value of the function_no-th shape function at the point_no-th quadrature point.
Definition: fe_values.hh:331
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:563
arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
Definition: fe_values.hh:455
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const unsigned int point_no) const
Return the gradient of the function_no-th shape function at the point_no-th quadrature point.
Definition: fe_values.hh:346
Discontinuous Lagrangean finite element on dim dimensional simplex.
Definition: fe_p.hh:108
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
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
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.
shared_ptr< FiniteElement< dim > > fe_
Finite element for the solution of the advection-diffusion equation.
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
FEValues< 3 > fe_values_
FEValues of object (of P disc finite element type)
static constexpr const char * name()
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
InitProjectionAssemblyDG(EqFields *eq_fields, EqData *eq_data)
Constructor.
~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:155
vector< PetscScalar > local_mass_balance_vector_
Same as previous.
Definition: assembly_dg.hh:158
PatchFEValues< 3 > fe_values_
FEValues of object (of P disc finite element type)
Definition: assembly_dg.hh:153
TransportDG< Model >::EqFields EqFields
Definition: assembly_dg.hh:39
void patch_reinit(PatchElementsList patch_elements) override
Reinit PatchFEValues objects (all computed elements in one step).
Definition: assembly_dg.hh:75
FieldSet used_fields_
Sub field set contains fields used in calculation.
Definition: assembly_dg.hh:150
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
Definition: assembly_dg.hh:156
void end() override
Implements AssemblyBase::end.
Definition: assembly_dg.hh:137
~MassAssemblyDG()
Destructor.
Definition: assembly_dg.hh:54
shared_ptr< FiniteElement< dim > > fe_
Finite element for the solution of the advection-diffusion equation.
Definition: assembly_dg.hh:143
EqData * eq_data_
Definition: assembly_dg.hh:147
void begin() override
Implements AssemblyBase::begin.
Definition: assembly_dg.hh:131
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:57
unsigned int ndofs_
Number of dofs.
Definition: assembly_dg.hh:152
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
Definition: assembly_dg.hh:82
TransportDG< Model >::EqData EqData
Definition: assembly_dg.hh:40
EqFields * eq_fields_
Data objects shared with TransportDG.
Definition: assembly_dg.hh:146
vector< PetscScalar > local_retardation_balance_vector_
Auxiliary vector for assemble mass matrix.
Definition: assembly_dg.hh:157
MassAssemblyDG(EqFields *eq_fields, EqData *eq_data)
Constructor.
Definition: assembly_dg.hh:45
void reinit(PatchElementsList patch_elements)
Reinit data.
double shape_value(const unsigned int function_no, const unsigned int point_no) const
Return the value of the function_no-th shape function at the point_no-th quadrature point.
Definition: fe_values.hh:531
void get_cell(const unsigned int patch_cell_idx)
Set element that is selected for processing. Element is given by index on patch.
Definition: fe_values.hh:515
double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
Definition: fe_values.hh:622
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:619
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:830
TransportDG< Model >::EqData EqData
Definition: assembly_dg.hh:723
TransportDG< Model >::EqFields EqFields
Definition: assembly_dg.hh:722
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
Definition: assembly_dg.hh:757
~SourcesAssemblyDG()
Destructor.
Definition: assembly_dg.hh:737
shared_ptr< FiniteElement< dim > > fe_
Finite element for the solution of the advection-diffusion equation.
Definition: assembly_dg.hh:817
vector< PetscScalar > local_source_balance_vector_
Auxiliary vector for set_sources method.
Definition: assembly_dg.hh:831
void end() override
Implements AssemblyBase::end.
Definition: assembly_dg.hh:810
FEValues< 3 > fe_values_
FEValues of object (of P disc finite element type)
Definition: assembly_dg.hh:827
static constexpr const char * name()
Definition: assembly_dg.hh:725
EqFields * eq_fields_
Data objects shared with TransportDG.
Definition: assembly_dg.hh:820
vector< PetscScalar > local_source_balance_rhs_
Auxiliary vector for set_sources method.
Definition: assembly_dg.hh:832
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
Definition: assembly_dg.hh:740
FieldSet used_fields_
Sub field set contains fields used in calculation.
Definition: assembly_dg.hh:824
unsigned int ndofs_
Number of dofs.
Definition: assembly_dg.hh:826
vector< LongIdx > dof_indices_
Vector of global DOF indices.
Definition: assembly_dg.hh:829
void begin() override
Implements AssemblyBase::begin.
Definition: assembly_dg.hh:804
SourcesAssemblyDG(EqFields *eq_fields, EqData *eq_data)
Constructor.
Definition: assembly_dg.hh:728
TransportDG< Model >::EqData EqData
Definition: assembly_dg.hh:249
void edge_integral(RangeConvert< DHEdgeSide, DHCellSide > edge_side_range)
Assembles the fluxes between sides of elements of the same dimension.
Definition: assembly_dg.hh:443
StiffnessAssemblyDG(EqFields *eq_fields, EqData *eq_data)
Constructor.
Definition: assembly_dg.hh:254
FieldSet used_fields_
Sub field set contains fields used in calculation.
Definition: assembly_dg.hh:690
unsigned int ndofs_
Number of dofs.
Definition: assembly_dg.hh:692
vector< vector< LongIdx > > side_dof_indices_
Vector of vectors of side DOF indices.
Definition: assembly_dg.hh:701
shared_ptr< FiniteElement< dim-1 > > fe_low_
Finite element for the solution of the advection-diffusion equation (dim-1).
Definition: assembly_dg.hh:683
shared_ptr< FiniteElement< dim > > fe_
Finite element for the solution of the advection-diffusion equation.
Definition: assembly_dg.hh:682
void boundary_side_integral(DHCellSide cell_side)
Assembles the fluxes on the boundary.
Definition: assembly_dg.hh:364
vector< double * > jumps
Auxiliary storage for jumps of shape functions.
Definition: assembly_dg.hh:707
vector< double * > waverages
Auxiliary storage for weighted averages of shape functions.
Definition: assembly_dg.hh:706
vector< FEValues< 3 > * > fv_sb_
Auxiliary vector, holds FEValues objects for assemble element-side.
Definition: assembly_dg.hh:698
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
Assembles the fluxes between elements of different dimensions.
Definition: assembly_dg.hh:606
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
Definition: assembly_dg.hh:275
vector< LongIdx > side_dof_indices_vb_
Vector of side DOF indices (assemble element-side fluxex)
Definition: assembly_dg.hh:702
FEValues< 3 > fe_values_side_
FEValues of object (of P disc finite element type)
Definition: assembly_dg.hh:696
vector< FEValues< 3 > > fe_values_vec_
Vector of FEValues of object (of P disc finite element types)
Definition: assembly_dg.hh:697
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assembles the cell (volume) integral into the stiffness matrix.
Definition: assembly_dg.hh:325
TransportDG< Model >::EqFields EqFields
Definition: assembly_dg.hh:248
vector< double * > averages
Auxiliary storage for averages of shape functions.
Definition: assembly_dg.hh:705
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
Definition: assembly_dg.hh:703
~StiffnessAssemblyDG()
Destructor.
Definition: assembly_dg.hh:268
static constexpr const char * name()
Definition: assembly_dg.hh:251
vector< LongIdx > dof_indices_
Vector of global DOF indices.
Definition: assembly_dg.hh:700
FEValues< 3 > fe_values_vb_
FEValues of dim-1 object (of P disc finite element type)
Definition: assembly_dg.hh:695
FEValues< 3 > fe_values_
FEValues of object (of P disc finite element type)
Definition: assembly_dg.hh:694
EqFields * eq_fields_
Data objects shared with TransportDG.
Definition: assembly_dg.hh:686
unsigned int qsize_lower_dim_
Size of quadrature of dim-1.
Definition: assembly_dg.hh:693
const TimeStep & step(int index=-1) const
unsigned int index() const
std::shared_ptr< DOFHandlerMultiDim > dh_
Object for distribution of dofs.
std::vector< VectorMPI > output_vec
Vector of solution data.
LinSys ** ls_dt
Linear algebra system for the time derivative (actually it is used only for handling the matrix struc...
int dg_variant
DG variant ((non-)symmetric/incomplete.
TimeGovernor * time_
std::shared_ptr< Balance > balance_
unsigned int dg_order
Polynomial order of finite elements.
LinSys ** ls
Linear algebra system for the transport equation.
std::vector< Vec > ret_vec
Auxiliary vectors for calculation of sources in balance due to retardation (e.g. sorption).
unsigned int max_edg_sides
Maximal number of edge sides (evaluate from dim 1,2,3)
MultiField< 3, FieldValue< 3 >::Scalar > dg_penalty
Penalty enforcing inter-element continuity of solution (for each substance).
MultiField< 3, FieldValue< 3 >::Scalar > fracture_sigma
Transition parameter for diffusive transfer on fractures (for each substance).
Definitions of basic Lagrangean finite elements with polynomial shape functions.
Class FEValues calculates finite element data on the actual cells such as shape function values,...
@ coupling
@ boundary
@ bulk
@ edge
#define DebugOut()
Macro defining 'debug' record of log.
Definition: logger.hh:284
unsigned int IntDim
A dimension index type.
Definition: mixed.hh:19
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.