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.
feMat.h
1 //
2 // Created by milinda on 10/30/18.
3 //
9 #ifndef DENDRO_5_0_FEMAT_H
10 #define DENDRO_5_0_FEMAT_H
11 
12 #include "oda.h"
13 #include "point.h"
14 #ifdef BUILD_WITH_PETSC
15 #include "petscdmda.h"
16 #endif
17 
18 class feMat {
19 
20 protected:
21 
24 
27 
30 
33 
34 
35 #ifdef BUILD_WITH_PETSC
36 
37  DM m_uiPETSC_DA;
38 #endif
39 
40 public:
45  {
46  m_uiOctDA=da;
47  }
48 
51  {
52 
53  }
54 
60  virtual void matVec(const VECType * in, VECType* out,double scale=1.0)=0;
61 
63  inline void setProblemDimensions(const Point& pt_min, const Point& pt_max)
64  {
65  m_uiPtMin=pt_min;
66  m_uiPtMax=pt_max;
67  }
68 
69 
70 #ifdef BUILD_WITH_PETSC
71 // all PETSC function should go here.
72 
78  virtual void matVec(const Vec& in, Vec& out,double scale=1.0)=0;
79 
86  virtual bool getAssembledMatrix(Mat *J, MatType mtype) = 0;
87 
88 #endif
89 
90 
91 };
92 
93 #endif //DENDRO_5_0_FEMAT_H
Abstract class for stiffness matrix computation.
Definition: feMat.h:18
void setProblemDimensions(const Point &pt_min, const Point &pt_max)
set the problem dimension
Definition: feMat.h:63
ot::DAType m_uiDaType
: type of the DA
Definition: feMat.h:26
A point class.
Definition: point.h:36
Point m_uiPtMax
problem domain max point
Definition: feMat.h:32
ot::DA * m_uiOctDA
: pointer to OCT DA
Definition: feMat.h:23
feMat(ot::DA *da)
: feMat constructor
Definition: feMat.h:44
Definition: oda.h:125
DAType
: DA type OCT: uses Dendro DA and PETSC uses petsc based DA
Definition: oda.h:105
virtual void matVec(const VECType *in, VECType *out, double scale=1.0)=0
Computes the LHS of the weak formulation, normally the stifness matrix times a given vector...
~feMat()
deconstructor
Definition: feMat.h:50
The class that manages the octree mesh that support FEM computations. Note that this file is a refact...
Point m_uiPtMin
problem domain min point
Definition: feMat.h:29