Opm::NewtonIterationBlackoilInterleavedImpl< np, ScalarT > Class Template Reference

This class solves the fully implicit black-oil system by solving the reduced system (after eliminating well variables) as a block-structured matrix (one block for all cell variables) for a fixed number of cell variables np . More...

Inheritance diagram for Opm::NewtonIterationBlackoilInterleavedImpl< np, ScalarT >:
Opm::NewtonIterationBlackoilInterface

List of all members.

Public Types

typedef
NewtonIterationBlackoilInterface::SolutionVector 
SolutionVector
 Return type for linearSolve(). A simple, non-ad vector type.

Public Member Functions

 NewtonIterationBlackoilInterleavedImpl (const NewtonIterationBlackoilInterleavedParameters &param, const boost::any &parallelInformation_arg=boost::any())
 Construct a system solver.
int iterations () const
 Solve the system of linear equations Ax = b, with A being the combined derivative matrix of the residual and b being the residual itself.
const boost::any & parallelInformation () const
void formInterleavedSystem (const std::vector< LinearisedBlackoilResidual::ADB > &eqs, Mat &istlA) const
SolutionVector computeNewtonIncrement (const LinearisedBlackoilResidual &residual) const
 Solve the linear system Ax = b, with A being the combined derivative matrix of the residual and b being the residual itself.

Protected Attributes

ISTLSolverType istlSolver_
NewtonIterationBlackoilInterleavedParameters parameters_

Detailed Description

template<int np, class ScalarT = double>
class Opm::NewtonIterationBlackoilInterleavedImpl< np, ScalarT >

This class solves the fully implicit black-oil system by solving the reduced system (after eliminating well variables) as a block-structured matrix (one block for all cell variables) for a fixed number of cell variables np .


Constructor & Destructor Documentation

template<int np, class ScalarT = double>
Opm::NewtonIterationBlackoilInterleavedImpl< np, ScalarT >::NewtonIterationBlackoilInterleavedImpl ( const NewtonIterationBlackoilInterleavedParameters param,
const boost::any &  parallelInformation_arg = boost::any() 
) [inline]

Construct a system solver.

Parameters:
[in] param parameters controlling the behaviour of the linear solvers
[in] parallelInformation In the case of a parallel run with dune-istl the information about the parallelization.

Member Function Documentation

template<int np, class ScalarT = double>
SolutionVector Opm::NewtonIterationBlackoilInterleavedImpl< np, ScalarT >::computeNewtonIncrement ( const LinearisedBlackoilResidual residual  )  const [inline, virtual]

Solve the linear system Ax = b, with A being the combined derivative matrix of the residual and b being the residual itself.

Parameters:
[in] residual residual object containing A and b.
Returns:
the solution x

Implements Opm::NewtonIterationBlackoilInterface.

template<int np, class ScalarT = double>
void Opm::NewtonIterationBlackoilInterleavedImpl< np, ScalarT >::formInterleavedSystem ( const std::vector< LinearisedBlackoilResidual::ADB > &  eqs,
Mat &  istlA 
) const [inline]

Go through all jacobians, and insert in correct spot

The straight forward way to do this would be to run through each element in the output matrix, and set all block entries by gathering from all "input matrices" (derivatives).

A faster alternative is to instead run through each "input matrix" and insert its elements in the correct spot in the output matrix.

template<int np, class ScalarT = double>
int Opm::NewtonIterationBlackoilInterleavedImpl< np, ScalarT >::iterations (  )  const [inline, virtual]

Solve the system of linear equations Ax = b, with A being the combined derivative matrix of the residual and b being the residual itself.

Parameters:
[in] residual residual object containing A and b.
Returns:
the solution x

Implements Opm::NewtonIterationBlackoilInterface.

template<int np, class ScalarT = double>
const boost::any& Opm::NewtonIterationBlackoilInterleavedImpl< np, ScalarT >::parallelInformation (  )  const [inline, virtual]

Get the information about the parallelization of the grid.

Implements Opm::NewtonIterationBlackoilInterface.


The documentation for this class was generated from the following file:

Generated on 26 Mar 2018 by  doxygen 1.6.1