Moertel  Development
 All Classes Namespaces Files Functions Enumerations Friends Pages
Protected Member Functions | Protected Attributes | List of all members
MOERTEL::Node Class Reference

A class to construct a single node. More...

#include <mrtr_pnode.H>

Inheritance diagram for MOERTEL::Node:
Inheritance graph
[legend]

Public Member Functions

 Node (int Id, const double *x, int ndof, const int *dof, bool isonboundary, int out)
 Constructor. More...
 
 Node (int out)
 Constructor. More...
 
 Node (const MOERTEL::Node &old)
 Copy-Constructor. More...
 
virtual ~Node ()
 Destructor.
 
int OutLevel ()
 Return the level of output written to stdout ( 0 - 10 )
 
int Id () const
 Return unique and positive id of this node.
 
bool Print () const
 Print this node to stdout.
 
void Reset ()
 Reset the internal state of this node and clear integrated values. More...
 
const double * XCoords () const
 Returns the view (a vector of length 3) to the current coordinates of this node.
 
bool SetX (double *x)
 Update the coordinates of this Node. More...
 
const double * Normal () const
 Return pointer to vector of length 3 of normal of this node.
 
bool SetN (double *n)
 Store a normal in this node. More...
 
int Ndof () const
 Return the number of degrees of freedom on this node.
 
int Nlmdof () const
 Return the number of Lagrange mutlipliers on this node. More...
 
const int * Dof () const
 Return view of degrees of freedom on this node. More...
 
const int * LMDof () const
 Return view of the Lagrange multipliers on this node. More...
 
int Nseg () const
 Return the number of segments adjacent to this node. More...
 
int * SegmentIds ()
 Returns a view to the vector of segment ids of segments adjacent to this node. More...
 
MOERTEL::Segment ** Segments ()
 Returns a view to the vector of pointers to segments adjacent to this node. More...
 
bool AddSegment (int sid)
 Adds a segment id to the list of segments adjacent to this node. More...
 
double * Pack (int *size)
 Packs most information stored in this node to a double vector for communication with MPI. More...
 
bool UnPack (double *pack)
 Unpacks information stored in a double vector to an instance of Node. More...
 
bool BuildAveragedNormal ()
 Build averaged normal from adjacent segments at this node.
 
bool GetPtrstoSegments (MOERTEL::Interface &interface)
 Construct vector of pointers to segments from adjacent segments id list.
 
bool SetProjectedNode (MOERTEL::ProjectedNode *pnode)
 Store a pointer to a projected node. More...
 
Teuchos::RCP
< MOERTEL::ProjectedNode > * 
GetProjectedNode (int &length)
 Returns a view of all projected nodes this node owns. More...
 
Teuchos::RCP
< MOERTEL::ProjectedNode
GetProjectedNode ()
 Returns a view of the projection of this node. More...
 
bool SetLagrangeMultiplierId (int LMId)
 Add a Lagrange multiplier id to this node's list.
 
void AddDValue (double val, int col)
 Add a value to the 'D' map of this node. More...
 
void AddMValue (double val, int col)
 Add a value to the 'M' map of this node. More...
 
void AddMmodValue (int row, double val, int col)
 Add a value to the 'M' map of this node. More...
 
Teuchos::RCP< std::map< int,
double > > 
GetD ()
 Get view of the 'D' map of this node.
 
Teuchos::RCP< std::map< int,
double > > 
GetM ()
 Get view of the 'M' map of this node.
 
Teuchos::RCP< std::vector
< std::map< int, double > > > 
GetMmod ()
 Get view of the 'Mmod' map of this node.
 
void SetCorner ()
 Set the flag indicating that node is corner node of 1D interface.
 
bool IsCorner () const
 Query the flag indicating that node is corner node of 1D interface.
 
bool IsOnBoundary () const
 Query the flag indicating that node is on the boundary of 2D interface.
 
int NSupportSet () const
 Get number of nodes that support this boundary node. More...
 
void AddSupportedByNode (MOERTEL::Node *suppnode)
 Add an internal node to the map of nodes supporting this node. More...
 
std::map< int, MOERTEL::Node * > & GetSupportedByNode ()
 Get the map of nodes supporting this node. More...
 
double Gap ()
 Return distance between projection and source node. More...
 
void SetGap (double gap)
 Sets the distance between projection and source node. More...
 

Protected Member Functions

Node operator= (const Node &old)
 

Protected Attributes

int Id_
 
double x_ [3]
 
double n_ [3]
 
int outputlevel_
 
bool iscorner_
 
bool isonboundary_
 
std::map< int, MOERTEL::Node * > supportedby_
 
std::vector< int > dof_
 
std::vector< int > LMdof_
 
std::vector< int > seg_
 
std::vector< MOERTEL::Segment * > segptr_
 
Teuchos::RCP< std::map< int,
double > > 
Drow_
 
Teuchos::RCP< std::map< int,
double > > 
Mrow_
 
Teuchos::RCP< std::vector
< std::map< int, double > > > 
Mmodrow_
 
std::vector< Teuchos::RCP
< MOERTEL::ProjectedNode > > 
pnode_
 
double gap_
 

Detailed Description

A class to construct a single node.

A class to construct a projection of a node on some segment.

A class to construct a single node

Date
Last update do Doxygen: 20-March-06

This class defines a node on one side of a (non-)conforming interface.

The MOERTEL::Node class supports the ostream& operator <<

Author
Glen Hansen (gahan.nosp@m.se@s.nosp@m.andia.nosp@m..gov)
Date
Last update do Doxygen: 10-July-2010

Constructor & Destructor Documentation

MOERTEL::Node::Node ( int  Id,
const double *  x,
int  ndof,
const int *  dof,
bool  isonboundary,
int  out 
)
explicit

Constructor.

Constructs an instance of this class.
Note that this is not a collective call as nodes shall only have one owning process.

Parameters
Id: A unique positive node id. Does not need to be continous among nodes
x: ptr to vector of length 3 holding coordinates of node. In the 2D problem case x[2] has to be 0.0
ndof: Number of degrees of freedom on this node, e.g. 2 in 2D and 3 in 3D
dof: Pointer to vector of length ndof holding degrees of freedom on this node
isonboundary: True if node is on the boundary of a 2D interface (3D problem), or a 1D interface (2D problem). Note that the number of lmdofs for nodes on the boundary are set to zero.
out: Level of output information written to stdout ( 0 - 10 )
MOERTEL::Node::Node ( int  out)
explicit

Constructor.

Constructs an empty instance of this class.
Note that this is not a collective call as nodes shall only have one owning process.
This constructor is used internally together with Pack(int* size) and UnPack(double* pack)
to allow for communication of nodes among processes

Parameters
out: Level of output information written to stdout ( 0 - 10 )
MOERTEL::Node::Node ( const MOERTEL::Node old)

Copy-Constructor.

Makes a deep copy of an existing node

Parameters
old: Node to be copied

References Id().

Member Function Documentation

void MOERTEL::Node::AddDValue ( double  val,
int  col 
)

Add a value to the 'D' map of this node.

The 'D' map is later assembled to the D matrix. Note that the 'D' map in this Node is scalar and the column id is a node id while the row map of the matrix D might have several Lagrange mutlipliers per node.

val : Value to be added col : nodal column id of the value added

void MOERTEL::Node::AddMmodValue ( int  row,
double  val,
int  col 
)

Add a value to the 'M' map of this node.

The 'M' map is later assembled to the M matrix. Note that Mmod here is a vector

row : Local LM id to add to (rowwise) val : Value to be added col : nodal column id of the value added

void MOERTEL::Node::AddMValue ( double  val,
int  col 
)

Add a value to the 'M' map of this node.

The 'M' map is later assembled to the M matrix. Note that the 'M' map in this Node is scalar and the column id is a node id while the row map of the matrix M might have several Lagrange mutlipliers per node.

val : Value to be added col : nodal column id of the value added

Referenced by MOERTEL::Integrator::Assemble().

bool MOERTEL::Node::AddSegment ( int  sid)

Adds a segment id to the list of segments adjacent to this node.

This method checks whether segment already exists in the adjacency

void MOERTEL::Node::AddSupportedByNode ( MOERTEL::Node suppnode)
inline

Add an internal node to the map of nodes supporting this node.

If this Node is on the boundary, its shape functions are added to the support of topologically close internal nodes. This method adds an internal node to the map of nodes supporting this node

References Id().

const int* MOERTEL::Node::Dof ( ) const
inline

Return view of degrees of freedom on this node.

Returns a view of the vector of length Ndof()

Referenced by MOERTEL::Integrator::Assemble(), and MOERTEL::Integrator::Assemble_2D_Mod().

double MOERTEL::Node::Gap ( )
inline

Return distance between projection and source node.

Returns the gap distance between the source node and projected image on the target segment. This distance is along the projection vector, and will be negative if the projected image is on the "wrong" side of the interface that contains the source node.

Teuchos::RCP< MOERTEL::ProjectedNode > * MOERTEL::Node::GetProjectedNode ( int &  length)

Returns a view of all projected nodes this node owns.

Returns a vector of RefCountPtr<MOERTEL::ProjectedNode> to all projected nodes this node owns. When using a continous normal field as projection method, length of this vector is 1. When using an orthogonal projection (1D interfaces only) a node can have more then 1 projection

length : returns length of RefCountPtr<MOERTEL::ProjectedNode>*

Teuchos::RCP< MOERTEL::ProjectedNode > MOERTEL::Node::GetProjectedNode ( )

Returns a view of the projection of this node.

Returns a RefCountPtr<MOERTEL::ProjectedNode> to the projected node. RefCountPtr<MOERTEL::ProjectedNode> is Teuchos::null if this Node does not have a projection onto an opposing segment.

References MOERTEL::length().

std::map<int,MOERTEL::Node*>& MOERTEL::Node::GetSupportedByNode ( )
inline

Get the map of nodes supporting this node.

If this Node is on the boundary, its shape functions are added to the support of topologically close internal nodes. This method returns the map of nodes supporting this node

Referenced by MOERTEL::Integrator::Assemble(), and MOERTEL::Integrator::Assemble_2D_Mod().

const int* MOERTEL::Node::LMDof ( ) const
inline

Return view of the Lagrange multipliers on this node.

Returns a view of the vector of length Nlmdof()

Referenced by MOERTEL::Integrator::Assemble().

int MOERTEL::Node::Nlmdof ( ) const
inline

Return the number of Lagrange mutlipliers on this node.

Returns the number of Lagrange mutlipliers on this node which is normally equal to the number of degrees of freedom on this node on the slave side of an interface and equal to zero on the master side of an interface

Referenced by MOERTEL::Integrator::Assemble().

int MOERTEL::Node::Nseg ( ) const
inline

Return the number of segments adjacent to this node.

Returns the number of MOERTEL::Segment adjacent to this node.

Warning
This information is computed in MOERTEL::Interface::Complete() and does not exist before
int MOERTEL::Node::NSupportSet ( ) const
inline

Get number of nodes that support this boundary node.

If this Node is on the boundary, its shape functions are added to the support of topologically close internal nodes. This method returns # of internal nodes supporting this node. It returns 0 for all internal nodes.

double * MOERTEL::Node::Pack ( int *  size)

Packs most information stored in this node to a double vector for communication with MPI.

This method stores most but not all information stored in a Node in a double vector so it can be send easily using MPI. The method Unpack(double* pack) is then used to unpack the information into an empty Node instance on the receiving side. It is used to provide processes that are member of an interface-local Epetra_Comm a redundant map of nodes on that interface.

Parameters
size: Returns the length of the double vector in size
Returns
pointer to double vector containing node information

References MOERTEL::ReportError().

void MOERTEL::Node::Reset ( )

Reset the internal state of this node and clear integrated values.

This function resets the node and prepares it to be re-integrated. In nonlinear problems, the interfaces (and nodes on the interface) may move relative to each other, each Newton iteration. Calling Reset() clears the state of each node in preparation for resetting the nodal coordinates, and re-integration each Newton iteration.

int* MOERTEL::Node::SegmentIds ( )
inline

Returns a view to the vector of segment ids of segments adjacent to this node.

Returns a view to a vector of length Nseg()

Warning
This information is computed in MOERTEL::Interface::Complete() and does not exist before
MOERTEL::Segment** MOERTEL::Node::Segments ( )
inline

Returns a view to the vector of pointers to segments adjacent to this node.

Returns a view to a vector of length Nseg()

Warning
This information is computed in MOERTEL::Interface::Complete() and does not exist before
void MOERTEL::Node::SetGap ( double  gap)
inline

Sets the distance between projection and source node.

Sets the gap distance between the source node and projected image on the target segment. This distance is along the projection vector, and will be negative if the projected image is on the "wrong" side of the interface that contains the source node.

bool MOERTEL::Node::SetN ( double *  n)
inline

Store a normal in this node.

This class makes a deep copy of the vector n

Parameters
n: pointer to vector of length 3 holding an outward normal at this node to be stored (deep copy) in this instance
bool MOERTEL::Node::SetProjectedNode ( MOERTEL::ProjectedNode pnode)

Store a pointer to a projected node.

Stores a vector and takes ownership of a MOERTEL::ProjectedNode which is this Node projected onto an opposing interface surface

bool MOERTEL::Node::SetX ( double *  x)
inline

Update the coordinates of this Node.

Parameters
x: a vector of length 3 holding the new coordinates of this node
bool MOERTEL::Node::UnPack ( double *  pack)

Unpacks information stored in a double vector to an instance of Node.

Unpacks a double vector that was created using Pack(int* size) and communicated
into an empty instance of Node. The length of the vector is stored in the vector itself and consistency of information is checked.

Parameters
pack: Returns true upon success

References MOERTEL::ReportError().


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