FEOrthoMeshDealWrapper Class Reference

#include <feorthomeshdealwrapper.h>

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

List of all members.

Public Member Functions

 FEOrthoMeshDealWrapper ()
 ~FEOrthoMeshDealWrapper ()
void wrap (dealii::FiniteElement< 3 > &fe)
virtual unsigned n_fields ()
virtual Index n_dofs_per_face ()
virtual const VecIndexface_to_cell_local_map (Index faceId)
virtual Index n_local_dofs ()
virtual Index local_dof_field (unsigned local_dof)
virtual double shape_value (unsigned i, const VecDouble &p)
virtual void shape_grad (unsigned i, const VecDouble &pt, VecDouble &vGrad)
virtual const ArrayOfVecDoubleget_support_points ()

Protected Attributes

dealii::FiniteElement< 3 > * _dealFE
std::vector< VecIndex_face_to_cell_map

Private Attributes

ArrayOfVecDouble m_support_points

Detailed Description

FEOrthoMeshDealWrapper It is a wrapper for finite

Definition at line 10 of file feorthomeshdealwrapper.h.


Constructor & Destructor Documentation

FEOrthoMeshDealWrapper::FEOrthoMeshDealWrapper (  )  [inline]

Definition at line 19 of file feorthomeshdealwrapper.h.

00020     :_dealFE(NULL)
00021   {
00022   }

FEOrthoMeshDealWrapper::~FEOrthoMeshDealWrapper (  )  [inline]

Definition at line 23 of file feorthomeshdealwrapper.h.

00023 {}


Member Function Documentation

virtual const VecIndex& FEOrthoMeshDealWrapper::face_to_cell_local_map ( Index  faceId  )  [inline, virtual]

Implements FiniteElementInterface.

Definition at line 54 of file feorthomeshdealwrapper.h.

00055   {
00056     return _face_to_cell_map[faceId];
00057   }

virtual const ArrayOfVecDouble& FEOrthoMeshDealWrapper::get_support_points (  )  [inline, virtual]

Get the support points for each local dof on the reference element. So the support point of the nth local dof is return as a 3d point on the nth position of the ArrayOfVecDouble.

Definition of support points: In many finite elements, each degree of freedom, ak represents the value of the function a0*N0(xi) + a1*N1(xi) .... an*Nn(xi) at a specific point xi==xk. This xk is the support point of the dof ak.

In the vector returned from this function the position i of this vector has the support point of the local dof i in the reference domain.

Implements FiniteElementInterface.

Definition at line 85 of file feorthomeshdealwrapper.h.

00086   {
00087     return m_support_points;
00088   }

virtual Index FEOrthoMeshDealWrapper::local_dof_field ( unsigned  local_dof  )  [inline, virtual]

Taking a local dof in the reference domain returns the number of the field where this dof belongs to.

Implements FiniteElementInterface.

Definition at line 68 of file feorthomeshdealwrapper.h.

00069   {
00070     return _dealFE->system_to_component_index(local_dof).first;
00071   }

virtual Index FEOrthoMeshDealWrapper::n_dofs_per_face (  )  [inline, virtual]

Implements FiniteElementInterface.

Definition at line 49 of file feorthomeshdealwrapper.h.

00050   {
00051     return _dealFE->n_dofs_per_face();
00052   }

virtual unsigned FEOrthoMeshDealWrapper::n_fields (  )  [inline, virtual]

This function returns the number of fields (components) of the finite element.

Returns:

Implements FiniteElementInterface.

Definition at line 45 of file feorthomeshdealwrapper.h.

00045 {return _dealFE->n_components();}

virtual Index FEOrthoMeshDealWrapper::n_local_dofs (  )  [inline, virtual]

Get number of local dofs This is the dimension of the basis functions in the element level. For hat functions in 1D for example, it should return 2

Implements FiniteElementInterface.

Definition at line 61 of file feorthomeshdealwrapper.h.

00062   {
00063     return _dealFE->n_dofs_per_cell();
00064   }

virtual void FEOrthoMeshDealWrapper::shape_grad ( unsigned  i,
const VecDouble pt,
VecDouble vGrad 
) [inline, virtual]

Return the derivative of the shape function i in the reference element at point xi

Pre: xi is in [0,1]

Implements FiniteElementInterface.

Definition at line 79 of file feorthomeshdealwrapper.h.

00080   {
00081     Point3D p(pt(0),pt(1),pt(2));
00082     _dealFE->shape_grad(i,p).unroll(vGrad);
00083   }

virtual double FEOrthoMeshDealWrapper::shape_value ( unsigned  i,
const VecDouble p 
) [inline, virtual]

Get the value of the ith shape function (also called basis function) at position xi in the reference domain. Parameters: i = The shape function number xi= The point in the reference domain Pre: xi is in [0,1]

Implements FiniteElementInterface.

Definition at line 73 of file feorthomeshdealwrapper.h.

00074   {
00075     Point3D pt(p(0),p(1),p(2));
00076     return _dealFE->shape_value(i,pt);
00077   }

void FEOrthoMeshDealWrapper::wrap ( dealii::FiniteElement< 3 > &  fe  )  [inline]

Definition at line 26 of file feorthomeshdealwrapper.h.

00027   {
00028     _dealFE = &fe;
00029     //Build the maps between the local dofs of the finite element
00030     //with the local dofs of the faces
00031     _face_to_cell_map.resize(6);
00032 
00033     for (unsigned i=0;i<6;i++) //for all faces
00034     {
00035       //initialize the map of that face
00036       _face_to_cell_map[i].resize(fe.n_dofs_per_face());
00037 
00038       //for each dof j on the face
00039       for (unsigned j=0;j<fe.n_dofs_per_face();j++)
00040         _face_to_cell_map[i][j]=fe.face_to_cell_index(j,i);
00041     }
00042 
00043     m_support_points=_dealFE->get_unit_support_points();
00044   }


Member Data Documentation

dealii::FiniteElement<3>* FEOrthoMeshDealWrapper::_dealFE [protected]

Definition at line 15 of file feorthomeshdealwrapper.h.

Definition at line 16 of file feorthomeshdealwrapper.h.

Definition at line 13 of file feorthomeshdealwrapper.h.


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