Main Page   Namespace List   Class Hierarchy   Compound List   File List   Namespace Members   Compound Members   File Members   Related Pages   Examples  

BasisFamilyBase.h

Go to the documentation of this file.
00001 #ifndef BASISFAMILYBASE_H
00002 #define BASISFAMILYBASE_H
00003 
00004 #include "SundanceDefs.h"
00005 #include "Cell.h"
00006 #include "CellType.h"
00007 #include "TSFArray.h"
00008 #include "Point.h"
00009 #include "DenseSerialVector.h"
00010 #include "XMLObject.h"
00011 #include "MultiIndex.h"
00012 #include "TemporaryMapper.h"
00013 
00014 namespace Sundance
00015 {
00016 
00017   using namespace TSF;
00018   using std::string;
00019 
00020   using std::ostream;
00021 
00022   class Expr;
00023 
00024   /** \ingroup LowLevelFE
00025    * Base class for basis functions.
00026    */
00027 
00028   class BasisFamilyBase
00029     {
00030     public:
00031       /** Construct a basis family of a given order. */
00032       BasisFamilyBase(const string& name, int order);
00033       /** TUVD */
00034       virtual ~BasisFamilyBase() {;}
00035 
00036       /** return the order of the basis family */
00037       int order() const {return order_;}
00038 
00039       /** return the number of DOFs associated with this basis on
00040        * the given cell type */
00041       int nNodes(const CellTopologyCode& cellType) const ;
00042 
00043       /** return an identifying name for the basis family */
00044       const string& name() const {return name_;}
00045 
00046       /** return a hash code */
00047       int hashCode() const {return hashCode_;}
00048 
00049       /** return true if the subtype is a
00050        * vector basis (e.g., Raviart-Thomas) */
00051       virtual bool isVectorBasis() const {return false;}
00052 
00053       /** return true if the subtype is a scalar basis */
00054       virtual bool isScalarBasis() const {return false;}
00055 
00056       /** evaluate on a set of points given in reference coordinates */
00057       virtual void refEval(const CellTopologyCode& cellType,
00058                            const TSFArray<Point>& pts,
00059                            const MultiIndex& deriv,
00060                            TSFArray<DenseSerialVector>& result) const ;
00061 
00062       /**
00063        * assign DOF numbers for a set of function IDs on a specified cell
00064        */
00065       virtual void setNodeNumbers(const Cell& cell,
00066                                   const TSFArray<int>& funcID,
00067                                   TemporaryMapper& tmpMap) const = 0 ;
00068 
00069       /** get the DOF numbers for a set of function IDs on a specified
00070        * cell
00071        */
00072       virtual void getNodeNumbers(const Cell& cell,
00073                                   const TSFArray<int>& funcID,
00074                                   TemporaryMapper& tmpMap,
00075                                   TSFArray<TSFArray<int> >& indices) const = 0 ;
00076 
00077       /**
00078        * find the DOF list for one of the facets of a specified cell
00079        */
00080       virtual void getFacetDOFs(const Cell& c, int facetDim, int facetIndex,
00081                                 const TSFArray<int>& maximalDOFs,
00082                                 TSFArray<int>& facetDOFs) const = 0 ;
00083 
00084       /**
00085        * Find the points at the nodes. Some points may have been processed already
00086        * and are marked as inactive with the bool array. Inactive points are not
00087        * returned.
00088        */
00089       virtual void getNodalPoints(const Cell& cell,
00090                                   const TSFArray<bool>& activeNodes,
00091                                   TSFArray<Point>& values) const ;
00092 
00093       /**
00094        * Write to XML
00095        */
00096       virtual XMLObject toXML() const ;
00097 
00098     protected:
00099       /** check that a point's dimension is appropriate for
00100        * a given cell type */
00101       void pointCheck(const Point& pt, int expectedDim,
00102                       const string& cellType) const ;
00103 
00104       /** get number of nodes for a point cell */
00105       virtual int nNodesOnPoint() const ;
00106 
00107       /** get number of nodes for a line cell */
00108       virtual int nNodesOnLine() const ;
00109 
00110       /** get number of nodes for a triangle cell */
00111       virtual int nNodesOnTriangle() const ;
00112 
00113       /** get number of nodes for a quad cell  */
00114       virtual int nNodesOnQuad() const ;
00115 
00116       /** get number of nodes for a tet cell */
00117       virtual int nNodesOnTet() const ;
00118 
00119       /** get number of nodes for a hex cell */
00120       virtual int nNodesOnBrick() const ;
00121 
00122       /** evaluate on a point cell */
00123       virtual void evalOnPoint(const Point& pt,
00124                                const MultiIndex& deriv,
00125                                DenseSerialVector& result) const ;
00126       /** evaluate on a line cell  */
00127       virtual void evalOnLine(const Point& pt,
00128                               const MultiIndex& deriv,
00129                               DenseSerialVector& result) const ;
00130 
00131       /** evaluate on a triangle cell  */
00132       virtual void evalOnTriangle(const Point& pt,
00133                                   const MultiIndex& deriv,
00134                                   DenseSerialVector& result) const ;
00135 
00136       /** evaluate on a quad cell  */
00137       virtual void evalOnQuad(const Point& pt,
00138                               const MultiIndex& deriv,
00139                               DenseSerialVector& result) const ;
00140 
00141       /** evaluate on a tet cell  */
00142       virtual void evalOnTet(const Point& pt,
00143                              const MultiIndex& deriv,
00144                              DenseSerialVector& result) const ;
00145       /** evaluate on a brick cell  */
00146       virtual void evalOnBrick(const Point& pt,
00147                                const MultiIndex& deriv,
00148                                DenseSerialVector& result) const ;
00149 
00150       /** order of the basis function */
00151       int order_;
00152       /** an identifying name for the basis function type */
00153       string name_;
00154       /** hash code */
00155       int hashCode_;
00156     };
00157 }
00158 
00159 #endif

Contact:
Kevin Long (krlong@ca.sandia.gov)


Documentation generated by