28 #include "petscviewer.h" 29 #include "petscerror.h" 77 const
it::Selection &
DarcyLMH::get_mh_mortar_selection() {
79 .
add_value(NoMortar,
"None",
"No Mortar method is applied.")
80 .
add_value(MortarP0,
"P0",
"Mortar space: P0 on elements of lower dimension.")
81 .
add_value(MortarP1,
"P1",
"Mortar space: P1 on intersections, using non-conforming pressures.")
91 "Boundary piezometric head for BC types: dirichlet, robin, and river." )
93 "Boundary switch piezometric head for BC types: seepage, river." )
95 "Initial condition for the pressure given as the piezometric head." )
97 return field_descriptor;
104 "Linear solver for MH problem.")
106 "Residual tolerance.")
108 "Minimum number of iterations (linear solutions) to use.\nThis is usefull if the convergence criteria " 109 "does not characterize your goal well enough so it converges prematurely, possibly even without a single linear solution." 110 "If greater then 'max_it' the value is set to 'max_it'.")
112 "Maximum number of iterations (linear solutions) of the non-linear solver.")
114 "If a stagnation of the nonlinear solver is detected the solver stops. " 115 "A divergence is reported by default, forcing the end of the simulation. By setting this flag to 'true', the solver " 116 "ends with convergence success on stagnation, but it reports warning about it.")
121 return it::Record(
"Flow_Darcy_LMH",
"Lumped Mixed-Hybrid solver for saturated Darcy flow.")
125 "Vector of the gravity force. Dimensionless.")
127 "Input data for Darcy flow model.")
129 "Non-linear solver for MH problem.")
131 "Output stream settings.\n Specify file format, precision etc.")
134 IT::Default(
"{ \"fields\": [ \"pressure_p0\", \"velocity_p0\" ] }"),
135 "Specification of output fields and output times.")
137 "Output settings specific to Darcy flow model.\n" 138 "Includes raw output and some experimental functionality.")
140 "Settings for computing mass balance.")
142 "Method for coupling Darcy flow between dimensions on incompatible meshes. [Experimental]" )
148 Input::register_class< DarcyLMH, Mesh &, const Input::Record >(
"Flow_Darcy_LMH") +
193 MessageOut() <<
"Duplicate key 'time', time in flow equation is already initialized from parent class!";
200 data_ = make_shared<EqData>();
203 data_->is_linear=
true;
232 gravity_array.copy_to(gvec);
235 data_->gravity_vec_ =
data_->gravity_.subvec(0,2);
237 data_->bc_pressure.add_factory(
239 (
data_->gravity_,
"bc_piezo_head") );
240 data_->bc_switch_pressure.add_factory(
242 (
data_->gravity_,
"bc_switch_piezo_head") );
243 data_->init_pressure.add_factory(
245 (
data_->gravity_,
"init_piezo_head") );
253 MessageOut() <<
"Missing the key 'time', obligatory for the transient problems." << endl;
264 std::shared_ptr< FiniteElement<1> > fe1_rt = std::make_shared<FE_RT0_disc<1>>();
265 std::shared_ptr< FiniteElement<2> > fe2_rt = std::make_shared<FE_RT0_disc<2>>();
266 std::shared_ptr< FiniteElement<3> > fe3_rt = std::make_shared<FE_RT0_disc<3>>();
267 std::shared_ptr< FiniteElement<0> > fe0_disc = std::make_shared<FE_P_disc<0>>(0);
268 std::shared_ptr< FiniteElement<1> > fe1_disc = std::make_shared<FE_P_disc<1>>(0);
269 std::shared_ptr< FiniteElement<2> > fe2_disc = std::make_shared<FE_P_disc<2>>(0);
270 std::shared_ptr< FiniteElement<3> > fe3_disc = std::make_shared<FE_P_disc<3>>(0);
271 std::shared_ptr< FiniteElement<0> > fe0_cr = std::make_shared<FE_CR<0>>();
272 std::shared_ptr< FiniteElement<1> > fe1_cr = std::make_shared<FE_CR<1>>();
273 std::shared_ptr< FiniteElement<2> > fe2_cr = std::make_shared<FE_CR<2>>();
274 std::shared_ptr< FiniteElement<3> > fe3_cr = std::make_shared<FE_CR<3>>();
276 FESystem<0> fe0_sys( {fe0_disc, fe0_disc, fe0_cr} );
282 std::shared_ptr<DiscreteSpace> ds = std::make_shared<EqualOrderDiscreteSpace>(
mesh_, fe_sys);
283 data_->dh_ = std::make_shared<DOFHandlerMultiDim>(*mesh_);
284 data_->dh_->distribute_dofs(ds);
291 uint rt_component = 0;
292 auto ele_flux_ptr = create_field_fe<3, FieldValue<3>::VectorFixed>(
data_->dh_, rt_component);
293 data_->full_solution = ele_flux_ptr->vec();
296 auto ele_velocity_ptr = std::make_shared< FieldDivide<3, FieldValue<3>::VectorFixed> >(ele_flux_ptr,
data_->cross_section);
299 uint p_ele_component = 0;
300 auto ele_pressure_ptr = create_field_fe<3, FieldValue<3>::Scalar>(
data_->dh_, p_ele_component, &
data_->full_solution);
303 uint p_edge_component = 1;
304 auto edge_pressure_ptr = create_field_fe<3, FieldValue<3>::Scalar>(
data_->dh_, p_edge_component, &
data_->full_solution);
307 arma::vec4 gravity = (-1) *
data_->gravity_;
313 uint p_element_component = 1;
314 data_->dh_p_ = std::make_shared<SubDOFHandlerMultiDim>(
data_->dh_,p_element_component);
318 uint p_edge_component = 2;
319 data_->dh_cr_ = std::make_shared<SubDOFHandlerMultiDim>(
data_->dh_,p_edge_component);
324 std::shared_ptr<DiscreteSpace> ds_cr_disc = std::make_shared<EqualOrderDiscreteSpace>(
mesh_, fe_cr_disc);
325 data_->dh_cr_disc_ = std::make_shared<DOFHandlerMultiDim>(*mesh_);
326 data_->dh_cr_disc_->distribute_dofs(ds_cr_disc);
333 data_->p_edge_solution =
data_->dh_cr_->create_vector();
334 data_->p_edge_solution_previous =
data_->dh_cr_->create_vector();
335 data_->p_edge_solution_previous_time =
data_->dh_cr_->create_vector();
344 .val<Input::AbstractRecord>(
"linear_solver");
356 data_->water_balance_idx =
balance_->add_quantity(
"water_volume");
365 data_->multidim_assembler = AssemblyBase::create< AssemblyLMH >(
data_);
397 DebugOut().fmt(
"Read initial condition\n");
401 LocDofVec p_indices = dh_cell.cell_with_other_dh(
data_->dh_p_.get()).get_loc_dof_indices();
403 LocDofVec l_indices = dh_cell.cell_with_other_dh(
data_->dh_cr_.get()).get_loc_dof_indices();
407 double init_value =
data_->init_pressure.value(ele.
centre(),ele);
408 unsigned int p_idx =
data_->dh_p_->parent_indices()[p_indices[0]];
409 data_->full_solution[p_idx] = init_value;
411 for (
unsigned int i=0; i<ele->
n_sides(); i++) {
413 unsigned int l_idx =
data_->dh_cr_->parent_indices()[l_indices[i]];
414 data_->full_solution[l_idx] += init_value/n_sides_of_edge;
416 data_->p_edge_solution[l_indices[i]] += init_value/n_sides_of_edge;
420 data_->full_solution.ghost_to_local_begin();
421 data_->full_solution.ghost_to_local_end();
423 data_->p_edge_solution.ghost_to_local_begin();
424 data_->p_edge_solution.ghost_to_local_end();
425 data_->p_edge_solution_previous_time.copy_from(
data_->p_edge_solution);
447 data_->p_edge_solution.zero_entries();
449 if (
data_->use_steady_assembly_) {
468 data_->full_solution.copy_from(temp);
490 data_->full_solution.local_to_ghost_begin();
491 data_->full_solution.local_to_ghost_end();
499 bool jump_time =
data_->storativity.is_jump_time();
500 if (! zero_time_term_from_left) {
505 data_->use_steady_assembly_ =
false;
511 WarningOut() <<
"Output of solution discontinuous in time not supported yet.\n";
520 if (! zero_time_term_from_left && ! jump_time && output)
527 if (zero_time_term_from_right) {
529 data_->use_steady_assembly_ =
true;
534 }
else if (! zero_time_term_from_left && jump_time) {
535 WarningOut() <<
"Discontinuous time term not supported yet.\n";
546 return (
data_->storativity.input_list_size() == 0);
559 MessageOut().fmt(
"[nonlinear solver] norm of initial residual: {}\n", residual_norm);
562 int is_linear_common;
567 this->
max_n_it_ = nl_solver_rec.
val<
unsigned int>(
"max_it");
568 this->
min_n_it_ = nl_solver_rec.
val<
unsigned int>(
"min_it");
569 if (this->min_n_it_ > this->max_n_it_) this->min_n_it_ = this->
max_n_it_;
571 if (! is_linear_common) {
577 while (nonlinear_iteration_ < this->min_n_it_ ||
578 (residual_norm > this->tolerance_ && nonlinear_iteration_ < this->max_n_it_ )) {
580 convergence_history.push_back(residual_norm);
584 if (convergence_history.size() >= 5 &&
585 convergence_history[ convergence_history.size() - 1]/convergence_history[ convergence_history.size() - 2] > 0.9 &&
586 convergence_history[ convergence_history.size() - 1]/convergence_history[ convergence_history.size() - 5] > 0.8) {
592 THROW(ExcSolverDiverge() << EI_Reason(
"Stagnation."));
596 if (! is_linear_common){
597 data_->p_edge_solution_previous.copy_from(
data_->p_edge_solution);
598 data_->p_edge_solution_previous.local_to_ghost_begin();
599 data_->p_edge_solution_previous.local_to_ghost_end();
603 MessageOut().fmt(
"[schur solver] lin. it: {}, reason: {}, residual: {}\n",
609 if (is_linear_common){
612 MessageOut().fmt(
"[nonlinear solver] lin. it: {}, reason: {}, residual: {}\n",
619 VecAXPBY(
data_->p_edge_solution.petsc_vec(), (1-alpha), alpha,
data_->p_edge_solution_previous.petsc_vec());
626 MessageOut().fmt(
"[nonlinear solver] it: {} lin. it: {}, reason: {}, residual: {}\n",
646 data_->p_edge_solution_previous_time.copy_from(
data_->p_edge_solution);
647 data_->p_edge_solution_previous_time.local_to_ghost_begin();
648 data_->p_edge_solution_previous_time.local_to_ghost_end();
662 balance_->calculate_cumulative(
data_->water_balance_idx,
data_->full_solution.petsc_vec());
663 balance_->calculate_instant(
data_->water_balance_idx,
data_->full_solution.petsc_vec());
670 return data_->lin_sys_schur->get_solution_precision();
682 START_TIMER(
"DarcyLMH::assembly_steady_mh_matrix");
695 unsigned int dim = dh_cell.dim();
696 assembler[dim-1]->assemble(dh_cell);
712 double zeros[100000];
713 for(
int i=0; i<100000; i++) zeros[i] = 0.0;
716 tmp_rows.reserve(200);
719 dofs.reserve(
data_->dh_cr_->max_elem_dofs());
720 dofs_ngh.reserve(
data_->dh_cr_->max_elem_dofs());
726 const uint ndofs = dh_cell.n_dofs();
727 dofs.resize(dh_cell.n_dofs());
728 dh_cell.get_dof_indices(dofs);
730 int* dofs_ptr = dofs.
data();
737 for (
DHCellSide neighb_side : dh_cell.neighb_sides() ) {
745 const uint ndofs_ngh = dh_neighb_cell.
n_dofs();
746 dofs_ngh.resize(ndofs_ngh);
750 tmp_rows.push_back(dofs_ngh[neighb_side.side().side_idx()]);
760 for(
auto &isec : isec_list ) {
764 const uint ndofs_slave = dh_cell_slave.
n_dofs();
765 dofs_ngh.resize(ndofs_slave);
769 for(
unsigned int i_side=0; i_side < dh_cell_slave.
elm()->
n_sides(); i_side++) {
770 tmp_rows.push_back( dofs_ngh[i_side] );
910 data_->lin_sys_schur = std::make_shared<LinSys_PETSC>( &(*
data_->dh_cr_->distr()) );
983 data_->full_solution.zero_entries();
984 data_->p_edge_solution.zero_entries();
988 xprintf(
Err,
"Unknown solver type. Internal error.\n");
999 START_TIMER(
"DarcyFlowMH::reconstruct_solution_from_schur");
1001 data_->full_solution.zero_entries();
1002 data_->p_edge_solution.local_to_ghost_begin();
1003 data_->p_edge_solution.local_to_ghost_end();
1010 unsigned int dim = dh_cell.dim();
1011 assembler[dim-1]->assemble_reconstruct(dh_cell);
1014 data_->full_solution.local_to_ghost_begin();
1015 data_->full_solution.local_to_ghost_end();
1023 START_TIMER(
"DarcyFlowMH::assembly_linear_system");
1026 data_->p_edge_solution.local_to_ghost_begin();
1027 data_->p_edge_solution.local_to_ghost_end();
1029 data_->is_linear=
true;
1061 std::string output_file;
1064 double d_max = std::numeric_limits<double>::max();
1065 double h1 = d_max, h2 = d_max, h3 = d_max;
1066 double he2 = d_max, he3 = d_max;
1069 case 1: h1 = std::min(h1,ele.measure());
break;
1070 case 2: h2 = std::min(h2,ele.measure());
break;
1071 case 3: h3 = std::min(h3,ele.measure());
break;
1074 for (
unsigned int j=0; j<ele->n_sides(); j++) {
1076 case 2: he2 = std::min(he2, ele.side(j)->measure());
break;
1077 case 3: he3 = std::min(he3, ele.side(j)->measure());
break;
1081 if(h1 == d_max) h1 = 0;
1082 if(h2 == d_max) h2 = 0;
1083 if(h3 == d_max) h3 = 0;
1084 if(he2 == d_max) he2 = 0;
1085 if(he3 == d_max) he3 = 0;
1088 file = fopen(output_file.c_str(),
"a");
1089 fprintf(file,
"nA = %d;\n",
data_->dh_cr_disc_->distr()->size());
1090 fprintf(file,
"nB = %d;\n",
data_->dh_->mesh()->get_el_ds()->size());
1091 fprintf(file,
"nBF = %d;\n",
data_->dh_cr_->distr()->size());
1092 fprintf(file,
"h1 = %e;\nh2 = %e;\nh3 = %e;\n", h1, h2, h3);
1093 fprintf(file,
"he2 = %e;\nhe3 = %e;\n", he2, he3);
1099 PetscViewerASCIIOpen(PETSC_COMM_WORLD, output_file.c_str(), &viewer);
1100 PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
1101 MatView( *const_cast<Mat*>(
lin_sys_schur().get_matrix()), viewer);
1102 VecView( *const_cast<Vec*>(
lin_sys_schur().get_rhs()), viewer);
1103 VecView( *const_cast<Vec*>(&(
lin_sys_schur().get_solution())), viewer);
1104 VecView( *const_cast<Vec*>(&(
data_->full_solution.petsc_vec())), viewer);
1303 if(
time_ !=
nullptr)
1311 unsigned int i, n_dofs, min, max;
1315 n_dofs = dh_cell.get_dof_indices(dof_indices);
1317 for (i=min; i<max; ++i) dof_vec.push_back(dof_indices[i]);
unsigned int n_sides() const
Returns number of sides aligned with the edge.
virtual bool zero_time_term(bool time_global=false)
static const Input::Type::Record & get_input_type()
Main balance input record type.
Output class for darcy_flow_mh model.
RegionSet get_region_set(const std::string &set_name) const
arma::Col< IntIdx > LocDofVec
void set_symmetric(bool flag=true)
Classes with algorithms for computation of intersections of meshes.
static const int registrar
Registrar of class to factory.
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
virtual void set_from_input(const Input::Record in_rec)
MixedMeshIntersections & mixed_intersections()
static const Input::Type::Record & get_input_type()
The specification of output stream.
Common base for intersection object.
unsigned int get_dof_indices(std::vector< LongIdx > &indices) const
Fill vector of the global indices of dofs associated to the cell.
Edge edge() const
Returns pointer to the edge connected to the side.
static const Input::Type::Record & get_input_type()
virtual void initialize_specific()
#define MessageOut()
Macro defining 'message' record of log.
virtual void start_add_assembly()
void assembly_mh_matrix(MultidimAssembly &assembler)
virtual PetscErrorCode mat_zero_entries()
Wrappers for linear systems based on MPIAIJ and MATIS format.
bool is_end() const
Returns true if the actual time is greater than or equal to the end time.
void next_time()
Proceed to the next time according to current estimated time step.
void zero_time_step() override
void create_linear_system(Input::AbstractRecord rec)
void solve_time_step(bool output=true)
Solve the problem without moving to next time and without output.
Lumped mixed-hybrid model of linear Darcy flow, possibly unsteady.
friend class DarcyFlowMHOutput
std::vector< int > get_component_indices_vec(unsigned int component) const
Get vector of all DOF indices of given component (0..side, 1..element, 2..edge)
virtual void start_allocation()
Cell accessor allow iterate over DOF handler cells.
virtual void finish_assembly()=0
std::vector< std::vector< ILpair > > element_intersections_
static const Input::Type::Instance & get_input_type(FieldSet &eq_data, const std::string &equation_name)
SideIter side(const unsigned int loc_index)
const RegionDB & region_db() const
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
const TimeStep & step(int index=-1) const
static const Input::Type::Selection & get_mh_mortar_selection()
Selection for enum MortarMethod.
Mixed-hybrid model of linear Darcy flow, possibly unsteady.
static const std::string field_descriptor_record_description(const string &record_name)
virtual void initial_condition_postprocess()
Basic time management functionality for unsteady (and steady) solvers (class Equation).
unsigned int bulk_ele_idx() const
Returns index of bulk element.
unsigned int n_dofs() const
Return number of dofs on given cell.
Basic time management class.
virtual void set_tolerances(double r_tol, double a_tol, unsigned int max_it)=0
DarcyFlowMHOutput * output_object
void view(const char *name="") const
static Input::Type::Record & record_template()
Template Record with common keys for derived equations.
Assembly explicit Schur complement for the given linear system. Provides method for resolution of the...
void initialize() override
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
void print_matlab_matrix(string matlab_file)
Print darcy flow matrix in matlab format into a file.
void copy_from(VectorMPI &other)
virtual void postprocess()
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
void read_initial_condition()
static const Input::Type::Instance & get_input_type_specific()
Compound finite element on dim dimensional simplex.
void allocate_mh_matrix()
FMT_FUNC int fprintf(std::ostream &os, CStringRef format, ArgList args)
Definitions of basic Lagrangean finite elements with polynomial shape functions.
unsigned int n_sides() const
virtual void output_data() override
Write computed fields.
#define START_TIMER(tag)
Starts a timer with specified tag.
static Input::Type::Abstract & get_input_type()
virtual Range< ElementAccessor< 3 > > elements_range() const
Returns range of bulk elements.
virtual ~DarcyLMH() override
unsigned int n_sides() const
std::shared_ptr< EqData > data_
unsigned int nonlinear_iteration_
virtual PetscErrorCode rhs_zero_entries()
virtual void accept_time_step()
postprocess velocity field (add sources)
void set_solution(Vec sol_vec)
Dedicated class for storing path to input and output files.
Support classes for parallel programing.
#define MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm)
int set_upper_constraint(double upper, std::string message)
Sets upper constraint for the next time step estimating.
void solve_nonlinear()
Solve method common to zero_time_step and update solution.
void update_solution() override
virtual unsigned int n_elements(bool boundary=false) const
Returns count of boundary or bulk elements.
LinSys & lin_sys_schur()
Getter for the linear system of the 2. Schur complement.
Input::Record input_record_
std::shared_ptr< Balance > balance_
virtual void assembly_linear_system()
#define WarningOut()
Macro defining 'warning' record of log.
virtual SolveInfo solve()=0
virtual double solution_precision() const
#define END_TIMER(tag)
Ends a timer with specified tag.
#define OLD_ASSERT_EQUAL(a, b)
void set_matrix_changed()
static const Input::Type::Record & get_input_type()
void dofs_range(unsigned int n_dofs, unsigned int &min, unsigned int &max, unsigned int component)
Helper method fills range (min and max) of given component.
void set_positive_definite(bool flag=true)
static const Input::Type::Record & type_field_descriptor()
static Input::Type::Abstract & get_input_type()
virtual void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)=0
unsigned int n_edges() const
Class for representation SI units of Fields.
void reconstruct_solution_from_schur(MultidimAssembly &assembler)
MortarMethod
Type of experimental Mortar-like method for non-compatible 1d-2d interaction.
unsigned int n_neighs_vb() const
Return number of neighbours.
#define DebugOut()
Macro defining 'debug' record of log.
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
virtual double compute_residual()=0
#define THROW(whole_exception_expr)
Wrapper for throw. Saves the throwing point.
Implementation of range helper class.
Side accessor allows to iterate over sides of DOF handler cell.
#define FLOW123D_FORCE_LINK_IN_CHILD(x)
DarcyLMH(Mesh &mesh, const Input::Record in_rec, TimeGovernor *tm=nullptr)
CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL.
void output()
Calculate values for output.
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Mixed-hybrid model of linear Darcy flow, possibly unsteady.