NewtonIterationBlackoilCPR.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2014 SINTEF ICT, Applied Mathematics.
3  Copyright 2015 IRIS AS
4  Copyright 2015 Dr. Blatt - HPC-Simulation-Software & Services
5  Copyright 2015 NTNU
6  Copyright 2015 Statoil AS
7 
8  This file is part of the Open Porous Media project (OPM).
9 
10  OPM is free software: you can redistribute it and/or modify
11  it under the terms of the GNU General Public License as published by
12  the Free Software Foundation, either version 3 of the License, or
13  (at your option) any later version.
14 
15  OPM is distributed in the hope that it will be useful,
16  but WITHOUT ANY WARRANTY; without even the implied warranty of
17  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18  GNU General Public License for more details.
19 
20  You should have received a copy of the GNU General Public License
21  along with OPM. If not, see <http://www.gnu.org/licenses/>.
22 */
23 
24 #ifndef OPM_NEWTONITERATIONBLACKOILCPR_HEADER_INCLUDED
25 #define OPM_NEWTONITERATIONBLACKOILCPR_HEADER_INCLUDED
26 
30 #include <opm/core/utility/parameters/ParameterGroup.hpp>
31 #include <opm/core/linalg/LinearSolverInterface.hpp>
32 #include <dune/istl/scalarproducts.hh>
33 #include <dune/istl/operators.hh>
34 #include <dune/istl/bvector.hh>
35 #include <memory>
36 
37 namespace Opm
38 {
39 
47  {
48  typedef Dune::FieldVector<double, 1 > VectorBlockType;
49  typedef Dune::FieldMatrix<double, 1, 1> MatrixBlockType;
50  typedef Dune::BCRSMatrix <MatrixBlockType> Mat;
51  typedef Dune::BlockVector<VectorBlockType> Vector;
52 
53  public:
54 
66  NewtonIterationBlackoilCPR(const parameter::ParameterGroup& param,
67  const boost::any& parallelInformation=boost::any());
68 
75 
77  virtual int iterations () const { return iterations_; }
78 
80  virtual const boost::any& parallelInformation() const;
81 
82  private:
83 
87  template<int category=Dune::SolverCategory::sequential, class O, class P>
88  void constructPreconditionerAndSolve(O& opA, DuneMatrix& istlAe,
89  Vector& x, Vector& istlb,
90  const P& parallelInformation_arg,
91  const P& parallelInformationAe,
92  Dune::InverseOperatorResult& result) const
93  {
94  typedef Dune::ScalarProductChooser<Vector,P,category> ScalarProductChooser;
95  std::unique_ptr<typename ScalarProductChooser::ScalarProduct>
96  sp(ScalarProductChooser::construct(parallelInformation_arg));
97  // Construct preconditioner.
98  // typedef Dune::SeqILU0<Mat,Vector,Vector> Preconditioner;
99  typedef Opm::CPRPreconditioner<Mat,Vector,Vector,P> Preconditioner;
100  parallelInformation_arg.copyOwnerToAll(istlb, istlb);
101  Preconditioner precond(cpr_param_, opA.getmat(), istlAe, parallelInformation_arg,
102  parallelInformationAe);
103 
104  // TODO: Revise when linear solvers interface opm-core is done
105  // Construct linear solver.
106  // GMRes solver
107  if ( newton_use_gmres_ ) {
108  Dune::RestartedGMResSolver<Vector> linsolve(opA, *sp, precond,
109  linear_solver_reduction_, linear_solver_restart_, linear_solver_maxiter_, linear_solver_verbosity_);
110  // Solve system.
111  linsolve.apply(x, istlb, result);
112  }
113  else { // BiCGstab solver
114  Dune::BiCGSTABSolver<Vector> linsolve(opA, *sp, precond,
115  linear_solver_reduction_, linear_solver_maxiter_, linear_solver_verbosity_);
116  // Solve system.
117  linsolve.apply(x, istlb, result);
118  }
119  }
120 
121  CPRParameter cpr_param_;
122 
123  mutable int iterations_;
124  boost::any parallelInformation_;
125 
126  const bool newton_use_gmres_;
127  const double linear_solver_reduction_;
128  const int linear_solver_maxiter_;
129  const int linear_solver_restart_;
130  const int linear_solver_verbosity_;
131  };
132 
133 } // namespace Opm
134 
135 #endif // OPM_NEWTONITERATIONBLACKOILCPR_HEADER_INCLUDED
NewtonIterationBlackoilCPR(const parameter::ParameterGroup &param, const boost::any &parallelInformation=boost::any())
Definition: LinearisedBlackoilResidual.hpp:47
LinearisedBlackoilResidual::ADB::V SolutionVector
Return type for linearSolve(). A simple, non-ad vector type.
Definition: NewtonIterationBlackoilInterface.hpp:35
virtual SolutionVector computeNewtonIncrement(const LinearisedBlackoilResidual &residual) const
Definition: AdditionalObjectDeleter.hpp:22
virtual int iterations() const
Definition: NewtonIterationBlackoilCPR.hpp:77
Definition: DuneMatrix.hpp:50
Definition: NewtonIterationBlackoilCPR.hpp:46
Interface class for (linear) solvers for the fully implicit black-oil system.
Definition: NewtonIterationBlackoilInterface.hpp:31
virtual const boost::any & parallelInformation() const
Get the information about the parallelization of the grid.
CPR preconditioner.
Definition: CPRPreconditioner.hpp:324