8 #ifndef SRC_ASSEMBLY_LMH_HH_ 9 #define SRC_ASSEMBLY_LMH_HH_ 46 velocity_interpolation_quad_(dim, 0),
54 auto fe = ad_->dh_->ds()->fe(ad_->dh_->own_range().begin()->elm())[
Dim<dim>{}];
56 loc_side_dofs = fe_system->
fe_dofs(0);
57 loc_ele_dof = fe_system->
fe_dofs(1)[0];
58 loc_edge_dofs = fe_system->
fe_dofs(2);
61 .error(
"Mortar methods are not supported in Lumped Mixed-Hybrid Method.");
70 schur_offset_ = loc_edge_dofs[0];
71 reconstructed_solution_.zeros(schur_offset_);
96 arma::vec schur_solution = ad_->p_edge_solution.get_subvec(loc_schur_.row_dofs);
98 loc_system_.reconstruct_solution_schur(schur_offset_, schur_solution, reconstructed_solution_);
103 ad_->full_solution.set_subvec(loc_system_.row_dofs.head(schur_offset_), reconstructed_solution_);
104 ad_->full_solution.set_subvec(loc_system_.row_dofs.tail(loc_schur_.row_dofs.n_elem), schur_solution);
110 save_local_system_ =
false;
111 bc_fluxes_reconstruted =
false;
115 loc_system_.compute_schur_complement(schur_offset_, loc_schur_,
true);
119 loc_schur_.eliminate_solution();
120 ad_->lin_sys_schur->set_local_system(loc_schur_, ad_->dh_cr_->get_local_to_global_map());
156 if(bc_fluxes_reconstruted)
160 auto ls = ad_->seepage_bc_systems.find(dh_cell.
elm_idx());
161 if (ls != ad_->seepage_bc_systems.end())
163 arma::vec schur_solution = ad_->p_edge_solution.get_subvec(loc_schur_.row_dofs);
165 ls->second.reconstruct_solution_schur(schur_offset_, schur_solution, reconstructed_solution_);
169 bc_fluxes_reconstruted =
true;
179 if(save_local_system_)
180 ad_->seepage_bc_systems[dh_cell.
elm_idx()] = loc_system_;
214 const unsigned int p =
size()+i;
216 const unsigned int t = dh_neighb_cell.
n_dofs() - dh_neighb_cr_cell.
n_dofs() + neighb_side.side().side_idx();
220 const unsigned int tt = dh_cr_cell.
n_dofs()+i;
222 [neighb_side.side().side_idx()];
226 loc_system_.reset(dofs, dofs);
227 loc_schur_.reset(dofs_schur, dofs_schur);
234 dirichlet_edge.resize(ele->
n_sides());
236 unsigned int sidx = dh_side.side_idx();
237 dirichlet_edge[sidx] = 0;
240 if (dh_side.side().is_boundary()) {
241 double cross_section = ad_->cross_section.value(ele.
centre(), ele);
244 ad_->balance->add_flux_values(ad_->water_balance_idx, dh_side,
245 {loc_system_.row_dofs[loc_side_dofs[sidx]]},
250 loc_system_.add_value(loc_side_dofs[sidx], loc_edge_dofs[sidx], 1.0);
251 loc_system_.add_value(loc_edge_dofs[sidx], loc_side_dofs[sidx], 1.0);
258 const unsigned int sidx = side.
side_idx();
259 const unsigned int side_row = loc_side_dofs[sidx];
260 const unsigned int edge_row = loc_edge_dofs[sidx];
268 double bc_pressure = ad_->bc_pressure.value(b_ele.
centre(), b_ele);
269 loc_schur_.set_solution(sidx, bc_pressure);
270 dirichlet_edge[sidx] = 1;
274 double bc_flux = -ad_->bc_flux.value(b_ele.
centre(), b_ele);
275 double bc_pressure = ad_->bc_pressure.value(b_ele.
centre(), b_ele);
276 double bc_sigma = ad_->bc_robin_sigma.value(b_ele.
centre(), b_ele);
278 dirichlet_edge[sidx] = 2;
279 loc_system_.add_value(edge_row, edge_row,
280 -b_ele.
measure() * bc_sigma * cross_section,
281 (bc_flux - bc_sigma * bc_pressure) * b_ele.
measure() * cross_section);
284 ad_->is_linear=
false;
286 char & switch_dirichlet = ad_->bc_switch_dirichlet[b_ele.
idx()];
287 double bc_pressure = ad_->bc_switch_pressure.value(b_ele.
centre(), b_ele);
288 double bc_flux = -ad_->bc_flux.value(b_ele.
centre(), b_ele);
289 double side_flux = bc_flux * b_ele.
measure() * cross_section;
292 if(use_dirichlet_switch){
293 if (switch_dirichlet) {
300 double solution_flux = reconstructed_solution_[side_row];
302 if ( solution_flux < side_flux) {
317 double solution_head = ad_->p_edge_solution[loc_schur_.row_dofs[sidx]];
319 if ( solution_head > bc_pressure) {
326 save_local_system_ = (bool) switch_dirichlet;
330 if (switch_dirichlet || ad_->force_no_neumann_bc ) {
332 loc_schur_.set_solution(sidx, bc_pressure);
333 dirichlet_edge[sidx] = 1;
336 loc_system_.add_value(edge_row, side_flux);
340 ad_->is_linear=
false;
342 double bc_pressure = ad_->bc_pressure.value(b_ele.
centre(), b_ele);
343 double bc_switch_pressure = ad_->bc_switch_pressure.value(b_ele.
centre(), b_ele);
344 double bc_flux = -ad_->bc_flux.value(b_ele.
centre(), b_ele);
345 double bc_sigma = ad_->bc_robin_sigma.value(b_ele.
centre(), b_ele);
347 double solution_head = ad_->p_edge_solution[loc_schur_.row_dofs[sidx]];
350 if (solution_head > bc_switch_pressure || ad_->force_no_neumann_bc) {
353 loc_system_.add_value(edge_row, edge_row,
354 -b_ele.
measure() * bc_sigma * cross_section,
355 b_ele.
measure() * cross_section * (bc_flux - bc_sigma * bc_pressure) );
359 double bc_total_flux = bc_flux + bc_sigma*(bc_switch_pressure - bc_pressure);
361 loc_system_.add_value(edge_row, bc_total_flux * b_ele.
measure() * cross_section);
373 double cs = ad_->cross_section.value(ele.
centre(), ele);
374 double conduct = ad_->conductivity.value(ele.
centre(), ele);
375 double scale = 1 / cs /conduct;
383 auto ele = dh_cell.
elm();
390 for (
unsigned int k=0; k<qsize; k++)
391 for (
unsigned int i=0; i<ndofs; i++){
393 arma::dot(gravity_vec,velocity.value(i,k))
395 loc_system_.add_value(i, rhs_val);
397 for (
unsigned int j=0; j<ndofs; j++){
399 arma::dot(velocity.value(i,k),
400 (ad_->anisotropy.value(ele.centre(), ele)).i()
401 * velocity.value(j,k))
404 loc_system_.add_value(i, j, mat_val);
435 for(
unsigned int side = 0; side < loc_side_dofs.size(); side++){
436 loc_system_.add_value(loc_ele_dof, loc_side_dofs[side], -1.0);
437 loc_system_.add_value(loc_side_dofs[side], loc_ele_dof, -1.0);
452 double alpha = 1.0 / ele->
n_sides();
453 double cross_section = ad_->cross_section.value(ele.
centre(), ele);
454 double coef = alpha * ele.
measure() * cross_section;
456 double source = ad_->water_source_density.value(ele.
centre(), ele)
457 + ad_->extra_source.value(ele.
centre(), ele);
458 double source_term = coef * source;
461 double storativity = 0.0;
462 double time_term_diag = 0.0, time_term = 0.0, time_term_rhs = 0.0;
464 if(! ad_->use_steady_assembly_)
466 storativity = ad_->storativity.value(ele.
centre(), ele)
467 + ad_->extra_storativity.value(ele.
centre(), ele);
468 time_term = coef * storativity;
471 for (
unsigned int i=0; i<ele->
n_sides(); i++)
473 if(! ad_->use_steady_assembly_)
475 time_term_diag = time_term / ad_->time_step_;
476 time_term_rhs = time_term_diag * ad_->p_edge_solution_previous_time[loc_schur_.row_dofs[i]];
478 ad_->balance->add_mass_values(ad_->water_balance_idx, dh_cell,
479 {loc_system_.row_dofs[loc_edge_dofs[i]]}, {time_term}, 0);
485 ad_->balance->add_mass_values(ad_->water_balance_idx, dh_cell,
486 {loc_system_.row_dofs[loc_edge_dofs[i]]}, {0}, 0);
489 this->loc_system_.add_value(loc_edge_dofs[i], loc_edge_dofs[i],
491 -source_term - time_term_rhs);
493 ad_->balance->add_source_values(ad_->water_balance_idx, ele.
region().
bulk_idx(),
494 {loc_system_.row_dofs[loc_edge_dofs[i]]}, {0},{source_term});
512 unsigned int p =
size()+i;
515 ngh_values_.fe_side_values_.reinit(neighb_side.side());
516 nv = ngh_values_.fe_side_values_.normal_vector(0);
518 double value = ad_->sigma.value( ele.
centre(), ele) *
519 2*ad_->conductivity.value( ele.
centre(), ele) *
520 arma::dot(ad_->anisotropy.value( ele.
centre(), ele)*nv, nv) *
521 ad_->cross_section.value( neighb_side.centre(), ele_higher ) *
522 ad_->cross_section.value( neighb_side.centre(), ele_higher ) /
523 ad_->cross_section.value( ele.
centre(), ele ) *
524 neighb_side.measure();
526 loc_system_.add_value(loc_ele_dof, loc_ele_dof, -value);
527 loc_system_.add_value(loc_ele_dof, p, value);
528 loc_system_.add_value(p,loc_ele_dof, value);
529 loc_system_.add_value(p,p, -value);
546 double edge_scale = ele.
measure()
547 * ad_->cross_section.value(ele.
centre(), ele)
550 double edge_source_term = edge_scale *
551 ( ad_->water_source_density.value(ele.
centre(), ele)
552 + ad_->extra_source.value(ele.
centre(), ele));
558 double edge_scale,
double edge_source_term)
562 double storativity = ad_->storativity.value(ele.
centre(), ele)
563 + ad_->extra_storativity.value(ele.
centre(), ele);
564 double new_pressure, old_pressure, time_term = 0.0;
566 for (
unsigned int i=0; i<ele->
n_sides(); i++) {
568 if( ! ad_->use_steady_assembly_)
570 new_pressure = ad_->p_edge_solution[loc_schur_.row_dofs[i]];
571 old_pressure = ad_->p_edge_solution_previous_time[loc_schur_.row_dofs[i]];
572 time_term = edge_scale * storativity / ad_->time_step_ * (new_pressure - old_pressure);
574 solution[loc_side_dofs[i]] += edge_source_term - time_term;
586 QGauss velocity_interpolation_quad_;
590 AssemblyDataPtrLMH ad_;
602 unsigned int loc_ele_dof;
607 unsigned int schur_offset_;
614 bool save_local_system_;
617 bool bc_fluxes_reconstruted;
virtual void assemble_sides(const DHCellAccessor &dh_cell)
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
void assemble_local_system(const DHCellAccessor &dh_cell, bool use_dirichlet_switch)
Range< DHCellSide > side_range() const
Returns range of cell sides.
static unsigned int size()
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
const FEValuesViews::Vector< spacedim > & vector_view(unsigned int i) const
Accessor to vector values of multicomponent FE.
unsigned int elm_idx() const
Return serial idx to element of loc_ele_idx_.
arma::Col< IntIdx > LocDofVec
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
void save_local_system(const DHCellAccessor &dh_cell)
unsigned int dim() const
Return dimension of element appropriate to cell.
Lumped mixed-hybrid model of linear Darcy flow, possibly unsteady.
unsigned int n_dofs() const
Returns the number of shape functions.
Cell accessor allow iterate over DOF handler cells.
Class FEValues calculates finite element data on the actual cells such as shape function values...
Class FESystem for compound finite elements.
LocDofVec get_loc_dof_indices() const
Returns the local indices of dofs associated to the cell on the local process.
unsigned int side_idx() const
unsigned int n_dofs() const
Return number of dofs on given cell.
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
static constexpr bool value
Symmetric Gauss-Legendre quadrature formulae on simplices.
Assembly explicit Schur complement for the given linear system. Provides method for resolution of the...
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
void load_local_system(const DHCellAccessor &dh_cell)
std::shared_ptr< DarcyLMH::EqData > AssemblyDataPtrLMH
virtual void assemble_source_term(const DHCellAccessor &dh_cell)
Transformed quadrature points.
void assemble_side_bc(const DHCellSide &side, double cross_section, bool use_dirichlet_switch)
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
Compound finite element on dim dimensional simplex.
void assemble(const DHCellAccessor &dh_cell) override
Raviart-Thomas element of order 0.
Definitions of basic Lagrangean finite elements with polynomial shape functions.
std::vector< unsigned int > fe_dofs(unsigned int fe_index)
Return dof indices belonging to given sub-FE.
unsigned int n_sides() const
RangeConvert< DHNeighbSide, DHCellSide > neighb_sides() const
Returns range of neighbour cell of lower dimension corresponding to cell of higher dimension...
void initialize(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags _flags)
Initialize structures and calculates cell-independent data.
void assembly_dim_connections(const DHCellAccessor &dh_cell)
void assemble_element(FMT_UNUSED const DHCellAccessor &dh_cell)
double measure() const
Computes the measure of the element.
Shape function gradients.
const DHCellAccessor & cell() const
Return DHCellAccessor appropriate to the side.
virtual void postprocess_velocity_specific(const DHCellAccessor &dh_cell, arma::vec &solution, double edge_scale, double edge_source_term)
void fix_velocity(const DHCellAccessor &) override
void assemble_bc(const DHCellAccessor &dh_cell, bool use_dirichlet_switch)
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
Definitions of particular quadrature rules on simplices.
Calculates finite element data on the actual cell.
void set_dofs(const DHCellAccessor &dh_cell)
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...
double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point...
unsigned int n_neighs_vb() const
Return number of neighbours.
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
void postprocess_velocity(const DHCellAccessor &dh_cell, arma::vec &solution)
void assemble_reconstruct(const DHCellAccessor &dh_cell) override
unsigned int n_points() const
Returns the number of quadrature points.
Definitions of Raviart-Thomas finite elements.
Side accessor allows to iterate over sides of DOF handler cell.
void update_water_content(const DHCellAccessor &) override
Updates water content in Richards.
ElementAccessor< 3 > element_accessor()
#define FEAL_ASSERT(expr)
Definition of assert for debug and release mode.
Transformed quadrature weights.
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
void assemble_sides_scale(const DHCellAccessor &dh_cell, double scale)