A class describing a quad grid. More...

#include <boundarygrid.hh>

Classes

struct  BoundedPredicate
 Predicate for checking if a vertex falls within a quads bounding box. More...
 
class  Quad
 A class describing a linear, quadrilateral element. More...
 
class  Vertex
 A class describing a 2D vertex. More...
 
struct  VertexLess
 Predicate for sorting vertices. More...
 

Public Types

typedef Dune::FieldVector
< double, 3 > 
GlobalCoordinate
 A coordinate on the underlying 3D grid. More...
 
typedef Dune::FieldVector
< double, 2 > 
FaceCoord
 A coordinate on the quad grid. More...
 

Public Member Functions

 BoundaryGrid ()
 Default (empty) constructor. More...
 
virtual ~BoundaryGrid ()
 Default (empty) destructor. More...
 
void add (const Quad &quad)
 Add a quad to the grid. More...
 
void addToColumn (size_t col, const Quad &quad)
 
Quadoperator[] (int index)
 Obtain a reference to a quad. More...
 
const Quadoperator[] (int index) const
 Obtain a const reference to a quad. More...
 
const QuadgetQuad (int col, int index) const
 
size_t size () const
 Obtain the number of quads in the grid. More...
 
size_t colSize (int i) const
 
size_t totalNodes () const
 Return the total number of nodes on the grid when known. More...
 
bool isFixed (int node) const
 Check if a given node is marked as fixed. More...
 
bool find (Vertex &res, const Vertex &coord) const
 Locate the position of a vertex on the grid. More...
 

Static Public Member Functions

static BoundaryGrid uniform (const FaceCoord &min, const FaceCoord &max, int k1, int k2, bool dc=false)
 Establish an uniform quad grid. More...
 
static void extract (FaceCoord &res, const GlobalCoordinate &coord, int dir)
 Helper function for extracting given 2D coordinates from a 3D coordinate. More...
 
static void extract (Vertex &res, const Dune::FieldVector< double, 3 > &coord, int dir)
 Helper function for extracting given 2D coordinates from a 3D vector. More...
 

Protected Member Functions

bool bilinearSolve (double epsilon, double order, const double *A, const double *B, std::vector< double > &X, std::vector< double > &Y) const
 Solves a bi-linear set of equations in x and y. A1 * x*y + A2 * x + A3 * y = A4 B1 * x*y + B2 * x + B3 * y = B4. More...
 
bool cubicSolve (double eps, double A, double B, double C, double D, std::vector< double > &X) const
 Solves the cubic equation A*x^3+B*x^2+C*x+D=0. More...
 
int Q4inv (FaceCoord &res, const Quad &q, const FaceCoord &coord, double epsZero, double epsOut) const
 Find the local coordinates of a given point within a given quad. More...
 

Protected Attributes

std::vector< Quadgrid
 Our quadrilateral elements. More...
 
std::vector< std::vector< Quad > > colGrids
 
std::vector< bool > fixNodes
 Whether or not a given node is marked as fixed. More...
 
size_t nodes
 Total number of nodes on grid. More...
 

Friends

std::ostream & operator<< (std::ostream &os, const BoundaryGrid &g)
 Print to a stream. More...
 

Detailed Description

A class describing a quad grid.

Member Typedef Documentation

typedef Dune::FieldVector<double,2> Opm::Elasticity::BoundaryGrid::FaceCoord

A coordinate on the quad grid.

typedef Dune::FieldVector<double,3> Opm::Elasticity::BoundaryGrid::GlobalCoordinate

A coordinate on the underlying 3D grid.

Constructor & Destructor Documentation

Opm::Elasticity::BoundaryGrid::BoundaryGrid ( )
inline

Default (empty) constructor.

virtual Opm::Elasticity::BoundaryGrid::~BoundaryGrid ( )
inlinevirtual

Default (empty) destructor.

Member Function Documentation

void Opm::Elasticity::BoundaryGrid::add ( const Quad quad)

Add a quad to the grid.

Parameters
[in]quadThe quad to add

Referenced by Opm::Elasticity::IMPL_FUNC().

void Opm::Elasticity::BoundaryGrid::addToColumn ( size_t  col,
const Quad quad 
)
inline

References colGrids.

Referenced by Opm::Elasticity::IMPL_FUNC().

bool Opm::Elasticity::BoundaryGrid::bilinearSolve ( double  epsilon,
double  order,
const double *  A,
const double *  B,
std::vector< double > &  X,
std::vector< double > &  Y 
) const
protected

Solves a bi-linear set of equations in x and y. A1 * x*y + A2 * x + A3 * y = A4 B1 * x*y + B2 * x + B3 * y = B4.

Parameters
[in]epsilonThe tolerance for equality checks with zero
[in]orderThe expected order of the solution (used for unique checks)
[in]AThe coefficients of the first equation
[in]BThe coefficients of the second equation
[out]XThe first component of the solutions
[out]YThe second component of the solutions
size_t Opm::Elasticity::BoundaryGrid::colSize ( int  i) const
inline

References colGrids.

Referenced by Opm::Elasticity::IMPL_FUNC().

bool Opm::Elasticity::BoundaryGrid::cubicSolve ( double  eps,
double  A,
double  B,
double  C,
double  D,
std::vector< double > &  X 
) const
protected

Solves the cubic equation A*x^3+B*x^2+C*x+D=0.

Parameters
[in]epsThe tolerance for equality checks with zero
[in]AEquation coefficient
[in]BEquation coefficient
[in]CEquation coefficient
[in]DEquation coefficient
[out]XThe obtained solutions
static void Opm::Elasticity::BoundaryGrid::extract ( FaceCoord res,
const GlobalCoordinate coord,
int  dir 
)
static

Helper function for extracting given 2D coordinates from a 3D coordinate.

Parameters
[in]coordThe 3D coordinates of the vertex
[in]dirThe direction to ignore
[out]resThe resulting coordinates

Referenced by Opm::Elasticity::HexGeometry< 2, cdim, GridImp >::HexGeometry(), and Opm::Elasticity::IMPL_FUNC().

static void Opm::Elasticity::BoundaryGrid::extract ( Vertex res,
const Dune::FieldVector< double, 3 > &  coord,
int  dir 
)
static

Helper function for extracting given 2D coordinates from a 3D vector.

Parameters
[in]coordThe 3D coordinates of the vertex
[in]dirThe direction to ignore
[out]resThe resulting coordinates
bool Opm::Elasticity::BoundaryGrid::find ( Vertex res,
const Vertex coord 
) const

Locate the position of a vertex on the grid.

Parameters
[in]coordThe coordinate of the vertex
[out]resThe resulting coordinates

Referenced by Opm::Elasticity::IMPL_FUNC().

const Quad& Opm::Elasticity::BoundaryGrid::getQuad ( int  col,
int  index 
) const
inline

References colGrids.

Referenced by Opm::Elasticity::IMPL_FUNC().

bool Opm::Elasticity::BoundaryGrid::isFixed ( int  node) const
inline

Check if a given node is marked as fixed.

Parameters
[in]nodeThe requested node
Returns
Whether or not the node is marked as fixed

References fixNodes.

Quad& Opm::Elasticity::BoundaryGrid::operator[] ( int  index)
inline

Obtain a reference to a quad.

Parameters
[in]indexThe index of the requested quad
Returns
A reference to the requested quad

References grid.

const Quad& Opm::Elasticity::BoundaryGrid::operator[] ( int  index) const
inline

Obtain a const reference to a quad.

Parameters
[in]indexThe index of the requested quad
Returns
A const reference to the requested quad

References grid.

int Opm::Elasticity::BoundaryGrid::Q4inv ( FaceCoord res,
const Quad q,
const FaceCoord coord,
double  epsZero,
double  epsOut 
) const
protected

Find the local coordinates of a given point within a given quad.

Parameters
[in]qThe quad to search within
[in]coordThe coordinates to search for
[in]epsZeroThe tolerance for equality checks with zero
[in]epsOutThe tolerance check for outside checks
[out]resThe obtained result
size_t Opm::Elasticity::BoundaryGrid::size ( ) const
inline

Obtain the number of quads in the grid.

References grid.

Referenced by Opm::Elasticity::IMPL_FUNC(), and Opm::Elasticity::renumber().

size_t Opm::Elasticity::BoundaryGrid::totalNodes ( ) const
inline

Return the total number of nodes on the grid when known.

See also
uniform

References nodes.

static BoundaryGrid Opm::Elasticity::BoundaryGrid::uniform ( const FaceCoord min,
const FaceCoord max,
int  k1,
int  k2,
bool  dc = false 
)
static

Establish an uniform quad grid.

Parameters
[in]minLower left corner
[in]maxUpper right corner
[in]k1Number of elements in the first direction
[in]k2Number of elements in the second direction
[in]dcIf true, order quads according to dune conventions
Returns
A quad grid spanning the given area. Nodes are numbered

Friends And Related Function Documentation

std::ostream& operator<< ( std::ostream &  os,
const BoundaryGrid g 
)
friend

Print to a stream.

Member Data Documentation

std::vector<std::vector<Quad> > Opm::Elasticity::BoundaryGrid::colGrids
protected

Referenced by addToColumn(), colSize(), and getQuad().

std::vector<bool> Opm::Elasticity::BoundaryGrid::fixNodes
protected

Whether or not a given node is marked as fixed.

Referenced by isFixed().

std::vector<Quad> Opm::Elasticity::BoundaryGrid::grid
protected

Our quadrilateral elements.

Referenced by operator[](), and size().

size_t Opm::Elasticity::BoundaryGrid::nodes
protected

Total number of nodes on grid.

Referenced by totalNodes().


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