DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_IntrepidBasisFactory.cpp
Go to the documentation of this file.
1 //---------------------------------------------------------------------------//
2 /*
3  Copyright (c) 2012, Stuart R. Slattery
4  All rights reserved.
5 
6  Redistribution and use in source and binary forms, with or without
7  modification, are permitted provided that the following conditions are
8  met:
9 
10  *: Redistributions of source code must retain the above copyright
11  notice, this list of conditions and the following disclaimer.
12 
13  *: Redistributions in binary form must reproduce the above copyright
14  notice, this list of conditions and the following disclaimer in the
15  documentation and/or other materials provided with the distribution.
16 
17  *: Neither the name of the University of Wisconsin - Madison nor the
18  names of its contributors may be used to endorse or promote products
19  derived from this software without specific prior written permission.
20 
21  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
22  "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
23  LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
24  A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
25  HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
26  SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
27  LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
28  DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
29  THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
30  (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
31  OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
32 */
33 //---------------------------------------------------------------------------//
39 //---------------------------------------------------------------------------//
40 
42 #include "DTK_DBC.hpp"
43 
44 #include <Teuchos_RCP.hpp>
45 
46 #include <Shards_BasicTopologies.hpp>
47 
48 #include "Intrepid_HGRAD_HEX_C1_FEM.hpp"
49 #include "Intrepid_HGRAD_HEX_C2_FEM.hpp"
50 #include "Intrepid_HGRAD_LINE_C1_FEM.hpp"
51 #include "Intrepid_HGRAD_PYR_C1_FEM.hpp"
52 #include "Intrepid_HGRAD_QUAD_C1_FEM.hpp"
53 #include "Intrepid_HGRAD_QUAD_C2_FEM.hpp"
54 #include "Intrepid_HGRAD_TET_C1_FEM.hpp"
55 #include "Intrepid_HGRAD_TET_C2_FEM.hpp"
56 #include "Intrepid_HGRAD_TRI_C1_FEM.hpp"
57 #include "Intrepid_HGRAD_TRI_C2_FEM.hpp"
58 #include "Intrepid_HGRAD_WEDGE_C1_FEM.hpp"
59 #include "Intrepid_HGRAD_WEDGE_C2_FEM.hpp"
60 #include <Intrepid_Basis.hpp>
61 
62 namespace DataTransferKit
63 {
64 //---------------------------------------------------------------------------//
69 //---------------------------------------------------------------------------//
70 Teuchos::RCP<Intrepid::Basis<double, Intrepid::FieldContainer<double>>>
71 IntrepidBasisFactory::create( const shards::CellTopology &cell_topo )
72 {
73 
74  Teuchos::RCP<Intrepid::Basis<double, Intrepid::FieldContainer<double>>>
75  basis;
76 
77  switch ( cell_topo.getKey() )
78  {
79 
80  case shards::Line<2>::key:
81  basis = Teuchos::rcp( new Intrepid::Basis_HGRAD_LINE_C1_FEM<
82  double, Intrepid::FieldContainer<double>>() );
83  break;
84 
85  case shards::Triangle<3>::key:
86  basis = Teuchos::rcp( new Intrepid::Basis_HGRAD_TRI_C1_FEM<
87  double, Intrepid::FieldContainer<double>>() );
88  break;
89 
90  case shards::Triangle<6>::key:
91  basis = Teuchos::rcp( new Intrepid::Basis_HGRAD_TRI_C2_FEM<
92  double, Intrepid::FieldContainer<double>>() );
93  break;
94 
95  case shards::Quadrilateral<4>::key:
96  basis = Teuchos::rcp( new Intrepid::Basis_HGRAD_QUAD_C1_FEM<
97  double, Intrepid::FieldContainer<double>>() );
98  break;
99 
100  case shards::Quadrilateral<9>::key:
101  basis = Teuchos::rcp( new Intrepid::Basis_HGRAD_QUAD_C2_FEM<
102  double, Intrepid::FieldContainer<double>>() );
103  break;
104 
105  case shards::Tetrahedron<4>::key:
106  basis = Teuchos::rcp( new Intrepid::Basis_HGRAD_TET_C1_FEM<
107  double, Intrepid::FieldContainer<double>>() );
108  break;
109 
110  case shards::Tetrahedron<10>::key:
111  basis = Teuchos::rcp( new Intrepid::Basis_HGRAD_TET_C2_FEM<
112  double, Intrepid::FieldContainer<double>>() );
113  break;
114 
115  case shards::Hexahedron<8>::key:
116  basis = Teuchos::rcp( new Intrepid::Basis_HGRAD_HEX_C1_FEM<
117  double, Intrepid::FieldContainer<double>>() );
118  break;
119 
120  case shards::Hexahedron<27>::key:
121  basis = Teuchos::rcp( new Intrepid::Basis_HGRAD_HEX_C2_FEM<
122  double, Intrepid::FieldContainer<double>>() );
123  break;
124 
125  case shards::Wedge<6>::key:
126  basis = Teuchos::rcp( new Intrepid::Basis_HGRAD_WEDGE_C1_FEM<
127  double, Intrepid::FieldContainer<double>>() );
128  break;
129 
130  case shards::Wedge<18>::key:
131  basis = Teuchos::rcp( new Intrepid::Basis_HGRAD_WEDGE_C2_FEM<
132  double, Intrepid::FieldContainer<double>>() );
133  break;
134 
135  case shards::Pyramid<5>::key:
136  case shards::Pyramid<13>::key:
137  case shards::Pyramid<14>::key:
138  basis = Teuchos::rcp( new Intrepid::Basis_HGRAD_PYR_C1_FEM<
139  double, Intrepid::FieldContainer<double>>() );
140  break;
141 
142  default:
143  bool topology_supported = false;
144  DTK_INSIST( topology_supported );
145  break;
146  }
147 
148  return basis;
149 }
150 
151 //---------------------------------------------------------------------------//
152 
153 } // end namespace DataTransferKit
154 
155 //---------------------------------------------------------------------------//
156 // end DTK_IntrepidBasisFactory.cpp
157 //---------------------------------------------------------------------------//
Assertions and Design-by-Contract for error handling.
Factory for Intrepid basis functions.
DTK_BasicEntitySet.cpp.