opm-upscaling
Loading...
Searching...
No Matches
Opm::Elasticity::ElasticityUpscale< GridType, PC > Class Template Reference

The main driver class. More...

#include <elasticity_upscale.hpp>

Public Types

typedef GridType::LeafGridView::ctype ctype
 A basic number.
typedef Dune::FieldVector< double, dimNodeValue
 A vectorial node value.
typedef GridType::LeafGridView::template Codim< 1 >::Geometry::GlobalCoordinate GlobalCoordinate
 A global coordinate.
typedef GridType::LeafGridView::IndexSet LeafIndexSet
 A set of indices.
typedef GridType::LeafGridView::template Codim< 0 >::Iterator LeafIterator
 An iterator over grid cells.
typedef PC::type PCType
 Our preconditioner type.
typedef std::shared_ptr< typename PC::type > PCPtr
 A pointer to our preconditioner.

Public Member Functions

 ElasticityUpscale (const GridType &gv_, ctype tol_, ctype Escale_, const std::string &file, const std::string &rocklist, bool verbose_)
 Main constructor.
void findBoundaries (double *min, double *max)
 Find boundary coordinates.
void addMPC (Direction dir, int slavenode, const BoundaryGrid::Vertex &m)
 Add a MPC equation.
void periodicBCs (const double *min, const double *max)
 Establish periodic boundaries using the MPC approach.
void periodicBCsMortar (const double *min, const double *max, int n1, int n2, int p1, int p2)
 Establish periodic boundaries using the mortar approach.
void fixCorners (const double *min, const double *max)
 Fix corner nodes.
void assemble (int loadcase, bool matrix)
 Assemble (optionally) stiffness matrix A and load vector.
template<int comp>
void averageStress (Dune::FieldVector< ctype, comp > &sigma, const Vector &u, int loadcase)
 Calculate the average stress vector for the given loadcase.
void solve (int loadcase)
 Solve Au = b for u.
void setupSolvers (const LinSolParams &params)

Public Attributes

ASMHandler< GridTypeA
 The linear operator.
Vector u [6]
 The solution vectors.
Vector b [6]
 The load vectors.
std::vector< double > volumeFractions
 Vector holding the volume fractions for materials (grouped by SATNUM).
bool bySat
 Are volume fractions grouped by SATNUM?
double upscaledRho
 Upscaled density.

Static Public Attributes

static const int dim = GridType::dimension
 Dimension of our grid.

Detailed Description

template<class GridType, class PC>
class Opm::Elasticity::ElasticityUpscale< GridType, PC >

The main driver class.

Constructor & Destructor Documentation

◆ ElasticityUpscale()

template<class GridType, class PC>
Opm::Elasticity::ElasticityUpscale< GridType, PC >::ElasticityUpscale ( const GridType & gv_,
ctype tol_,
ctype Escale_,
const std::string & file,
const std::string & rocklist,
bool verbose_ )
inline

Main constructor.

Parameters
[in]gv_The grid to operate on
[in]tol_The tolerance to use when deciding whether or not a coordinate falls on a plane/line/point.
See also
tol
Parameters
[in]Escale_A scale value for E-moduluses to avoid numerical issues
[in]fileThe eclipse grid file
[in]rocklistIf not blank, file is a rocklist
[in]verbose_If true, give verbose output

Member Function Documentation

◆ addMPC()

template<class GridType, class PC>
void Opm::Elasticity::ElasticityUpscale< GridType, PC >::addMPC ( Direction dir,
int slavenode,
const BoundaryGrid::Vertex & m )

Add a MPC equation.

Parameters
[in]dirThe direction of the MPC
[in]slavenodeThe slave node index
[in]mThe vertices on the master grid

◆ assemble()

template<class GridType, class PC>
void Opm::Elasticity::ElasticityUpscale< GridType, PC >::assemble ( int loadcase,
bool matrix )

Assemble (optionally) stiffness matrix A and load vector.

Parameters
[in]loadcaseThe strain load case. Set to -1 to skip
[in]matrixWhether or not to assemble the matrix

◆ averageStress()

template<class GridType, class PC>
template<int comp>
void Opm::Elasticity::ElasticityUpscale< GridType, PC >::averageStress ( Dune::FieldVector< ctype, comp > & sigma,
const Vector & u,
int loadcase )

Calculate the average stress vector for the given loadcase.

Parameters
[out]sigmaThe stress vector
[in]uThe displacement vector
[in]loadcaseThe strain load case considered

◆ findBoundaries()

template<class GridType, class PC>
void Opm::Elasticity::ElasticityUpscale< GridType, PC >::findBoundaries ( double * min,
double * max )

Find boundary coordinates.

Parameters
[out]minThe miminum coordinates of the grid
[out]maxThe maximum coordinates of the grid

◆ fixCorners()

template<class GridType, class PC>
void Opm::Elasticity::ElasticityUpscale< GridType, PC >::fixCorners ( const double * min,
const double * max )

Fix corner nodes.

Parameters
[in]minThe minimum coordinates on the grid
[in]maxThe maximum coordinates on the grid

◆ periodicBCs()

template<class GridType, class PC>
void Opm::Elasticity::ElasticityUpscale< GridType, PC >::periodicBCs ( const double * min,
const double * max )

Establish periodic boundaries using the MPC approach.

Parameters
[in]minThe minimum coordinates of the grid
[in]maxThe maximum coordinates of the grid

◆ periodicBCsMortar()

template<class GridType, class PC>
void Opm::Elasticity::ElasticityUpscale< GridType, PC >::periodicBCsMortar ( const double * min,
const double * max,
int n1,
int n2,
int p1,
int p2 )

Establish periodic boundaries using the mortar approach.

Parameters
[in]minThe minimum coordinates of the grid
[in]maxThe maximum coordinates of the grid
[in]n1The number of elements on the lambda grid in the X direction
[in]n2The number of elements on the lambda grid in the Y direction
[in]p1The order of multipliers in the X direction
[in]p2The order of multipliers in the Y direction

◆ setupSolvers()

template<class GridType, class PC>
void Opm::Elasticity::ElasticityUpscale< GridType, PC >::setupSolvers ( const LinSolParams & params)
Parameters
[in]paramsThe linear solver parameters

◆ solve()

template<class GridType, class PC>
void Opm::Elasticity::ElasticityUpscale< GridType, PC >::solve ( int loadcase)

Solve Au = b for u.

Parameters
[in]loadcaseThe load case to solve

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