A light weight node. More...
#include <mrtr_point.H>
Public Member Functions | |
Point (const int id, const double *xi, int out) | |
Constructor. More... | |
virtual | ~Point () |
Destructor. | |
int | OutLevel () |
Return the level of output written to stdout ( 0 - 10 ) | |
void | Print () const |
Print this node to stdout. | |
int | Id () |
Return id of this point. | |
const double * | Xi () |
Return view of segment local coordinates of this point (2D) | |
const double * | XCoords () |
Return view of global coordinates of this point (3D) More... | |
Teuchos::RCP< MOERTEL::Node > | Node () |
Return view of Node. More... | |
bool | SetXi (const double *xi) |
Set segment local coordinates of this point (2D) in a segment. | |
bool | SetNode (MOERTEL::Node *node) |
Set a Node to this point. More... | |
void | StoreFunctionValues (int place, double *val, int valdim) |
Store finite element function values at the Point 's coordinate Xi() More... | |
std::vector< double > * | FunctionValues () |
Return view of function values stored in this Point. More... | |
A light weight node.
A light weight version of a node
This class defines a point on a segment. It is a light weight version of a node. It is used in the integration of 2D interfaces where the mortar side is imprinted to the slave side. The overlap between a mortar and a slave segment leads to a polygon defined by points on the slave segment. The polygon is then discretized by triangle finite elements (eventually adding more points) to perform the integration on the polygon region. A point might therefore become a node of the polygon discretization and therefore has capabilities to store a Node class.
The MOERTEL::Point class supports the ostream& operator <<
MOERTEL::Point::Point | ( | const int | id, |
const double * | xi, | ||
int | out | ||
) |
Constructor.
Constructs an instance of this class.
Note that this is not a collective call as points shall only have one owning process.
id | : A unique positive point id. |
xi | : Coordinates of point in a segment (2D) |
out | : Level of output information written to stdout ( 0 - 10 ) |
|
inline |
Return view of function values stored in this Point.
Returns a view of the function values that were stored in this Point using StoreFunctionValues
|
inline |
|
inline |
void MOERTEL::Point::StoreFunctionValues | ( | int | place, |
double * | val, | ||
int | valdim | ||
) |
Store finite element function values at the Point 's coordinate Xi()
place | : Place in internal data structure where to store function values. place=0 is used to store trace space function values of the slave segment. place=1 is used to store Lagrange multiplier space function values of the slave segment. place=2 is used to store trace space function values of the master segment. |
val | : Vector of length valdim holding function values |
valdim | : Dimension of val |
|
inline |