Panzer  Version of the Day
Panzer_IntrepidBasisFactory.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_INTREPID_BASIS_FACTORY_H
44 #define PANZER_INTREPID_BASIS_FACTORY_H
45 
46 #include <sstream>
47 #include <string>
48 #include <map>
49 #include "Teuchos_RCP.hpp"
50 #include "Intrepid2_Basis.hpp"
51 
52 #include "Shards_CellTopology.hpp"
53 
54 #include "Intrepid2_HGRAD_QUAD_C1_FEM.hpp"
55 #include "Intrepid2_HGRAD_QUAD_C2_FEM.hpp"
56 
57 #include "Intrepid2_HGRAD_HEX_C1_FEM.hpp"
58 #include "Intrepid2_HGRAD_HEX_C2_FEM.hpp"
59 
60 #include "Intrepid2_HGRAD_TET_C1_FEM.hpp"
61 #include "Intrepid2_HGRAD_TET_C2_FEM.hpp"
62 
63 #include "Intrepid2_HGRAD_TRI_C1_FEM.hpp"
64 #include "Intrepid2_HGRAD_TRI_C2_FEM.hpp"
65 
66 #include "Intrepid2_HGRAD_LINE_C1_FEM.hpp"
67 
68 #include "Intrepid2_HCURL_TRI_I1_FEM.hpp"
69 #include "Intrepid2_HCURL_TET_I1_FEM.hpp"
70 #include "Intrepid2_HCURL_QUAD_I1_FEM.hpp"
71 #include "Intrepid2_HCURL_HEX_I1_FEM.hpp"
72 
73 #include "Intrepid2_HDIV_TRI_I1_FEM.hpp"
74 #include "Intrepid2_HDIV_QUAD_I1_FEM.hpp"
75 #include "Intrepid2_HDIV_TET_I1_FEM.hpp"
76 #include "Intrepid2_HDIV_HEX_I1_FEM.hpp"
77 
79 
80 namespace panzer {
81 
82 
96  template <typename ScalarT, typename ArrayT>
98  createIntrepid2Basis(const std::string basis_type, int basis_order,
99  const shards::CellTopology & cell_topology)
100  {
101  // Shards supports extended topologies so the names have a "size"
102  // associated with the number of nodes. We prune the size to
103  // avoid combinatorial explosion of checks.
104  std::string cell_topology_type = cell_topology.getName();
105  std::size_t end_position = 0;
106  end_position = cell_topology_type.find("_");
107  std::string cell_type = cell_topology_type.substr(0,end_position);
108 
110 
111  if ( (basis_type == "Const") && (basis_order == 0) )
113 
114  else if ( (basis_type == "HGrad") && (cell_type == "Hexahedron") && (basis_order == 1) )
115  basis = Teuchos::rcp( new Intrepid2::Basis_HGRAD_HEX_C1_FEM<ScalarT,ArrayT> );
116 
117  else if ( (basis_type == "HGrad") && (cell_type == "Hexahedron") && (basis_order == 2) )
118  basis = Teuchos::rcp( new Intrepid2::Basis_HGRAD_HEX_C2_FEM<ScalarT,ArrayT> );
119 
120  else if ( (basis_type == "HCurl") && (cell_type == "Hexahedron") && (basis_order == 1) )
121  basis = Teuchos::rcp( new Intrepid2::Basis_HCURL_HEX_I1_FEM<ScalarT,ArrayT> );
122 
123  else if ( (basis_type == "HDiv") && (cell_type == "Hexahedron") && (basis_order == 1) )
124  basis = Teuchos::rcp( new Intrepid2::Basis_HDIV_HEX_I1_FEM<ScalarT,ArrayT> );
125 
126  else if ( (basis_type == "HGrad") && (cell_type == "Tetrahedron") && (basis_order == 1) )
127  basis = Teuchos::rcp( new Intrepid2::Basis_HGRAD_TET_C1_FEM<ScalarT,ArrayT> );
128 
129  else if ( (basis_type == "HGrad") && (cell_type == "Tetrahedron") && (basis_order == 2) )
130  basis = Teuchos::rcp( new Intrepid2::Basis_HGRAD_TET_C1_FEM<ScalarT,ArrayT> );
131 
132  else if ( (basis_type == "HCurl") && (cell_type == "Tetrahedron") && (basis_order == 1) )
133  basis = Teuchos::rcp( new Intrepid2::Basis_HCURL_TET_I1_FEM<ScalarT,ArrayT> );
134 
135  else if ( (basis_type == "HDiv") && (cell_type == "Tetrahedron") && (basis_order == 1) )
136  { basis = Teuchos::rcp( new Intrepid2::Basis_HDIV_TET_I1_FEM<ScalarT,ArrayT> ); }
137 
138  else if ( (basis_type == "HGrad") && (cell_type == "Quadrilateral") && (basis_order == 1) )
139  basis = Teuchos::rcp( new Intrepid2::Basis_HGRAD_QUAD_C1_FEM<ScalarT,ArrayT> );
140 
141  else if ( (basis_type == "HGrad") && (cell_type == "Quadrilateral") && (basis_order == 2) )
142  basis = Teuchos::rcp( new Intrepid2::Basis_HGRAD_QUAD_C2_FEM<ScalarT,ArrayT> );
143 
144  else if ( (basis_type == "HCurl") && (cell_type == "Quadrilateral") && (basis_order == 1) )
145  basis = Teuchos::rcp( new Intrepid2::Basis_HCURL_QUAD_I1_FEM<ScalarT,ArrayT> );
146 
147  else if ( (basis_type == "HDiv") && (cell_type == "Quadrilateral") && (basis_order == 1) )
148  basis = Teuchos::rcp( new Intrepid2::Basis_HDIV_QUAD_I1_FEM<ScalarT,ArrayT> );
149 
150  else if ( (basis_type == "HGrad") && (cell_type == "Triangle") && (basis_order == 1) )
151  basis = Teuchos::rcp( new Intrepid2::Basis_HGRAD_TRI_C1_FEM<ScalarT,ArrayT> );
152 
153  else if ( (basis_type == "HGrad") && (cell_type == "Triangle") && (basis_order == 2) )
154  basis = Teuchos::rcp( new Intrepid2::Basis_HGRAD_TRI_C2_FEM<ScalarT,ArrayT> );
155 
156  else if ( (basis_type == "HCurl") && (cell_type == "Triangle") && (basis_order == 1) )
157  basis = Teuchos::rcp( new Intrepid2::Basis_HCURL_TRI_I1_FEM<ScalarT,ArrayT> );
158 
159  else if ( (basis_type == "HDiv") && (cell_type == "Triangle") && (basis_order == 1) )
160  basis = Teuchos::rcp( new Intrepid2::Basis_HDIV_TRI_I1_FEM<ScalarT,ArrayT> );
161 
162  else if ( (basis_type == "HGrad") && (cell_type == "Line") && (basis_order == 1) )
163  basis = Teuchos::rcp( new Intrepid2::Basis_HGRAD_LINE_C1_FEM<ScalarT,ArrayT> );
164 
166  "Failed to create the requestedbasis with basis_type=\"" << basis_type <<
167  "\", basis_order=\"" << basis_order << "\", and cell_type=\"" << cell_type << "\"!\n");
168 
169  TEUCHOS_TEST_FOR_EXCEPTION(cell_topology!=basis->getBaseCellTopology(),
170  std::runtime_error,
171  "Failed to create basis. Intrepid2 basis topology does not match mesh cell topology!");
172 
173  return basis;
174  }
175 
189  template <typename ScalarT, typename ArrayT>
191  createIntrepid2Basis(const std::string basis_type, int basis_order,
192  const Teuchos::RCP<const shards::CellTopology> & cell_topology)
193  {
194  return createIntrepid2Basis<ScalarT,ArrayT>(basis_type,basis_order,*cell_topology);
195  }
196 
197 }
198 
199 
200 #endif
bool is_null(const std::shared_ptr< T > &p)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< Intrepid2::Basis< ScalarT, ArrayT > > createIntrepid2Basis(const std::string basis_type, int basis_order, const shards::CellTopology &cell_topology)
Creates an Intrepid2::Basis object given the basis, order and cell topology.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< const panzer::PureBasis > basis
Interpolates basis DOF values to IP DOF values.