00001 #ifndef _MY_OrthoMeshEqVar_
00002 #define _MY_OrthoMeshEqVar_
00003 #include "globals.h"
00004 #include <base/quadrature_lib.h>
00005 #include <fe/fe.h>
00006 #include "orthomesh.h"
00007
00008
00009 typedef std::vector<FEFaceValues<3>* > FaceValuesVec;
00010
00014 class DealFEWrapper
00015 {
00016 private:
00017
00018
00019
00020 FEValues<3> *m_fev;
00021 Triangulation<3> *m_pTria;
00022 DoFHandler<3> *m_pDof;
00023 DoFHandler<DIM>::active_cell_iterator *m_pDoFCell;
00024 OrthoMesh::Cell_It *m_pCell;
00025 FaceDirection3D m_faceDir;
00026 FaceValuesVec m_facesValues;
00027 const FiniteElement<3> *m_fe;
00028 OrthoMesh *m_pMesh;
00029 protected:
00030
00031
00032
00033 public:
00034 DealFEWrapper();
00035 void initFEValues(OrthoMesh &mesh,const FiniteElement<3> &fe, const Quadrature<3> &quadrature, const Quadrature<2> &face_quadrature,const UpdateFlags update_flags,const UpdateFlags face_update_flags);
00036 void clear();
00037 OrthoMesh& getMesh(){assert(m_pMesh);return *m_pMesh;}
00038 const FiniteElement<3>* getFE(){return m_fe;}
00039 void setCell(OrthoMesh::Cell_It &cell);
00040 void setFace(OrthoMesh::Face_It &face);
00041
00042 FaceDirection3D getFaceAtBoundaryDir(){return m_faceDir;}
00043
00044 OrthoMesh::Cell_It& getCell(){assert(m_pCell);return *m_pCell;}
00046 ~DealFEWrapper();
00047 FEValues<3>& getFEValues(){assert(m_fev);return *m_fev;}
00048
00049
00050 unsigned n_dofs_per_face();
00051 unsigned n_dofs_per_cell();
00052 unsigned n_quadrature_points();
00053 unsigned n_face_quadrature_points();
00054 unsigned n_dofs();
00056
00057 Index local_dof_component(Index cDof){return m_fe->system_to_component_index(cDof).first;}
00058 Index local_face_to_local_cell_map(unsigned face_no, unsigned local_face_dof);
00059 void printDealFEInfo();
00060
00061
00062
00063 double shape_value(unsigned i,unsigned qPoint) {assert(m_fev);return m_fev->shape_value(i,qPoint);}
00064 const Tensor<1,DIM>& shape_grad(unsigned i,unsigned qPoint) {assert(m_fev);return m_fev->shape_grad(i,qPoint);}
00065 double shape_value(unsigned i,unsigned component,unsigned qPoint) {assert(m_fev);return m_fev->shape_value_component(i,qPoint,component);}
00066 Tensor<1,DIM> shape_grad (unsigned i,unsigned component,unsigned qPoint) {assert(m_fev);return m_fev->shape_grad_component(i,qPoint,component);}
00067 double JxW(unsigned qPoint) {assert(m_fev);return m_fev->JxW(qPoint);}
00068 Point3D get_qpoint(unsigned qPoint);
00069
00070 double JxWFace(unsigned qPoint) {return m_facesValues[m_faceDir]->JxW(qPoint);}
00071 const Tensor<1,DIM>& normal_vector(unsigned face_no,unsigned qPoint) {return m_facesValues[face_no]->normal_vector(qPoint);}
00072 const Tensor<1,DIM>& shape_grad_face (unsigned face_no,unsigned i,unsigned qPoint) {return m_facesValues[face_no]->shape_grad(i,qPoint);}
00073 double shape_value_face(unsigned face_no,unsigned i,unsigned qPoint) {return m_facesValues[face_no]->shape_value(i,qPoint);}
00074 Point3D get_face_qpoint(unsigned qPoint);
00075
00076
00077 void changedCellEvent(){}
00078 };
00079
00080
00081 inline unsigned DealFEWrapper::n_dofs_per_face()
00082 {
00083 assert(m_fe);
00084 return m_fe->n_dofs_per_face();
00085 }
00086
00087 inline unsigned DealFEWrapper::n_dofs_per_cell()
00088 {
00089 assert(m_fe);
00090 return m_fe->n_dofs_per_cell();
00091 }
00092
00093
00094 inline unsigned DealFEWrapper::n_dofs()
00095 {
00096 assert(m_fe);
00097 return m_fe->n_dofs_per_face()*getMesh().numFaces() +
00098 m_fe->n_dofs_per_vertex()*getMesh().numVertices() +
00099 m_fe->n_dofs_per_hex()*getMesh().numCells();
00100 }
00101
00102 inline unsigned DealFEWrapper::n_quadrature_points()
00103 {
00104 assert(m_fe);
00105 return m_fev->n_quadrature_points;
00106 }
00107
00108
00109 inline unsigned DealFEWrapper::n_face_quadrature_points()
00110 {
00111 assert(m_fe);
00112 assert(m_facesValues.size() != 0);
00113 return m_facesValues[0]->n_quadrature_points;
00114 }
00115
00121 inline Index DealFEWrapper::local_face_to_local_cell_map(unsigned face_no, unsigned local_face_dof)
00122 {
00123 assert(m_fe);
00124 return m_fe->face_to_cell_index(local_face_dof,face_no);
00125 }
00126
00127 #endif