00001 #ifndef __FE_MATRIX_H_
00002 #define __FE_MATRIX_H_
00003
00004 #include <string>
00005 #include "feMat.h"
00006 #include "timeInfo.h"
00007
00008 namespace ot {
00009 namespace fem {
00010
00011 template <typename T>
00012 class feMatrix : public feMat {
00013 public:
00014
00015 enum stdElemType {
00016 ST_0,ST_1,ST_2,ST_3,ST_4,ST_5,ST_6,ST_7
00017 };
00018
00019 enum exhaustiveElemType {
00020
00021
00022 ET_N = 0,
00023 ET_Y = 2,
00024 ET_X = 1,
00025 ET_XY = 3,
00026 ET_Z = 8,
00027 ET_ZY = 10,
00028 ET_ZX = 9,
00029 ET_ZXY = 11,
00030 ET_XY_XY = 7,
00031 ET_XY_ZXY = 15,
00032 ET_YZ_ZY = 42,
00033 ET_YZ_ZXY = 43,
00034 ET_YZ_XY_ZXY = 47,
00035 ET_ZX_ZX = 25,
00036 ET_ZX_ZXY = 27,
00037 ET_ZX_XY_ZXY = 31,
00038 ET_ZX_YZ_ZXY = 59,
00039 ET_ZX_YZ_XY_ZXY = 63
00040 };
00041
00042 feMatrix();
00043 feMatrix(daType da);
00044 ~feMatrix();
00045
00046
00047 void setStencil(void* stencil);
00048
00049 void setName(std::string name);
00050
00064 virtual bool MatVec(Vec _in, Vec _out, double scale=1.0);
00065
00066 virtual bool MatGetDiagonal(Vec _diag, double scale=1.0);
00067
00081 inline bool ElementalMatVec(int i, int j, int k, PetscScalar ***in, PetscScalar ***out, double scale) {
00082 return asLeaf().ElementalMatVec(i,j,k,in,out,scale);
00083 }
00084
00085 inline bool ElementalMatGetDiagonal(int i, int j, int k, PetscScalar ***diag, double scale) {
00086 return asLeaf().ElementalMatGetDiagonal(i,j,k,diag,scale);
00087 }
00088
00102 inline bool ElementalMatVec(unsigned int idx, PetscScalar *in, PetscScalar *out, double scale) {
00103 return asLeaf().ElementalMatVec(idx, in, out, scale);
00104 }
00105
00106 inline bool ElementalMatGetDiagonal(unsigned int idx, PetscScalar *diag, double scale) {
00107 return asLeaf().ElementalMatGetDiagonal(idx, diag, scale);
00108 }
00109
00110
00111
00117 T& asLeaf() { return static_cast<T&>(*this);}
00118
00119 bool initStencils() {
00120 return asLeaf().initStencils();
00121 }
00122
00123 bool preMatVec() {
00124 return asLeaf().preMatVec();
00125 }
00126
00127 bool postMatVec() {
00128 return asLeaf().postMatVec();
00129 }
00130
00131 void setDof(unsigned int dof) { m_uiDof = dof; }
00132 unsigned int getDof() { return m_uiDof; }
00133
00134 void setTimeInfo(timeInfo *t) { m_time =t; }
00135 timeInfo* getTimeInfo() { return m_time; }
00136
00137 void initOctLut();
00138
00139
00140 inline PetscErrorCode alignElementAndVertices(ot::DA * da,
00141 stdElemType & sType, unsigned int* indices);
00142 inline PetscErrorCode mapVtxAndFlagsToOrientation(int childNum,
00143 unsigned int* indices, unsigned char & mask);
00144 inline PetscErrorCode reOrderIndices(unsigned char eType,
00145 unsigned int* indices);
00146
00147 protected:
00148 void * m_stencil;
00149
00150 std::string m_strMatrixType;
00151
00152 timeInfo *m_time;
00153
00154 unsigned int m_uiDof;
00155
00156
00157 unsigned char ** m_ucpLut;
00158 };
00159
00160 #include "feMatrix.txx"
00161
00162 }
00163 }
00164 #endif