Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_STK_QuadraticToLinearMeshFactory.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Panzer: A partial differential equation assembly
5 // engine for strongly coupled complex multiphysics systems
6 // Copyright (2011) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Roger P. Pawlowski (rppawlo@sandia.gov) and
39 // Eric C. Cyr (eccyr@sandia.gov)
40 // ***********************************************************************
41 // @HEADER
42 
43 #ifndef __Panzer_STK_QuadraticToLinearMeshFactory_hpp__
44 #define __Panzer_STK_QuadraticToLinearMeshFactory_hpp__
45 
47 #include <Panzer_STK_Interface.hpp>
48 #include <vector>
49 #include <string>
50 
51 namespace panzer_stk {
52 
53 class STK_Interface;
54 
60 public:
61 
62  QuadraticToLinearMeshFactory(const std::string& quadMeshFileName,
63  stk::ParallelMachine mpi_comm = MPI_COMM_WORLD, // CHECK: ALLOW MPI_COMM_WORLD
64  const bool print_debug = false);
65 
67  const bool print_debug = false);
68 
69  Teuchos::RCP<STK_Interface> buildMesh(stk::ParallelMachine parallelMach) const;
70  virtual Teuchos::RCP<STK_Interface> buildUncommitedMesh(stk::ParallelMachine parallelMach) const;
71  virtual void completeMeshConstruction(STK_Interface & mesh,stk::ParallelMachine parallelMach) const;
72 
75 
78 
79 protected:
80 
81  void buildMetaData(stk::ParallelMachine parallelMach,STK_Interface & mesh) const;
82  void buildElements(stk::ParallelMachine parallelMach,STK_Interface & mesh) const;
83 
84  void addSideSets(STK_Interface & mesh) const;
85  void addNodeSets(STK_Interface & mesh) const;
86  void addEdgeBlocks(STK_Interface & mesh) const;
87  void copyCellFieldData(STK_Interface & mesh) const;
88 
94  void getOutputTopology();
95 
97 
98  mutable unsigned int machRank_, machSize_;
99 
101 
104  mutable stk::mesh::EntityId offset_;
105 
107 
108  std::string edgeBlockName_;
109 
110  const CellTopologyData * outputTopoData_;
111 
112  unsigned int nDim_;
113  unsigned int nNodes_;
114 
116  std::vector<shards::CellTopology> supportedInputTopos_ = {
117  shards::CellTopology(shards::getCellTopologyData<shards::Quadrilateral<8>>()),
118  shards::CellTopology(shards::getCellTopologyData<shards::Triangle<6>>()),
119  shards::CellTopology(shards::getCellTopologyData<shards::Tetrahedron<10>>()),
120  shards::CellTopology(shards::getCellTopologyData<shards::Hexahedron<20>>())
121  };
122 
125  std::map<const std::string,const CellTopologyData *> outputTopoMap_ = {
126  {shards::getCellTopologyData<shards::Quadrilateral<8>>()->name,
127  shards::getCellTopologyData<shards::Quadrilateral<4>>()},
128  {shards::getCellTopologyData<shards::Triangle<6>>()->name,
129  shards::getCellTopologyData<shards::Triangle<3>>()},
130  {shards::getCellTopologyData<shards::Tetrahedron<10>>()->name,
131  shards::getCellTopologyData<shards::Tetrahedron<4>>()},
132  {shards::getCellTopologyData<shards::Hexahedron<20>>()->name,
133  shards::getCellTopologyData<shards::Hexahedron<8>>()}
134  };
135 };
136 
137 }
138 
139 #endif
QuadraticToLinearMeshFactory(const std::string &quadMeshFileName, stk::ParallelMachine mpi_comm=MPI_COMM_WORLD, const bool print_debug=false)
std::map< const std::string, const CellTopologyData * > outputTopoMap_
void buildElements(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
Teuchos::RCP< STK_Interface > buildMesh(stk::ParallelMachine parallelMach) const
Build the mesh object.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Derived from ParameterListAcceptor.
void buildMetaData(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &paramList)
Derived from ParameterListAcceptor.
bool offsetGIDs_
If true, offset mesh GIDs to exercise 32-bit limits.
virtual void completeMeshConstruction(STK_Interface &mesh, stk::ParallelMachine parallelMach) const
virtual Teuchos::RCP< STK_Interface > buildUncommitedMesh(stk::ParallelMachine parallelMach) const
std::vector< shards::CellTopology > supportedInputTopos_
Nodes in one element of the linear basis.