Panzer
Version of the Day
Toggle main menu visibility
Loading...
Searching...
No Matches
adapters-stk
src
stk_interface
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
46
#include <
Panzer_STK_MeshFactory.hpp
>
47
#include <
Panzer_STK_Interface.hpp
>
48
#include <vector>
49
#include <string>
50
51
namespace
panzer_stk
{
52
53
class
STK_Interface
;
54
59
class
QuadraticToLinearMeshFactory
:
public
STK_MeshFactory
{
60
public
:
61
62
QuadraticToLinearMeshFactory
(
const
std::string& quadMeshFileName,
63
stk::ParallelMachine mpi_comm = MPI_COMM_WORLD,
64
const
bool
print_debug =
false
);
65
66
QuadraticToLinearMeshFactory
(
const
Teuchos::RCP<panzer_stk::STK_Interface>& quadMesh,
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
74
void
setParameterList
(
const
Teuchos::RCP<Teuchos::ParameterList> & paramList);
75
77
Teuchos::RCP<const Teuchos::ParameterList>
getValidParameters
()
const
;
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
96
Teuchos::RCP<panzer_stk::STK_Interface>
quadMesh_
;
97
98
mutable
unsigned
int
machRank_
,
machSize_
;
99
100
bool
createEdgeBlocks_
;
101
103
bool
offsetGIDs_
;
104
mutable
stk::mesh::EntityId
offset_
;
105
106
bool
print_debug_
;
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
Panzer_STK_Interface.hpp
Panzer_STK_MeshFactory.hpp
panzer_stk::QuadraticToLinearMeshFactory::offsetGIDs_
bool offsetGIDs_
If true, offset mesh GIDs to exercise 32-bit limits.
Definition
Panzer_STK_QuadraticToLinearMeshFactory.hpp:103
panzer_stk::QuadraticToLinearMeshFactory::nNodes_
unsigned int nNodes_
Dimension of the mesh.
Definition
Panzer_STK_QuadraticToLinearMeshFactory.hpp:113
panzer_stk::QuadraticToLinearMeshFactory::offset_
stk::mesh::EntityId offset_
Definition
Panzer_STK_QuadraticToLinearMeshFactory.hpp:104
panzer_stk::QuadraticToLinearMeshFactory::completeMeshConstruction
virtual void completeMeshConstruction(STK_Interface &mesh, stk::ParallelMachine parallelMach) const
Definition
Panzer_STK_QuadraticToLinearMeshFactory.cpp:167
panzer_stk::QuadraticToLinearMeshFactory::buildMetaData
void buildMetaData(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
Definition
Panzer_STK_QuadraticToLinearMeshFactory.cpp:247
panzer_stk::QuadraticToLinearMeshFactory::print_debug_
bool print_debug_
Definition
Panzer_STK_QuadraticToLinearMeshFactory.hpp:106
panzer_stk::QuadraticToLinearMeshFactory::outputTopoData_
const CellTopologyData * outputTopoData_
Definition
Panzer_STK_QuadraticToLinearMeshFactory.hpp:110
panzer_stk::QuadraticToLinearMeshFactory::outputTopoMap_
std::map< const std::string, const CellTopologyData * > outputTopoMap_
Definition
Panzer_STK_QuadraticToLinearMeshFactory.hpp:125
panzer_stk::QuadraticToLinearMeshFactory::nDim_
unsigned int nDim_
Output mesh topology data.
Definition
Panzer_STK_QuadraticToLinearMeshFactory.hpp:112
panzer_stk::QuadraticToLinearMeshFactory::QuadraticToLinearMeshFactory
QuadraticToLinearMeshFactory(const std::string &quadMeshFileName, stk::ParallelMachine mpi_comm=MPI_COMM_WORLD, const bool print_debug=false)
Definition
Panzer_STK_QuadraticToLinearMeshFactory.cpp:57
panzer_stk::QuadraticToLinearMeshFactory::getValidParameters
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Derived from ParameterListAcceptor.
Definition
Panzer_STK_QuadraticToLinearMeshFactory.cpp:222
panzer_stk::QuadraticToLinearMeshFactory::setParameterList
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > ¶mList)
Derived from ParameterListAcceptor.
Definition
Panzer_STK_QuadraticToLinearMeshFactory.cpp:207
panzer_stk::QuadraticToLinearMeshFactory::machSize_
unsigned int machSize_
Definition
Panzer_STK_QuadraticToLinearMeshFactory.hpp:98
panzer_stk::QuadraticToLinearMeshFactory::buildMesh
Teuchos::RCP< STK_Interface > buildMesh(stk::ParallelMachine parallelMach) const
Build the mesh object.
Definition
Panzer_STK_QuadraticToLinearMeshFactory.cpp:94
panzer_stk::QuadraticToLinearMeshFactory::buildElements
void buildElements(stk::ParallelMachine parallelMach, STK_Interface &mesh) const
Definition
Panzer_STK_QuadraticToLinearMeshFactory.cpp:333
panzer_stk::QuadraticToLinearMeshFactory::createEdgeBlocks_
bool createEdgeBlocks_
Definition
Panzer_STK_QuadraticToLinearMeshFactory.hpp:100
panzer_stk::QuadraticToLinearMeshFactory::buildUncommitedMesh
virtual Teuchos::RCP< STK_Interface > buildUncommitedMesh(stk::ParallelMachine parallelMach) const
Definition
Panzer_STK_QuadraticToLinearMeshFactory.cpp:149
panzer_stk::QuadraticToLinearMeshFactory::addEdgeBlocks
void addEdgeBlocks(STK_Interface &mesh) const
Definition
Panzer_STK_QuadraticToLinearMeshFactory.cpp:436
panzer_stk::QuadraticToLinearMeshFactory::getOutputTopology
void getOutputTopology()
Definition
Panzer_STK_QuadraticToLinearMeshFactory.cpp:117
panzer_stk::QuadraticToLinearMeshFactory::addNodeSets
void addNodeSets(STK_Interface &mesh) const
Definition
Panzer_STK_QuadraticToLinearMeshFactory.cpp:433
panzer_stk::QuadraticToLinearMeshFactory::machRank_
unsigned int machRank_
Second order mesh.
Definition
Panzer_STK_QuadraticToLinearMeshFactory.hpp:98
panzer_stk::QuadraticToLinearMeshFactory::quadMesh_
Teuchos::RCP< panzer_stk::STK_Interface > quadMesh_
Definition
Panzer_STK_QuadraticToLinearMeshFactory.hpp:96
panzer_stk::QuadraticToLinearMeshFactory::addSideSets
void addSideSets(STK_Interface &mesh) const
Definition
Panzer_STK_QuadraticToLinearMeshFactory.cpp:404
panzer_stk::QuadraticToLinearMeshFactory::edgeBlockName_
std::string edgeBlockName_
Definition
Panzer_STK_QuadraticToLinearMeshFactory.hpp:108
panzer_stk::QuadraticToLinearMeshFactory::supportedInputTopos_
std::vector< shards::CellTopology > supportedInputTopos_
Nodes in one element of the linear basis.
Definition
Panzer_STK_QuadraticToLinearMeshFactory.hpp:116
panzer_stk::QuadraticToLinearMeshFactory::copyCellFieldData
void copyCellFieldData(STK_Interface &mesh) const
Definition
Panzer_STK_QuadraticToLinearMeshFactory.cpp:454
panzer_stk::STK_Interface
Definition
Panzer_STK_Interface.hpp:104
panzer_stk::STK_MeshFactory::STK_MeshFactory
STK_MeshFactory()
Definition
Panzer_STK_MeshFactory.hpp:62
panzer_stk
Definition
Panzer_STK_GatherFields_decl.hpp:58
Generated by
1.17.0