DealFEWrapper Class Reference

#include <dealfewrapper.h>

Inheritance diagram for DealFEWrapper:
Inheritance graph
[legend]
Collaboration diagram for DealFEWrapper:
Collaboration graph
[legend]

List of all members.

Public Member Functions

 DealFEWrapper ()
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)
void clear ()
OrthoMeshgetMesh ()
const FiniteElement< 3 > * getFE ()
void setCell (OrthoMesh::Cell_It &cell)
void setFace (OrthoMesh::Face_It &face)
FaceDirection3D getFaceAtBoundaryDir ()
OrthoMesh::Cell_ItgetCell ()
 ~DealFEWrapper ()
FEValues< 3 > & getFEValues ()
unsigned n_dofs_per_face ()
unsigned n_dofs_per_cell ()
unsigned n_quadrature_points ()
unsigned n_face_quadrature_points ()
unsigned n_dofs ()
Index local_dof_component (Index cDof)
Index local_face_to_local_cell_map (unsigned face_no, unsigned local_face_dof)
void printDealFEInfo ()
double shape_value (unsigned i, unsigned qPoint)
const Tensor< 1, DIM > & shape_grad (unsigned i, unsigned qPoint)
double shape_value (unsigned i, unsigned component, unsigned qPoint)
Tensor< 1, DIM > shape_grad (unsigned i, unsigned component, unsigned qPoint)
double JxW (unsigned qPoint)
Point3D get_qpoint (unsigned qPoint)
double JxWFace (unsigned qPoint)
const Tensor< 1, DIM > & normal_vector (unsigned face_no, unsigned qPoint)
const Tensor< 1, DIM > & shape_grad_face (unsigned face_no, unsigned i, unsigned qPoint)
double shape_value_face (unsigned face_no, unsigned i, unsigned qPoint)
Point3D get_face_qpoint (unsigned qPoint)
void changedCellEvent ()

Private Attributes

FEValues< 3 > * m_fev
Triangulation< 3 > * m_pTria
DoFHandler< 3 > * m_pDof
DoFHandler< DIM >
::active_cell_iterator * 
m_pDoFCell
OrthoMesh::Cell_Itm_pCell
FaceDirection3D m_faceDir
FaceValuesVec m_facesValues
const FiniteElement< 3 > * m_fe
OrthoMeshm_pMesh

Detailed Description

This class implements a wrapper class of the finite elements in Deal library

Definition at line 14 of file dealfewrapper.h.


Constructor & Destructor Documentation

DealFEWrapper::DealFEWrapper (  ) 

Definition at line 40 of file dealfewrapper.cpp.

00041 {
00042   m_fev=NULL;
00043   m_pTria=NULL;
00044   m_pDof=NULL;
00045   m_fe=NULL;
00046   m_pMesh=NULL;
00047 }

DealFEWrapper::~DealFEWrapper (  ) 

Definition at line 100 of file dealfewrapper.cpp.

00101 {
00102   if (m_pDof)
00103   {
00104     this->clear();
00105   }
00106 }


Member Function Documentation

void DealFEWrapper::changedCellEvent (  )  [inline]

Reimplemented in PoissonEqVar.

Definition at line 77 of file dealfewrapper.h.

00077 {}

void DealFEWrapper::clear (  ) 

Definition at line 83 of file dealfewrapper.cpp.

00084 {
00085   assert(m_pDof);
00086   m_pDof->clear();
00087   delete m_pDoFCell;
00088   delete m_pDof;
00089   delete m_pTria;
00090   delete m_fev;
00091   for (unsigned i=0;i<6;i++)
00092     delete m_facesValues[i];
00093 
00094   m_pDoFCell=NULL;
00095   m_pDof=NULL;
00096   m_pTria=NULL;
00097   m_fev=NULL;
00098 }

Point3D DealFEWrapper::get_face_qpoint ( unsigned  qPoint  ) 

Definition at line 123 of file dealfewrapper.cpp.

00124 {
00125   assert(m_facesValues[m_faceDir]);
00126   Point3D pt = m_facesValues[m_faceDir]->get_quadrature_points()[qPoint];
00127   pt+=(*m_pCell)->vertex(VERTEX_000);
00128   return pt;
00129 }

Point3D DealFEWrapper::get_qpoint ( unsigned  qPoint  ) 

Definition at line 115 of file dealfewrapper.cpp.

00116 {
00117   assert(m_fev);
00118   Point3D pt = m_fev->get_quadrature_points()[qPoint];
00119   pt+=getCell()->vertex(VERTEX_000);
00120   return pt;
00121 }

OrthoMesh::Cell_It& DealFEWrapper::getCell (  )  [inline]

Return the current cell element

Definition at line 44 of file dealfewrapper.h.

FaceDirection3D DealFEWrapper::getFaceAtBoundaryDir (  )  [inline]

Definition at line 42 of file dealfewrapper.h.

00042 {return m_faceDir;}

const FiniteElement<3>* DealFEWrapper::getFE (  )  [inline]

Definition at line 38 of file dealfewrapper.h.

00038 {return m_fe;}

FEValues<3>& DealFEWrapper::getFEValues (  )  [inline]

Definition at line 47 of file dealfewrapper.h.

00047 {assert(m_fev);return *m_fev;}

OrthoMesh& DealFEWrapper::getMesh (  )  [inline]

Definition at line 37 of file dealfewrapper.h.

00037 {assert(m_pMesh);return *m_pMesh;}

void DealFEWrapper::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 
)

Definition at line 50 of file dealfewrapper.cpp.

00051 {
00052   m_pMesh=&mesh;
00053   m_fe=&fe;
00054   m_pTria= new Triangulation<3>();
00055   Point3D P1(0,0,0);
00056   GridGenerator::hyper_rectangle(*m_pTria,P1,mesh.getDX());
00057   m_pDof = new DoFHandler<3>(*m_pTria);
00058   m_pDof->distribute_dofs(fe);
00059   m_pDoFCell = new DoFHandler<DIM>::active_cell_iterator(m_pDof->begin_active());
00060   m_fev = new FEValues<3>(fe,quadrature,update_flags);
00061   m_fev->reinit(*m_pDoFCell);
00062 
00063   
00064   //INIT VECTOR OF FACE VALUES, One for each face
00065   m_facesValues.resize(6,NULL);
00066   m_facesValues[UP_FACE]    =new FEFaceValues<3>(fe,face_quadrature,face_update_flags);
00067   m_facesValues[BOTTOM_FACE]=new FEFaceValues<3>(fe,face_quadrature,face_update_flags);
00068   m_facesValues[LEFT_FACE]  =new FEFaceValues<3>(fe,face_quadrature,face_update_flags);
00069   m_facesValues[RIGHT_FACE] =new FEFaceValues<3>(fe,face_quadrature,face_update_flags);
00070   m_facesValues[BACK_FACE]  =new FEFaceValues<3>(fe,face_quadrature,face_update_flags);
00071   m_facesValues[FRONT_FACE] =new FEFaceValues<3>(fe,face_quadrature,face_update_flags);
00072 
00073   //Set the FEFaceValues
00074   m_facesValues[LEFT_FACE  ]->reinit(*m_pDoFCell,LEFT_FACE  );
00075   m_facesValues[RIGHT_FACE ]->reinit(*m_pDoFCell,RIGHT_FACE );
00076   m_facesValues[BOTTOM_FACE]->reinit(*m_pDoFCell,BOTTOM_FACE);
00077   m_facesValues[UP_FACE    ]->reinit(*m_pDoFCell,UP_FACE    );
00078   m_facesValues[FRONT_FACE ]->reinit(*m_pDoFCell,FRONT_FACE );
00079   m_facesValues[BACK_FACE  ]->reinit(*m_pDoFCell,BACK_FACE) ;
00080 }

double DealFEWrapper::JxW ( unsigned  qPoint  )  [inline]

Definition at line 67 of file dealfewrapper.h.

00067 {assert(m_fev);return m_fev->JxW(qPoint);}

double DealFEWrapper::JxWFace ( unsigned  qPoint  )  [inline]

Definition at line 70 of file dealfewrapper.h.

00070 {return m_facesValues[m_faceDir]->JxW(qPoint);}

Index DealFEWrapper::local_dof_component ( Index  cDof  )  [inline]

Given the local cell dof, this function returns its component

Reimplemented in PoissonEqVar.

Definition at line 57 of file dealfewrapper.h.

Index DealFEWrapper::local_face_to_local_cell_map ( unsigned  face_no,
unsigned  local_face_dof 
) [inline]

Give the local dof of a face, the method returns the local dof indice with respect to the cell

Parameters:
face The local index of the face (LEFT_FACE, RIGHT_FACE ..... SO ON)
local The local dof number with respect to the face.
Returns:
Local DOF with respect to the cell

Reimplemented in PoissonEqVar.

Definition at line 121 of file dealfewrapper.h.

00122 {
00123   assert(m_fe); //Assert m_fe is initialized
00124   return m_fe->face_to_cell_index(local_face_dof,face_no);
00125 }

unsigned DealFEWrapper::n_dofs (  )  [inline]

Number of dofs in the domain

Reimplemented in PoissonEqVar.

Definition at line 94 of file dealfewrapper.h.

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 }

unsigned DealFEWrapper::n_dofs_per_cell (  )  [inline]

Reimplemented in PoissonEqVar.

Definition at line 87 of file dealfewrapper.h.

00088 {
00089   assert(m_fe);
00090   return m_fe->n_dofs_per_cell();
00091 }

unsigned DealFEWrapper::n_dofs_per_face (  )  [inline]

Number of dofs per face

Reimplemented in PoissonEqVar.

Definition at line 81 of file dealfewrapper.h.

00082 {
00083   assert(m_fe);
00084   return m_fe->n_dofs_per_face();
00085 }

unsigned DealFEWrapper::n_face_quadrature_points (  )  [inline]

M Number of quadrature points per cell Number of quadrature points at faces

Reimplemented in PoissonEqVar.

Definition at line 109 of file dealfewrapper.h.

00110 {
00111   assert(m_fe);
00112   assert(m_facesValues.size() != 0);
00113   return m_facesValues[0]->n_quadrature_points;
00114 }

unsigned DealFEWrapper::n_quadrature_points (  )  [inline]

M Number of dofs per cell (8 for each component)

Reimplemented in PoissonEqVar.

Definition at line 102 of file dealfewrapper.h.

00103 {
00104   assert(m_fe);
00105   return m_fev->n_quadrature_points;
00106 }

const Tensor<1,DIM>& DealFEWrapper::normal_vector ( unsigned  face_no,
unsigned  qPoint 
) [inline]

Definition at line 71 of file dealfewrapper.h.

00071 {return m_facesValues[face_no]->normal_vector(qPoint);}

void DealFEWrapper::printDealFEInfo (  ) 

Definition at line 142 of file dealfewrapper.cpp.

00143 {
00144   assert(m_fe);
00145   printf("====Printing Finete Element====\n");
00146   printf("DOFS PER CELL %d\n",m_fe->n_dofs_per_cell());
00147   printf("DOFS PER FACE %d\n",m_fe->n_dofs_per_face());
00148   printf("DOFS PER HEX %d\n",m_fe->n_dofs_per_hex());
00149   printf("COMPONENTS %d\n",m_fe->n_components());
00150 
00151   printf("LOCAL FACE LOCAL CELL MAPPING\n\n");
00152   for (unsigned face=0;face < GeometryInfo<3>::faces_per_cell;face++)
00153   {
00154     printf("FACE %d\n",face);
00155     for (unsigned i=0;i<m_fe->n_dofs_per_face();i++)
00156     {
00157       printf("\tDOF %d = %d\n",i,local_face_to_local_cell_map(face,i));
00158     }
00159   }
00160     printf("DOF TO COMPONENT MAPPING\n\n");
00161   for (unsigned dof=0;dof<n_dofs_per_cell();dof++)
00162   {
00163     printf("\tDof %d is Component %d\n",dof,local_dof_component(dof));
00164   }
00165     
00166 }

void DealFEWrapper::setCell ( OrthoMesh::Cell_It cell  ) 

Reimplemented in IsotropicElasticEqVar, and PoissonEqVar.

Definition at line 107 of file dealfewrapper.cpp.

00108 {
00109   m_pCell=&cell;
00110 }

void DealFEWrapper::setFace ( OrthoMesh::Face_It face  ) 

Definition at line 133 of file dealfewrapper.cpp.

00134 {
00135   m_faceDir = face->getFaceDir();//thids method only makes sense for faces at boundary
00136   
00137 }

Tensor<1,DIM> DealFEWrapper::shape_grad ( unsigned  i,
unsigned  component,
unsigned  qPoint 
) [inline]

Definition at line 66 of file dealfewrapper.h.

00066 {assert(m_fev);return m_fev->shape_grad_component(i,qPoint,component);}

const Tensor<1,DIM>& DealFEWrapper::shape_grad ( unsigned  i,
unsigned  qPoint 
) [inline]

Get the gradient of the shape function

Definition at line 64 of file dealfewrapper.h.

const Tensor<1,DIM>& DealFEWrapper::shape_grad_face ( unsigned  face_no,
unsigned  i,
unsigned  qPoint 
) [inline]

Definition at line 72 of file dealfewrapper.h.

00072 {return m_facesValues[face_no]->shape_grad(i,qPoint);}

double DealFEWrapper::shape_value ( unsigned  i,
unsigned  component,
unsigned  qPoint 
) [inline]

Definition at line 65 of file dealfewrapper.h.

00065 {assert(m_fev);return m_fev->shape_value_component(i,qPoint,component);}

double DealFEWrapper::shape_value ( unsigned  i,
unsigned  qPoint 
) [inline]

Get the shape value at the current quadrature point

Definition at line 63 of file dealfewrapper.h.

double DealFEWrapper::shape_value_face ( unsigned  face_no,
unsigned  i,
unsigned  qPoint 
) [inline]

Definition at line 73 of file dealfewrapper.h.

00073 {return m_facesValues[face_no]->shape_value(i,qPoint);}


Member Data Documentation

Definition at line 25 of file dealfewrapper.h.

Definition at line 26 of file dealfewrapper.h.

const FiniteElement<3>* DealFEWrapper::m_fe [private]

Definition at line 27 of file dealfewrapper.h.

FEValues<3>* DealFEWrapper::m_fev [private]

Definition at line 20 of file dealfewrapper.h.

Definition at line 24 of file dealfewrapper.h.

DoFHandler<3>* DealFEWrapper::m_pDof [private]

Definition at line 22 of file dealfewrapper.h.

DoFHandler<DIM>::active_cell_iterator* DealFEWrapper::m_pDoFCell [private]

Definition at line 23 of file dealfewrapper.h.

Definition at line 28 of file dealfewrapper.h.

Triangulation<3>* DealFEWrapper::m_pTria [private]

Definition at line 21 of file dealfewrapper.h.


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines
Generated on Sun Apr 8 23:13:00 2012 for CO2INJECTION by  doxygen 1.6.3