28 #include "petscviewer.h"
29 #include "petscerror.h"
82 .
add_value(
MortarP1,
"P1",
"Mortar space: P1 on intersections, using non-conforming pressures.")
90 "Homogeneous Neumann boundary condition\n(zero normal flux over the boundary).")
92 "Dirichlet boundary condition. "
93 "Specify the pressure head through the ``bc_pressure`` field "
94 "or the piezometric head through the ``bc_piezo_head`` field.")
95 .
add_value(
total_flux,
"total_flux",
"Flux boundary condition (combines Neumann and Robin type). "
96 "Water inflow equal to (($ \\delta_d(q_d^N + \\sigma_d (h_d^R - h_d) )$)). "
97 "Specify the water inflow by the ``bc_flux`` field, the transition coefficient by ``bc_robin_sigma`` "
98 "and the reference pressure head or piezometric head through ``bc_pressure`` or ``bc_piezo_head`` respectively.")
100 "Seepage face boundary condition. Pressure and inflow bounded from above. Boundary with potential seepage flow "
101 "is described by the pair of inequalities: "
102 "(($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. "
103 "Caution: setting (($q_d^N$)) strictly negative "
104 "may lead to an ill posed problem since a positive outflow is enforced. "
105 "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."
108 "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: "
109 "(( $ \\delta_d(q_d^N + \\sigma_d(H_d^D - H_d) )$)). For the water level under the bedrock, constant infiltration is used: "
110 "(( $ \\delta_d(q_d^N + \\sigma_d(H_d^D - H_d^S) )$)). Parameters: ``bc_pressure``, ``bc_switch_pressure``, "
111 " ``bc_sigma``, ``bc_flux``."
122 "Boundary piezometric head for BC types: dirichlet, robin, and river." )
124 "Boundary switch piezometric head for BC types: seepage, river." )
126 "Initial condition for the pressure given as the piezometric head." )
128 return field_descriptor;
135 "Linear solver for MH problem.")
137 "Residual tolerance.")
139 "Minimum number of iterations (linear solutions) to use.\nThis is usefull if the convergence criteria "
140 "does not characterize your goal well enough so it converges prematurely, possibly even without a single linear solution."
141 "If greater then 'max_it' the value is set to 'max_it'.")
143 "Maximum number of iterations (linear solutions) of the non-linear solver.")
145 "If a stagnation of the nonlinear solver is detected the solver stops. "
146 "A divergence is reported by default, forcing the end of the simulation. By setting this flag to 'true', the solver "
147 "ends with convergence success on stagnation, but it reports warning about it.")
152 return it::Record(
"Flow_Darcy_MH",
"Mixed-Hybrid solver for saturated Darcy flow.")
156 "Vector of the gravity force. Dimensionless.")
158 "Input data for Darcy flow model.")
160 "Non-linear solver for MH problem.")
162 "Output stream settings.\n Specify file format, precision etc.")
165 IT::Default(
"{ \"fields\": [ \"pressure_p0\", \"velocity_p0\" ] }"),
166 "Specification of output fields and output times.")
168 "Output settings specific to Darcy flow model.\n"
169 "Includes raw output and some experimental functionality.")
171 "Settings for computing mass balance.")
173 "Number of Schur complements to perform when solving MH system.")
175 "Method for coupling Darcy flow between dimensions on incompatible meshes. [Experimental]" )
181 Input::register_class< DarcyMH, Mesh &, const Input::Record >(
"Flow_Darcy_MH") +
190 *
this += field_ele_pressure.name(
"pressure_p0")
193 .description(
"Pressure solution - P0 interpolation.");
195 *
this += field_edge_pressure.name(
"pressure_edge")
198 .description(
"Pressure solution - Crouzeix-Raviart interpolation.");
200 *
this += field_ele_piezo_head.name(
"piezo_head_p0")
203 .description(
"Piezo head solution - P0 interpolation.");
205 *
this += field_ele_velocity.name(
"velocity_p0")
206 .units(
UnitSI().m().s(-1))
208 .description(
"Velocity solution - P0 interpolation.");
210 *
this += flux.name(
"flux")
211 .units(
UnitSI().m().s(-1))
213 .description(
"Darcy flow flux.");
215 *
this += anisotropy.name(
"anisotropy")
216 .description(
"Anisotropy of the conductivity tensor.")
217 .input_default(
"1.0")
220 *
this += cross_section.name(
"cross_section")
221 .description(
"Complement dimension parameter (cross section for 1D, thickness for 2D).")
222 .input_default(
"1.0")
223 .units(
UnitSI().m(3).md() );
225 *
this += conductivity.name(
"conductivity")
226 .description(
"Isotropic conductivity scalar.")
227 .input_default(
"1.0")
228 .units(
UnitSI().m().s(-1) )
231 *
this += sigma.name(
"sigma")
232 .description(
"Transition coefficient between dimensions.")
233 .input_default(
"1.0")
236 *
this += water_source_density.name(
"water_source_density")
237 .description(
"Water source density.")
238 .input_default(
"0.0")
241 *
this += bc_type.name(
"bc_type")
242 .description(
"Boundary condition type.")
243 .input_selection( get_bc_type_selection() )
244 .input_default(
"\"none\"")
248 .disable_where(bc_type, {
none, seepage} )
250 .description(
"Prescribed pressure value on the boundary. Used for all values of ``bc_type`` except ``none`` and ``seepage``. "
251 "See documentation of ``bc_type`` for exact meaning of ``bc_pressure`` in individual boundary condition types.")
252 .input_default(
"0.0")
256 .disable_where(bc_type, {
none, dirichlet} )
258 .description(
"Incoming water boundary flux. Used for bc_types : ``total_flux``, ``seepage``, ``river``.")
259 .input_default(
"0.0")
260 .units(
UnitSI().m().s(-1) );
262 *
this += bc_robin_sigma
263 .disable_where(bc_type, {
none, dirichlet, seepage} )
264 .name(
"bc_robin_sigma")
265 .description(
"Conductivity coefficient in the ``total_flux`` or the ``river`` boundary condition type.")
266 .input_default(
"0.0")
269 *
this += bc_switch_pressure
270 .disable_where(bc_type, {
none, dirichlet, total_flux} )
271 .name(
"bc_switch_pressure")
272 .description(
"Critical switch pressure for ``seepage`` and ``river`` boundary conditions.")
273 .input_default(
"0.0")
278 *
this += init_pressure.name(
"init_pressure")
279 .description(
"Initial condition for pressure in time dependent problems.")
280 .input_default(
"0.0")
283 *
this += storativity.name(
"storativity")
284 .description(
"Storativity (in time dependent problems).")
285 .input_default(
"0.0")
288 *
this += extra_storativity.name(
"extra_storativity")
289 .description(
"Storativity added from upstream equation.")
291 .input_default(
"0.0")
292 .flags( input_copy );
294 *
this += extra_source.name(
"extra_water_source_density")
295 .description(
"Water source density added from upstream equation.")
296 .input_default(
"0.0")
298 .flags( input_copy );
342 MessageOut() <<
"Duplicate key 'time', time in flow equation is already initialized from parent class!";
349 data_ = make_shared<EqData>();
352 data_->is_linear=
true;
382 gravity_array.copy_to(gvec);
385 data_->gravity_vec_ =
data_->gravity_.subvec(0,2);
387 data_->bc_pressure.add_factory(
389 (
data_->gravity_,
"bc_piezo_head") );
390 data_->bc_switch_pressure.add_factory(
392 (
data_->gravity_,
"bc_switch_piezo_head") );
393 data_->init_pressure.add_factory(
395 (
data_->gravity_,
"init_piezo_head") );
403 MessageOut() <<
"Missing the key 'time', obligatory for the transient problems." << endl;
416 std::shared_ptr< FiniteElement<1> > fe1_rt = std::make_shared<FE_RT0_disc<1>>();
417 std::shared_ptr< FiniteElement<2> > fe2_rt = std::make_shared<FE_RT0_disc<2>>();
418 std::shared_ptr< FiniteElement<3> > fe3_rt = std::make_shared<FE_RT0_disc<3>>();
419 std::shared_ptr< FiniteElement<0> > fe0_disc = std::make_shared<FE_P_disc<0>>(0);
420 std::shared_ptr< FiniteElement<1> > fe1_disc = std::make_shared<FE_P_disc<1>>(0);
421 std::shared_ptr< FiniteElement<2> > fe2_disc = std::make_shared<FE_P_disc<2>>(0);
422 std::shared_ptr< FiniteElement<3> > fe3_disc = std::make_shared<FE_P_disc<3>>(0);
423 std::shared_ptr< FiniteElement<0> > fe0_cr = std::make_shared<FE_CR<0>>();
424 std::shared_ptr< FiniteElement<1> > fe1_cr = std::make_shared<FE_CR<1>>();
425 std::shared_ptr< FiniteElement<2> > fe2_cr = std::make_shared<FE_CR<2>>();
426 std::shared_ptr< FiniteElement<3> > fe3_cr = std::make_shared<FE_CR<3>>();
428 FESystem<0> fe0_sys( {fe0_disc, fe0_disc, fe0_cr} );
434 std::shared_ptr<DiscreteSpace> ds = std::make_shared<EqualOrderDiscreteSpace>(
mesh_, fe_sys);
435 data_->dh_ = std::make_shared<DOFHandlerMultiDim>(*
mesh_);
436 data_->dh_->distribute_dofs(ds);
440 this->
data_->multidim_assembler = AssemblyBase::create< AssemblyMH >(
data_);
444 uint rt_component = 0;
445 data_->full_solution =
data_->dh_->create_vector();
446 auto ele_flux_ptr = create_field_fe<3, FieldValue<3>::VectorFixed>(
data_->dh_, &
data_->full_solution, rt_component);
447 data_->flux.set(ele_flux_ptr, 0.0);
449 auto ele_velocity_ptr = std::make_shared< FieldDivide<3, FieldValue<3>::VectorFixed> >(ele_flux_ptr,
data_->cross_section);
450 data_->field_ele_velocity.set(ele_velocity_ptr, 0.0);
452 uint p_ele_component = 1;
453 auto ele_pressure_ptr = create_field_fe<3, FieldValue<3>::Scalar>(
data_->dh_, &
data_->full_solution, p_ele_component);
454 data_->field_ele_pressure.set(ele_pressure_ptr, 0.0);
456 uint p_edge_component = 2;
457 auto edge_pressure_ptr = create_field_fe<3, FieldValue<3>::Scalar>(
data_->dh_, &
data_->full_solution, p_edge_component);
458 data_->field_edge_pressure.set(edge_pressure_ptr, 0.0);
460 arma::vec4 gravity = (-1) *
data_->gravity_;
462 data_->field_ele_piezo_head.set(ele_piezo_head_ptr, 0.0);
466 uint p_edge_component = 2;
467 data_->dh_cr_ = std::make_shared<SubDOFHandlerMultiDim>(
data_->dh_,p_edge_component);
472 std::shared_ptr<DiscreteSpace> ds_cr_disc = std::make_shared<EqualOrderDiscreteSpace>(
mesh_, fe_cr_disc);
473 data_->dh_cr_disc_ = std::make_shared<DOFHandlerMultiDim>(*
mesh_);
474 data_->dh_cr_disc_->distribute_dofs(ds_cr_disc);
484 .val<Input::AbstractRecord>(
"linear_solver");
494 VecCreateMPI(PETSC_COMM_WORLD,
data_->dh_->distr()->lsize(),PETSC_DETERMINE,&(
steady_diagonal));
503 data_->water_balance_idx =
balance_->add_quantity(
"water_volume");
533 if (zero_time_term_from_right) {
564 data_->full_solution.local_to_ghost_begin();
565 data_->full_solution.local_to_ghost_end();
575 bool jump_time =
data_->storativity.is_jump_time();
576 if (! zero_time_term_from_left) {
586 WarningOut() <<
"Output of solution discontinuous in time not supported yet.\n";
595 if (! zero_time_term_from_left && ! jump_time && output)
output_data();
601 if (zero_time_term_from_right) {
606 }
else if (! zero_time_term_from_left && jump_time) {
607 WarningOut() <<
"Discontinuous time term not supported yet.\n";
619 return (
data_->storativity.input_list_size() == 0);
632 MessageOut().fmt(
"[nonlinear solver] norm of initial residual: {}\n", residual_norm);
635 int is_linear_common;
640 this->
max_n_it_ = nl_solver_rec.
val<
unsigned int>(
"max_it");
641 this->
min_n_it_ = nl_solver_rec.
val<
unsigned int>(
"min_it");
642 if (this->min_n_it_ > this->max_n_it_) this->min_n_it_ = this->
max_n_it_;
644 if (! is_linear_common) {
652 while (nonlinear_iteration_ < this->
min_n_it_ ||
653 (residual_norm > this->tolerance_ && nonlinear_iteration_ < this->
max_n_it_ )) {
655 convergence_history.push_back(residual_norm);
658 if (convergence_history.size() >= 5 &&
659 convergence_history[ convergence_history.size() - 1]/convergence_history[ convergence_history.size() - 2] > 0.9 &&
660 convergence_history[ convergence_history.size() - 1]/convergence_history[ convergence_history.size() - 5] > 0.8) {
666 THROW(ExcSolverDiverge() << EI_Reason(
"Stagnation."));
670 if (! is_linear_common)
676 if (is_linear_common){
679 MessageOut().fmt(
"[nonlinear solver] lin. it: {}, reason: {}, residual: {}\n",
693 MessageOut().fmt(
"[nonlinear solver] it: {} lin. it: {}, reason: {}, residual: {}\n",
696 chkerr(VecDestroy(&save_solution));
720 if(
data_->mortar_method_ != MortarMethod::NoMortar){
721 auto multidim_assembler = AssemblyBase::create< AssemblyMH >(
data_);
723 unsigned int dim = dh_cell.dim();
724 multidim_assembler[dim-1]->fix_velocity(dh_cell);
791 unsigned int dim = dh_cell.dim();
792 assembler[dim-1]->assemble(dh_cell);
810 double zeros[100000];
811 for(
int i=0; i<100000; i++) zeros[i] = 0.0;
814 tmp_rows.reserve(200);
817 dofs.reserve(
data_->dh_->max_elem_dofs());
818 dofs_ngh.reserve(
data_->dh_->max_elem_dofs());
823 const uint ndofs = dh_cell.n_dofs();
825 dh_cell.get_dof_indices(dofs);
834 for (
DHCellSide neighb_side : dh_cell.neighb_sides() ) {
842 const uint ndofs_ngh = dh_neighb_cell.
n_dofs();
843 dofs_ngh.resize(ndofs_ngh);
849 const unsigned int t = dh_neighb_cell.
n_dofs() - (dh_neighb_cell.
dim()+1) + neighb_side.side().side_idx();
850 tmp_rows.push_back(dofs_ngh[t]);
864 for(
auto &isec : isec_list ) {
868 const uint ndofs_slave = dh_cell_slave.
n_dofs();
869 dofs_ngh.resize(ndofs_slave);
873 for(
unsigned int i_side=0; i_side < dh_cell_slave.
elm()->
n_sides(); i_side++) {
875 tmp_rows.push_back( dofs_ngh[(ndofs_slave+1)/2+i_side] );
885 edge_rows = dofs.
data() + nsides +1;
924 const uint ndofs = dh_cell.n_dofs();
925 global_dofs.resize(ndofs);
926 dh_cell.get_dof_indices(global_dofs);
928 double cs =
data_->cross_section.value(ele.
centre(), ele);
931 double source = ele.
measure() * cs *
932 (
data_->water_source_density.value(ele.
centre(), ele)
938 {dh_cell.get_loc_dof_indices()[ndofs/2]}, {0}, {source});
958 #ifdef FLOW123D_HAVE_BDDCML
959 WarningOut() <<
"For BDDC no Schur complements are used.";
972 THROW( ExcBddcmlNotSupported() );
973 #endif // FLOW123D_HAVE_BDDCML
978 WarningOut() <<
"Invalid number of Schur Complements. Using 2.";
991 ls->LinSys::set_from_input(in_rec);
1000 ISCreateGeneral(PETSC_COMM_SELF, side_dofs_vec.size(), &(side_dofs_vec[0]), PETSC_COPY_VALUES, &is);
1015 const PetscInt *b_indices;
1016 ISGetIndices(ls->
IsB, &b_indices);
1018 for(
uint i_b=0, i_bb=0; i_b < b_size && i_bb < elem_dofs_vec.size(); i_b++) {
1019 if (b_indices[i_b] == elem_dofs_vec[i_bb])
1020 elem_dofs_vec[i_bb++] = i_b + ds->
begin();
1022 ISRestoreIndices(ls->
IsB, &b_indices);
1025 ISCreateGeneral(PETSC_COMM_SELF, elem_dofs_vec.size(), &(elem_dofs_vec[0]), PETSC_COPY_VALUES, &is);
1053 THROW( ExcUnknownSolver() );
1065 START_TIMER(
"DarcyFlowMH_Steady::assembly_linear_system");
1067 data_->is_linear=
true;
1122 std::string output_file;
1134 PetscViewerASCIIOpen(PETSC_COMM_WORLD, output_file.c_str(), &viewer);
1135 PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
1147 double d_max = std::numeric_limits<double>::max();
1148 double h1 = d_max, h2 = d_max, h3 = d_max;
1149 double he2 = d_max, he3 = d_max;
1152 case 1: h1 = std::min(h1,ele.measure());
break;
1153 case 2: h2 = std::min(h2,ele.measure());
break;
1154 case 3: h3 = std::min(h3,ele.measure());
break;
1157 for (
unsigned int j=0; j<ele->n_sides(); j++) {
1159 case 2: he2 = std::min(he2, ele.side(j)->measure());
break;
1160 case 3: he3 = std::min(he3, ele.side(j)->measure());
break;
1164 if(h1 == d_max) h1 = 0;
1165 if(h2 == d_max) h2 = 0;
1166 if(h3 == d_max) h3 = 0;
1167 if(he2 == d_max) he2 = 0;
1168 if(he3 == d_max) he3 = 0;
1171 file = fopen(output_file.c_str(),
"a");
1172 fprintf(file,
"nA = %d;\n",
data_->dh_cr_disc_->distr()->size());
1173 fprintf(file,
"nB = %d;\n",
data_->dh_->mesh()->get_el_ds()->size());
1174 fprintf(file,
"nBF = %d;\n",
data_->dh_cr_->distr()->size());
1175 fprintf(file,
"h1 = %e;\nh2 = %e;\nh3 = %e;\n", h1, h2, h3);
1176 fprintf(file,
"he2 = %e;\nhe3 = %e;\n", he2, he3);
1194 START_TIMER(
"DarcyFlowMH_Steady::set_mesh_data_for_bddc");
1222 dh_cell.get_dof_indices(cell_dofs_global);
1224 inet.insert(inet.end(), cell_dofs_global.begin(), cell_dofs_global.end());
1225 uint n_inet = cell_dofs_global.size();
1228 uint dim = dh_cell.elm().dim();
1229 elDimMax = std::max( elDimMax, dim );
1230 elDimMin = std::min( elDimMin, dim );
1234 isegn.push_back( dh_cell.elm_idx());
1237 for (
unsigned int si=0; si<dh_cell.elm()->n_sides(); si++) {
1238 arma::vec3 coord = dh_cell.elm().side(si)->centre();
1241 localDofMap.insert( std::make_pair( cell_dofs_global[si], coord ) );
1243 localDofMap.insert( std::make_pair( cell_dofs_global[si+dim+2], coord ) );
1246 arma::vec3 elm_centre = dh_cell.elm().centre();
1247 localDofMap.insert( std::make_pair( cell_dofs_global[dim+1], elm_centre ) );
1251 for(
DHCellSide side : dh_cell.neighb_sides()) {
1252 uint neigh_dim = side.cell().elm().dim();
1253 side.cell().get_dof_indices(cell_dofs_global);
1254 int edge_row = cell_dofs_global[neigh_dim+2+side.side_idx()];
1255 localDofMap.insert( std::make_pair( edge_row, side.centre() ) );
1256 inet.push_back( edge_row );
1259 nnet.push_back(n_inet);
1264 double conduct =
data_->conductivity.value( elm_centre , dh_cell.elm() );
1265 auto aniso =
data_->anisotropy.value( elm_centre , dh_cell.elm() );
1269 for (
int i = 0; i < 3; i++) {
1270 coef = coef + aniso.at(i,i);
1273 coef = conduct*coef / 3;
1276 "Zero coefficient of hydrodynamic resistance %f . \n ", coef );
1277 element_permeability.push_back( 1. / coef );
1287 auto distr =
data_->dh_->distr();
1296 int numNodeSub = localDofMap.size();
1307 for ( ; itB != localDofMap.end(); ++itB ) {
1308 isngn[ind] = itB -> first;
1311 for (
int j = 0; j < 3; j++ ) {
1312 xyz[ j*numNodeSub + ind ] = coord[j];
1317 localDofMap.clear();
1325 Global2LocalMap_ global2LocalNodeMap;
1326 for (
unsigned ind = 0; ind < isngn.size(); ++ind ) {
1327 global2LocalNodeMap.insert( std::make_pair(
static_cast<unsigned>( isngn[ind] ), ind ) );
1332 for (
unsigned int iEle = 0; iEle < isegn.size(); iEle++ ) {
1333 int nne = nnet[ iEle ];
1334 for (
int ien = 0; ien < nne; ien++ ) {
1336 int indGlob = inet[indInet];
1338 Global2LocalMap_::iterator pos = global2LocalNodeMap.find( indGlob );
1339 OLD_ASSERT( pos != global2LocalNodeMap.end(),
1340 "Cannot remap node index %d to local indices. \n ", indGlob );
1341 int indLoc =
static_cast<int> ( pos -> second );
1344 inet[ indInet++ ] = indLoc;
1348 int numNodes =
size;
1349 int numDofsInt =
size;
1351 int meshDim = elDimMax;
1363 bddc_ls -> load_mesh( LinSys_BDDC::BDDCMatrixType::SPD_VIA_SYMMETRICGENERAL, spaceDim, numNodes, numDofsInt, inet, nnet, nndf, isegn, isngn, isngn, xyz, element_permeability, meshDim );
1385 if(
time_ !=
nullptr)
1438 const IntIdx p_ele_dof = dh_cell.get_loc_dof_indices()[dh_cell.n_dofs()/2];
1439 local_sol[p_ele_dof] =
data_->init_pressure.value(ele.
centre(),ele);
1450 PetscScalar *local_diagonal;
1458 dofs.reserve(
data_->dh_->max_elem_dofs());
1462 dofs.resize(dh_cell.n_dofs());
1463 dh_cell.get_dof_indices(dofs);
1466 const uint p_ele_dof = dh_cell.n_dofs() / 2;
1468 double diagonal_coeff =
data_->cross_section.value(ele.
centre(), ele)
1470 +
data_->extra_storativity.value(ele.
centre(), ele))
1472 local_diagonal[dh_cell.get_loc_dof_indices()[p_ele_dof]]= - diagonal_coeff /
time_->
dt();
1474 balance_->add_mass_values(
data_->water_balance_idx, dh_cell,
1475 {dh_cell.get_loc_dof_indices()[p_ele_dof]},
1476 {diagonal_coeff}, 0.0);
1478 VecRestoreArray(new_diagonal,& local_diagonal);
1479 MatDiagonalSet(*( schur0->get_matrix() ), new_diagonal, ADD_VALUES);
1481 schur0->set_matrix_changed();
1483 balance_->finish_mass_assembly(data_->water_balance_idx);
1490 if (scale_factor != 1.0) {
1511 void dofs_range(
unsigned int n_dofs,
unsigned int &min,
unsigned int &max,
unsigned int component) {
1515 }
else if (component==1) {
1527 unsigned int i, n_dofs, min, max;
1531 n_dofs = dh_cell.get_dof_indices(dof_indices);
1533 for (i=min; i<max; ++i) dof_vec.push_back(dof_indices[i]);