3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
sequential/cellcentered/pressure.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
19#ifndef DUMUX_FVPRESSURE_HH
20#define DUMUX_FVPRESSURE_HH
21
22// dumux environment
23#include <type_traits>
24#include <dumux/common/math.hh>
26#include <map>
32namespace Dumux
33{
35
47template<class TypeTag> class FVPressure
48{
49 //the model implementation
50 using Implementation = typename GET_PROP_TYPE(TypeTag, PressureModel);
51
52 using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
53 using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
54 using Problem = typename GET_PROP_TYPE(TypeTag, Problem);
55 using CellData = typename GET_PROP_TYPE(TypeTag, CellData);
56
57 // using declarations to abbreviate several dune classes...
58 using Element = typename GridView::Traits::template Codim<0>::Entity;
59 using Intersection = typename GridView::Intersection;
60
61 // the typenames used for the stiffness matrix and solution vector
62 using Matrix = typename GET_PROP_TYPE(TypeTag, PressureCoefficientMatrix);
63 using RHSVector = typename GET_PROP_TYPE(TypeTag, PressureRHSVector);
64 using PressureSolution = typename GET_PROP_TYPE(TypeTag, PressureSolutionVector);
65
66 using SolutionTypes = typename GET_PROP(TypeTag, SolutionTypes);
67 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
68
69 using Indices = typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
70
71protected:
72
77 using EntryType = Dune::FieldVector<Scalar, 2>;
78
80
85 enum
86 {
87 rhs = 1,
88 matrix = 0
89
90 };
91
92 enum
93 {
94 pressEqIdx = Indices::pressureEqIdx,
95 };
96
98 void initializeMatrix();
101
102
113 void assemble(bool first);
114
116 void solve();
117
119 PressureSolution& pressure()
120 { return pressure_;}
121
123 const PressureSolution& pressure() const
124 { return pressure_;}
125
127 //result in better convergence of the linear solver!
129 {
130 for (const auto& element : elements(problem_.gridView()))
131 {
132 PrimaryVariables initValues;
133 problem_.initial(initValues, element);
134
135 int eIdxGlobal = problem_.variables().index(element);
136
137 pressure_[eIdxGlobal] = initValues[pressEqIdx];
138 }
139 }
140
141public:
151 void getSource(EntryType& entry, const Element& element, const CellData& cellData, const bool first);
152
162 void getStorage(EntryType& entry, const Element& element, const CellData& cellData, const bool first);
163
173 void getFlux(EntryType& entry, const Intersection& intersection, const CellData& cellData, const bool first);
174
185 const Intersection& intersection, const CellData& cellData, const bool first);
186
193 const Scalar pressure(const int eIdxGlobal) const
194 { return pressure_[eIdxGlobal];}
195
197 const Matrix& globalMatrix()
198 {
199 return A_;
200 }
201
203 const RHSVector& rightHandSide()
204 {
205 return f_;
206 }
207
213 {
214 int size = problem_.gridView().size(0);//resize to make sure the final grid size (after the problem was completely built) is used!
215 A_.setSize(size, size);
216 A_.setBuildMode(Matrix::random);
217 f_.resize(size);
218 pressure_.resize(size);
220 asImp_().initializeMatrix();// initialize sparse matrix
221 }
222
227 void update()
228 {
229 asImp_().assemble(false); Dune::dinfo << "pressure calculation"<< std::endl;
230 solve();
231
232 return;
233 }
234
236 {
237 DUNE_THROW(Dune::NotImplemented,"Velocity calculation not implemented in pressure model!");
238 }
239
241 {} //empty function for the case the velocity is calculated in the transport model
242
250 void serializeEntity(std::ostream &outstream, const Element &element)
251 {
252 int eIdxGlobal = problem_.variables().index(element);
253 outstream << pressure_[eIdxGlobal][0];
254 }
255
263 void deserializeEntity(std::istream &instream, const Element &element)
264 {
265 int eIdxGlobal = problem_.variables().index(element);
266 instream >> pressure_[eIdxGlobal][0];
267 }
268
278 void setFixPressureAtIndex(Scalar pressure, int eIdxGlobal)
279 {
280 fixPressure_.insert(std::make_pair(eIdxGlobal, pressure));
281 }
282
289 void unsetFixPressureAtIndex(int eIdxGlobal)
290 {
291 fixPressure_.erase(eIdxGlobal);
292 }
293
295 {
296 fixPressure_.clear();
297 }
298
303 FVPressure(Problem& problem) :
304 problem_(problem)
305 {}
306
307private:
309 Implementation &asImp_()
310 { return *static_cast<Implementation *>(this);}
311
313 const Implementation &asImp_() const
314 { return *static_cast<const Implementation *>(this);}
315
316 Problem& problem_;
317
318 PressureSolution pressure_;
319
320 std::string solverName_;
321 std::string preconditionerName_;
322protected:
323 Matrix A_;
324 RHSVector f_;
325private:
326 std::map<int, Scalar> fixPressure_;
327};
328
330template<class TypeTag>
332{
333 initializeMatrixRowSize();
334 A_.endrowsizes();
335 initializeMatrixIndices();
336 A_.endindices();
337}
338
340template<class TypeTag>
342{
343 // determine matrix row sizes
344 for (const auto& element : elements(problem_.gridView()))
345 {
346 // cell index
347 int eIdxGlobalI = problem_.variables().index(element);
348
349 // initialize row size
350 int rowSize = 1;
351
352 // run through all intersections with neighbors
353 for (const auto& intersection : intersections(problem_.gridView(), element))
354 {
355 if (intersection.neighbor())
356 rowSize++;
357 }
358 A_.setrowsize(eIdxGlobalI, rowSize);
359 }
360}
361
363template<class TypeTag>
365{
366 // determine position of matrix entries
367 for (const auto& element : elements(problem_.gridView()))
368 {
369 // cell index
370 int eIdxGlobalI = problem_.variables().index(element);
371
372 // add diagonal index
373 A_.addindex(eIdxGlobalI, eIdxGlobalI);
374
375 // run through all intersections with neighbors
376 for (const auto& intersection : intersections(problem_.gridView(), element))
377 if (intersection.neighbor())
378 {
379 // access neighbor
380 int eIdxGlobalJ = problem_.variables().index(intersection.outside());
381
382 // add off diagonal index
383 A_.addindex(eIdxGlobalI, eIdxGlobalJ);
384 }
385 }
386}
387
388
399template<class TypeTag>
401{
402 // initialization: set matrix A_ to zero
403 A_ = 0;
404 f_ = 0;
405
406 for (const auto& element : elements(problem_.gridView()))
407 {
408 // get the global index of the cell
409 int eIdxGlobalI = problem_.variables().index(element);
410
411 // assemble interior element contributions
412 if (element.partitionType() == Dune::InteriorEntity)
413 {
414 // get the cell data
415 CellData& cellDataI = problem_.variables().cellData(eIdxGlobalI);
416
417 EntryType entries(0.);
418
419 /***** source term ***********/
420 asImp_().getSource(entries, element, cellDataI, first);
421 f_[eIdxGlobalI] += entries[rhs];
422
423 /***** flux term ***********/
424 // iterate over all faces of the cell
425 for (const auto& intersection : intersections(problem_.gridView(), element))
426 {
427 /************* handle interior face *****************/
428 if (intersection.neighbor())
429 {
430 auto elementNeighbor = intersection.outside();
431
432 int eIdxGlobalJ = problem_.variables().index(elementNeighbor);
433
434 // check for hanging nodes
435 bool haveSameLevel = (element.level() == elementNeighbor.level());
436 // calculate only from one side (except for hanging nodes), but add matrix entries for both sides
437 // the last condition is needed to properly assemble in the presence
438 // of ghost elements
440 && (eIdxGlobalI > eIdxGlobalJ) && haveSameLevel
441 && elementNeighbor.partitionType() == Dune::InteriorEntity)
442 continue;
443
444 entries = 0;
445 asImp_().getFlux(entries, intersection, cellDataI, first);
446
447 //set right hand side
448 f_[eIdxGlobalI] -= entries[rhs];
449
450 // set diagonal entry
451 A_[eIdxGlobalI][eIdxGlobalI] += entries[matrix];
452
453 // set off-diagonal entry
454 A_[eIdxGlobalI][eIdxGlobalJ] -= entries[matrix];
455
456 // The second condition is needed to not spoil the ghost element entries
458 && elementNeighbor.partitionType() == Dune::InteriorEntity)
459 {
460 f_[eIdxGlobalJ] += entries[rhs];
461 A_[eIdxGlobalJ][eIdxGlobalJ] += entries[matrix];
462 A_[eIdxGlobalJ][eIdxGlobalI] -= entries[matrix];
463 }
464
465 } // end neighbor
466
467 /************* boundary face ************************/
468 else
469 {
470 entries = 0;
471 asImp_().getFluxOnBoundary(entries, intersection, cellDataI, first);
472
473 //set right hand side
474 f_[eIdxGlobalI] += entries[rhs];
475 // set diagonal entry
476 A_[eIdxGlobalI][eIdxGlobalI] += entries[matrix];
477 }
478 } //end interfaces loop
479 // printmatrix(std::cout, A_, "global stiffness matrix", "row", 11, 3);
480
481 /***** storage term ***********/
482 entries = 0;
483 asImp_().getStorage(entries, element, cellDataI, first);
484 f_[eIdxGlobalI] += entries[rhs];
485 // set diagonal entry
486 A_[eIdxGlobalI][eIdxGlobalI] += entries[matrix];
487 }
488 // assemble overlap and ghost element contributions
489 else
490 {
491 A_[eIdxGlobalI] = 0.0;
492 A_[eIdxGlobalI][eIdxGlobalI] = 1.0;
493 f_[eIdxGlobalI] = pressure_[eIdxGlobalI];
494 }
495 } // end grid traversal
496// printmatrix(std::cout, A_, "global stiffness matrix after assempling", "row", 11,3);
497// printvector(std::cout, f_, "right hand side", "row", 10);
498}
499
500// forward declaration
501template<class GV, class T>
503
504namespace Detail {
505template<class T> struct isParallelAMGBackend : public std::false_type {};
506template<class GV, class T>
507struct isParallelAMGBackend<Dumux::ParallelAMGBackend<GV, T>> : public std::true_type {};
508} // end namespace Detail
509
510template<class Solver, class Problem>
511typename std::enable_if_t<!Detail::isParallelAMGBackend<Solver>::value, Solver>
512getSolver(const Problem& problem)
513{
514 return Solver();
515}
516
517template<class Solver, class Problem>
518typename std::enable_if_t<Detail::isParallelAMGBackend<Solver>::value, Solver>
519getSolver(const Problem& problem)
520{
521 return Solver(problem.gridView(), problem.model().dofMapper());
522}
523
525template<class TypeTag>
527{
528 using Solver = typename GET_PROP_TYPE(TypeTag, LinearSolver);
529
530 int verboseLevelSolver = getParam<int>("LinearSolver.Verbosity");
531
532 if (verboseLevelSolver)
533 std::cout << __FILE__ << ": solve for pressure" << std::endl;
534
535 //set a fixed pressure for a certain cell
536 if (fixPressure_.size() > 0)
537 {
538 for (auto it = fixPressure_.begin(); it != fixPressure_.end(); ++it)
539 {
540 A_[it->first] = 0;
541 A_[it->first][it->first] = 1;
542 f_[it->first] = it->second;
543 }
544 }
545
546// printmatrix(std::cout, A_, "global stiffness matrix", "row", 11, 3);
547// printvector(std::cout, f_, "right hand side", "row", 10, 1, 3);
548
549 auto solver = getSolver<Solver>(problem_);
550 solver.solve(A_, pressure_, f_);
551
552// printvector(std::cout, pressure_, "pressure", "row", 200, 1, 3);
553}
554
555} //end namespace Dumux
556#endif
#define GET_PROP_VALUE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:282
#define GET_PROP(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:281
#define GET_PROP_TYPE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:283
Define some often used mathematical functions.
Dune::IteratorRange< typename MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper >::Intersections::const_iterator > intersections(const MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper > &glue)
Range generator to iterate with range-based for loops over all intersections as follows: for (const a...
Definition: glue.hh:62
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
std::enable_if_t<!Detail::isParallelAMGBackend< Solver >::value, Solver > getSolver(const Problem &problem)
Definition: sequential/cellcentered/pressure.hh:512
Property tag PressureRHSVector
Type of the right hand side vector given to the linear solver.
Definition: sequential/pressureproperties.hh:60
Property tag VisitFacesOnlyOnce
Type of solution vector or pressure system.
Definition: sequential/pressureproperties.hh:62
Property tag PressureSolutionVector
Definition: sequential/pressureproperties.hh:61
Property tag PressureCoefficientMatrix
Gives maximum number of intersections of an element and neighboring elements.
Definition: porousmediumflow/sequential/properties.hh:74
Property tag PressureModel
The type of the discretization of a pressure model.
Definition: porousmediumflow/sequential/properties.hh:65
A linear solver based on the ISTL AMG preconditioner and the ISTL BiCGSTAB solver.
Definition: amgbackend.hh:81
Base class for linear solvers.
Definition: dumux/linear/solver.hh:37
The finite volume base class for the solution of a pressure equation.
Definition: sequential/cellcentered/pressure.hh:48
void getSource(EntryType &entry, const Element &element, const CellData &cellData, const bool first)
Function which calculates the source entry.
void initializeMatrixIndices()
Initialize the global matrix of the system of equations to solve.
Definition: sequential/cellcentered/pressure.hh:364
void assemble(bool first)
Function which assembles the system of equations to be solved.
Definition: sequential/cellcentered/pressure.hh:400
@ pressEqIdx
Definition: sequential/cellcentered/pressure.hh:94
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:303
void initialize()
Initialize pressure model.
Definition: sequential/cellcentered/pressure.hh:212
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:240
void deserializeEntity(std::istream &instream, const Element &element)
Function for deserialization of the pressure field.
Definition: sequential/cellcentered/pressure.hh:263
void resetFixPressureAtIndex()
Definition: sequential/cellcentered/pressure.hh:294
void unsetFixPressureAtIndex(int eIdxGlobal)
Reset the fixed pressure state.
Definition: sequential/cellcentered/pressure.hh:289
void initializePressure()
Initialization of the pressure solution vector: Initialization with meaningful values may.
Definition: sequential/cellcentered/pressure.hh:128
PressureSolution & pressure()
Returns the vector containing the pressure solution.
Definition: sequential/cellcentered/pressure.hh:119
const Matrix & globalMatrix()
Returns the global matrix of the last pressure solution step.
Definition: sequential/cellcentered/pressure.hh:197
void serializeEntity(std::ostream &outstream, const Element &element)
Function for serialization of the pressure field.
Definition: sequential/cellcentered/pressure.hh:250
const PressureSolution & pressure() const
Returns the vector containing the pressure solution.
Definition: sequential/cellcentered/pressure.hh:123
@ rhs
index for the right hand side entry
Definition: sequential/cellcentered/pressure.hh:87
@ matrix
index for the global matrix entry
Definition: sequential/cellcentered/pressure.hh:88
Matrix A_
Global stiffnes matrix (sparse matrix which is build by the initializeMatrix() function)
Definition: sequential/cellcentered/pressure.hh:323
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:203
void solve()
Solves the global system of equations to get the spatial distribution of the pressure.
Definition: sequential/cellcentered/pressure.hh:526
const Scalar pressure(const int eIdxGlobal) const
Public access function for the primary pressure variable.
Definition: sequential/cellcentered/pressure.hh:193
void initializeMatrixRowSize()
Initialize the global matrix of the system of equations to solve.
Definition: sequential/cellcentered/pressure.hh:341
void calculateVelocity()
Definition: sequential/cellcentered/pressure.hh:235
RHSVector f_
Right hand side vector.
Definition: sequential/cellcentered/pressure.hh:324
void setFixPressureAtIndex(Scalar pressure, int eIdxGlobal)
Set a pressure to be fixed at a certain cell.
Definition: sequential/cellcentered/pressure.hh:278
void update()
Pressure update.
Definition: sequential/cellcentered/pressure.hh:227
Dune::FieldVector< Scalar, 2 > EntryType
Definition: sequential/cellcentered/pressure.hh:77
void initializeMatrix()
Initialize the global matrix of the system of equations to solve.
Definition: sequential/cellcentered/pressure.hh:331
Definition: sequential/cellcentered/pressure.hh:505
Base file for properties related to sequential IMPET algorithms.