Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_PureBasis.cpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Panzer: A partial differential equation assembly
4// engine for strongly coupled complex multiphysics systems
5//
6// Copyright 2011 NTESS and the Panzer contributors.
7// SPDX-License-Identifier: BSD-3-Clause
8// *****************************************************************************
9// @HEADER
10
11#include "Panzer_PureBasis.hpp"
12#include "Panzer_Dimension.hpp"
13#include "Panzer_CellData.hpp"
16#include "Teuchos_Assert.hpp"
17#include "Phalanx_DataLayout_MDALayout.hpp"
18#include <sstream>
19
21PureBasis(const std::string & basis_type,
22 const int basis_order,
23 const int num_cells,
24 const Teuchos::RCP<const shards::CellTopology> & cell_topology) :
25 topology_(cell_topology),
26 num_cells_(num_cells)
27{
28 initialize(basis_type,basis_order);
29}
30
32PureBasis(const std::string & in_basis_type,
33 const int in_basis_order,
34 const CellData & in_cell_data) :
35 topology_(in_cell_data.getCellTopology()),
36 num_cells_(in_cell_data.numCells())
37{
38 initialize(in_basis_type,in_basis_order);
39}
40
42PureBasis(const panzer::BasisDescriptor & description,
43 const Teuchos::RCP<const shards::CellTopology> & cell_topology,
44 const int num_cells):
45 topology_(cell_topology),
46 num_cells_(num_cells)
47{
48 initialize(description.getType(), description.getOrder());
49}
50
51void panzer::PureBasis::initialize(const std::string & in_basis_type,const int in_basis_order)
52{
53 // Support for deprecated basis descriptions
54 std::string basis_type = in_basis_type;
55 int basis_order = in_basis_order;
56
57 if (basis_type=="Q1" || basis_type=="T1") {
58 basis_type = "HGrad";
59 basis_order = 1;
60 }
61 else if (basis_type == "Q2" || basis_type=="T2") {
62 basis_type = "HGrad";
63 basis_order = 2;
64 }
65 else if (basis_type == "TEdge1" || basis_type=="QEdge1") {
66 basis_type = "HCurl";
67 basis_order = 1;
68 }
69 else if(basis_type == "Const") {
70 basis_type = "Const";
71 basis_order = 0;
72 }
73 // End deprecated basis support
74
75 intrepid_basis_ = panzer::createIntrepid2Basis<PHX::Device::execution_space,double,double>(basis_type, basis_order, *topology_);
76
77 basis_type_ = basis_type;
78
79 std::ostringstream os;
80 os << basis_type_ << ":" << basis_order;
81 basis_name_ = os.str();
82
83 field_basis_name_ = "Basis: " + basis_name_;
84 field_basis_name_D1_ = "Grad Basis: " + basis_name_;
85 field_basis_name_D2_ = "D2 Basis: " + basis_name_;
86
87 if( basis_type_ == "HGrad")
88 element_space_ = HGRAD;
89 else if(basis_type_=="HCurl")
90 element_space_ = HCURL;
91 else if(basis_type_=="HDiv")
92 element_space_ = HDIV;
93 else if(basis_type_=="Const")
94 element_space_ = CONST;
95 else if(basis_type_=="HVol")
96 element_space_ = HVOL;
97 else { TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
98 "PureBasis::initializeIntrospection - Invalid basis name \""
99 << basis_type_ << "\""); }
100
101 switch(getElementSpace()) {
102 case CONST:
103 basis_rank_ = 0;
104 break;
105 case HVOL:
106 basis_rank_ = 0;
107 break;
108 case HGRAD:
109 basis_rank_ = 0;
110 break;
111 case HCURL:
112 basis_rank_ = 1;
113 break;
114 case HDIV:
115 basis_rank_ = 1;
116 break;
117 default:
118 TEUCHOS_ASSERT(false);
119 break;
120 };
121
122 using Teuchos::rcp;
123 using PHX::MDALayout;
124
125 cell_data = rcp(new MDALayout<Cell>(numCells()));
126
127 functional = rcp(new MDALayout<Cell,BASIS>(numCells(), cardinality()));
128
129 functional_grad = rcp(new MDALayout<Cell,BASIS,Dim>(numCells(),
130 cardinality(),
131 dimension()));
132
133 coordinates = rcp(new MDALayout<Cell,BASIS,Dim>(numCells(),
134 cardinality(),
135 dimension()));
136
137 functional_D2 = rcp(new MDALayout<Cell,BASIS,Dim,Dim>(numCells(),
138 cardinality(),
139 dimension(),
140 dimension()));
141
142 local_mat_layout = Teuchos::rcp(new PHX::MDALayout<panzer::Cell, panzer::BASIS, panzer::BASIS>(
143 this->numCells(), this->cardinality(), this->cardinality()));
144
145}
146
148{
149 return intrepid_basis_->getCardinality();
150}
151
153{
154 return num_cells_;
155}
156
158{
159 return topology_->getDimension();
160}
161
162std::string panzer::PureBasis::type() const
163{
164 return basis_type_;
165}
166
168{
169 return intrepid_basis_->getDegree();
170}
171
172std::string panzer::PureBasis::name() const
173{
174 return basis_name_;
175}
176
178{
179 return field_basis_name_;
180}
181
183{
184 return field_basis_name_D1_;
185}
186
188{
189 return field_basis_name_D2_;
190}
191
192Teuchos::RCP< Intrepid2::Basis<PHX::Device::execution_space,double,double> >
194{
195 return intrepid_basis_;
196}
197
198bool
200{
201 // typedef Kokkos::DynRankView<double,PHX::Device> Array;
202 // Teuchos::RCP<const Intrepid2::DofCoordsInterface<Array> > coord_interface
203 // = Teuchos::rcp_dynamic_cast<const Intrepid2::DofCoordsInterface<Array> >(getIntrepid2Basis());
204
205 // return !Teuchos::is_null(coord_interface);
206
207 return true;
208}
int getOrder() const
Get order of basis.
const std::string & getType() const
Get type of basis.
Data for determining cell topology and dimensionality.
std::string fieldName() const
void initialize(const std::string &basis_type, const int basis_order)
Initialize the basis object.
int numCells() const
Returns the number of cells in the data layouts.
bool supportsBasisCoordinates() const
Teuchos::RCP< Intrepid2::Basis< PHX::Device::execution_space, double, double > > getIntrepid2Basis() const
std::string type() const
Returns the basis type.
std::string fieldNameD1() const
PureBasis(const std::string &basis_type, const int basis_order, const CellData &cell_data)
int cardinality() const
Returns the number of basis coefficients.
std::string name() const
A unique key that is the combination of the basis type and basis order.
int order() const
Returns the polynomial order of the basis.
int dimension() const
Returns the dimension of the basis from the topology.
std::string fieldNameD2() const