A class describing a quad grid.
More...
#include <boundarygrid.hh>
|
typedef Dune::FieldVector< double, 3 > | GlobalCoordinate |
| A coordinate on the underlying 3D grid.
|
|
typedef Dune::FieldVector< double, 2 > | FaceCoord |
| A coordinate on the quad grid.
|
|
|
| BoundaryGrid () |
| Default (empty) constructor.
|
|
virtual | ~BoundaryGrid () |
| Default (empty) destructor.
|
|
void | add (const Quad &quad) |
| Add a quad to the grid.
|
|
void | addToColumn (size_t col, const Quad &quad) |
|
Quad & | operator[] (int index) |
| Obtain a reference to a quad.
|
|
const Quad & | operator[] (int index) const |
| Obtain a const reference to a quad.
|
|
const Quad & | getQuad (int col, int index) const |
|
size_t | size () const |
| Obtain the number of quads in the grid.
|
|
size_t | colSize (int i) const |
|
size_t | totalNodes () const |
| Return the total number of nodes on the grid when known.
|
|
bool | isFixed (int node) const |
| Check if a given node is marked as fixed.
|
|
bool | find (Vertex &res, const Vertex &coord) const |
| Locate the position of a vertex on the grid.
|
|
|
static BoundaryGrid | uniform (const FaceCoord &min, const FaceCoord &max, int k1, int k2, bool dc=false) |
| Establish an uniform quad grid.
|
|
static void | extract (FaceCoord &res, const GlobalCoordinate &coord, int dir) |
| Helper function for extracting given 2D coordinates from a 3D coordinate.
|
|
static void | extract (Vertex &res, const Dune::FieldVector< double, 3 > &coord, int dir) |
| Helper function for extracting given 2D coordinates from a 3D vector.
|
|
|
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.
|
|
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.
|
|
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.
|
|
|
std::vector< Quad > | grid |
| Our quadrilateral elements.
|
|
std::vector< std::vector< Quad > > | colGrids |
|
std::vector< bool > | fixNodes |
| Whether or not a given node is marked as fixed.
|
|
size_t | nodes |
| Total number of nodes on grid.
|
|
|
std::ostream & | operator<< (std::ostream &os, const BoundaryGrid &g) |
| Print to a stream.
|
|
A class describing a quad grid.
◆ add()
void Opm::Elasticity::BoundaryGrid::add |
( |
const Quad & |
quad | ) |
|
Add a quad to the grid.
- Parameters
-
◆ bilinearSolve()
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] | epsilon | The tolerance for equality checks with zero |
[in] | order | The expected order of the solution (used for unique checks) |
[in] | A | The coefficients of the first equation |
[in] | B | The coefficients of the second equation |
[out] | X | The first component of the solutions |
[out] | Y | The second component of the solutions |
◆ cubicSolve()
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] | eps | The tolerance for equality checks with zero |
[in] | A | Equation coefficient |
[in] | B | Equation coefficient |
[in] | C | Equation coefficient |
[in] | D | Equation coefficient |
[out] | X | The obtained solutions |
◆ extract() [1/2]
Helper function for extracting given 2D coordinates from a 3D coordinate.
- Parameters
-
[in] | coord | The 3D coordinates of the vertex |
[in] | dir | The direction to ignore |
[out] | res | The resulting coordinates |
◆ extract() [2/2]
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] | coord | The 3D coordinates of the vertex |
[in] | dir | The direction to ignore |
[out] | res | The resulting coordinates |
◆ find()
bool Opm::Elasticity::BoundaryGrid::find |
( |
Vertex & |
res, |
|
|
const Vertex & |
coord |
|
) |
| const |
Locate the position of a vertex on the grid.
- Parameters
-
[in] | coord | The coordinate of the vertex |
[out] | res | The resulting coordinates |
◆ isFixed()
bool Opm::Elasticity::BoundaryGrid::isFixed |
( |
int |
node | ) |
const |
|
inline |
Check if a given node is marked as fixed.
- Parameters
-
[in] | node | The requested node |
- Returns
- Whether or not the node is marked as fixed
◆ operator[]() [1/2]
Quad & Opm::Elasticity::BoundaryGrid::operator[] |
( |
int |
index | ) |
|
|
inline |
Obtain a reference to a quad.
- Parameters
-
[in] | index | The index of the requested quad |
- Returns
- A reference to the requested quad
◆ operator[]() [2/2]
const Quad & Opm::Elasticity::BoundaryGrid::operator[] |
( |
int |
index | ) |
const |
|
inline |
Obtain a const reference to a quad.
- Parameters
-
[in] | index | The index of the requested quad |
- Returns
- A const reference to the requested quad
◆ Q4inv()
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] | q | The quad to search within |
[in] | coord | The coordinates to search for |
[in] | epsZero | The tolerance for equality checks with zero |
[in] | epsOut | The tolerance check for outside checks |
[out] | res | The obtained result |
◆ totalNodes()
size_t Opm::Elasticity::BoundaryGrid::totalNodes |
( |
| ) |
const |
|
inline |
Return the total number of nodes on the grid when known.
- See also
- uniform
◆ uniform()
Establish an uniform quad grid.
- Parameters
-
[in] | min | Lower left corner |
[in] | max | Upper right corner |
[in] | k1 | Number of elements in the first direction |
[in] | k2 | Number of elements in the second direction |
[in] | dc | If true, order quads according to dune conventions |
- Returns
- A quad grid spanning the given area. Nodes are numbered
The documentation for this class was generated from the following files: