• Main Page
  • Related Pages
  • Namespaces
  • Classes
  • Files
  • File List
  • File Members

Albany_AbstractMeshStruct.hpp

Go to the documentation of this file.
00001 //*****************************************************************//
00002 //    Albany 2.0:  Copyright 2012 Sandia Corporation               //
00003 //    This Software is released under the BSD license detailed     //
00004 //    in the file "license.txt" in the top-level Albany directory  //
00005 //*****************************************************************//
00006 
00007 #ifndef ALBANY_ABSTRACTMESHSTRUCT_HPP
00008 #define ALBANY_ABSTRACTMESHSTRUCT_HPP
00009 
00010 #include "Teuchos_ParameterList.hpp"
00011 #include "Epetra_Comm.h"
00012 #include "Epetra_Map.h"
00013 #include "Albany_StateInfoStruct.hpp"
00014 #include "Albany_AbstractFieldContainer.hpp"
00015 
00016 #include "Shards_CellTopology.hpp"
00017 #include "Albany_Layouts.hpp"
00018 #include "Albany_ProblemUtils.hpp"
00019 #include "Intrepid_DefaultCubatureFactory.hpp"
00020 #include "Intrepid_FunctionSpaceTools.hpp"
00021 #include "Adapt_NodalDataBlock.hpp"
00022 
00023 
00024 namespace Albany {
00025 
00026 template <typename T>
00027 struct DynamicDataArray {
00028    typedef Teuchos::ArrayRCP<Teuchos::RCP<T> > type;
00029 };
00030 
00031 class CellSpecs {
00032 
00033   public:
00034 
00035     CellSpecs(const CellTopologyData &ctd, const int worksetSize, const int cubdegree, 
00036                   const int numdim, const int vecdim = -1, const int numface = 0, bool compositeTet = false) :
00037         cellTopologyData(ctd),
00038         cellType(shards::CellTopology (&ctd)),
00039         intrepidBasis(Albany::getIntrepidBasis(ctd, compositeTet)),
00040         cellCubature(cubFactory.create(cellType, cubdegree)),
00041         dl(worksetSize, cellType.getNodeCount(),
00042                           intrepidBasis->getCardinality(), cellCubature->getNumPoints(), numdim, vecdim, numface)
00043      { }
00044 
00045      unsigned int getNumVertices(){ return cellType.getNodeCount(); }
00046      unsigned int getNumNodes(){ return intrepidBasis->getCardinality(); }
00047      unsigned int getNumQPs(){ return cellCubature->getNumPoints(); }
00048 
00049    private:
00050 
00051      static Intrepid::DefaultCubatureFactory<RealType> cubFactory;
00052 
00053      const CellTopologyData &cellTopologyData; // Information about the topology of the elements contained in the workset
00054      const shards::CellTopology &cellType; // the topology of the elements contained in the workset
00055      const Teuchos::RCP<Intrepid::Basis<RealType, Intrepid::FieldContainer<RealType> > > intrepidBasis; // The basis
00056      const Teuchos::RCP<Intrepid::Cubature<RealType> > cellCubature; // The cubature of the cells in the workset
00057      // Make sure this appears after the above, as it depends on the above being initialized prior to
00058      // dl being initialized
00059      const Albany::Layouts dl; // the data layout for the elements in the workset
00060 
00061 };
00062 
00063 
00064 struct AbstractMeshStruct {
00065 
00066     virtual ~AbstractMeshStruct() {}
00067 
00068   public:
00069 
00071     enum msType { STK_MS, FMDB_VTK_MS, FMDB_EXODUS_MS };
00072 
00073     virtual void setFieldAndBulkData(
00074       const Teuchos::RCP<const Epetra_Comm>& comm,
00075       const Teuchos::RCP<Teuchos::ParameterList>& params,
00076       const unsigned int neq_,
00077       const AbstractFieldContainer::FieldContainerRequirements& req,
00078       const Teuchos::RCP<Albany::StateInfoStruct>& sis,
00079       const unsigned int worksetSize) = 0;
00080 
00081     virtual Teuchos::ArrayRCP<Teuchos::RCP<Albany::MeshSpecsStruct> >& getMeshSpecs() = 0;
00082 
00083     virtual const Albany::DynamicDataArray<Albany::CellSpecs>::type& getMeshDynamicData() const = 0;
00084 
00085     virtual msType meshSpecsType() = 0;
00086 
00087     Teuchos::RCP<Adapt::NodalDataBlock> nodal_data_block;
00088 
00089 };
00090 }
00091 
00092 #endif // ALBANY_ABSTRACTMESHSTRUCT_HPP

Generated on Wed Mar 26 2014 18:36:35 for Albany: a Trilinos-based PDE code by  doxygen 1.7.1