Moertel  Development
Public Types | Public Member Functions | Protected Attributes | List of all members
MOERTEL::Segment Class Referenceabstract

A virtual class to serve as base class for different types of interface segmentsconstruct a single interface. More...

#include <mrtr_segment.H>

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

Public Types

enum  SegmentType { seg_none, seg_Linear1D, seg_BiLinearQuad, seg_BiLinearTri }
 Type of segment. More...
 

Public Member Functions

 Segment (int id, int nnode, int *nodeId, int outlevel)
 Standard Constructor. More...
 
 Segment (int id, const std::vector< int > &nodeId, int outlevel)
 
 Segment (int outlevel)
 Empty Constructor. More...
 
 Segment (MOERTEL::Segment &old)
 Copy Constructor. More...
 
virtual ~Segment ()
 Destructor. More...
 
int OutLevel ()
 Return level of output to be generated by this class (0-10) More...
 
int Id () const
 Return unique id of this Segment. More...
 
int Nnode () const
 Return number of nodes attached to this Segment. More...
 
MOERTEL::Segment::SegmentType Type () const
 Return type of Segment. More...
 
const int * NodeIds () const
 Return view of node ids of nodes attached to this Segment. More...
 
MOERTEL::Node ** Nodes ()
 Return pointer to vector of length Nnode() of pointers to Nodes attached to this Segment. More...
 
int Nfunctions ()
 Return number of functions defined on this Segment. More...
 
MOERTEL::Function::FunctionType FunctionType (int id)
 Return FunctionType of a function with the Id id. More...
 
bool SetFunction (int id, MOERTEL::Function *func)
 Attach a function to this Segment. More...
 
bool EvaluateFunction (int id, const double *xi, double *val, int valdim, double *deriv)
 Evaluate a function with a certain id. More...
 
double * BuildNormalAtNode (int nid)
 Build normal at a node adjacent to this Segment. More...
 
bool GetPtrstoNodes (MOERTEL::Interface &interface)
 Get pointers to Nodes attached to this Segment from the Interface this Segment resides on. More...
 
bool GetPtrstoNodes (std::vector< MOERTEL::Node * > &nodes)
 Get pointers to Nodes attached to this Segment from a vector of Node pointers. More...
 
virtual bool Print () const
 Print this Segment. More...
 
int GetLocalNodeId (int nid)
 Get segment-local node id from global node id nid. More...
 
virtual MOERTEL::SegmentClone ()=0
 Deep copy the derived class and return pointer to it. More...
 
virtual int * Pack (int *size)=0
 Pack some data from this class to an int vector of length size so it can be communicated using MPI. More...
 
virtual bool UnPack (int *pack)=0
 Unpack some data an int vector and store data in this class. More...
 
virtual double * BuildNormal (double *xi)=0
 Build an outward normal at segment coordinates xi. More...
 
virtual double Area ()=0
 Compute and return the area of this Segment. More...
 
virtual double Metric (double *xi, double g[], double G[][3])=0
 Build the basis vectors and metric tensor at a given local coord in this segment. More...
 
virtual bool LocalCoordinatesOfNode (int lid, double *xi)=0
 Get local coords of a node attached to this segment with local node Id lid. More...
 

Protected Attributes

int Id_
 
int outputlevel_
 
SegmentType stype_
 
std::vector< int > nodeId_
 
std::vector< MOERTEL::Node * > nodeptr_
 
std::map< int, Teuchos::RCP< MOERTEL::Function > > functions_
 

Detailed Description

A virtual class to serve as base class for different types of interface segmentsconstruct a single interface.

A virtual class as a basis for different types of interface segments

Date
Last update to Doxygen: 15-Dec-05

This class serves as a (not pure) virtual base class to several types of interface segments.

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

Author
Glen Hansen (gahan.nosp@m.se@s.nosp@m.andia.nosp@m..gov)

Member Enumeration Documentation

◆ SegmentType

Type of segment.

Parameters
seg_none: default value
seg_Linear1D: linear 1D segment with 2 nodes
seg_BiLinearQuad: linear 2D triangle with 3 nodes
seg_BiLinearTri: linear 2D quadrilateral with 4 nodes

Constructor & Destructor Documentation

◆ Segment() [1/3]

MOERTEL::Segment::Segment ( int  id,
int  nnode,
int *  nodeId,
int  outlevel 
)

Standard Constructor.

Parameters
Id: A unique positive Segment id. Does not need to be continous among segments
nnode: Number of nodes this segment is attached to
nodeId: Pointer to vector length nnode holding unique positive node ids of nodes this segment is attached to
outlevel: Level of output to stdout to be generated by this class (0-10)

◆ Segment() [2/3]

MOERTEL::Segment::Segment ( int  outlevel)

Empty Constructor.

This constructor is used together with Pack and UnPack for communicating segments

Parameters
outlevel: Level of output to stdout to be generated by this class (0-10)

◆ Segment() [3/3]

MOERTEL::Segment::Segment ( MOERTEL::Segment old)

Copy Constructor.

Makes a deep copy of a Segment

References MOERTEL::ReportError().

◆ ~Segment()

MOERTEL::Segment::~Segment ( )
virtual

Destructor.

Member Function Documentation

◆ Area()

virtual double MOERTEL::Segment::Area ( )
pure virtual

◆ BuildNormal()

virtual double* MOERTEL::Segment::BuildNormal ( double *  xi)
pure virtual

Build an outward normal at segment coordinates xi.

Implemented in MOERTEL::Segment_BiLinearQuad, MOERTEL::Segment_BiLinearTri, and MOERTEL::Segment_Linear1D.

◆ BuildNormalAtNode()

double * MOERTEL::Segment::BuildNormalAtNode ( int  nid)

Build normal at a node adjacent to this Segment.

Parameters
nid: global unique node id

◆ Clone()

virtual MOERTEL::Segment* MOERTEL::Segment::Clone ( )
pure virtual

Deep copy the derived class and return pointer to it.

Implemented in MOERTEL::Segment_BiLinearQuad, MOERTEL::Segment_BiLinearTri, and MOERTEL::Segment_Linear1D.

Referenced by MOERTEL::Interface::AddSegment().

◆ EvaluateFunction()

bool MOERTEL::Segment::EvaluateFunction ( int  id,
const double *  xi,
double *  val,
int  valdim,
double *  deriv 
)

Evaluate a function with a certain id.

Will evaluate the function with Id id at a given local coordinate

Parameters
id(in): unique function id
xi(in): Segment local coordinates where to evaluate the function
val(out): Vector holding function values at xi on output. If NULL on input, function will not evaluate values.
valdim(in): length of val
deriv(out): Vector holding function derivatives at xi on output, should be of length 2*valdim in most cases. If NULL on input, function will not evaluate derivatives.

References MOERTEL::ReportError().

Referenced by MOERTEL::Integrator::Integrate(), and MOERTEL::Integrator::Integrate_2D_Mmod().

◆ FunctionType()

MOERTEL::Function::FunctionType MOERTEL::Segment::FunctionType ( int  id)

Return FunctionType of a function with the Id id.

Parameters
id: function id to lookup the type for

◆ GetLocalNodeId()

int MOERTEL::Segment::GetLocalNodeId ( int  nid)

Get segment-local node id from global node id nid.

References MOERTEL::ReportError().

Referenced by MOERTEL::Projector::ProjectNodetoSegment_Orthogonal_to_Slave().

◆ GetPtrstoNodes() [1/2]

bool MOERTEL::Segment::GetPtrstoNodes ( MOERTEL::Interface interface)

◆ GetPtrstoNodes() [2/2]

bool MOERTEL::Segment::GetPtrstoNodes ( std::vector< MOERTEL::Node * > &  nodes)

Get pointers to Nodes attached to this Segment from a vector of Node pointers.

References MOERTEL::ReportError().

◆ Id()

int MOERTEL::Segment::Id ( ) const
inline

◆ LocalCoordinatesOfNode()

virtual bool MOERTEL::Segment::LocalCoordinatesOfNode ( int  lid,
double *  xi 
)
pure virtual

Get local coords of a node attached to this segment with local node Id lid.

Implemented in MOERTEL::Segment_BiLinearQuad, MOERTEL::Segment_Linear1D, and MOERTEL::Segment_BiLinearTri.

Referenced by MOERTEL::Projector::ProjectNodetoSegment_Orthogonal_to_Slave().

◆ Metric()

virtual double MOERTEL::Segment::Metric ( double *  xi,
double  g[],
double  G[][3] 
)
pure virtual

◆ Nfunctions()

int MOERTEL::Segment::Nfunctions ( )
inline

Return number of functions defined on this Segment.

◆ Nnode()

int MOERTEL::Segment::Nnode ( ) const
inline

◆ NodeIds()

const int* MOERTEL::Segment::NodeIds ( ) const
inline

Return view of node ids of nodes attached to this Segment.

Referenced by MOERTEL::Interface::AddSegment().

◆ Nodes()

MOERTEL::Node** MOERTEL::Segment::Nodes ( )
inline

Return pointer to vector of length Nnode() of pointers to Nodes attached to this Segment.

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

◆ OutLevel()

int MOERTEL::Segment::OutLevel ( )
inline

Return level of output to be generated by this class (0-10)

Referenced by MOERTEL::Interface::AddSegment().

◆ Pack()

virtual int* MOERTEL::Segment::Pack ( int *  size)
pure virtual

Pack some data from this class to an int vector of length size so it can be communicated using MPI.

Implemented in MOERTEL::Segment_BiLinearQuad, MOERTEL::Segment_BiLinearTri, and MOERTEL::Segment_Linear1D.

◆ Print()

bool MOERTEL::Segment::Print ( ) const
virtual

Print this Segment.

◆ SetFunction()

bool MOERTEL::Segment::SetFunction ( int  id,
MOERTEL::Function func 
)

Attach a function to this Segment.

Will attach a function to this Segment under the function Id id. Segment will not take ownership of func but will store a deep copy

Parameters
id: unique function id to store function
func: Function to store in this Segment

References MOERTEL::Function::Clone().

◆ Type()

MOERTEL::Segment::SegmentType MOERTEL::Segment::Type ( ) const
inline

Return type of Segment.

Referenced by MOERTEL::Interface::AddSegment().

◆ UnPack()

virtual bool MOERTEL::Segment::UnPack ( int *  pack)
pure virtual

Unpack some data an int vector and store data in this class.

Implemented in MOERTEL::Segment_BiLinearQuad, MOERTEL::Segment_BiLinearTri, and MOERTEL::Segment_Linear1D.


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