SCIRun  5.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
HexTriquadraticLgn.h
Go to the documentation of this file.
1 //
2 // For more information, please see: http://software.sci.utah.edu
3 //
4 // The MIT License
5 //
6 // Copyright (c) 2009 Scientific Computing and Imaging Institute,
7 // University of Utah.
8 //
9 //
10 // Permission is hereby granted, free of charge, to any person obtaining a
11 // copy of this software and associated documentation files (the "Software"),
12 // to deal in the Software without restriction, including without limitation
13 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
14 // and/or sell copies of the Software, and to permit persons to whom the
15 // Software is furnished to do so, subject to the following conditions:
16 //
17 // The above copyright notice and this permission notice shall be included
18 // in all copies or substantial portions of the Software.
19 //
20 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
21 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
22 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
23 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
24 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
25 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
26 // DEALINGS IN THE SOFTWARE.
27 //
28 /// @file HexTriquadraticLgn.h
29 /// @author Martin Cole, Frank B. Sachse
30 /// @date Dec 3 2004
31 
32 #ifndef CORE_BASIS_HEXTRIQUADRATICLGN_H
33 #define CORE_BASIS_HEXTRIQUADRATICLGN_H 1
34 
36 
37 #include <Core/Basis/share.h>
38 
39 namespace SCIRun {
40 namespace Core {
41 namespace Basis {
42 
43 /// Class for describing unit geometry of HexTriquadraticLgn
45 public:
46  static double unit_vertices[20][3]; ///< Parametric coordinates of vertices of unit edge
47 
50 
51  static int number_of_vertices()
52  { return 20; } ///< return number of vertices
53  static int dofs()
54  { return 20; } ///< return degrees of freedom
55 };
56 
57 
58 /// Class for handling of element of type hexahedron with
59 /// triquadratic lagrangian interpolation
60 template <class T>
61 class HexTriquadraticLgn : public BasisAddNodes<T>,
62  public HexApprox,
63  public HexGaussian3<double>,
64  public HexSamplingSchemes,
66  public HexElementWeights
67 {
68 public:
69  typedef T value_type;
70 
72  virtual ~HexTriquadraticLgn() {}
73 
74  static int polynomial_order() { return 2; }
75 
76 
77  template<class VECTOR>
78  inline void get_weights(const VECTOR& coords, double *w) const
79  { get_quadratic_weights(coords,w); }
80 
81  template<class VECTOR>
82  inline void get_derivate_weights(const VECTOR& coords, double *w) const
83  { get_quadratic_derivate_weights(coords,w); }
84 
85  /// get value at parametric coordinate
86  template <class ElemData, class VECTOR>
87  T interpolate(const VECTOR &coords, const ElemData &cd) const
88  {
89  double w[20];
90  get_quadratic_weights(coords, w);
91 
92  return (T)(
93  w[0] * cd.node0() +
94  w[1] * cd.node1() +
95  w[2] * cd.node2() +
96  w[3] * cd.node3() +
97  w[4] * cd.node4() +
98  w[5] * cd.node5() +
99  w[6] * cd.node6() +
100  w[7] * cd.node7() +
101  w[8] * this->nodes_[cd.edge0_index()] +
102  w[9] * this->nodes_[cd.edge1_index()] +
103  w[10] * this->nodes_[cd.edge2_index()] +
104  w[11] * this->nodes_[cd.edge3_index()] +
105  w[12] * this->nodes_[cd.edge4_index()] +
106  w[13] * this->nodes_[cd.edge5_index()] +
107  w[14] * this->nodes_[cd.edge6_index()] +
108  w[15] * this->nodes_[cd.edge7_index()] +
109  w[16] * this->nodes_[cd.edge8_index()] +
110  w[17] * this->nodes_[cd.edge9_index()] +
111  w[18] * this->nodes_[cd.edge10_index()] +
112  w[19] * this->nodes_[cd.edge11_index()]);
113  }
114 
115  /// get first derivative at parametric coordinate
116  template <class ElemData, class VECTOR1, class VECTOR2>
117  void derivate(const VECTOR1 &coords, const ElemData &cd,
118  VECTOR2 &derivs) const
119  {
120  double w[60];
121  get_quadratic_derivate_weights(coords, w);
122 
123  derivs.resize(3);
124 
125  derivs[0]=static_cast<typename VECTOR2::value_type>(
126  w[0] * cd.node0() +
127  w[1] * cd.node1() +
128  w[2] * cd.node2() +
129  w[3] * cd.node3() +
130  w[4] * cd.node4() +
131  w[5] * cd.node5() +
132  w[6] * cd.node6() +
133  w[7] * cd.node7() +
134  w[8] * this->nodes_[cd.edge0_index()] +
135  w[9] * this->nodes_[cd.edge1_index()] +
136  w[10] * this->nodes_[cd.edge2_index()] +
137  w[11] * this->nodes_[cd.edge3_index()] +
138  w[12] * this->nodes_[cd.edge4_index()] +
139  w[13] * this->nodes_[cd.edge5_index()] +
140  w[14] * this->nodes_[cd.edge6_index()] +
141  w[15] * this->nodes_[cd.edge7_index()] +
142  w[16] * this->nodes_[cd.edge8_index()] +
143  w[17] * this->nodes_[cd.edge9_index()] +
144  w[18] * this->nodes_[cd.edge10_index()] +
145  w[19] * this->nodes_[cd.edge11_index()]);
146 
147  derivs[1]=static_cast<typename VECTOR2::value_type>(
148  w[20] * cd.node0() +
149  w[21] * cd.node1() +
150  w[22] * cd.node2() +
151  w[23] * cd.node3() +
152  w[24] * cd.node4() +
153  w[25] * cd.node5() +
154  w[26] * cd.node6() +
155  w[27] * cd.node7() +
156  w[28] * this->nodes_[cd.edge0_index()] +
157  w[29] * this->nodes_[cd.edge1_index()] +
158  w[30] * this->nodes_[cd.edge2_index()] +
159  w[31] * this->nodes_[cd.edge3_index()] +
160  w[32] * this->nodes_[cd.edge4_index()] +
161  w[33] * this->nodes_[cd.edge5_index()] +
162  w[34] * this->nodes_[cd.edge6_index()] +
163  w[35] * this->nodes_[cd.edge7_index()] +
164  w[36] * this->nodes_[cd.edge8_index()] +
165  w[37] * this->nodes_[cd.edge9_index()] +
166  w[38] * this->nodes_[cd.edge10_index()] +
167  w[39] * this->nodes_[cd.edge11_index()]);
168 
169  derivs[2]=static_cast<typename VECTOR2::value_type>(
170  w[40] * cd.node0() +
171  w[41] * cd.node1() +
172  w[42] * cd.node2() +
173  w[43] * cd.node3() +
174  w[44] * cd.node4() +
175  w[45] * cd.node5() +
176  w[46] * cd.node6() +
177  w[47] * cd.node7() +
178  w[48] * this->nodes_[cd.edge0_index()] +
179  w[49] * this->nodes_[cd.edge1_index()] +
180  w[50] * this->nodes_[cd.edge2_index()] +
181  w[51] * this->nodes_[cd.edge3_index()] +
182  w[52] * this->nodes_[cd.edge4_index()] +
183  w[53] * this->nodes_[cd.edge5_index()] +
184  w[54] * this->nodes_[cd.edge6_index()] +
185  w[55] * this->nodes_[cd.edge7_index()] +
186  w[56] * this->nodes_[cd.edge8_index()] +
187  w[57] * this->nodes_[cd.edge9_index()] +
188  w[58] * this->nodes_[cd.edge10_index()] +
189  w[59] * this->nodes_[cd.edge11_index()]);
190 
191  }
192 
193  /// get parametric coordinate for value within the element
194  template <class ElemData, class VECTOR>
195  bool get_coords(VECTOR &coords, const T& value,
196  const ElemData &cd) const
197  {
199  return CL.get_coords(this, coords, value, cd);
200  }
201 
202  /// get arc length for edge
203  template <class ElemData>
204  double get_arc_length(const unsigned edge, const ElemData &cd) const
205  {
206  return get_arc3d_length<CrvGaussian2<double> >(this, edge, cd);
207  }
208 
209  /// get area
210  template <class ElemData>
211  double get_area(const unsigned face, const ElemData &cd) const
212  {
213  return get_area3<QuadGaussian3<double> >(this, face, cd);
214  }
215 
216  /// get volume
217  template <class ElemData>
218  double get_volume(const ElemData & cd) const
219  {
220  return get_volume3(this, cd);
221  }
222 
223  static const std::string type_name(int n = -1);
224 
225  virtual void io (Piostream& str);
226 
227 };
228 
229 template <class T>
230 const std::string
232 {
233  ASSERT((n >= -1) && n <= 1);
234  if (n == -1)
235  {
236  static const std::string name = TypeNameGenerator::make_template_id(type_name(0), type_name(1));
237  return name;
238  }
239  else if (n == 0)
240  {
241  static const std::string nm("HexTriquadraticLgn");
242  return nm;
243  } else {
244  return find_type_name((T *)0);
245  }
246 }
247 
248 
249 
250 
251 }}
252 template <class T>
254 {
255  static TypeDescription* td = 0;
256  if(!td){
257  const TypeDescription *sub = get_type_description((T*)0);
259  (*subs)[0] = sub;
260  td = new TypeDescription("HexTriquadraticLgn", subs,
261  std::string(__FILE__),
262  "SCIRun",
264  }
265  return td;
266 }
267 
269 template <class T>
270 void
272 {
273  stream.begin_class(get_type_description(this)->get_name(),
275  Pio(stream, this->nodes_);
276  stream.end_class();
277 }
278 }
279 
280 #endif
double get_arc_length(const unsigned edge, const ElemData &cd) const
get arc length for edge
Definition: HexTriquadraticLgn.h:204
T value_type
Definition: HexTriquadraticLgn.h:69
T interpolate(const VECTOR &coords, const ElemData &cd) const
get value at parametric coordinate
Definition: HexTriquadraticLgn.h:87
void get_derivate_weights(const VECTOR &coords, double *w) const
Definition: HexTriquadraticLgn.h:82
Definition: HexElementsWeights.h:36
Definition: Persistent.h:89
static const std::string type_name(int n=-1)
Definition: HexTriquadraticLgn.h:231
static int number_of_vertices()
return number of vertices
Definition: HexTriquadraticLgn.h:51
void get_weights(const VECTOR &coords, double *w) const
Definition: HexTriquadraticLgn.h:78
double get_area(const unsigned face, const ElemData &cd) const
get area
Definition: HexTriquadraticLgn.h:211
HexTriquadraticLgnUnitElement()
Definition: HexTriquadraticLgn.h:48
Definition: TypeDescription.h:45
#define SCISHARE
Definition: share.h:39
std::vector< const TypeDescription * > td_vec
Definition: TypeDescription.h:56
std::vector< T > nodes_
Definition: Basis.h:159
Class for creating geometrical approximations of Hex meshes.
Definition: HexTrilinearLgn.h:112
#define ASSERT(condition)
Definition: Assert.h:110
virtual int begin_class(const std::string &name, int current_version)
Definition: Persistent.cc:143
void get_quadratic_weights(const VECTOR &coords, double *w) const
get weight factors at parametric coordinate
Definition: HexElementsWeights.h:87
const string find_type_name(float *)
Definition: TypeName.cc:63
const int HEXTRIQUADRATICLGN_VERSION
Definition: HexTriquadraticLgn.h:268
virtual ~HexTriquadraticLgn()
Definition: HexTriquadraticLgn.h:72
Class with weights and coordinates for 3rd order Gaussian integration.
Definition: HexTrilinearLgn.h:331
const char * name[]
Definition: BoostGraphExampleTests.cc:87
virtual void io(Piostream &str)
Definition: HexTriquadraticLgn.h:271
Class for describing interfaces to basis elements with additional nodes.
Definition: Basis.h:169
static int dofs()
return degrees of freedom
Definition: HexTriquadraticLgn.h:53
static int polynomial_order()
Definition: HexTriquadraticLgn.h:74
static const std::string make_template_id(const std::string &templateName, const std::string &templateParam)
Definition: TypeName.h:62
double get_volume(const ElemData &cd) const
get volume
Definition: HexTriquadraticLgn.h:218
void Pio(Piostream &stream, Array1< T > &array)
Definition: Array1.h:65
virtual ~HexTriquadraticLgnUnitElement()
Definition: HexTriquadraticLgn.h:49
virtual void end_class()
Definition: Persistent.cc:178
Class for describing unit geometry of HexTriquadraticLgn.
Definition: HexTriquadraticLgn.h:44
HexTriquadraticLgn()
Definition: HexTriquadraticLgn.h:71
Definition: HexTriquadraticLgn.h:61
void derivate(const VECTOR1 &coords, const ElemData &cd, VECTOR2 &derivs) const
get first derivative at parametric coordinate
Definition: HexTriquadraticLgn.h:117
Definition: HexSamplingSchemes.h:43
Definition: HexTrilinearLgn.h:197
int n
Definition: eab.py:9
bool get_coords(const ElemBasis *pEB, VECTOR &coords, const T &value, const ElemData &cd) const
find value in interpolation for given value
Definition: HexTrilinearLgn.h:206
void get_quadratic_derivate_weights(const VECTOR &coords, double *w) const
get weight factors of derivative at parametric coordinate
Definition: HexElementsWeights.h:114
double get_volume3(const ElemBasis *pEB, const ElemData &cd)
Definition: Locate.h:179
bool get_coords(VECTOR &coords, const T &value, const ElemData &cd) const
get parametric coordinate for value within the element
Definition: HexTriquadraticLgn.h:195
Class for describing unit geometry of HexTrilinearLgn.
Definition: HexTrilinearLgn.h:49
Definition: TypeDescription.h:49
const TypeDescription * get_type_description(Core::Basis::ConstantBasis< T > *)
Definition: Constant.h:209