#include <fedealwrapper.h>
Public Member Functions | |
FEDealWrapper () | |
void | wrap (dealii::FiniteElement< 3 > &fe) |
virtual unsigned | n_fields () |
~FEDealWrapper () | |
virtual Index | n_dofs_per_face () |
virtual const VecIndex & | face_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 ArrayOfVecDouble & | get_support_points () |
Protected Attributes | |
dealii::FiniteElement< 3 > * | _dealFE |
std::vector< VecIndex > | _face_to_cell_map |
Private Attributes | |
ArrayOfVecDouble | m_support_points |
This class creates a wrapper for the finite elements of the deal library. So when you call the methods of this class you are actually calling the finite element defined on that library.
Design Details: Generally the user create a derived class that instantiates a finite element and then passes this finite element to the base class FEDealWrapper. Note that if this finite element is created as a member of the derived class and we pass it to the FEDealWrapper we create a serious problem, because in the derived class the FEDEalWrapper constructor is called before the constructor of the finite element and then we are passing to FEDEalWrapper an initialized element. That will result in unpredictable errors. A definitely recipe for disaster. Thats why the construction of the FEDealWrapper is empty and defined as protected. To not allow programmers to commit such mistake.
So the FEDealWrapper does not have a constructor that receives the deal finite element. This element is set after the constructor with the method wrap() that assign a finite element to a pointer inside FEDealWrapper. This pointer is set as NULL in the default constructor, so in case the user forget to use the wrap() function, to set the finite element we will have an assert error or a segmentation fault in the case of release mode. At least the error will not be unpredictable. So last do the final message
BE SURE TO SET THIS CLASS WITH THE METHOD FEDealWrapper::wrap() BEFORE USE ANY OTHER METHODS OF THIS CLASS
Definition at line 38 of file fedealwrapper.h.
FEDealWrapper::FEDealWrapper | ( | ) | [inline] |
Definition at line 47 of file fedealwrapper.h.
00048 :_dealFE(NULL) 00049 { 00050 }
FEDealWrapper::~FEDealWrapper | ( | ) | [inline] |
Definition at line 75 of file fedealwrapper.h.
Implements FiniteElementInterface.
Definition at line 83 of file fedealwrapper.h.
00084 { 00085 return _face_to_cell_map[faceId]; 00086 }
virtual const ArrayOfVecDouble& FEDealWrapper::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 114 of file fedealwrapper.h.
00115 { 00116 return m_support_points; 00117 }
virtual Index FEDealWrapper::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 97 of file fedealwrapper.h.
00098 { 00099 return _dealFE->system_to_component_index(local_dof).first; 00100 }
virtual Index FEDealWrapper::n_dofs_per_face | ( | ) | [inline, virtual] |
Implements FiniteElementInterface.
Definition at line 78 of file fedealwrapper.h.
00079 { 00080 return _dealFE->n_dofs_per_face(); 00081 }
virtual unsigned FEDealWrapper::n_fields | ( | ) | [inline, virtual] |
This function returns the number of fields (components) of the finite element.
Implements FiniteElementInterface.
Definition at line 73 of file fedealwrapper.h.
00073 {return _dealFE->n_components();}
virtual Index FEDealWrapper::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 90 of file fedealwrapper.h.
00091 { 00092 return _dealFE->n_dofs_per_cell(); 00093 }
virtual void FEDealWrapper::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 108 of file fedealwrapper.h.
virtual double FEDealWrapper::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 102 of file fedealwrapper.h.
void FEDealWrapper::wrap | ( | dealii::FiniteElement< 3 > & | fe | ) | [inline] |
Definition at line 53 of file fedealwrapper.h.
00054 { 00055 _dealFE = &fe; 00056 //Build the maps between the local dofs of the finite element 00057 //with the local dofs of the faces 00058 _face_to_cell_map.resize(6); 00059 00060 for (unsigned i=0;i<6;i++) //for all faces 00061 { 00062 //initialize the map of that face 00063 _face_to_cell_map[i].resize(fe.n_dofs_per_face()); 00064 00065 //for each dof j on the face 00066 for (unsigned j=0;j<fe.n_dofs_per_face();j++) 00067 _face_to_cell_map[i][j]=fe.face_to_cell_index(j,i); 00068 } 00069 00070 m_support_points=_dealFE->get_unit_support_points(); 00071 00072 }
dealii::FiniteElement<3>* FEDealWrapper::_dealFE [protected] |
Definition at line 43 of file fedealwrapper.h.
std::vector<VecIndex > FEDealWrapper::_face_to_cell_map [protected] |
Definition at line 44 of file fedealwrapper.h.
Definition at line 41 of file fedealwrapper.h.