18 #ifndef ASSEMBLY_DG_HH_
19 #define ASSEMBLY_DG_HH_
35 template <
unsigned int dim,
class Model>
42 static constexpr
const char *
name() {
return "MassAssemblyDG"; }
77 ASSERT_EQ(cell.
dim(), dim).error(
"Dimension of element mismatch!");
81 for (
unsigned int sbi=0; sbi<
eq_data_->n_substances(); ++sbi)
84 for (
unsigned int i=0; i<
ndofs_; i++)
86 for (
unsigned int j=0; j<
ndofs_; j++)
89 for (
auto p : this->
bulk_points(element_patch_idx) )
97 for (
unsigned int i=0; i<
ndofs_; i++)
101 for (
auto p : this->
bulk_points(element_patch_idx) )
145 template <
template<
IntDim...>
class DimAssembly>
154 double h_max = 0, h_min = numeric_limits<double>::infinity();
155 for (
unsigned int i=0; i<e->
n_nodes(); i++)
156 for (
unsigned int j=i+1; j<e->
n_nodes(); j++)
158 double dist = arma::norm(*e.
node(i) - *e.
node(j));
159 h_max = max(h_max, dist);
160 h_min = min(h_min, dist);
179 delta += dot(diff_coef(p)*nv, nv);
182 return n == 0 ? 0 : (delta/n);
197 const double &diff_delta,
214 template <
class Po
intType>
217 double side_flux = 0;
219 side_flux += arma::dot(advection_coef(p), normal(p))*JxW(p);
228 template <
unsigned int dim,
class Model>
235 static constexpr
const char *
name() {
return "StiffnessAssemblyDG"; }
261 for (
auto a :
averages)
if (a !=
nullptr)
delete[] a;
262 for (
auto a :
waverages)
if (a !=
nullptr)
delete[] a;
263 for (
auto a :
jumps)
if (a !=
nullptr)
delete[] a;
294 for (
unsigned int s=0; s<2; s++)
305 ASSERT_EQ(cell.
dim(), dim).error(
"Dimension of element mismatch!");
306 if (!cell.
is_own())
return;
311 for (
unsigned int sbi=0; sbi<
eq_data_->n_substances(); sbi++)
313 for (
unsigned int i=0; i<
ndofs_; i++)
314 for (
unsigned int j=0; j<
ndofs_; j++)
317 for (
auto p : this->
bulk_points(element_patch_idx) )
319 for (
unsigned int i=0; i<
ndofs_; i++)
324 for (
unsigned int j=0; j<
ndofs_; j++)
338 ASSERT_EQ(cell_side.
dim(), dim).error(
"Dimension of element mismatch!");
348 for (
unsigned int sbi=0; sbi<
eq_data_->n_substances(); sbi++)
353 double transport_flux = side_flux/side.
measure();
359 unsigned int bc_type =
eq_fields_->bc_type[sbi](p_bdr);
368 transport_flux += gamma_l;
375 double flux_times_JxW;
390 flux_times_JxW = transport_flux*
JxW_side_(p);
392 for (
unsigned int i=0; i<
ndofs_; i++)
394 for (
unsigned int j=0; j<
ndofs_; j++)
416 ASSERT_EQ(edge_side_range.
begin()->element().dim(), dim).error(
"Dimension of element mismatch!");
419 double gamma_l, omega[2], transport_flux, delta[2], delta_sum;
420 double aniso1, aniso2;
421 double local_alpha=0.0;
425 auto dh_edge_cell =
eq_data_->
dh_->cell_accessor_from_element( edge_side.elem_idx() );
429 auto zero_edge_side = *edge_side_range.
begin();
430 auto p = *( this->
edge_points(zero_edge_side).begin() );
434 for (
unsigned int sbi=0; sbi<
eq_data_->n_substances(); sbi++)
437 double pflux = 0, nflux = 0;
443 pflux += fluxes[sid];
445 nflux += fluxes[sid];
456 for (
unsigned int i=0; i<
ndofs_; i++)
466 for(
DHCellSide edge_side1 : edge_side_range )
469 for(
DHCellSide edge_side2 : edge_side_range )
472 if (s2<=s1)
continue;
473 ASSERT(edge_side1.is_valid()).error(
"Invalid side of edge.");
475 auto p = *( this->
edge_points(edge_side1).begin() );
480 if (fluxes[s2] > 0 && fluxes[s1] < 0)
481 transport_flux = fluxes[s1]*fabs(fluxes[s2]/pflux);
482 else if (fluxes[s2] < 0 && fluxes[s1] > 0)
483 transport_flux = fluxes[s1]*fabs(fluxes[s2]/nflux);
487 gamma_l = 0.5*fabs(transport_flux);
493 auto p2 = p1.point_on(edge_side2);
494 delta[0] += dot(
eq_fields_->diffusion_coef[sbi](p1)*normal_vector,normal_vector);
495 delta[1] += dot(
eq_fields_->diffusion_coef[sbi](p2)*normal_vector,normal_vector);
501 delta_sum = delta[0] + delta[1];
504 if (fabs(delta_sum) > 0)
506 omega[0] = delta[1]/delta_sum;
507 omega[1] = delta[0]/delta_sum;
508 double h = edge_side1.diameter();
511 gamma_l += local_alpha/h*aniso1*aniso2*(delta[0]*delta[1]/delta_sum);
514 for (
int i=0; i<2; i++) omega[i] = 0;
517 int sd[2];
bool is_side_own[2];
518 sd[0] = s1; is_side_own[0] = edge_side1.cell().is_own();
519 sd[1] = s2; is_side_own[1] = edge_side2.cell().is_own();
525 auto p2 = p1.point_on(edge_side2);
526 for (
unsigned int i=0; i<
ndofs_; i++)
537 for (
int n=0; n<2; n++)
539 if (!is_side_own[n])
continue;
541 for (
int m=0; m<2; m++)
543 for (
unsigned int i=0; i<
ndofs_; i++)
544 for (
unsigned int j=0; j<
ndofs_; j++)
551 for (
unsigned int i=0; i<
ndofs_; i++)
553 for (
unsigned int j=0; j<
ndofs_; j++)
584 if (dim == 1)
return;
585 ASSERT_EQ(cell_lower_dim.
dim(), dim-1).error(
"Dimension of element mismatch!");
591 for(
unsigned int i=0; i<n_indices; ++i) {
598 for(
unsigned int i=0; i<n_indices; ++i) {
604 bool own_element_id[2];
605 own_element_id[0] = cell_lower_dim.
is_own();
606 own_element_id[1] = cell_higher_dim.
is_own();
608 for (
unsigned int sbi=0; sbi<
eq_data_->n_substances(); sbi++)
617 auto p_low = p_high.lower_dim(cell_lower_dim);
628 double transport_flux = arma::dot(
eq_fields_->advection_coef[sbi](p_high),
normal_(p_high));
631 uint is_high_i = conc_shape_i.is_high_dim();
632 if (!own_element_id[is_high_i])
continue;
633 uint i_mat_idx = conc_shape_i.join_idx();
634 double diff_shape_i = conc_shape_i(p_high) - conc_shape_i(p_low);
636 uint j_mat_idx = conc_shape_j.join_idx();
638 sigma * diff_shape_i * (conc_shape_j(p_high) - conc_shape_j(p_low))
639 + diff_shape_i * ( max(0.,transport_flux) * conc_shape_j(p_high) + min(0.,transport_flux) * conc_shape_j(p_low))
678 template <
template<
IntDim...>
class DimAssembly>
687 template <
unsigned int dim,
class Model>
694 static constexpr
const char *
name() {
return "SourcesAssemblyDG"; }
729 ASSERT_EQ(cell.
dim(), dim).error(
"Dimension of element mismatch!");
737 for (
unsigned int sbi=0; sbi<
eq_data_->n_substances(); sbi++)
743 for (
auto p : this->
bulk_points(element_patch_idx) )
747 for (
unsigned int i=0; i<
ndofs_; i++)
752 for (
unsigned int i=0; i<
ndofs_; i++)
754 for (
auto p : this->
bulk_points(element_patch_idx) )
798 template <
template<
IntDim...>
class DimAssembly>
807 template <
unsigned int dim,
class Model>
814 static constexpr
const char *
name() {
return "BdrConditionAssemblyDG"; }
860 for (
unsigned int sbi=0; sbi<
eq_data_->n_substances(); sbi++)
867 double transport_flux = side_flux/cell_side.
measure();
871 unsigned int bc_type =
eq_fields_->bc_type[sbi](p_bdr);
876 auto p_bdr = p.point_bdr(bc_elm);
877 double bc_term = -transport_flux*
eq_fields_->bc_dirichlet_value[sbi](p_bdr)*
JxW_(p);
878 for (
unsigned int i=0; i<
ndofs_; i++)
881 for (
unsigned int i=0; i<
ndofs_; i++)
887 double transport_flux = side_flux/cell_side.
measure();
897 auto p_bdr = p.point_bdr(bc_elm);
898 double bc_term = gamma_l*
eq_fields_->bc_dirichlet_value[sbi](p_bdr)*
JxW_(p);
900 for (
unsigned int i=0; i<
ndofs_; i++)
906 for (
unsigned int i=0; i<
ndofs_; i++)
914 for (
unsigned int i=0; i<
ndofs_; i++)
921 auto p_bdr = p.point_bdr(bc_elm);
924 for (
unsigned int i=0; i<
ndofs_; i++)
928 for (
unsigned int i=0; i<
ndofs_; i++)
931 auto p_bdr = p.point_bdr(bc_elm);
942 auto p_bdr = p.point_bdr(bc_elm);
945 for (
unsigned int i=0; i<
ndofs_; i++)
949 for (
unsigned int i=0; i<
ndofs_; i++)
952 auto p_bdr = p.point_bdr(bc_elm);
963 for (
unsigned int i=0; i<
ndofs_; i++)
1008 template <
template<
IntDim...>
class DimAssembly>
1017 template <
unsigned int dim,
class Model>
1024 static constexpr
const char *
name() {
return "InitProjectionAssemblyDG"; }
1056 ASSERT_EQ(cell.
dim(), dim).error(
"Dimension of element mismatch!");
1060 for (
unsigned int sbi=0; sbi<
eq_data_->n_substances(); sbi++)
1062 for (
unsigned int i=0; i<
ndofs_; i++)
1065 for (
unsigned int j=0; j<
ndofs_; j++)
1069 for (
auto p : this->
bulk_points(element_patch_idx) )
1073 for (
unsigned int i=0; i<
ndofs_; i++)
1075 for (
unsigned int j=0; j<
ndofs_; j++)
1102 template <
template<
IntDim...>
class DimAssembly>
1111 template <
unsigned int dim,
class Model>
1118 static constexpr
const char *
name() {
return "InitConditionAssemblyDG"; }
1127 for(
unsigned int i = 0; i<RefElement<dim>::n_nodes; i++)
1146 ASSERT_EQ(cell.
dim(), dim).error(
"Dimension of element mismatch!");
1150 for (
unsigned int sbi=0; sbi<
eq_data_->n_substances(); sbi++)
1153 for (
auto p : this->
bulk_points(element_patch_idx) )
1155 double val =
eq_fields_->init_condition[sbi](p);
1172 template <
template<
IntDim...>
class DimAssembly>
double elem_anisotropy(ElementAccessor< 3 > e)
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.
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.
double advective_flux(Field< 3, FieldValue< 3 >::VectorFixed > &advection_coef, Range< PointType > pts, ElQ< Scalar > &JxW, ElQ< Vector > normal)
Computes advective flux.
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
JoinValues< dim > join_values()
Return JoinValues object.
PatchFEValues< 3 > * fe_values_
Common FEValues object over all dimensions.
BulkValues< dim > bulk_values()
Return BulkValues object.
unsigned int n_dofs()
Return BulkValues object.
SideValues< dim > side_values()
Return SideValues object.
Range< BulkPoint > bulk_points(unsigned int element_patch_idx) const
Return BulkPoint range of appropriate dimension.
Range< EdgePoint > edge_points(const DHCellSide &cell_side) const
Return EdgePoint range of appropriate dimension.
Quadrature * quad_
Quadrature used in assembling methods.
std::string print_update_flags(UpdateFlags u) const
Print update flags to string format.
Quadrature * quad_low_
Quadrature used in assembling methods (dim-1).
Range< CouplingPoint > coupling_points(const DHCellSide &cell_side) const
Return CouplingPoint range of appropriate dimension.
int active_integrals_
Holds mask of active integrals.
Range< BoundaryPoint > boundary_points(const DHCellSide &cell_side) const
Return BoundaryPoint range of appropriate dimension.
ElementCacheMap * element_cache_map_
ElementCacheMap shared with GenericAssembly object.
static constexpr const char * name()
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.
BdrConditionAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
TransportDG< Model >::EqData EqData
TransportDG< Model >::EqFields EqFields
void begin() override
Implements AssemblyBase::begin.
~BdrConditionAssemblyDG()
Destructor.
PetscScalar local_flux_balance_rhs_
Auxiliary variable for set_boundary_conditions method.
FeQ< Scalar > conc_shape_
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
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.
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.
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
NodeAccessor< 3 > node(unsigned int ni) const
unsigned int idx() const
We need this method after replacing Region by RegionIdx, and movinf RegionDB instance into particular...
Directing class of FieldValueCache.
unsigned int n_nodes() const
Container for various descendants of FieldCommonBase.
Class template representing a field with values dependent on: point, element, and region.
Generic class of assemblation.
FieldSet used_fields_
Sub field set contains fields used in calculation.
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
static constexpr const char * name()
~InitConditionAssemblyDG()
Destructor.
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
InitConditionAssemblyDG(EqFields *eq_fields, EqData *eq_data)
Constructor.
EqFields * eq_fields_
Data objects shared with TransportDG.
TransportDG< Model >::EqData EqData
TransportDG< Model >::EqFields EqFields
EqFields * eq_fields_
Data objects shared with TransportDG.
vector< PetscScalar > local_rhs_
Auxiliary vector for set_sources method.
vector< LongIdx > dof_indices_
Vector of global DOF indices.
TransportDG< Model >::EqData EqData
FieldSet used_fields_
Sub field set contains fields used in calculation.
unsigned int ndofs_
Number of dofs.
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
FeQ< Scalar > init_shape_
InitProjectionAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
static constexpr const char * name()
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
~InitProjectionAssemblyDG()
Destructor.
TransportDG< Model >::EqFields EqFields
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
virtual void rhs_set_values(int nrow, int *rows, double *vals)=0
virtual void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)=0
void set_values(int nrow, int *rows, int ncol, int *cols, PetscScalar *mat_vals, PetscScalar *rhs_vals)
Set values in the system matrix and values in the right-hand side vector on corresponding rows.
static constexpr double almost_zero
vector< LongIdx > dof_indices_
Vector of global DOF indices.
vector< PetscScalar > local_mass_balance_vector_
Same as previous.
TransportDG< Model >::EqFields EqFields
MassAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
FieldSet used_fields_
Sub field set contains fields used in calculation.
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
void end() override
Implements AssemblyBase::end.
~MassAssemblyDG()
Destructor.
void begin() override
Implements AssemblyBase::begin.
static constexpr const char * name()
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
FeQ< Scalar > conc_shape_
unsigned int ndofs_
Number of dofs.
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
TransportDG< Model >::EqData EqData
EqFields * eq_fields_
Data objects shared with TransportDG.
vector< PetscScalar > local_retardation_balance_vector_
Auxiliary vector for assemble mass matrix.
Base class for quadrature rules on simplices in arbitrary dimensions.
unsigned int size() const
Returns number of quadrature points.
Armor::Array< double >::ArrayMatSet set(uint i)
double weight(unsigned int i) const
Returns the ith weight.
IterConvert< ObjectIn, ObjectOut > begin()
Iterator to begin item of range.
Iter< Object > begin()
Iterator to begin item of range.
static LocalPoint node_coords(unsigned int nid)
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
double measure() const
Calculate metrics of the side.
ElementAccessor< 3 > element() const
Returns iterator to the element of the side.
double diameter() const
Calculate the side diameter.
vector< PetscScalar > local_rhs_
Auxiliary vector for set_sources method.
TransportDG< Model >::EqData EqData
TransportDG< Model >::EqFields EqFields
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
~SourcesAssemblyDG()
Destructor.
SourcesAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
vector< PetscScalar > local_source_balance_vector_
Auxiliary vector for set_sources method.
void end() override
Implements AssemblyBase::end.
static constexpr const char * name()
EqFields * eq_fields_
Data objects shared with TransportDG.
vector< PetscScalar > local_source_balance_rhs_
Auxiliary vector for set_sources method.
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
FieldSet used_fields_
Sub field set contains fields used in calculation.
unsigned int ndofs_
Number of dofs.
vector< LongIdx > dof_indices_
Vector of global DOF indices.
void begin() override
Implements AssemblyBase::begin.
FeQ< Scalar > conc_shape_
TransportDG< Model >::EqData EqData
FeQ< Scalar > conc_shape_
void edge_integral(RangeConvert< DHEdgeSide, DHCellSide > edge_side_range)
Assembles the fluxes between sides of elements of the same dimension.
FeQ< Vector > conc_grad_sidw_
FieldSet used_fields_
Sub field set contains fields used in calculation.
unsigned int ndofs_
Number of dofs.
vector< vector< LongIdx > > side_dof_indices_
Vector of vectors of side DOF indices.
Range< JoinShapeAccessor< Scalar > > conc_join_shape_
StiffnessAssemblyDG(EqFields *eq_fields, EqData *eq_data, PatchFEValues< 3 > *fe_values)
Constructor.
void boundary_side_integral(DHCellSide cell_side)
Assembles the fluxes on the boundary.
vector< double * > jumps
Auxiliary storage for jumps of shape functions.
vector< double * > waverages
Auxiliary storage for weighted averages of shape functions.
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
Assembles the fluxes between elements of different dimensions.
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
vector< LongIdx > side_dof_indices_vb_
Vector of side DOF indices (assemble element-side fluxex)
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assembles the cell (volume) integral into the stiffness matrix.
TransportDG< Model >::EqFields EqFields
vector< double * > averages
Auxiliary storage for averages of shape functions.
vector< PetscScalar > local_matrix_
Auxiliary vector for assemble methods.
~StiffnessAssemblyDG()
Destructor.
static constexpr const char * name()
vector< LongIdx > dof_indices_
Vector of global DOF indices.
FeQ< Scalar > conc_shape_side_
EqFields * eq_fields_
Data objects shared with TransportDG.
unsigned int qsize_lower_dim_
Size of quadrature of dim-1.
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.
std::shared_ptr< Balance > balance_
LinSys ** ls
Linear algebra system for the transport equation.
std::vector< Vec > ret_vec
Auxiliary vectors for calculation of sources in balance due to retardation (e.g. sorption).
unsigned int max_edg_sides
Maximal number of edge sides (evaluate from dim 1,2,3)
MultiField< 3, FieldValue< 3 >::Scalar > dg_penalty
Penalty enforcing inter-element continuity of solution (for each substance).
MultiField< 3, FieldValue< 3 >::Scalar > fracture_sigma
Transition parameter for diffusive transfer on fractures (for each substance).
Definitions of basic Lagrangean finite elements with polynomial shape functions.
#define DebugOut()
Macro defining 'debug' record of log.
unsigned int IntDim
A dimension index type.
Class FEValues calculates finite element data on the actual cells such as shape function values,...
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.
@ update_values
Shape function values.
@ 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.
@ update_quadrature_points
Transformed quadrature points.