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

Albany_TmplSTKMeshStruct.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_TMPLSTKMESHSTRUCT_HPP
00008 #define ALBANY_TMPLSTKMESHSTRUCT_HPP
00009 
00010 #include "Albany_GenericSTKMeshStruct.hpp"
00011 
00012 namespace Albany {
00013 
00018 
00019 template<unsigned Dim>
00020 struct albany_stk_mesh_traits { };
00021 
00023 template<unsigned Dim, class traits = albany_stk_mesh_traits<Dim> >
00024 struct EBSpecsStruct {
00025 
00026     EBSpecsStruct(){}
00027 
00029     typedef traits traits_type;
00030 
00032     void Initialize(unsigned int nnelems[], double blLength[]);
00033 
00035     void Initialize(int i, const Teuchos::RCP<Teuchos::ParameterList>& params);
00036 
00038     // Note that elemNo is the logical lower left corner of the element being queried
00039     bool inEB(const std::vector<int>& elemNo){ for(std::size_t i = 0; i < elemNo.size(); i++){
00040           if(elemNo[i] < min[i]) return false;
00041           if(elemNo[i] >= max[i]) return false;
00042         }
00043         return true;
00044     }
00045 
00047     int numElems(int dim){ return (max[dim] - min[dim]);}
00048 
00050 //    void calcElemSizes(std::vector<std::vector<double> > &h){ 
00051     void calcElemSizes(std::vector<double> h[]){ 
00052 //        for(std::size_t i = 0; i < h.size(); i++)
00053         for(unsigned i = 0; i < Dim; i++)
00054           for(unsigned j = min[i]; j < max[i]; j++)
00055             h[i][j] = blLength[i] / double(max[i] - min[i]);
00056     }
00057 
00058     std::string name;      // Name of element block
00059     int min[traits_type::size];       // Minimum logical coordinate of the block 
00060     int max[traits_type::size];       // Maximum logical coordinate of the block 
00061     double blLength[traits_type::size];      
00062 
00063 };
00064 
00066 
00067 template<unsigned Dim, class traits = albany_stk_mesh_traits<Dim> >
00068 
00069   class TmplSTKMeshStruct : public GenericSTKMeshStruct {
00070 
00071     public:
00072 
00074     typedef traits traits_type;
00076     typedef typename traits_type::default_element_type default_element_type;
00077     typedef typename traits_type::optional_element_type optional_element_type;
00078     typedef typename traits_type::default_element_side_type default_element_side_type;
00079 
00081     TmplSTKMeshStruct(const Teuchos::RCP<Teuchos::ParameterList>& params, 
00082                       const Teuchos::RCP<Teuchos::ParameterList>& adaptParams, 
00083                       const Teuchos::RCP<const Epetra_Comm>& comm);
00084 
00085     ~TmplSTKMeshStruct() {};
00086 
00088     void setFieldAndBulkData(
00089                   const Teuchos::RCP<const Epetra_Comm>& comm,
00090                   const Teuchos::RCP<Teuchos::ParameterList>& params,
00091                   const unsigned int neq_,
00092                   const AbstractFieldContainer::FieldContainerRequirements& req,
00093                   const Teuchos::RCP<Albany::StateInfoStruct>& sis,
00094                   const unsigned int worksetSize);
00095 
00097     bool hasRestartSolution() const {return false; }
00098 
00100     double restartDataTime() const {return -1.0; }
00101 
00102 
00103     private:
00104 
00106     void buildMesh(const Teuchos::RCP<const Epetra_Comm>& comm);
00107 
00108  
00110     Teuchos::RCP<const Teuchos::ParameterList>
00111       getValidDiscretizationParameters() const;
00112 
00114     void DeclareParts(std::vector<EBSpecsStruct<Dim, traits>  > ebStructArray, 
00115          std::vector<std::string> ssNames,
00116          std::vector<std::string> nsNames);
00117 
00118     unsigned int nelem[traits_type::size];
00119     double scale[traits_type::size];
00120     unsigned int numEB;
00121     std::vector<double> x[traits_type::size];
00122 //    std::vector<std::vector<double> > x;
00123     Teuchos::RCP<Epetra_Map> elem_map;
00124     std::vector<EBSpecsStruct<Dim, traits> > EBSpecs;
00125 
00126     bool periodic_x, periodic_y, periodic_z;
00127     bool triangles; // Defaults to false, meaning quad elements
00128 
00129   };
00130 
00131 // Explicit template definitions in support of the above
00132 
00133   template <>
00134   struct albany_stk_mesh_traits<0> { 
00135 
00136     enum { size = 1 }; // stk wants one dimension
00137     typedef shards::Particle default_element_type;
00138     typedef shards::Particle optional_element_type;
00139     typedef shards::Particle default_element_side_type; // No sides in 0D
00140 
00141   };
00142 
00143   template <>
00144   struct albany_stk_mesh_traits<1> { 
00145 
00146     enum { size = 1 };
00147     typedef shards::Line<2> default_element_type;
00148     typedef shards::Line<2> optional_element_type;
00149     typedef shards::Particle default_element_side_type; // No sides in 1D
00150 
00151   };
00152 
00153   template <>
00154   struct albany_stk_mesh_traits<2> { 
00155 
00156     enum { size = 2 };
00157     typedef shards::Quadrilateral<4> default_element_type;
00158     typedef shards::Triangle<3> optional_element_type;
00159     typedef shards::Line<2> default_element_side_type;
00160 
00161   };
00162 
00163   template <>
00164   struct albany_stk_mesh_traits<3> { 
00165 
00166     enum { size = 3 };
00167     typedef shards::Hexahedron<8> default_element_type;
00168     typedef shards::Hexahedron<8> optional_element_type;
00169     typedef shards::Quadrilateral<4> default_element_side_type;
00170 
00171   };
00172 
00173 // Now, the explicit template function declarations (templated on dimension)
00174 
00175   template<> int  EBSpecsStruct<0>::numElems(int i);
00176   template<> void EBSpecsStruct<0>::calcElemSizes(std::vector<double> h[]);
00177 
00178   template<> void EBSpecsStruct<0>::Initialize(unsigned int nelems[], double blLen[]);
00179   template<> void EBSpecsStruct<0>::Initialize(int i, const Teuchos::RCP<Teuchos::ParameterList>& params);
00180   template<> void EBSpecsStruct<1>::Initialize(int i, const Teuchos::RCP<Teuchos::ParameterList>& params);
00181   template<> void EBSpecsStruct<2>::Initialize(int i, const Teuchos::RCP<Teuchos::ParameterList>& params);
00182   template<> void EBSpecsStruct<3>::Initialize(int i, const Teuchos::RCP<Teuchos::ParameterList>& params);
00183 
00184   template<> void TmplSTKMeshStruct<0>::buildMesh(const Teuchos::RCP<const Epetra_Comm>& comm);
00185   template<> void TmplSTKMeshStruct<1>::buildMesh(const Teuchos::RCP<const Epetra_Comm>& comm);
00186   template<> void TmplSTKMeshStruct<2>::buildMesh(const Teuchos::RCP<const Epetra_Comm>& comm);
00187   template<> void TmplSTKMeshStruct<3>::buildMesh(const Teuchos::RCP<const Epetra_Comm>& comm);
00188 
00189   template<> void TmplSTKMeshStruct<0, albany_stk_mesh_traits<0> >::setFieldAndBulkData(
00190                   const Teuchos::RCP<const Epetra_Comm>& comm,
00191                   const Teuchos::RCP<Teuchos::ParameterList>& params,
00192                   const unsigned int neq_,
00193                   const AbstractFieldContainer::FieldContainerRequirements& req,
00194                   const Teuchos::RCP<Albany::StateInfoStruct>& sis,
00195                   const unsigned int worksetSize);
00196 
00197   template<> Teuchos::RCP<const Teuchos::ParameterList> TmplSTKMeshStruct<0>::getValidDiscretizationParameters() const;
00198   template<> Teuchos::RCP<const Teuchos::ParameterList> TmplSTKMeshStruct<1>::getValidDiscretizationParameters() const;
00199   template<> Teuchos::RCP<const Teuchos::ParameterList> TmplSTKMeshStruct<2>::getValidDiscretizationParameters() const;
00200   template<> Teuchos::RCP<const Teuchos::ParameterList> TmplSTKMeshStruct<3>::getValidDiscretizationParameters() const;
00201 
00202 } // namespace Albany
00203 
00204 // Define macro for explicit template instantiation
00205 #define TMPLSTK_INSTANTIATE_TEMPLATE_CLASS_0D(name) \
00206   template class name<0>;
00207 #define TMPLSTK_INSTANTIATE_TEMPLATE_CLASS_1D(name) \
00208   template class name<1>;
00209 #define TMPLSTK_INSTANTIATE_TEMPLATE_CLASS_2D(name) \
00210   template class name<2>;
00211 #define TMPLSTK_INSTANTIATE_TEMPLATE_CLASS_3D(name) \
00212   template class name<3>;
00213 
00214 #define TMPLSTK_INSTANTIATE_TEMPLATE_CLASS(name) \
00215   TMPLSTK_INSTANTIATE_TEMPLATE_CLASS_0D(name) \
00216   TMPLSTK_INSTANTIATE_TEMPLATE_CLASS_1D(name) \
00217   TMPLSTK_INSTANTIATE_TEMPLATE_CLASS_2D(name) \
00218   TMPLSTK_INSTANTIATE_TEMPLATE_CLASS_3D(name)
00219 
00220 #endif // ALBANY_TMPLSTKMESHSTRUCT_HPP

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