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

Cell.h

Go to the documentation of this file.
00001 #ifndef CELL_H
00002 #define CELL_H
00003 
00004 #include "SundanceDefs.h"
00005 
00006 #include "TSFArray.h"
00007 #include "TSFSmartPtr.h"
00008 #include "CellBase.h"
00009 #include "FacetSetBase.h"
00010 
00011 
00012 namespace Sundance
00013 {
00014 
00015   using namespace TSF;
00016   using std::string;
00017 
00018   using std::ostream;
00019 
00020   class CellBase;
00021 
00022   class Mesh;
00023   class CellMap;
00024   class MeshData;
00025 
00026   /** \ingroup UserLevelGeometry
00027       Class Cell represents a cell in an unstructured mesh.
00028 
00029 
00030 
00031   */
00032 
00033 
00034   class Cell
00035     {
00036     public:
00037       /** empty ctor */
00038       Cell() : ptr_(0) {;}
00039       /** construct, assuming control of a ptr to a CellBase subtype */
00040       Cell(CellBase* ptr) : ptr_(ptr) {;}
00041 
00042       /** identify self as a zero-cell */
00043       bool isZeroCell() const {return ptr_->isZeroCell();}
00044       /** identify self as a maximal cell */
00045       bool isMaximalCell() const {return ptr_->isMaximalCell();}
00046       /** identify self as an affine cell */
00047       bool isAffine() const {return ptr_->isAffine();}
00048 
00049       /** return the rank of the processor that owns this cell */
00050       int ownerProcID() const {return ptr_->ownerProcID();}
00051       /** return the local index that the current processor uses to identify
00052        * this cell */
00053       int localIndex() const {return ptr_->localIndex();}
00054       /** return the global index of this cell */
00055       int globalIndex() const {return ptr_->globalIndex();}
00056       /** set the global index of this cell */
00057       void setGlobalIndex(int globalIndex) {ptr_->setGlobalIndex(globalIndex);}
00058       /** Inform a cell that it has a parent */
00059       void addParent(int parentIndex, int myFacetIndex)
00060         {ptr_->addParent(parentIndex, myFacetIndex);}
00061 
00062       /** return the spatial dimension of this cell */
00063       int dim() const {return ptr_->dim();}
00064       /** return the spatial dimension of the mesh in which this cell lives */
00065       int meshDim() const ;
00066       /** return an ID code for the mesh in which this cell lives */
00067       int meshID() const ;
00068 
00069 
00070       // topological methods
00071       /** return the number of nodes defining this cell */
00072       int numNodes() const {return ptr_->numNodes();}
00073       /** return the number of vertices */
00074       int numVertices() const {return ptr_->numVertices();}
00075       /** return the number of facets of a given dimension */
00076       int numFacets(int dim) const {return ptr_->numFacets(dim);}
00077       /** return the number of cofacets of a given dimension */
00078       int numCofacets(int dim) const {return ptr_->numCofacets(dim);}
00079       /** get a list of the cell's node numbers */
00080       void getNodes(TSFArray<int>& nodes) const {ptr_->getNodes(nodes);}
00081       /** get a list of the cell's facet indices */
00082       void getFacetIndices(TSFArray<TSFArray<int> >& facetIndices) const
00083         {ptr_->getFacetIndices(facetIndices);}
00084       /** determine if this cell is topologically equivalent to the
00085        * given list of nodes */
00086       bool topologicallyEquals(const TSFArray<int>& nodes) const
00087         {return ptr_->topologicallyEquals(nodes);}
00088       /** return the (facetNumber)-th facet cell of dimension (facetDim). */
00089       const Cell& facet(int facetDim, int facetNumber) const
00090         {return ptr_->facet(facetDim, facetNumber);}
00091       /** return the (cofacetNumber)-th cofacet cell of dimension (cofacetDim). */
00092       const Cell& cofacet(int cofacetDim, int cofacetNumber) const
00093         {return ptr_->cofacet(cofacetDim, cofacetNumber);}
00094 
00095       /** */
00096       int numParents() const {return ptr_->numParents();}
00097       /** return the index of the i-th maximal cell that owns this cell */
00098       int parentIndex(int i) const {return ptr_->parentIndex(i);}
00099       /** return the facet index by which the parent refers to this cell */
00100       int myFacetIndex(int i) const {return ptr_->myFacetIndex(i);}
00101       /** return the i-th maximal parent cell */
00102       const Cell& maximalParent(int i) const ;
00103       /** return the list of the cell's node numbers (MaximalCells only) */
00104       const TSFArray<int>& nodes() const {return ptr_->nodes();}
00105       /** return the node number of this cell (ZeroCells only) */
00106       int node() const {return ptr_->node();}
00107       //  const TSFArray<TSFArray<int> >& facets() const {return ptr_->facets();}
00108       const FacetSetBase* facets() const {return ptr_->facets();}
00109       /** return list of indices of cofacets */
00110       const TSFArray<TSFArray<int> >& cofacetIndices() const
00111         {return ptr_->cofacetIndices();}
00112       /** identify self as boundary cell or not */
00113       bool isBoundaryCell() const {return ptr_->isBoundaryCell();}
00114 
00115       /** return cell's label */
00116       const string& label() const {return ptr_->label();}
00117       /** return my label's index in the label list */
00118       int labelIndex() const {return ptr_->labelIndex();}
00119       /** set my label */
00120       void setLabel(const string& newLabel) {ptr_->setLabel(newLabel);}
00121 
00122       /** type id */
00123       const string& typeName() const {return ptr_->typeName();}
00124       /** topology code */
00125       CellTopologyCode topologyCode() const {return ptr_->topologyCode();}
00126 
00127       /** return the reference cell used to compute geometric quantities
00128        * for this cell */
00129       const ReferenceCell& refCell() const {return ptr_->refCell();}
00130 
00131 
00132       /** return the point coordinates of the i-th node */
00133       const Point& point(int i) const {return ptr_->point(i);}
00134 
00135       /** */
00136       Point position() const ;
00137 
00138       /** compute the diameter of this cell */
00139       double diameter() const ;
00140 
00141       /** compute the outward normal of this cell */
00142       Point normal() const ;
00143 
00144       /** get jacobian at an array of points */
00145       void jacobian(const TSFArray<Point>& refpts,
00146                     TSFArray<CellJacobian>& J) const ;
00147 
00148       /** get determinant of the Jacobian at an array of points */
00149       void detJacobian(const TSFArray<Point>& refpts,
00150                        TSFArray<double>& detJ) const ;
00151 
00152 
00153       /** get jacobian (specialized to affine cells) */
00154       void constantJacobian(CellJacobian& J) const ;
00155 
00156       /** get determinant of jacobian (specialized to affine cells) */
00157       void constantDetJacobian(double& detJ) const ;
00158 
00159       /** do pullback map (physical  to reference coordinates)  */
00160       void pullBack(const TSFArray<Point>& phys,
00161                     TSFArray<Point>& ref) const ;
00162       /** do push-forward map (reference to physical coordinates) */
00163       void pushFwd(const TSFArray<Point>& ref,
00164                    TSFArray<Point>& phys) const ;
00165 
00166       /** do push-forward map (reference to physical coordinates) */
00167       void pushFwd(int numPts, const double* ref, double* phys) const ;
00168 
00169       /** do pullback map (physical to reference coordinates) */
00170       void pullBack(int numPts, const double* phys, double* ref) const ;
00171 
00172       /** register a new facet with this cell */
00173       void registerFacet(int dim, int facetNum, int facetCellIndex)
00174         {ptr_->registerFacet(dim, facetNum, facetCellIndex);}
00175       /** register a new cofacet with this cell */
00176       void registerCofacet(int dim, int cofacetCellIndex)
00177         {ptr_->registerCofacet(dim, cofacetCellIndex);}
00178 
00179       /** print to stream */
00180       ostream& print(ostream& os) const {return ptr_->print(os);}
00181       /** write in a standard text form */
00182       void textWrite(ostream& os) const {ptr_->textWrite(os);}
00183       /** write to a string */
00184       string toString() const {return ptr_->toString();}
00185 
00186       /** return number of bytes required by this cell */
00187       int byteCount() const {return ptr_->byteCount();}
00188       TSFArray<int> wastage() const {return ptr_->wastage();}
00189 
00190     private:
00191       TSFSmartPtr<CellBase> ptr_;
00192     };
00193 
00194   inline ostream& operator<<(ostream& os, const Cell& cell)
00195     {
00196       return cell.print(os);
00197     }
00198 
00199 
00200 
00201 }
00202 
00203 namespace TSF
00204 {
00205   using namespace std;
00206 
00207   inline string toString(const Sundance::Cell& cell) {return cell.toString();}
00208 }
00209 #endif
00210 
00211 
00212 
00213 
00214 
00215 
00216 
00217 
00218 
00219 

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


Documentation generated by