Flow123d
JB_transport-9331eee
|
Go to the documentation of this file.
18 #ifndef ASSEMBLY_LMH_HH_
19 #define ASSEMBLY_LMH_HH_
45 template <
unsigned int dim>
52 static constexpr
const char *
name() {
return "ReadInitCondAssemblyLMH"; }
73 ASSERT_EQ(cell.
dim(), dim).error(
"Dimension of element mismatch!");
79 auto p = *( this->
bulk_points(element_patch_idx).begin() );
82 for (
unsigned int i=0; i<cell.
elm()->
n_sides(); i++) {
110 template <
template<
IntDim...>
class DimAssembly>
115 template <
unsigned int dim>
124 static constexpr
const char *
name() {
return "MHMatrixAssemblyLMH"; }
153 fe_ = std::make_shared< FE_P_disc<dim> >(0);
182 ASSERT_EQ(cell.
dim(), dim).error(
"Dimension of element mismatch!");
185 auto p = *( this->
bulk_points(element_patch_idx).begin() );
201 ASSERT_EQ(cell_side.
dim(), dim).error(
"Dimension of element mismatch!");
224 if (dim == 1)
return;
225 ASSERT_EQ(cell_lower_dim.
dim(), dim-1).error(
"Dimension of element mismatch!");
227 unsigned int neigh_idx =
ngh_idx(cell_lower_dim, neighb_side);
228 unsigned int loc_dof_higher = (2*(cell_lower_dim.
dim()+1) + 1) + neigh_idx;
233 auto p_low = p_high.lower_dim(cell_lower_dim);
296 unsigned int i_constr = 0;
344 schur_solution, this->reconstructed_solution_);
349 this->eq_data_->schur_offset_[dh_cell.dim()-1]), this->reconstructed_solution_);
351 this->loc_schur_.row_dofs.n_elem), schur_solution);
363 auto ele = cell.
elm();
377 arma::dot( velocity.value(i,k),
447 for (
unsigned int i=0; i<n_sides; i++)
479 cross_section_ = eq_fields_->cross_section(p_side);
483 side_row_ = eq_data_->loc_side_dofs[dim-1][sidx_];
484 edge_row_ = eq_data_->loc_edge_dofs[dim-1][sidx_];
490 char & switch_dirichlet = eq_data_->bc_switch_dirichlet[b_ele.
idx()];
491 if (switch_dirichlet) {
497 auto reconstr_solution = this->load_local_system(cell_side.
cell());
498 double side_flux = -eq_fields_->bc_flux(p_bdr) * b_ele.
measure() * cross_section_;
500 if ( reconstr_solution[side_row_] < side_flux) {
502 switch_dirichlet = 0;
517 double solution_head = eq_data_->p_edge_solution.get(loc_dof_vec[sidx_]);
518 double bc_pressure = eq_fields_->bc_switch_pressure(p_bdr);
520 if ( solution_head > bc_pressure) {
535 eq_data_->loc_constraint_.emplace_back( bulk_local_idx_, sidx_, eq_fields_->bc_pressure(p_bdr) );
539 eq_data_->loc_system_[bulk_local_idx_].add_value(edge_row_, edge_row_,
540 -b_ele.
measure() * eq_fields_->bc_robin_sigma(p_bdr) * cross_section_,
541 (-eq_fields_->bc_flux(p_bdr) - eq_fields_->bc_robin_sigma(p_bdr) * eq_fields_->bc_pressure(p_bdr)) * b_ele.
measure() * cross_section_);
543 eq_data_->is_linear=
false;
545 double bc_pressure = eq_fields_->bc_switch_pressure(p_bdr);
546 double side_flux = -eq_fields_->bc_flux(p_bdr) * b_ele.
measure() * cross_section_;
548 eq_data_->save_local_system_[bulk_local_idx_] = (bool) eq_data_->bc_switch_dirichlet[b_ele.
idx()];
552 if (eq_data_->save_local_system_[bulk_local_idx_] || eq_data_->force_no_neumann_bc ) {
554 eq_data_->loc_constraint_.emplace_back(bulk_local_idx_, sidx_, bc_pressure);
557 eq_data_->loc_system_[bulk_local_idx_].add_value(edge_row_, side_flux);
561 eq_data_->is_linear=
false;
563 double bc_pressure = eq_fields_->bc_pressure(p_bdr);
564 double bc_switch_pressure = eq_fields_->bc_switch_pressure(p_bdr);
565 double bc_flux = -eq_fields_->bc_flux(p_bdr);
566 double bc_sigma = eq_fields_->bc_robin_sigma(p_bdr);
570 double solution_head = eq_data_->p_edge_solution.get(loc_dof_vec[sidx_]);
573 if (solution_head > bc_switch_pressure || eq_data_->force_no_neumann_bc) {
576 eq_data_->loc_system_[bulk_local_idx_].add_value(edge_row_, edge_row_,
577 -b_ele.
measure() * bc_sigma * cross_section_,
578 b_ele.
measure() * cross_section_ * (bc_flux - bc_sigma * bc_pressure) );
582 double bc_total_flux = bc_flux + bc_sigma*(bc_switch_pressure - bc_pressure);
584 eq_data_->loc_system_[bulk_local_idx_].add_value(edge_row_, bc_total_flux * b_ele.
measure() * cross_section_);
588 THROW( ExcBCNotSupported() );
591 eq_data_->balance_->add_flux_values(eq_data_->water_balance_idx, cell_side,
592 {eq_data_->loc_system_[bulk_local_idx_].row_dofs[eq_data_->loc_side_dofs[dim-1][sidx_]]},
598 unsigned int size, loc_size, elm_dim;;
603 size = (elm_dim+1) + 1 + (elm_dim+1);
608 dofs.head(dh_cell.n_dofs()) = dh_cell.get_loc_dof_indices();
615 for (
DHCellSide neighb_side : dh_cell.neighb_sides() ) {
626 const unsigned int p = size+i;
628 const unsigned int t = dh_neighb_cell.
n_dofs() - dh_neighb_cr_cell.
n_dofs() + neighb_side.side().side_idx();
633 eq_data_->loc_system_[dh_cell.local_idx()].reset(dofs, dofs);
636 for (
DHCellSide dh_side : dh_cell.side_range()) {
637 unsigned int sidx = dh_side.side_idx();
639 eq_data_->loc_system_[dh_cell.local_idx()].add_value(eq_data_->loc_side_dofs[elm_dim-1][sidx], eq_data_->loc_edge_dofs[elm_dim-1][sidx], 1.0);
640 eq_data_->loc_system_[dh_cell.local_idx()].add_value(eq_data_->loc_edge_dofs[elm_dim-1][sidx], eq_data_->loc_side_dofs[elm_dim-1][sidx], 1.0);
668 const unsigned int tt = dh_cr_cell.
n_dofs()+i;
675 loc_schur_.reset(dofs_schur, dofs_schur);
681 for (
uint n_i=0; n_i<dh_cell.
elm()->n_neighs_vb(); ++n_i) {
683 if ( (side->elem_idx() == neighb_side.
elem_idx()) && (side->side_idx() == neighb_side.
side_idx()) )
return n_i;
698 auto ls = eq_data_->seepage_bc_systems.find(dh_cell.
elm_idx());
699 if (ls != eq_data_->seepage_bc_systems.end())
703 arma::vec schur_solution = eq_data_->p_edge_solution.get_subvec(loc_dof_vec);
705 ls->second.reconstruct_solution_schur(eq_data_->schur_offset_[dim-1], schur_solution, reconstructed_solution_);
707 unsigned int pos_in_cache = this->element_cache_map_->position_in_cache(dh_cell.
elm_idx());
708 auto p = *( this->bulk_points(pos_in_cache).begin() );
709 postprocess_velocity_darcy(dh_cell, p, reconstructed_solution_);
711 eq_data_->bc_fluxes_reconstruted[bulk_local_idx_] =
true;
714 return reconstructed_solution_;
726 * eq_fields_->cross_section(p)
729 edge_source_term_ = edge_scale_ *
730 ( eq_fields_->water_source_density(p)
731 + eq_fields_->extra_source(p));
738 postprocess_velocity(dh_cell, p);
744 for (
unsigned int i=0; i<dh_cell.
elm()->n_sides(); i++) {
746 if( ! eq_data_->use_steady_assembly_)
748 double new_pressure = eq_data_->p_edge_solution.get(loc_dof_vec[i]);
749 double old_pressure = eq_data_->p_edge_solution_previous_time.get(loc_dof_vec[i]);
750 time_term_ = edge_scale_ * (eq_fields_->storativity(p) + eq_fields_->extra_storativity(p))
751 / eq_data_->time_step_ * (new_pressure - old_pressure);
753 solution[eq_data_->loc_side_dofs[dim-1][i]] += edge_source_term_ - time_term_;
770 shared_ptr<FiniteElement<dim>>
fe_;
780 unsigned int sidx_, side_row_, edge_row_;
788 template <
template<
IntDim...>
class DimAssembly>
793 template <
unsigned int dim>
800 static constexpr
const char *
name() {
return "ReconstructSchurAssemblyLMH"; }
809 ASSERT_EQ(cell.
dim(), dim).error(
"Dimension of element mismatch!");
812 auto p = *( this->
bulk_points(element_patch_idx).begin() );
unsigned int n_neighs_vb() const
Return number of neighbours.
void precompute_boundary_side(DHCellSide &cell_side, BoundaryPoint &p_side, BulkPoint &p_bdr)
Precompute values used in boundary side integral on given DHCellSide.
std::vector< bool > save_local_system_
Flag for saving the local system. Currently used only in case of seepage BC.
arma::Col< IntIdx > LocDofVec
LocDofVec l_indices_
Vector of pre-computed local DOF indices.
Raviart-Thomas element of order 0.
Edge edge() const
Returns pointer to the edge connected to the side.
void postprocess_velocity_darcy(const DHCellAccessor &dh_cell, BulkPoint &p, arma::vec &solution)
Postprocess velocity during loading of local system and after calculating of cell integral.
void begin() override
Implements AssemblyBase::begin.
BCField< 3, FieldValue< 3 >::Scalar > bc_flux
DarcyLMH::EqFields EqFields
EqFields * eq_fields_
Data objects shared with Flow equation.
FE_RT0< dim > fe_rt_
Assembly volume integrals.
Field< 3, FieldValue< 3 >::Scalar > conductivity
void boundary_side_integral_in(DHCellSide cell_side, const ElementAccessor< 3 > &b_ele, BulkPoint &p_bdr)
Definitions of Raviart-Thomas finite elements.
void reset()
Reset data members.
const FEValuesViews::Vector< spacedim > & vector_view(unsigned int i) const
Accessor to vector values of multicomponent FE.
EqFields * eq_fields_
Data objects shared with DarcyFlow.
unsigned int n_points() const
Returns the number of quadrature points.
Side side() const
Return Side of given cell and side_idx.
std::shared_ptr< DOFHandlerMultiDim > dh_
full DOF handler represents DOFs of sides, elements and edges
void end() override
Implements AssemblyBase::end.
VectorMPI p_edge_solution
Quadrature * quad_low_
Quadrature used in assembling methods (dim-1).
arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
FieldSet used_fields_
Sub field set contains fields used in calculation.
double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
void ghost_to_local_begin()
ghost_to_local_{begin,end} updates the local values by adding ghost values from neighbouring processo...
unsigned int n_sides() const
Returns number of sides aligned with the edge.
std::map< LongIdx, LocalSystem > seepage_bc_systems
Directing class of FieldValueCache.
Class FEValues calculates finite element data on the actual cells such as shape function values,...
void set_loc_schur(const DHCellAccessor dh_cr_cell)
Precompute loc_schur data member of given cell.
bool is_own() const
Return true if accessor represents own element (false for ghost element)
virtual ~ReadInitCondAssemblyLMH()
Destructor.
Field< 3, FieldValue< 3 >::Tensor > anisotropy
void end() override
Implements AssemblyBase::end.
@ update_values
Shape function values.
Base point accessor class.
void set_solution(uint loc_dof, double solution, double diag=0.0)
Set the position and value of known solution. E.g. Dirichlet boundary condition.
void postprocess_velocity(const DHCellAccessor &dh_cell, BulkPoint &p)
std::array< unsigned int, 3 > schur_offset_
Index offset in the local system for the Schur complement (of dim = 1,2,3).
void initialize(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags _flags)
Initialize structures and calculates cell-independent data.
#define THROW(whole_exception_expr)
Wrapper for throw. Saves the throwing point.
void dimjoin_intergral(DHCellAccessor cell_lower_dim, DHCellSide neighb_side)
DarcyLMH::EqFields::BC_Type type_
Type of boundary condition.
@ update_quadrature_points
Transformed quadrature points.
std::shared_ptr< LinSys > lin_sys_schur
void boundary_side_integral(DHCellSide cell_side)
unsigned int dim() const
Return dimension of element appropriate to cell.
@ update_normal_vectors
Normal vectors.
unsigned int dim() const
Return dimension of element appropriate to the side.
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Assemble integral over element.
FieldSet used_fields_
Sub field set contains fields used in calculation.
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
void end_mh_matrix()
Common code of end method of MH matrix assembly (Darcy and Richards)
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
Lumped mixed-hybrid model of linear Darcy flow, possibly unsteady.
const Element * element() const
void begin_reconstruct_schur()
Common code of begin method of Reconstruct Schur assembly (Darcy and Richards)
Class FESystem for compound finite elements.
void asm_source_term_darcy(const DHCellAccessor &cell, BulkPoint &p)
Part of cell_integral method, specialized in Darcy equation.
Field< 3, FieldValue< 3 >::Scalar > cross_section
BCField< 3, FieldValue< 3 >::Scalar > bc_switch_pressure
Definitions of basic Lagrangean finite elements with polynomial shape functions.
Side accessor allows to iterate over sides of DOF handler cell.
BCField< 3, FieldValue< 3 >::Scalar > bc_robin_sigma
ReadInitCondAssemblyLMH(EqFields *eq_fields, EqData *eq_data)
Constructor.
unsigned int n_dofs() const
Returns the number of shape functions.
std::array< std::vector< unsigned int >, 3 > loc_edge_dofs
DECLARE_EXCEPTION(ExcBCNotSupported,<< "BC type not supported.\n")
DHCellAccessor cell_with_other_dh(const DOFHandlerMultiDim *dh) const
Create new accessor with same local idx and given DOF handler. Actual and given DOF handler must be c...
const DHCellAccessor & cell() const
Return DHCellAccessor appropriate to the side.
RangeConvert< DHNeighbSide, DHCellSide > neighb_sides() const
Returns range of neighbour cell of lower dimension corresponding to cell of higher dimension.
static constexpr const char * name()
#define ASSERT_PERMANENT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
std::vector< arma::vec > postprocess_solution_
Definitions of particular quadrature rules on simplices.
void boundary_side_integral(FMT_UNUSED DHCellSide cell_side)
Assembles between boundary element and corresponding side on bulk element.
unsigned int local_idx() const
Return local index to element (index of DOF handler).
void begin() override
Implements AssemblyBase::begin.
void end() override
Implements AssemblyBase::end.
Assembly explicit Schur complement for the given linear system. Provides method for resolution of the...
std::vector< LocalSystem > loc_system_
Following data members are stored in vectors, one item for every cell.
double cross_section_
Precomputed cross-section value.
void dimjoin_intergral(FMT_UNUSED DHCellAccessor cell_lower_dim, FMT_UNUSED DHCellSide neighb_side)
Assembles the fluxes between elements of different dimensions.
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
static constexpr const char * name()
shared_ptr< FiniteElement< dim > > fe_
Finite element for the solution of the advection-diffusion equation.
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
double source_term_
Variables used in compute lumped source.
arma::vec get_subvec(const LocDofVec &loc_indices)
BCField< 3, FieldValue< 3 >::Scalar > bc_pressure
DarcyLMH::EqFields EqFields
double time_term_rhs_
Variables used in compute time term (unsteady)
void ghost_to_local_end()
ghost_to_local_{begin,end} updates the local values by adding ghost values from neighbouring processo...
void asm_element()
Part of cell_integral method, common in all descendants.
void local_to_ghost_begin()
local_to_ghost_{begin,end} updates the ghost values on neighbouring processors from local values
Container for various descendants of FieldCommonBase.
unsigned int n_dofs() const
Return number of dofs on given cell.
ElementAccessor< 3 > element_accessor()
unsigned int elm_idx() const
Return serial idx to element of loc_ele_idx_.
unsigned int nonlinear_iteration_
void add(unsigned int pos, double val)
Add value to item on given position.
Field< 3, FieldValue< 3 >::Scalar > extra_source
Externally added storativity.
Symmetric Gauss-Legendre quadrature formulae on simplices.
bool use_steady_assembly_
DarcyLMH::EqFields EqFields
unsigned int elem_idx() const
Range< BulkPoint > bulk_points(unsigned int element_patch_idx) const
Return BulkPoint range of appropriate dimension.
arma::vec reconstructed_solution_
Vector for reconstructed solution (velocity and pressure on element) from Schur complement.
Field< 3, FieldValue< 3 >::Scalar > sigma
MHMatrixAssemblyLMH(EqFields *eq_fields, EqData *eq_data)
Constructor.
FEValues< 3 > fe_values_side_
FEValues of object (of P disc finite element type)
unsigned int n_sides() const
Point accessor allow iterate over quadrature points of given side defined in local element coordinate...
std::array< std::vector< unsigned int >, 3 > loc_side_dofs
void eliminate_solution()
Cell accessor allow iterate over DOF handler cells.
void set_subvec(const LocDofVec &loc_indices, const arma::vec &values)
void begin_mh_matrix()
Common code of begin method of MH matrix assembly (Darcy and Richards)
virtual ~MHMatrixAssemblyLMH()
Destructor.
std::vector< unsigned int > fe_dofs(unsigned int fe_index)
Return dof indices belonging to given sub-FE.
unsigned int ngh_idx(DHCellAccessor &dh_cell, DHCellSide &neighb_side)
Temporary method find neighbour index in higher-dim cell.
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
Integral over element.
bool force_no_neumann_bc
auxiliary flag for switchting Dirichlet like BC
std::shared_ptr< SubDOFHandlerMultiDim > dh_cr_
DOF handler represents DOFs of edges.
int active_integrals_
Holds mask of active integrals.
arma::vec load_local_system(const DHCellAccessor &dh_cell)
arma::vec3 nv_
Normal vector.
Field< 3, FieldValue< 3 >::Scalar > init_pressure
@ update_JxW_values
Transformed quadrature weights.
double edge_source_term_
Precomputed values in postprocess_velocity.
void local_to_ghost_end()
local_to_ghost_{begin,end} updates the ghost values on neighbouring processors from local values
void use_dirichlet_switch(DHCellSide &cell_side, const ElementAccessor< 3 > &b_ele, BulkPoint &p_bdr)
Update BC switch dirichlet in MH matrix assembly if BC type is seepage.
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
unsigned int bulk_local_idx_
Local idx of bulk element.
Field< 3, FieldValue< 3 >::Scalar > water_source_density
unsigned int idx() const
We need this method after replacing Region by RegionIdx, and movinf RegionDB instance into particular...
double ngh_value_
Precomputed ngh value.
double get(unsigned int pos) const
Return value on given position.
std::vector< LocalConstraint > loc_constraint_
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
std::shared_ptr< Balance > balance_
Shared Balance object.
Field< 3, FieldValue< 3 >::Scalar > storativity
static constexpr const char * name()
BCField< 3, FieldValue< 3 >::Enum > bc_type
Range< CouplingPoint > coupling_points(const DHCellSide &cell_side) const
Return CouplingPoint range of appropriate dimension.
unsigned int side_idx() const
void set_dofs()
Precompute loc_system and loc_schur data members.
void initialize(ElementCacheMap *element_cache_map)
Initialize auxiliary vectors and other data members.
SideIter side(const unsigned int loc_index)
Field< 3, FieldValue< 3 >::Scalar > extra_storativity
std::array< unsigned int, 3 > loc_ele_dof
void end_reconstruct_schur()
Common code of end method of Reconstruct Schur assembly (Darcy and Richards)
Generic class of assemblation.
Compound finite element on dim dimensional simplex.
ReconstructSchurAssemblyLMH(EqFields *eq_fields, EqData *eq_data)
Constructor.
Range< BoundaryPoint > boundary_points(const DHCellSide &cell_side) const
Return BoundaryPoint range of appropriate dimension.
void cell_integral(DHCellAccessor cell, unsigned int element_patch_idx)
ElementCacheMap * element_cache_map_
ElementCacheMap shared with GenericAssembly object.
double measure() const
Computes the measure of the element.
VectorMPI p_edge_solution_previous_time
void asm_sides(const DHCellAccessor &cell, BulkPoint &p, double conductivity)
Part of cell_integral method, common in all descendants.
unsigned int IntDim
A dimension index type.