Home | Doxygen Documentation | Tutorials | Developer Tools (restricted)

hp3D/hexahedron.hh
Go to the documentation of this file.
00001 /* hexahedral element for 3D hp FEM
00002  */
00003 
00004 #ifndef hexahedron_hh
00005 #define hexahedron_hh
00006 
00007 #include <iostream>
00008 #include <memory>
00009 #include "basics/typedefs.hh"
00010 #include "basics/vectorsMatrices.hh"
00011 #include "geometry/cell3D.hh"
00012 #include "integration/quadrature.hh"
00013 #include "integration/karniadakis.hh"
00014 #include "space/tmatrix.hh"
00015 #include "space/hpMethod.hh"
00016 #include "element.hh"
00017 #include "hexahedronGraphics.hh"
00018 
00019 // 0: Gauss Lobatto (good for graphics)
00020 // 4: Gauss Jacobi (highest order)
00021 #ifdef INTEGRATION_RULE
00022 #undef INTEGRATION_RULE
00023 #endif
00024 #define INTEGRATION_RULE 4
00025 
00026 namespace hp3D {
00027 
00028   using concepts::Real;
00029 
00030   // ************************************************************ Hexahedron **
00031 
00036   class Hexahedron : public Element<Real> {
00037   public:
00045     Hexahedron(concepts::Hexahedron3d& cell, const ushort* p, 
00046          concepts::TColumn<Real>* T0, concepts::TColumn<Real>* T1);
00047     virtual ~Hexahedron();
00048 
00049     virtual const concepts::Hexahedron& support() const {
00050       return cell_.connector(); }
00051     virtual concepts::Real3d vertex(uint i) const { return cell_.vertex(i); }
00052     virtual const concepts::Hexahedron3d& cell() const { return cell_; }
00053 
00055     inline const concepts::Karniadakis<1,0>* shpfctX() const {
00056       return shpfctX_.get(); }
00058     inline const concepts::Karniadakis<1,0>* shpfctY() const {
00059       return shpfctY_.get(); }
00061     inline const concepts::Karniadakis<1,0>* shpfctZ() const {
00062       return shpfctZ_.get(); }
00063 
00065     inline const concepts::Karniadakis<1,1>* shpfctDX() const {
00066       return shpfctDX_.get(); }
00068     inline const concepts::Karniadakis<1,1>* shpfctDY() const {
00069       return shpfctDY_.get(); }
00071     inline const concepts::Karniadakis<1,1>* shpfctDZ() const {
00072       return shpfctDZ_.get(); }
00073 
00075     inline const concepts::Quadrature<INTEGRATION_RULE>* integrationX() const {
00076       return intX_.get(); }
00078     inline const concepts::Quadrature<INTEGRATION_RULE>* integrationY() const {
00079       return intY_.get(); }
00081     inline const concepts::Quadrature<INTEGRATION_RULE>* integrationZ() const {
00082       return intZ_.get(); }
00083 
00086     inline concepts::Real3d chi(const Real x, const Real y,
00087         const Real z) const {
00088       return cell_.chi(x, y, z); }
00089 
00091     inline concepts::MapReal3d jacobian(const Real x, const Real y,
00092           const Real z) const {
00093       concepts::MapReal3d jac = cell_.jacobian(x, y, z);
00094       jac *= 0.5;
00095       return jac; }
00096 
00098     inline concepts::MapReal3d jacobianInverse(const Real x, const Real y, 
00099                  const Real z) const {
00100       return jacobian(x, y, z).inverse(); }
00101 
00103     inline Real jacobianDeterminant(const Real x, const Real y,
00104             const Real z) const {
00105       return jacobian(x, y, z).determinant(); }
00106 
00114     void setStrategy(const concepts::Hex3dSubdivision* strategy = 0)
00115       throw(concepts::StrategyChange) { cell_.setStrategy(strategy); }
00116 
00120     const concepts::Hex3dSubdivision* getStrategy() const {
00121       return cell_.getStrategy(); }
00122 
00123     virtual const concepts::ElementGraphics<Real>* graphics() const;
00124 
00125     virtual bool operator<(const Element<Real>& elm) const;
00126 
00132     bool hasSameMatrix(const Hexahedron& elm) const;
00133 
00135     void edgeP(const uint i, const concepts::AdaptiveControlP<1>& p) {
00136       conceptsAssert(i < 12, concepts::Assertion());
00137       edges_[i] = &p; }
00139     const concepts::AdaptiveControlP<1>& edgeP(const uint i) const {
00140       conceptsAssert(i < 12, concepts::Assertion());
00141       conceptsAssert(edges_[i], concepts::Assertion());
00142       return *(edges_[i]); }
00144     void faceP(const uint i, const concepts::AdaptiveControlP<2>& p) {
00145       conceptsAssert(i < 6, concepts::Assertion());
00146       faces_[i] = &p; }
00148     const concepts::AdaptiveControlP<2>& faceP(const uint i) const {
00149       conceptsAssert(i < 6, concepts::Assertion());
00150       conceptsAssert(faces_[i], concepts::Assertion());
00151       return *(faces_[i]); }
00152   protected:
00153     virtual std::ostream& info(std::ostream& os) const;
00154   private:
00156     concepts::Hexahedron3d& cell_;
00158     std::auto_ptr<concepts::Karniadakis<1,0> > shpfctX_, shpfctY_, shpfctZ_;
00160     std::auto_ptr<concepts::Karniadakis<1,1> > shpfctDX_, shpfctDY_, shpfctDZ_;
00162     std::auto_ptr<concepts::Quadrature<INTEGRATION_RULE> > intX_, intY_, intZ_;
00163 
00165     const concepts::AdaptiveControlP<1>* edges_[12];
00167     const concepts::AdaptiveControlP<2>* faces_[6];
00168 
00169     static std::auto_ptr<HexahedronGraphics> graphics_;
00170   };
00171 
00172 } // namespace hp3D
00173 
00174 #endif // hexahedron_hh

Home | Doxygen Documentation | Tutorials | Developer Tools (restricted)