Dendro  5.01
Dendro in Greek language means tree. The Dendro library is a large scale (262K cores on ORNL's Titan) distributed memory adaptive octree framework. The main goal of Dendro is to perform large scale multiphysics simulations efficeiently in mordern supercomputers. Dendro consists of efficient parallel data structures and algorithms to perform variational ( finite element) methods and finite difference mthods on 2:1 balanced arbitary adaptive octrees which enables the users to perform simulations raning from black holes (binary black hole mergers) to blood flow in human body, where applications ranging from relativity, astrophysics to biomedical engineering.
Public Member Functions | Protected Attributes | List of all members
feMatrix< T > Class Template Referenceabstract
Inheritance diagram for feMatrix< T >:
Inheritance graph
[legend]
Collaboration diagram for feMatrix< T >:
Collaboration graph
[legend]

Public Member Functions

 feMatrix (ot::DA *da, unsigned int dof=1)
 constructs an FEM stiffness matrix class. More...
 
virtual void matVec (const VECType *in, VECType *out, double scale=1.0)
 Computes the LHS of the weak formulation, normally the stifness matrix times a given vector. More...
 
virtual void elementalMatVec (const VECType *in, VECType *out, double *coords=NULL, double scale=1.0)=0
 Computes the elemental matvec. More...
 
T & asLeaf ()
 static cast to the leaf node of the inheritance
 
bool preMatVec (const VECType *in, VECType *out, double scale=1.0)
 executed just before the matVec loop in matvec function More...
 
bool postMatVec (const VECType *in, VECType *out, double scale=1.0)
 executed just after the matVec loop in matvec function More...
 
bool preMat ()
 executed before the matrix assembly
 
bool postMat ()
 executed after the matrix assembly
 
void getElementalMatrix (unsigned int eleID, std::vector< ot::MatRecord > &records, double *coords)
 Compute the elemental Matrix. More...
 
- Public Member Functions inherited from feMat
 feMat (ot::DA *da)
 : feMat constructor More...
 
 ~feMat ()
 deconstructor
 
void setProblemDimensions (const Point &pt_min, const Point &pt_max)
 set the problem dimension
 

Protected Attributes

unsigned int m_uiDof
 number of dof
 
VECType * m_uiEleVecIn
 element nodal vec in
 
VECType * m_uiEleVecOut
 
double * m_uiEleCoords
 
- Protected Attributes inherited from feMat
ot::DAm_uiOctDA
 : pointer to OCT DA
 
ot::DAType m_uiDaType
 : type of the DA
 
Point m_uiPtMin
 problem domain min point
 
Point m_uiPtMax
 problem domain max point
 

Constructor & Destructor Documentation

◆ feMatrix()

template<typename T >
feMatrix< T >::feMatrix ( ot::DA da,
unsigned int  dof = 1 
)

constructs an FEM stiffness matrix class.

Parameters
[in]daoctree DA

Member Function Documentation

◆ elementalMatVec()

template<typename T>
virtual void feMatrix< T >::elementalMatVec ( const VECType *  in,
VECType *  out,
double *  coords = NULL,
double  scale = 1.0 
)
pure virtual

Computes the elemental matvec.

Parameters
[in]ininput vector u
[out]outoutput vector Ku
[in]defaultparameter scale vector by scale*Ku

Implemented in HeatEq::HeatMat.

◆ getElementalMatrix()

template<typename T>
void feMatrix< T >::getElementalMatrix ( unsigned int  eleID,
std::vector< ot::MatRecord > &  records,
double *  coords 
)
inline

Compute the elemental Matrix.

Parameters
[in]eleIDelement ID
[in]coords: elemental coordinates
[out]recordsrecords corresponding to the elemental matrix.

◆ matVec()

template<typename T >
void feMatrix< T >::matVec ( const VECType *  in,
VECType *  out,
double  scale = 1.0 
)
virtual

Computes the LHS of the weak formulation, normally the stifness matrix times a given vector.

Parameters
[in]ininput vector u
[out]outoutput vector Ku
[in]defaultparameter scale vector by scale*Ku

Implements feMat.

◆ postMatVec()

template<typename T>
bool feMatrix< T >::postMatVec ( const VECType *  in,
VECType *  out,
double  scale = 1.0 
)
inline

executed just after the matVec loop in matvec function

Parameters
[in]in: input Vector
[out]outoutput vector
[in]scalescalaing factror

◆ preMatVec()

template<typename T>
bool feMatrix< T >::preMatVec ( const VECType *  in,
VECType *  out,
double  scale = 1.0 
)
inline

executed just before the matVec loop in matvec function

Parameters
[in]in: input Vector
[out]outoutput vector
[in]scalescalaing factror

Member Data Documentation

◆ m_uiEleCoords

template<typename T>
double* feMatrix< T >::m_uiEleCoords
protected

elemental coordinates


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