28 #include "petscviewer.h" 29 #include "petscerror.h" 77 const
it::Selection &
DarcyMH::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.")
89 "Homogeneous Neumann boundary condition\n(zero normal flux over the boundary).")
91 "Dirichlet boundary condition. " 92 "Specify the pressure head through the ``bc_pressure`` field " 93 "or the piezometric head through the ``bc_piezo_head`` field.")
94 .
add_value(
total_flux,
"total_flux",
"Flux boundary condition (combines Neumann and Robin type). " 95 "Water inflow equal to (($ \\delta_d(q_d^N + \\sigma_d (h_d^R - h_d) )$)). " 96 "Specify the water inflow by the ``bc_flux`` field, the transition coefficient by ``bc_robin_sigma`` " 97 "and the reference pressure head or piezometric head through ``bc_pressure`` or ``bc_piezo_head`` respectively.")
99 "Seepage face boundary condition. Pressure and inflow bounded from above. Boundary with potential seepage flow " 100 "is described by the pair of inequalities: " 101 "(($h_d \\le h_d^D$)) and (($ -\\boldsymbol q_d\\cdot\\boldsymbol n \\le \\delta q_d^N$)), where the equality holds in at least one of them. " 102 "Caution: setting (($q_d^N$)) strictly negative " 103 "may lead to an ill posed problem since a positive outflow is enforced. " 104 "Parameters (($h_d^D$)) and (($q_d^N$)) are given by the fields ``bc_switch_pressure`` (or ``bc_switch_piezo_head``) and ``bc_flux`` respectively." 107 "River boundary condition. For the water level above the bedrock, (($H_d > H_d^S$)), the Robin boundary condition is used with the inflow given by: " 108 "(( $ \\delta_d(q_d^N + \\sigma_d(H_d^D - H_d) )$)). For the water level under the bedrock, constant infiltration is used: " 109 "(( $ \\delta_d(q_d^N + \\sigma_d(H_d^D - H_d^S) )$)). Parameters: ``bc_pressure``, ``bc_switch_pressure``, " 110 " ``bc_sigma``, ``bc_flux``." 121 "Boundary piezometric head for BC types: dirichlet, robin, and river." )
123 "Boundary switch piezometric head for BC types: seepage, river." )
125 "Initial condition for the pressure given as the piezometric head." )
127 return field_descriptor;
134 "Linear solver for MH problem.")
136 "Residual tolerance.")
138 "Minimum number of iterations (linear solutions) to use.\nThis is usefull if the convergence criteria " 139 "does not characterize your goal well enough so it converges prematurely, possibly even without a single linear solution." 140 "If greater then 'max_it' the value is set to 'max_it'.")
142 "Maximum number of iterations (linear solutions) of the non-linear solver.")
144 "If a stagnation of the nonlinear solver is detected the solver stops. " 145 "A divergence is reported by default, forcing the end of the simulation. By setting this flag to 'true', the solver " 146 "ends with convergence success on stagnation, but it reports warning about it.")
151 return it::Record(
"Flow_Darcy_MH",
"Mixed-Hybrid solver for saturated Darcy flow.")
155 "Vector of the gravity force. Dimensionless.")
157 "Input data for Darcy flow model.")
159 "Non-linear solver for MH problem.")
161 "Output stream settings.\n Specify file format, precision etc.")
164 IT::Default(
"{ \"fields\": [ \"pressure_p0\", \"velocity_p0\" ] }"),
165 "Specification of output fields and output times.")
167 "Output settings specific to Darcy flow model.\n" 168 "Includes raw output and some experimental functionality.")
170 "Settings for computing mass balance.")
172 "Number of Schur complements to perform when solving MH system.")
174 "Method for coupling Darcy flow between dimensions on incompatible meshes. [Experimental]" )
180 Input::register_class< DarcyMH, Mesh &, const Input::Record >(
"Flow_Darcy_MH") +
192 .
description(
"Pressure solution - P0 interpolation.");
197 .
description(
"Pressure solution - Crouzeix-Raviart interpolation.");
202 .
description(
"Piezo head solution - P0 interpolation.");
207 .
description(
"Velocity solution - P0 interpolation.");
215 .
description(
"Anisotropy of the conductivity tensor.")
220 .
description(
"Complement dimension parameter (cross section for 1D, thickness for 2D).")
231 .
description(
"Transition coefficient between dimensions.")
249 .description(
"Prescribed pressure value on the boundary. Used for all values of ``bc_type`` except ``none`` and ``seepage``. " 250 "See documentation of ``bc_type`` for exact meaning of ``bc_pressure`` in individual boundary condition types.")
251 .input_default(
"0.0")
257 .description(
"Incoming water boundary flux. Used for bc_types : ``total_flux``, ``seepage``, ``river``.")
258 .input_default(
"0.0")
259 .units(
UnitSI().m().s(-1) );
263 .name(
"bc_robin_sigma")
264 .description(
"Conductivity coefficient in the ``total_flux`` or the ``river`` boundary condition type.")
265 .input_default(
"0.0")
270 .name(
"bc_switch_pressure")
271 .description(
"Critical switch pressure for ``seepage`` and ``river`` boundary conditions.")
272 .input_default(
"0.0")
278 .
description(
"Initial condition for pressure in time dependent problems.")
283 .
description(
"Storativity (in time dependent problems).")
288 .
description(
"Storativity added from upstream equation.")
294 .
description(
"Water source density added from upstream equation.")
341 MessageOut() <<
"Duplicate key 'time', time in flow equation is already initialized from parent class!";
348 data_ = make_shared<EqData>();
351 data_->is_linear=
true;
381 gravity_array.copy_to(gvec);
384 data_->gravity_vec_ =
data_->gravity_.subvec(0,2);
386 data_->bc_pressure.add_factory(
388 (
data_->gravity_,
"bc_piezo_head") );
389 data_->bc_switch_pressure.add_factory(
391 (
data_->gravity_,
"bc_switch_piezo_head") );
392 data_->init_pressure.add_factory(
394 (
data_->gravity_,
"init_piezo_head") );
402 MessageOut() <<
"Missing the key 'time', obligatory for the transient problems." << endl;
415 std::shared_ptr< FiniteElement<1> > fe1_rt = std::make_shared<FE_RT0_disc<1>>();
416 std::shared_ptr< FiniteElement<2> > fe2_rt = std::make_shared<FE_RT0_disc<2>>();
417 std::shared_ptr< FiniteElement<3> > fe3_rt = std::make_shared<FE_RT0_disc<3>>();
418 std::shared_ptr< FiniteElement<0> > fe0_disc = std::make_shared<FE_P_disc<0>>(0);
419 std::shared_ptr< FiniteElement<1> > fe1_disc = std::make_shared<FE_P_disc<1>>(0);
420 std::shared_ptr< FiniteElement<2> > fe2_disc = std::make_shared<FE_P_disc<2>>(0);
421 std::shared_ptr< FiniteElement<3> > fe3_disc = std::make_shared<FE_P_disc<3>>(0);
422 std::shared_ptr< FiniteElement<0> > fe0_cr = std::make_shared<FE_CR<0>>();
423 std::shared_ptr< FiniteElement<1> > fe1_cr = std::make_shared<FE_CR<1>>();
424 std::shared_ptr< FiniteElement<2> > fe2_cr = std::make_shared<FE_CR<2>>();
425 std::shared_ptr< FiniteElement<3> > fe3_cr = std::make_shared<FE_CR<3>>();
427 FESystem<0> fe0_sys( {fe0_disc, fe0_disc, fe0_cr} );
433 std::shared_ptr<DiscreteSpace>
ds = std::make_shared<EqualOrderDiscreteSpace>(
mesh_, fe_sys);
434 data_->dh_ = std::make_shared<DOFHandlerMultiDim>(*mesh_);
435 data_->dh_->distribute_dofs(ds);
439 this->
data_->multidim_assembler = AssemblyBase::create< AssemblyMH >(
data_);
443 uint rt_component = 0;
444 auto ele_flux_ptr = create_field_fe<3, FieldValue<3>::VectorFixed>(
data_->dh_, rt_component);
445 data_->full_solution = ele_flux_ptr->vec();
448 auto ele_velocity_ptr = std::make_shared< FieldDivide<3, FieldValue<3>::VectorFixed> >(ele_flux_ptr,
data_->cross_section);
451 uint p_ele_component = 0;
452 auto ele_pressure_ptr = create_field_fe<3, FieldValue<3>::Scalar>(
data_->dh_, p_ele_component, &
data_->full_solution);
455 uint p_edge_component = 1;
456 auto edge_pressure_ptr = create_field_fe<3, FieldValue<3>::Scalar>(
data_->dh_, p_edge_component, &
data_->full_solution);
459 arma::vec4 gravity = (-1) *
data_->gravity_;
465 uint p_edge_component = 2;
466 data_->dh_cr_ = std::make_shared<SubDOFHandlerMultiDim>(
data_->dh_,p_edge_component);
471 std::shared_ptr<DiscreteSpace> ds_cr_disc = std::make_shared<EqualOrderDiscreteSpace>(
mesh_, fe_cr_disc);
472 data_->dh_cr_disc_ = std::make_shared<DOFHandlerMultiDim>(*mesh_);
473 data_->dh_cr_disc_->distribute_dofs(ds_cr_disc);
483 .val<Input::AbstractRecord>(
"linear_solver");
493 VecCreateMPI(PETSC_COMM_WORLD,
data_->dh_->distr()->lsize(),PETSC_DETERMINE,&(
steady_diagonal));
502 data_->water_balance_idx =
balance_->add_quantity(
"water_volume");
532 if (zero_time_term_from_right) {
563 data_->full_solution.local_to_ghost_begin();
564 data_->full_solution.local_to_ghost_end();
574 bool jump_time =
data_->storativity.is_jump_time();
575 if (! zero_time_term_from_left) {
585 WarningOut() <<
"Output of solution discontinuous in time not supported yet.\n";
594 if (! zero_time_term_from_left && ! jump_time && output)
output_data();
600 if (zero_time_term_from_right) {
605 }
else if (! zero_time_term_from_left && jump_time) {
606 WarningOut() <<
"Discontinuous time term not supported yet.\n";
618 return (
data_->storativity.input_list_size() == 0);
631 MessageOut().fmt(
"[nonlinear solver] norm of initial residual: {}\n", residual_norm);
634 int is_linear_common;
639 this->
max_n_it_ = nl_solver_rec.
val<
unsigned int>(
"max_it");
640 this->
min_n_it_ = nl_solver_rec.
val<
unsigned int>(
"min_it");
641 if (this->min_n_it_ > this->max_n_it_) this->min_n_it_ = this->
max_n_it_;
643 if (! is_linear_common) {
651 while (nonlinear_iteration_ < this->min_n_it_ ||
652 (residual_norm > this->tolerance_ && nonlinear_iteration_ < this->max_n_it_ )) {
654 convergence_history.push_back(residual_norm);
657 if (convergence_history.size() >= 5 &&
658 convergence_history[ convergence_history.size() - 1]/convergence_history[ convergence_history.size() - 2] > 0.9 &&
659 convergence_history[ convergence_history.size() - 1]/convergence_history[ convergence_history.size() - 5] > 0.8) {
665 THROW(ExcSolverDiverge() << EI_Reason(
"Stagnation."));
669 if (! is_linear_common)
675 if (is_linear_common){
678 MessageOut().fmt(
"[nonlinear solver] lin. it: {}, reason: {}, residual: {}\n",
679 si.n_iterations, si.converged_reason, residual_norm);
692 MessageOut().fmt(
"[nonlinear solver] it: {} lin. it: {}, reason: {}, residual: {}\n",
695 chkerr(VecDestroy(&save_solution));
719 if(
data_->mortar_method_ != MortarMethod::NoMortar){
720 auto multidim_assembler = AssemblyBase::create< AssemblyMH >(
data_);
722 unsigned int dim = dh_cell.dim();
723 multidim_assembler[dim-1]->fix_velocity(dh_cell);
790 unsigned int dim = dh_cell.dim();
791 assembler[dim-1]->assemble(dh_cell);
809 double zeros[100000];
810 for(
int i=0; i<100000; i++) zeros[i] = 0.0;
813 tmp_rows.reserve(200);
816 dofs.reserve(
data_->dh_->max_elem_dofs());
817 dofs_ngh.reserve(
data_->dh_->max_elem_dofs());
822 const uint ndofs = dh_cell.n_dofs();
824 dh_cell.get_dof_indices(dofs);
833 for (
DHCellSide neighb_side : dh_cell.neighb_sides() ) {
841 const uint ndofs_ngh = dh_neighb_cell.
n_dofs();
842 dofs_ngh.resize(ndofs_ngh);
848 const unsigned int t = dh_neighb_cell.
n_dofs() - (dh_neighb_cell.
dim()+1) + neighb_side.side().side_idx();
849 tmp_rows.push_back(dofs_ngh[t]);
863 for(
auto &isec : isec_list ) {
867 const uint ndofs_slave = dh_cell_slave.
n_dofs();
868 dofs_ngh.resize(ndofs_slave);
872 for(
unsigned int i_side=0; i_side < dh_cell_slave.
elm()->
n_sides(); i_side++) {
874 tmp_rows.push_back( dofs_ngh[(ndofs_slave+1)/2+i_side] );
884 edge_rows = dofs.
data() + nsides +1;
923 const uint ndofs = dh_cell.n_dofs();
924 global_dofs.resize(ndofs);
925 dh_cell.get_dof_indices(global_dofs);
927 double cs =
data_->cross_section.value(ele.
centre(), ele);
930 double source = ele.
measure() * cs *
931 (
data_->water_source_density.value(ele.
centre(), ele)
937 {dh_cell.get_loc_dof_indices()[ndofs/2]}, {0}, {source});
957 #ifdef FLOW123D_HAVE_BDDCML 958 WarningOut() <<
"For BDDC no Schur complements are used.";
971 xprintf(
Err,
"Flow123d was not build with BDDCML support.\n");
972 #endif // FLOW123D_HAVE_BDDCML 977 WarningOut() <<
"Invalid number of Schur Complements. Using 2.";
990 ls->LinSys::set_from_input(in_rec);
999 ISCreateGeneral(PETSC_COMM_SELF, side_dofs_vec.size(), &(side_dofs_vec[0]), PETSC_COPY_VALUES, &is);
1014 const PetscInt *b_indices;
1015 ISGetIndices(ls->
IsB, &b_indices);
1017 for(
uint i_b=0, i_bb=0; i_b < b_size && i_bb < elem_dofs_vec.size(); i_b++) {
1018 if (b_indices[i_b] == elem_dofs_vec[i_bb])
1019 elem_dofs_vec[i_bb++] = i_b + ds->
begin();
1021 ISRestoreIndices(ls->
IsB, &b_indices);
1024 ISCreateGeneral(PETSC_COMM_SELF, elem_dofs_vec.size(), &(elem_dofs_vec[0]), PETSC_COPY_VALUES, &is);
1052 xprintf(
Err,
"Unknown solver type. Internal error.\n");
1064 START_TIMER(
"DarcyFlowMH_Steady::assembly_linear_system");
1066 data_->is_linear=
true;
1120 std::string output_file;
1132 PetscViewerASCIIOpen(PETSC_COMM_WORLD, output_file.c_str(), &viewer);
1133 PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
1145 double d_max = std::numeric_limits<double>::max();
1146 double h1 = d_max, h2 = d_max, h3 = d_max;
1147 double he2 = d_max, he3 = d_max;
1150 case 1: h1 = std::min(h1,ele.measure());
break;
1151 case 2: h2 = std::min(h2,ele.measure());
break;
1152 case 3: h3 = std::min(h3,ele.measure());
break;
1155 for (
unsigned int j=0; j<ele->n_sides(); j++) {
1157 case 2: he2 = std::min(he2, ele.side(j)->measure());
break;
1158 case 3: he3 = std::min(he3, ele.side(j)->measure());
break;
1162 if(h1 == d_max) h1 = 0;
1163 if(h2 == d_max) h2 = 0;
1164 if(h3 == d_max) h3 = 0;
1165 if(he2 == d_max) he2 = 0;
1166 if(he3 == d_max) he3 = 0;
1169 file = fopen(output_file.c_str(),
"a");
1170 fprintf(file,
"nA = %d;\n",
data_->dh_cr_disc_->distr()->size());
1171 fprintf(file,
"nB = %d;\n",
data_->dh_->mesh()->get_el_ds()->size());
1172 fprintf(file,
"nBF = %d;\n",
data_->dh_cr_->distr()->size());
1173 fprintf(file,
"h1 = %e;\nh2 = %e;\nh3 = %e;\n", h1, h2, h3);
1174 fprintf(file,
"he2 = %e;\nhe3 = %e;\n", he2, he3);
1192 START_TIMER(
"DarcyFlowMH_Steady::set_mesh_data_for_bddc");
1220 dh_cell.get_dof_indices(cell_dofs_global);
1222 inet.insert(inet.end(), cell_dofs_global.begin(), cell_dofs_global.end());
1223 uint n_inet = cell_dofs_global.size();
1226 uint dim = dh_cell.elm().dim();
1227 elDimMax = std::max( elDimMax, dim );
1228 elDimMin = std::min( elDimMin, dim );
1232 isegn.push_back( dh_cell.elm_idx());
1235 for (
unsigned int si=0; si<dh_cell.elm()->n_sides(); si++) {
1236 arma::vec3 coord = dh_cell.elm().side(si)->centre();
1239 localDofMap.insert( std::make_pair( cell_dofs_global[si], coord ) );
1241 localDofMap.insert( std::make_pair( cell_dofs_global[si+dim+2], coord ) );
1244 arma::vec3 elm_centre = dh_cell.elm().centre();
1245 localDofMap.insert( std::make_pair( cell_dofs_global[dim+1], elm_centre ) );
1249 for(
DHCellSide side : dh_cell.neighb_sides()) {
1250 uint neigh_dim = side.cell().elm().dim();
1251 side.cell().get_dof_indices(cell_dofs_global);
1252 int edge_row = cell_dofs_global[neigh_dim+2+side.side_idx()];
1253 localDofMap.insert( std::make_pair( edge_row, side.centre() ) );
1254 inet.push_back( edge_row );
1257 nnet.push_back(n_inet);
1262 double conduct =
data_->conductivity.value( elm_centre , dh_cell.elm() );
1263 auto aniso =
data_->anisotropy.value( elm_centre , dh_cell.elm() );
1267 for (
int i = 0; i < 3; i++) {
1268 coef = coef + aniso.at(i,i);
1271 coef = conduct*coef / 3;
1274 "Zero coefficient of hydrodynamic resistance %f . \n ", coef );
1275 element_permeability.push_back( 1. / coef );
1285 auto distr =
data_->dh_->distr();
1294 int numNodeSub = localDofMap.size();
1305 for ( ; itB != localDofMap.end(); ++itB ) {
1306 isngn[ind] = itB -> first;
1309 for (
int j = 0; j < 3; j++ ) {
1310 xyz[ j*numNodeSub + ind ] = coord[j];
1315 localDofMap.clear();
1323 Global2LocalMap_ global2LocalNodeMap;
1324 for (
unsigned ind = 0; ind < isngn.size(); ++ind ) {
1325 global2LocalNodeMap.insert( std::make_pair( static_cast<unsigned>( isngn[ind] ), ind ) );
1330 for (
unsigned int iEle = 0; iEle < isegn.size(); iEle++ ) {
1331 int nne = nnet[ iEle ];
1332 for (
int ien = 0; ien < nne; ien++ ) {
1334 int indGlob = inet[indInet];
1336 Global2LocalMap_::iterator pos = global2LocalNodeMap.find( indGlob );
1337 OLD_ASSERT( pos != global2LocalNodeMap.end(),
1338 "Cannot remap node index %d to local indices. \n ", indGlob );
1339 int indLoc =
static_cast<int> ( pos -> second );
1342 inet[ indInet++ ] = indLoc;
1346 int numNodes =
size;
1347 int numDofsInt =
size;
1349 int meshDim = elDimMax;
1361 bddc_ls -> load_mesh( LinSys_BDDC::BDDCMatrixType::SPD_VIA_SYMMETRICGENERAL, spaceDim, numNodes, numDofsInt, inet, nnet, nndf, isegn, isngn, isngn, xyz, element_permeability, meshDim );
1383 if(
time_ !=
nullptr)
1436 const IntIdx p_ele_dof = dh_cell.get_loc_dof_indices()[dh_cell.n_dofs()/2];
1437 local_sol[p_ele_dof] =
data_->init_pressure.value(ele.
centre(),ele);
1448 PetscScalar *local_diagonal;
1456 dofs.reserve(
data_->dh_->max_elem_dofs());
1460 dofs.resize(dh_cell.n_dofs());
1461 dh_cell.get_dof_indices(dofs);
1464 const uint p_ele_dof = dh_cell.n_dofs() / 2;
1466 double diagonal_coeff =
data_->cross_section.value(ele.
centre(), ele)
1468 +
data_->extra_storativity.value(ele.
centre(), ele))
1470 local_diagonal[dh_cell.get_loc_dof_indices()[p_ele_dof]]= - diagonal_coeff /
time_->
dt();
1472 balance_->add_mass_values(
data_->water_balance_idx, dh_cell,
1473 {dh_cell.get_loc_dof_indices()[p_ele_dof]},
1474 {diagonal_coeff}, 0.0);
1488 if (scale_factor != 1.0) {
1509 void dofs_range(
unsigned int n_dofs,
unsigned int &min,
unsigned int &max,
unsigned int component) {
1513 }
else if (component==1) {
1525 unsigned int i, n_dofs, min, max;
1529 n_dofs = dh_cell.get_dof_indices(dof_indices);
1531 for (i=min; i<max; ++i) dof_vec.push_back(dof_indices[i]);
Field< 3, FieldValue< 3 >::VectorFixed > field_ele_velocity
Field< 3, FieldValue< 3 >::Scalar > extra_source
Externally added storativity.
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)
static const Input::Type::Record & get_input_type()
Main balance input record type.
Output class for darcy_flow_mh model.
Field< 3, FieldValue< 3 >::Scalar > init_pressure
RegionSet get_region_set(const std::string &set_name) const
virtual void initialize_specific()
void set_complement(LinSys_PETSC *ls)
Set complement LinSys object.
SchurComplement SchurComplement
void assembly_mh_matrix(MultidimAssembly &assembler)
void set_symmetric(bool flag=true)
Classes with algorithms for computation of intersections of meshes.
FieldCommon & input_selection(Input::Type::Selection element_selection)
Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construc...
virtual ~DarcyMH() override
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.
void set_from_input(const Input::Record in_rec) override
virtual void assembly_source_term()
Source term is implemented differently in LMH version.
friend class DarcyFlowMHOutput
#define MessageOut()
Macro defining 'message' record of log.
virtual void start_add_assembly()
virtual void output_data() override
Write computed fields.
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.
static const int registrar
Registrar of class to factory.
Field< 3, FieldValue< 3 >::Scalar > extra_storativity
Field< 3, FieldValue< 3 >::Scalar > sigma
unsigned int dim() const
Return dimension of element appropriate to cell.
static const Input::Type::Selection & get_mh_mortar_selection()
Selection for enum MortarMethod.
void initialize() override
void set_negative_definite(bool flag=true)
static const Input::Type::Record & type_field_descriptor()
virtual void start_allocation()
Cell accessor allow iterate over DOF handler cells.
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
virtual void finish_assembly()=0
virtual double get_solution_precision()=0
std::vector< std::vector< ILpair > > element_intersections_
static const Input::Type::Instance & get_input_type(FieldSet &eq_data, const std::string &equation_name)
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 std::string field_descriptor_record_description(const string &record_name)
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.
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
void solve_nonlinear()
Solve method common to zero_time_step and update solution.
static const Input::Type::Record & get_input_type()
std::shared_ptr< EqData > data_
Basic time management class.
std::shared_ptr< DiscreteSpace > ds
virtual void set_tolerances(double r_tol, double a_tol, unsigned int max_it)=0
Field< 3, FieldValue< 3 >::Scalar > field_ele_pressure
Externally added water source.
void view(const char *name="") const
static Input::Type::Record & record_template()
Template Record with common keys for derived equations.
unsigned int nonlinear_iteration_
void set_mesh_data_for_bddc(LinSys_BDDC *bddc_ls)
Assembly explicit Schur complement for the given linear system. Provides method for resolution of the...
BCField< 3, FieldValue< 3 >::Scalar > bc_flux
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
void update_solution() override
void solve_time_step(bool output=true)
Solve the problem without moving to next time and without output.
void set_from_input(const Input::Record in_rec) override
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
BCField< 3, FieldValue< 3 >::Scalar > bc_pressure
virtual void postprocess()
postprocess velocity field (add sources)
static constexpr Mask equation_result
Match result fields. These are never given by input or copy of input.
Field< 3, FieldValue< 3 >::Scalar > water_source_density
static const Input::Type::Instance & get_input_type_specific()
double * get_solution_array()
BCField< 3, FieldValue< 3 >::Scalar > bc_switch_pressure
Field< 3, FieldValue< 3 >::Scalar > storativity
static constexpr Mask input_copy
Compound finite element on dim dimensional simplex.
std::shared_ptr< Balance > balance_
FieldCommon & input_default(const string &input_default)
const Vec & get_solution()
FMT_FUNC int fprintf(std::ostream &os, CStringRef format, ArgList args)
Definitions of basic Lagrangean finite elements with polynomial shape functions.
unsigned int begin(int proc) const
get starting local index
unsigned int n_sides() const
bool is_changed_dt() const
#define START_TIMER(tag)
Starts a timer with specified tag.
void set_from_input(const Input::Record in_rec) override
static Input::Type::Abstract & get_input_type()
virtual Range< ElementAccessor< 3 > > elements_range() const
Returns range of bulk elements.
bool use_steady_assembly_
double measure() const
Computes the measure of the element.
unsigned int n_sides() const
void allocate_mh_matrix()
static const Input::Type::Selection & get_bc_type_selection()
Return a Selection corresponding to enum BC_Type.
virtual void assembly_linear_system()
virtual const Vec * get_rhs()
virtual PetscErrorCode rhs_zero_entries()
void set_solution(Vec sol_vec)
Field< 3, FieldValue< 3 >::VectorFixed > flux
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.
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Field< 3, FieldValue< 3 >::Scalar > field_ele_piezo_head
FieldCommon & description(const string &description)
void create_linear_system(Input::AbstractRecord rec)
virtual unsigned int n_elements(bool boundary=false) const
Returns count of boundary or bulk elements.
Distribution * make_complement_distribution()
get distribution of complement object if complement is defined
Input::Record input_record_
auto disable_where(const Field< spacedim, typename FieldValue< spacedim >::Enum > &control_field, const vector< FieldEnum > &value_list) -> Field &
Field< 3, FieldValue< 3 >::TensorFixed > anisotropy
DarcyMH(Mesh &mesh, const Input::Record in_rec, TimeGovernor *tm=nullptr)
CREATE AND FILL GLOBAL MH MATRIX OF THE WATER MODEL.
unsigned int bulk_idx() const
Returns index of the region in the bulk set.
#define WarningOut()
Macro defining 'warning' record of log.
FieldCommon & name(const string &name)
virtual SolveInfo solve()=0
#define END_TIMER(tag)
Ends a timer with specified tag.
#define OLD_ASSERT_EQUAL(a, b)
EqData()
Creation of all fields.
void set_matrix_changed()
static const Input::Type::Record & get_input_type()
void print_matlab_matrix(string matlab_file)
Print darcy flow matrix in matlab format into a file.
mixed-hybrid model of linear Darcy flow, possibly unsteady.
unsigned int index() const
Field< 3, FieldValue< 3 >::Scalar > field_edge_pressure
void zero_time_step() override
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)
virtual void read_initial_condition()
FieldCommon & set_limits(double min, double max=std::numeric_limits< double >::max())
Field< 3, FieldValue< 3 >::Scalar > cross_section
static Input::Type::Abstract & get_input_type()
virtual void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)=0
BCField< 3, FieldValue< 3 >::Enum > bc_type
DarcyFlowMHOutput * output_object
Field< 3, FieldValue< 3 >::Scalar > conductivity
FieldCommon & flags(FieldFlag::Flags::Mask mask)
unsigned int n_edges() const
void prepare_new_time_step()
Class for representation SI units of Fields.
virtual const Mat * get_matrix()
MortarMethod
Type of experimental Mortar-like method for non-compatible 1d-2d interaction.
virtual double solution_precision() const
BCField< 3, FieldValue< 3 >::Scalar > bc_robin_sigma
unsigned int n_neighs_vb() const
Return number of neighbours.
virtual void setup_time_term()
static UnitSI & dimensionless()
Returns dimensionless unit.
#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.
MortarMethod mortar_method_
#define FLOW123D_FORCE_LINK_IN_CHILD(x)
void rhs_set_value(int row, double val)
Input::Type::Record make_field_descriptor_type(const std::string &equation_name) const
static const Input::Type::Record & get_input_type()
Solver based on Multilevel BDDC - using corresponding class of OpenFTL package.
virtual bool zero_time_term(bool time_global=false)
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.