19#ifndef DUMUX_FVPRESSURE_HH
20#define DUMUX_FVPRESSURE_HH
59 using Element =
typename GridView::Traits::template Codim<0>::Entity;
60 using Intersection =
typename GridView::Intersection;
68 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
131 for (
const auto& element : elements(problem_.gridView()))
133 PrimaryVariables initValues;
134 problem_.initial(initValues, element);
136 int eIdxGlobal = problem_.variables().index(element);
138 pressure_[eIdxGlobal] = initValues[
pressEqIdx];
152 void getSource(
EntryType& entry,
const Element& element,
const CellData& cellData,
const bool first);
174 void getFlux(
EntryType& entry,
const Intersection& intersection,
const CellData& cellData,
const bool first);
186 const Intersection& intersection,
const CellData& cellData,
const bool first);
195 {
return pressure_[eIdxGlobal];}
215 int size = problem_.gridView().size(0);
216 A_.setSize(size, size);
217 A_.setBuildMode(Matrix::random);
219 pressure_.resize(size);
221 asImp_().initializeMatrix();
230 asImp_().assemble(
false); Dune::dinfo <<
"pressure calculation"<< std::endl;
238 DUNE_THROW(Dune::NotImplemented,
"Velocity calculation not implemented in pressure model!");
253 int eIdxGlobal = problem_.variables().index(element);
254 outstream << pressure_[eIdxGlobal][0];
266 int eIdxGlobal = problem_.variables().index(element);
267 instream >> pressure_[eIdxGlobal][0];
281 fixPressure_.insert(std::make_pair(eIdxGlobal,
pressure));
292 fixPressure_.erase(eIdxGlobal);
297 fixPressure_.clear();
310 Implementation &asImp_()
311 {
return *
static_cast<Implementation *
>(
this);}
314 const Implementation &asImp_()
const
315 {
return *
static_cast<const Implementation *
>(
this);}
319 PressureSolution pressure_;
321 std::string solverName_;
322 std::string preconditionerName_;
327 std::map<int, Scalar> fixPressure_;
331template<
class TypeTag>
334 initializeMatrixRowSize();
336 initializeMatrixIndices();
341template<
class TypeTag>
345 for (
const auto& element : elements(problem_.gridView()))
348 int eIdxGlobalI = problem_.variables().index(element);
354 for (
const auto& intersection : intersections(problem_.gridView(), element))
356 if (intersection.neighbor())
359 A_.setrowsize(eIdxGlobalI, rowSize);
364template<
class TypeTag>
368 for (
const auto& element : elements(problem_.gridView()))
371 int eIdxGlobalI = problem_.variables().index(element);
374 A_.addindex(eIdxGlobalI, eIdxGlobalI);
377 for (
const auto& intersection : intersections(problem_.gridView(), element))
378 if (intersection.neighbor())
381 int eIdxGlobalJ = problem_.variables().index(intersection.outside());
384 A_.addindex(eIdxGlobalI, eIdxGlobalJ);
400template<
class TypeTag>
407 for (
const auto& element : elements(problem_.gridView()))
410 int eIdxGlobalI = problem_.variables().index(element);
413 if (element.partitionType() == Dune::InteriorEntity)
416 CellData& cellDataI = problem_.variables().cellData(eIdxGlobalI);
421 asImp_().getSource(entries, element, cellDataI, first);
422 f_[eIdxGlobalI] += entries[rhs];
426 for (
const auto& intersection : intersections(problem_.gridView(), element))
429 if (intersection.neighbor())
431 auto elementNeighbor = intersection.outside();
433 int eIdxGlobalJ = problem_.variables().index(elementNeighbor);
436 bool haveSameLevel = (element.level() == elementNeighbor.level());
440 if (getPropValue<TypeTag, Properties::VisitFacesOnlyOnce>()
441 && (eIdxGlobalI > eIdxGlobalJ) && haveSameLevel
442 && elementNeighbor.partitionType() == Dune::InteriorEntity)
446 asImp_().getFlux(entries, intersection, cellDataI, first);
449 f_[eIdxGlobalI] -= entries[rhs];
452 A_[eIdxGlobalI][eIdxGlobalI] += entries[matrix];
455 A_[eIdxGlobalI][eIdxGlobalJ] -= entries[matrix];
458 if (getPropValue<TypeTag, Properties::VisitFacesOnlyOnce>()
459 && elementNeighbor.partitionType() == Dune::InteriorEntity)
461 f_[eIdxGlobalJ] += entries[rhs];
462 A_[eIdxGlobalJ][eIdxGlobalJ] += entries[matrix];
463 A_[eIdxGlobalJ][eIdxGlobalI] -= entries[matrix];
472 asImp_().getFluxOnBoundary(entries, intersection, cellDataI, first);
475 f_[eIdxGlobalI] += entries[rhs];
477 A_[eIdxGlobalI][eIdxGlobalI] += entries[matrix];
484 asImp_().getStorage(entries, element, cellDataI, first);
485 f_[eIdxGlobalI] += entries[rhs];
487 A_[eIdxGlobalI][eIdxGlobalI] += entries[matrix];
492 A_[eIdxGlobalI] = 0.0;
493 A_[eIdxGlobalI][eIdxGlobalI] = 1.0;
494 f_[eIdxGlobalI] = pressure_[eIdxGlobalI];
511template<
class Solver,
class Problem>
512typename std::enable_if_t<!Detail::isParallelAMGBackend<Solver>::value, Solver>
518template<
class Solver,
class Problem>
519typename std::enable_if_t<Detail::isParallelAMGBackend<Solver>::value, Solver>
522 return Solver(problem.gridView(), problem.model().dofMapper());
526template<
class TypeTag>
531 int verboseLevelSolver = getParam<int>(
"LinearSolver.Verbosity", 0);
533 if (verboseLevelSolver)
534 std::cout << __FILE__ <<
": solve for pressure" << std::endl;
537 if (fixPressure_.size() > 0)
539 for (
auto it = fixPressure_.begin(); it != fixPressure_.end(); ++it)
542 A_[it->first][it->first] = 1;
543 f_[it->first] = it->second;
550 auto solver = getSolver<Solver>(problem_);
551 bool converged = solver.solve(A_, pressure_, f_);
Some exceptions thrown in DuMux
Define some often used mathematical functions.
std::enable_if_t<!Detail::isParallelAMGBackend< Solver >::value, Solver > getSolver(const Problem &problem)
Definition: sequential/cellcentered/pressure.hh:513
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type GetProp
get the type of a property (equivalent to old macro GET_PROP(...))
Definition: propertysystem.hh:140
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
A linear solver based on the ISTL AMG preconditioner and the ISTL BiCGSTAB solver.
Definition: amgbackend.hh:50
The finite volume base class for the solution of a pressure equation.
Definition: sequential/cellcentered/pressure.hh:49
void getSource(EntryType &entry, const Element &element, const CellData &cellData, const bool first)
Function which calculates the source entry.
@ pressEqIdx
Definition: sequential/cellcentered/pressure.hh:95
void initializeMatrixIndices()
Initialize the global matrix of the system of equations to solve.
Definition: sequential/cellcentered/pressure.hh:365
void assemble(bool first)
Function which assembles the system of equations to be solved.
Definition: sequential/cellcentered/pressure.hh:401
void getStorage(EntryType &entry, const Element &element, const CellData &cellData, const bool first)
Function which calculates the storage entry.
FVPressure(Problem &problem)
Constructs a FVPressure object.
Definition: sequential/cellcentered/pressure.hh:304
void initialize()
Initialize pressure model.
Definition: sequential/cellcentered/pressure.hh:213
void getFlux(EntryType &entry, const Intersection &intersection, const CellData &cellData, const bool first)
Function which calculates the flux entry.
void updateVelocity()
Definition: sequential/cellcentered/pressure.hh:241
void deserializeEntity(std::istream &instream, const Element &element)
Function for deserialization of the pressure field.
Definition: sequential/cellcentered/pressure.hh:264
void resetFixPressureAtIndex()
Definition: sequential/cellcentered/pressure.hh:295
void unsetFixPressureAtIndex(int eIdxGlobal)
Reset the fixed pressure state.
Definition: sequential/cellcentered/pressure.hh:290
void initializePressure()
Initialization of the pressure solution vector: Initialization with meaningful values may.
Definition: sequential/cellcentered/pressure.hh:129
PressureSolution & pressure()
Returns the vector containing the pressure solution.
Definition: sequential/cellcentered/pressure.hh:120
const Matrix & globalMatrix()
Returns the global matrix of the last pressure solution step.
Definition: sequential/cellcentered/pressure.hh:198
void serializeEntity(std::ostream &outstream, const Element &element)
Function for serialization of the pressure field.
Definition: sequential/cellcentered/pressure.hh:251
const PressureSolution & pressure() const
Returns the vector containing the pressure solution.
Definition: sequential/cellcentered/pressure.hh:124
@ rhs
index for the right hand side entry
Definition: sequential/cellcentered/pressure.hh:88
@ matrix
index for the global matrix entry
Definition: sequential/cellcentered/pressure.hh:89
Matrix A_
Global stiffnes matrix (sparse matrix which is build by the initializeMatrix() function)
Definition: sequential/cellcentered/pressure.hh:324
void getFluxOnBoundary(EntryType &entry, const Intersection &intersection, const CellData &cellData, const bool first)
Function which calculates the boundary flux entry.
const RHSVector & rightHandSide()
Returns the right hand side vector of the last pressure solution step.
Definition: sequential/cellcentered/pressure.hh:204
void solve()
Solves the global system of equations to get the spatial distribution of the pressure.
Definition: sequential/cellcentered/pressure.hh:527
const Scalar pressure(const int eIdxGlobal) const
Public access function for the primary pressure variable.
Definition: sequential/cellcentered/pressure.hh:194
void initializeMatrixRowSize()
Initialize the global matrix of the system of equations to solve.
Definition: sequential/cellcentered/pressure.hh:342
void calculateVelocity()
Definition: sequential/cellcentered/pressure.hh:236
RHSVector f_
Right hand side vector.
Definition: sequential/cellcentered/pressure.hh:325
void setFixPressureAtIndex(Scalar pressure, int eIdxGlobal)
Set a pressure to be fixed at a certain cell.
Definition: sequential/cellcentered/pressure.hh:279
void update()
Pressure update.
Definition: sequential/cellcentered/pressure.hh:228
Dune::FieldVector< Scalar, 2 > EntryType
Definition: sequential/cellcentered/pressure.hh:78
void initializeMatrix()
Initialize the global matrix of the system of equations to solve.
Definition: sequential/cellcentered/pressure.hh:332
Definition: sequential/cellcentered/pressure.hh:506
Base file for properties related to sequential IMPET algorithms.