Flow123d  JS_constraints-e651b99
assembly_convection.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_convection.hh
15  * @brief
16  */
17 
18 #ifndef ASSEMBLY_CONVECTION_HH_
19 #define ASSEMBLY_CONVECTION_HH_
20 
23 #include "transport/transport.h"
24 #include "fem/patch_fe_values.hh"
25 #include "fem/patch_op_impl.hh"
27 #include "coupling/balance.hh"
28 #include "fem/element_cache_map.hh"
29 
30 
31 /**
32  * Auxiliary container class for Finite element and related objects of given dimension.
33  */
34 template <unsigned int dim, class TEqData>
36 {
37 public:
38  typedef typename TEqData::EqFields EqFields;
39  typedef TEqData EqData;
40 
41  static constexpr const char * name() { return "Convection_Mass_Assembly"; }
42 
43  /// Constructor.
45  : AssemblyBasePatch<dim>(0, asm_internals), eq_fields_(eq_data->eq_fields_.get()), eq_data_(eq_data),
46  mass_integral_( this->create_bulk_integral(this->quad_) ) {
47  this->used_fields_ += eq_fields_->cross_section;
48  this->used_fields_ += eq_fields_->water_content;
49  }
50 
51  /// Destructor.
53 
54  /// Initialize auxiliary vectors and other data members
55  void initialize() {}
56 
57 
58  /// Assemble integral over element
59  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
60  {
61  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
62 
63  ElementAccessor<3> elm = cell.elm();
64  // we have currently zero order P_Disc FE
65  ASSERT(cell.get_loc_dof_indices().size() == 1);
66  IntIdx local_p0_dof = cell.get_loc_dof_indices()[0];
67 
68  auto p = *( mass_integral_->points(element_patch_idx).begin() );
69  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); ++sbi)
70  eq_data_->balance_->add_mass_values(eq_data_->subst_idx[sbi], cell, {local_p0_dof},
71  {eq_fields_->cross_section(p)*eq_fields_->water_content(p)*elm.measure()}, 0);
72 
73  VecSetValue(eq_data_->mass_diag, eq_data_->dh_->get_local_to_global_map()[local_p0_dof],
74  eq_fields_->cross_section(p)*eq_fields_->water_content(p), INSERT_VALUES);
75  }
76 
77  /// Implements @p AssemblyBase::begin.
78  void begin() override
79  {
80  VecZeroEntries(eq_data_->mass_diag);
81  eq_data_->balance_->start_mass_assembly(eq_data_->subst_idx);
82  }
83 
84  /// Implements @p AssemblyBase::end.
85  void end() override
86  {
87  eq_data_->balance_->finish_mass_assembly(eq_data_->subst_idx);
88 
89  VecAssemblyBegin(eq_data_->mass_diag);
90  VecAssemblyEnd(eq_data_->mass_diag);
91 
92  eq_data_->is_mass_diag_changed = true;
93  }
94 
95  private:
96  shared_ptr<FiniteElement<dim>> fe_; ///< Finite element for the solution of the advection-diffusion equation.
97 
98  /// Data objects shared with TransportDG
101 
102  /// Sub field set contains fields used in calculation.
104 
105  /// Bulk integral of assembly class
106  std::shared_ptr<BulkIntegralAcc<dim>> mass_integral_;
107 
108  template < template<IntDim...> class DimAssembly>
109  friend class GenericAssembly;
110 
111 };
112 
113 
114 
115 /**
116  * Auxiliary container class for Finite element and related objects of given dimension.
117  */
118 template <unsigned int dim, class TEqData>
120 {
121 public:
122  typedef typename TEqData::EqFields EqFields;
123  typedef TEqData EqData;
124 
125  static constexpr const char * name() { return "Convection_InitCond_Assembly"; }
126 
127  /// Constructor.
129  : AssemblyBasePatch<dim>(0, asm_internals), eq_fields_(eq_data->eq_fields_.get()), eq_data_(eq_data),
130  bulk_integral_( this->create_bulk_integral(this->quad_) ) {
131  this->used_fields_ += eq_fields_->init_conc;
132  }
133 
134  /// Destructor.
136 
137  /// Initialize auxiliary vectors and other data members
138  void initialize() {
139  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++) {
140  vecs_.push_back(eq_fields_->conc_mobile_fe[sbi]->vec());
141  }
142  }
143 
144 
145  /// Assemble integral over element
146  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
147  {
148  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
149 
150  LongIdx index = cell.local_idx();
151  auto p = *( bulk_integral_->points(element_patch_idx).begin() );
152 
153  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++) {
154  vecs_[sbi].set( index, eq_fields_->init_conc[sbi](p) );
155  }
156  }
157 
158 private:
159  //shared_ptr<FiniteElement<dim>> fe_; ///< Finite element for the solution of the advection-diffusion equation.
160 
161  /// Data objects shared with TransportDG
164 
165  std::vector<VectorMPI> vecs_; ///< Set of data vectors of conc_mobile_fe objects.
166 
167  /// Sub field set contains fields used in calculation.
169 
170  /// Bulk integral of assembly class
171  std::shared_ptr<BulkIntegralAcc<dim>> bulk_integral_;
172 
173  template < template<IntDim...> class DimAssembly>
174  friend class GenericAssembly;
175 
176 };
177 
178 
179 /**
180  * Auxiliary container class for Finite element and related objects of given dimension.
181  *
182  * Assembles concentration sources for each substance and set boundary conditions.
183  * note: the source of concentration is multiplied by time interval (gives the mass, not the flow like before)
184  */
185 template <unsigned int dim, class TEqData>
187 {
188 public:
189  typedef typename TEqData::EqFields EqFields;
190  typedef TEqData EqData;
191 
192  static constexpr const char * name() { return "Convection_ConcSourcesBdr_Assembly"; }
193 
194  /// Constructor.
196  : AssemblyBasePatch<dim>(0, asm_internals), eq_fields_(eq_data->eq_fields_.get()), eq_data_(eq_data),
199  JxW_bdr_( bdr_integral_->JxW() ),
200  normal_bdr_( bdr_integral_->normal_vector() ) {
201  this->used_fields_ += eq_fields_->cross_section;
202  this->used_fields_ += eq_fields_->sources_sigma;
203  this->used_fields_ += eq_fields_->sources_density;
204  this->used_fields_ += eq_fields_->sources_conc;
205  this->used_fields_ += eq_fields_->flow_flux;
206  this->used_fields_ += eq_fields_->bc_conc;
207  }
208 
209  /// Destructor.
211 
212  /// Initialize auxiliary vectors and other data members
213  void initialize() {}
214 
215 
216  /// Assemble integral over element
217  inline void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
218  {
219  ASSERT_EQ(cell.dim(), dim).error("Dimension of element mismatch!");
220  if (!eq_data_->sources_changed_) return;
221 
222  // read for all substances
223  double max_cfl=0;
224  double source, diag;
225  ElementAccessor<3> elm = cell.elm();
226  // we have currently zero order P_Disc FE
227  ASSERT(cell.get_loc_dof_indices().size() == 1);
228  IntIdx local_p0_dof = cell.get_loc_dof_indices()[0];
229 
230  auto p = *( bulk_integral_->points(element_patch_idx).begin() );
231  for (unsigned int sbi = 0; sbi < eq_data_->n_substances(); sbi++)
232  {
233  source = eq_fields_->cross_section(p) * (eq_fields_->sources_density[sbi](p)
234  + eq_fields_->sources_sigma[sbi](p) * eq_fields_->sources_conc[sbi](p));
235  // addition to RHS
236  eq_data_->corr_vec[sbi].set(local_p0_dof, source);
237  // addition to diagonal of the transport matrix
238  diag = eq_fields_->sources_sigma[sbi](p) * eq_fields_->cross_section(p);
239  eq_data_->tm_diag[sbi].set(local_p0_dof, - diag);
240 
241  // compute maximal cfl condition over all substances
242  max_cfl = std::max(max_cfl, fabs(diag));
243 
244  eq_data_->balance_->add_source_values(sbi, elm.region().bulk_idx(), {local_p0_dof},
245  {- eq_fields_->sources_sigma[sbi](p) * elm.measure() * eq_fields_->cross_section(p)},
246  {source * elm.measure()});
247  }
248 
249  eq_data_->cfl_source_.set(local_p0_dof, max_cfl);
250  }
251 
252  /// Assembles the fluxes on the boundary.
253  inline void boundary_side_integral(DHCellSide cell_side)
254  {
255  ASSERT_EQ(cell_side.dim(), dim).error("Dimension of element mismatch!");
256  if (!cell_side.cell().is_own()) return;
257 
258  // we have currently zero order P_Disc FE
259  ASSERT(cell_side.cell().get_loc_dof_indices().size() == 1);
260  IntIdx local_p0_dof = cell_side.cell().get_loc_dof_indices()[0];
261  LongIdx glob_p0_dof = eq_data_->dh_->get_local_to_global_map()[local_p0_dof];
262 
263  unsigned int sbi;
264  auto p_side = *( bdr_integral_->points(cell_side).begin() );
265  auto p_bdr = p_side.point_bdr(cell_side.cond().element_accessor() );
266  double flux = eq_fields_->side_flux(p_side, normal_bdr_(p_side), JxW_bdr_(p_side));
267  if (flux < 0.0) {
268  double aij = -(flux / cell_side.element().measure() );
269 
270  for (sbi=0; sbi<eq_data_->n_substances(); sbi++)
271  {
272  double value = eq_fields_->bc_conc[sbi](p_bdr);
273 
274  VecSetValue(eq_data_->bcvcorr[sbi], glob_p0_dof, value * aij, ADD_VALUES);
275 
276  // CAUTION: It seems that PETSc possibly optimize allocated space during assembly.
277  // So we have to add also values that may be non-zero in future due to changing velocity field.
278  eq_data_->balance_->add_flux_values(eq_data_->subst_idx[sbi], cell_side,
279  {local_p0_dof}, {0.0}, flux*value);
280  }
281  } else {
282  for (sbi=0; sbi<eq_data_->n_substances(); sbi++)
283  VecSetValue(eq_data_->bcvcorr[sbi], glob_p0_dof, 0, ADD_VALUES);
284 
285  for (sbi=0; sbi<eq_data_->n_substances(); sbi++)
286  eq_data_->balance_->add_flux_values(eq_data_->subst_idx[sbi], cell_side,
287  {local_p0_dof}, {flux}, 0.0);
288  }
289  }
290 
291  /// Implements @p AssemblyBase::begin.
292  void begin() override
293  {
294  eq_data_->sources_changed_ = ( (eq_fields_->sources_density.changed() )
295  || (eq_fields_->sources_conc.changed() )
296  || (eq_fields_->sources_sigma.changed() )
297  || (eq_fields_->cross_section.changed()));
298 
299  //TODO: would it be possible to check the change in data for chosen substance? (may be in multifields?)
300 
301  if (eq_data_->sources_changed_) eq_data_->balance_->start_source_assembly(eq_data_->subst_idx);
302  // Assembly bcvcorr vector
303  for(unsigned int sbi=0; sbi < eq_data_->n_substances(); sbi++) VecZeroEntries(eq_data_->bcvcorr[sbi]);
304 
305  eq_data_->balance_->start_flux_assembly(eq_data_->subst_idx);
306  }
307 
308  /// Implements @p AssemblyBase::end.
309  void end() override
310  {
311  eq_data_->balance_->finish_flux_assembly(eq_data_->subst_idx);
312  if (eq_data_->sources_changed_) eq_data_->balance_->finish_source_assembly(eq_data_->subst_idx);
313 
314  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++) VecAssemblyBegin(eq_data_->bcvcorr[sbi]);
315  for (unsigned int sbi=0; sbi<eq_data_->n_substances(); sbi++) VecAssemblyEnd(eq_data_->bcvcorr[sbi]);
316 
317  // we are calling set_boundary_conditions() after next_time() and
318  // we are using data from t() before, so we need to set corresponding bc time
319  eq_data_->transport_bc_time = eq_data_->time_->last_t();
320  }
321 
322  private:
323  /// Data objects shared with ConvectionTransport
326 
327  /// Sub field set contains fields used in calculation.
329 
330  std::shared_ptr<BulkIntegralAcc<dim>> bulk_integral_; ///< Bulk integral of assembly class
331  std::shared_ptr<BoundaryIntegralAcc<dim>> bdr_integral_; ///< Boundary integral of assembly class
332 
335 
336  template < template<IntDim...> class DimAssembly>
337  friend class GenericAssembly;
338 
339 };
340 
341 
342 /**
343  * Auxiliary container class for Finite element and related objects of given dimension.
344  *
345  *
346  * Assembly convection term part of the matrix and boundary matrix for application of boundary conditions.
347  *
348  * Discretization of the convection term use explicit time scheme and finite volumes with full upwinding.
349  * We count on with exchange between dimensions and mixing on edges where more then two elements connect (can happen for 2D and 1D elements in
350  * 3D embedding space)
351  *
352  * In order to get multiplication matrix for explicit transport one have to scale the convection part by the acctual time step and
353  * add time term, i. e. unit matrix (see. transport_matrix_step_mpi)
354  *
355  * Updates CFL time step constrain.
356  */
357 template <unsigned int dim, class TEqData>
359 {
360 public:
361  typedef typename TEqData::EqFields EqFields;
362  typedef TEqData EqData;
363 
364  static constexpr const char * name() { return "Convection_MatrixMpi_Assembly"; }
365 
366  /// Constructor.
368  : AssemblyBasePatch<dim>(0, asm_internals), eq_fields_(eq_data->eq_fields_.get()), eq_data_(eq_data),
371  JxW_side_( edge_integral_->JxW() ),
372  JxW_join_( coupling_integral_->JxW() ),
373  normal_side_( edge_integral_->normal_vector() ),
374  normal_join_( coupling_integral_->normal_vector() ) {
375  this->used_fields_ += eq_fields_->flow_flux;
376  }
377 
378  /// Destructor.
380 
381  /// Initialize auxiliary vectors and other data members
382  void initialize() {
383  side_dofs_.resize(eq_data_->max_edg_sides);
384  side_flux_.resize(eq_data_->max_edg_sides);
385  elm_meassures_.resize(eq_data_->max_edg_sides);
386  all_elem_dofs_.resize(eq_data_->max_edg_sides-1);
387  row_values_.resize(eq_data_->max_edg_sides-1);
388  dof_indices_i_.resize(1);
389  dof_indices_j_.resize(1);
390  }
391 
392  /// Assembles the fluxes between sides of elements of the same dimension.
394  ASSERT_EQ(edge_side_range.begin()->element().dim(), dim).error("Dimension of element mismatch!");
395 
396  unsigned int sid=0, s1, s2, i_col;
397  edg_flux = 0.0;
398  for( DHCellSide edge_side : edge_side_range )
399  {
400  edge_side.cell().get_dof_indices(dof_indices_i_);
401  side_dofs_[sid] = dof_indices_i_[0];
402  elm_meassures_[sid] = edge_side.element().measure();
403  auto p = *( edge_integral_->points(edge_side).begin() );
404  side_flux_[sid] = eq_fields_->side_flux(p, normal_side_(p), JxW_side_(p));
405  if (side_flux_[sid] > 0.0) {
406  eq_data_->cfl_flow_.add_global(side_dofs_[sid], -(side_flux_[sid] / edge_side.element().measure()) );
407  edg_flux += side_flux_[sid];
408  }
409  ++sid;
410  }
411 
412  unsigned int n_edge_sides = edge_side_range.begin()->n_edge_sides();
413  if (n_edge_sides<2) return; // following loop has no effect if edge contains only 1 side
414  for( s1=0; s1<n_edge_sides; s1++ )
415  {
416  for( s2=0, i_col=0; s2<n_edge_sides; s2++ )
417  {
418  if (s2==s1) continue;
419 
420  if (side_flux_[s2] > 0.0 && side_flux_[s1] < 0.0)
421  row_values_[i_col] = -(side_flux_[s1] * side_flux_[s2] / ( edg_flux * elm_meassures_[s1] ) );
422  else row_values_[i_col] = 0;
423  all_elem_dofs_[i_col++] = side_dofs_[s2];
424  }
425  MatSetValues(eq_data_->tm, 1, &(side_dofs_[s1]), n_edge_sides-1, &(all_elem_dofs_[0]), &(row_values_[0]), INSERT_VALUES);
426  }
427  }
428 
429 
430  /// Assembles the fluxes between elements of different dimensions.
431  inline void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side) {
432  if (dim == 3) return;
433  ASSERT_EQ(cell_lower_dim.dim(), dim).error("Dimension of element mismatch!");
434 
435  auto p_high = *( coupling_integral_->points(neighb_side).begin() );
436 
437  cell_lower_dim.get_dof_indices(dof_indices_i_);
438  neighb_side.cell().get_dof_indices(dof_indices_j_);
439  flux = eq_fields_->side_flux(p_high, normal_join_(p_high), JxW_join_(p_high));
440 
441  // volume source - out-flow from higher dimension
442  if (flux > 0.0) aij = flux / cell_lower_dim.elm().measure();
443  else aij=0;
444  MatSetValue(eq_data_->tm, dof_indices_i_[0], dof_indices_j_[0], aij, INSERT_VALUES);
445  // out flow from higher dim. already accounted
446 
447  // volume drain - in-flow to higher dimension
448  if (flux < 0.0) {
449  eq_data_->cfl_flow_.add_global( dof_indices_i_[0], (flux / cell_lower_dim.elm().measure()) ); // diagonal drain
450  aij = (-flux) / neighb_side.element().measure();
451  } else aij=0;
452  MatSetValue(eq_data_->tm, dof_indices_j_[0], dof_indices_i_[0], aij, INSERT_VALUES);
453  }
454 
455 
456  /// Implements @p AssemblyBase::begin.
457  void begin() override
458  {
459  MatZeroEntries(eq_data_->tm);
460  eq_data_->cfl_flow_.zero_entries();
461  }
462 
463  /// Implements @p AssemblyBase::end.
464  void end() override
465  {
466  for ( DHCellAccessor dh_cell : eq_data_->dh_->own_range() ) {
467  dh_cell.get_dof_indices(dof_indices_i_);
468  MatSetValue(eq_data_->tm, dof_indices_i_[0], dof_indices_i_[0], eq_data_->cfl_flow_.get(dh_cell.local_idx()), INSERT_VALUES);
469 
470  eq_data_->cfl_flow_.set(dh_cell.local_idx(), fabs(eq_data_->cfl_flow_.get(dh_cell.local_idx())) );
471  }
472 
473  MatAssemblyBegin(eq_data_->tm, MAT_FINAL_ASSEMBLY);
474  MatAssemblyEnd(eq_data_->tm, MAT_FINAL_ASSEMBLY);
475  VecAssemblyBegin(eq_data_->cfl_flow_.petsc_vec());
476  VecAssemblyEnd(eq_data_->cfl_flow_.petsc_vec());
477 
478  eq_data_->is_convection_matrix_scaled = false;
479  eq_data_->transport_matrix_time = eq_data_->time_->t();
480  }
481 
482 private:
483  /// Data objects shared with ConvectionTransport
486 
487  /// Sub field set contains fields used in calculation.
489 
490  vector<LongIdx> dof_indices_i_, dof_indices_j_; ///< Global DOF indices.
491 
497  double aij;
498  double edg_flux, flux;
499 
500  std::shared_ptr<EdgeIntegralAcc<dim>> edge_integral_; ///< Edge integral of assembly class
501  std::shared_ptr<CouplingIntegralAcc<dim>> coupling_integral_; ///< Coupling integral of assembly class
502 
507 
508  template < template<IntDim...> class DimAssembly>
509  friend class GenericAssembly;
510 
511 };
512 
513 
514 #endif /* ASSEMBLY_CONVECTION_HH_ */
#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)
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.
std::shared_ptr< EdgeIntegralAcc< dim > > create_edge_integral(Quadrature *quad)
std::shared_ptr< BoundaryIntegralAcc< dim > > create_boundary_integral(Quadrature *quad)
ElementAccessor< 3 > element_accessor()
ConcSourcesBdrAssemblyConvection(EqData *eq_data, AssemblyInternals *asm_internals)
Constructor.
std::shared_ptr< BulkIntegralAcc< dim > > bulk_integral_
Bulk integral of assembly class.
static constexpr const char * name()
void initialize()
Initialize auxiliary vectors and other data members.
FieldSet used_fields_
Sub field set contains fields used in calculation.
void boundary_side_integral(DHCellSide cell_side)
Assembles the fluxes on the boundary.
std::shared_ptr< BoundaryIntegralAcc< dim > > bdr_integral_
Boundary integral of assembly class.
void end() override
Implements AssemblyBase::end.
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
EqFields * eq_fields_
Data objects shared with ConvectionTransport.
void begin() override
Implements AssemblyBase::begin.
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_.
unsigned int local_idx() const
Return local index to element (index of DOF handler).
Side accessor allows to iterate over sides of DOF handler cell.
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.
ElementAccessor< 3 > element() const
double measure() const
Computes the measure of the element.
Region region() const
Definition: accessors.hh:198
Container for various descendants of FieldCommonBase.
Definition: field_set.hh:160
Generic class of assemblation.
FieldSet used_fields_
Sub field set contains fields used in calculation.
void initialize()
Initialize auxiliary vectors and other data members.
std::vector< VectorMPI > vecs_
Set of data vectors of conc_mobile_fe objects.
std::shared_ptr< BulkIntegralAcc< dim > > bulk_integral_
Bulk integral of assembly class.
EqFields * eq_fields_
Data objects shared with TransportDG.
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
static constexpr const char * name()
InitCondAssemblyConvection(EqData *eq_data, AssemblyInternals *asm_internals)
Constructor.
MassAssemblyConvection(EqData *eq_data, AssemblyInternals *asm_internals)
Constructor.
FieldSet used_fields_
Sub field set contains fields used in calculation.
~MassAssemblyConvection()
Destructor.
std::shared_ptr< BulkIntegralAcc< dim > > mass_integral_
Bulk integral of assembly class.
shared_ptr< FiniteElement< dim > > fe_
Finite element for the solution of the advection-diffusion equation.
TEqData::EqFields EqFields
void begin() override
Implements AssemblyBase::begin.
static constexpr const char * name()
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
void initialize()
Initialize auxiliary vectors and other data members.
void end() override
Implements AssemblyBase::end.
EqFields * eq_fields_
Data objects shared with TransportDG.
std::vector< double > elm_meassures_
void initialize()
Initialize auxiliary vectors and other data members.
std::vector< LongIdx > all_elem_dofs_
std::shared_ptr< EdgeIntegralAcc< dim > > edge_integral_
Edge integral of assembly class.
void end() override
Implements AssemblyBase::end.
static constexpr const char * name()
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
Assembles the fluxes between elements of different dimensions.
EqFields * eq_fields_
Data objects shared with ConvectionTransport.
FieldSet used_fields_
Sub field set contains fields used in calculation.
std::vector< double > side_flux_
std::shared_ptr< CouplingIntegralAcc< dim > > coupling_integral_
Coupling integral of assembly class.
std::vector< LongIdx > side_dofs_
void begin() override
Implements AssemblyBase::begin.
vector< LongIdx > dof_indices_j_
Global DOF indices.
std::vector< double > row_values_
MatrixMpiAssemblyConvection(EqData *eq_data, AssemblyInternals *asm_internals)
Constructor.
void edge_integral(RangeConvert< DHEdgeSide, DHCellSide > edge_side_range)
Assembles the fluxes between sides of elements of the same dimension.
IterConvert< ObjectIn, ObjectOut > begin()
Iterator to begin item of range.
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
Definition: region.hh:90
int IntIdx
Definition: index_types.hh:25
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: index_types.hh:24
static constexpr bool value
Definition: json.hpp:87
unsigned int IntDim
A dimension index type.
Definition: mixed.hh:19
Class FEValues calculates finite element data on the actual cells such as shape function values,...
Store finite element reinit functions.
Definitions of particular quadrature rules on simplices.
Holds common data shared between GenericAssemblz and Assembly<dim> classes.