BEM++  2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
linear_scalar_shapeset_barycentric.hpp
1 // Copyright (C) 2011-2013 by the BEM++ Authors
2 //
3 // Permission is hereby granted, free of charge, to any person obtaining a copy
4 // of this software and associated documentation files (the "Software"), to deal
5 // in the Software without restriction, including without limitation the rights
6 // to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
7 // copies of the Software, and to permit persons to whom the Software is
8 // furnished to do so, subject to the following conditions:
9 //
10 // The above copyright notice and this permission notice shall be included in
11 // all copies or substantial portions of the Software.
12 //
13 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
14 // IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
15 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
16 // AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
17 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
18 // OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
19 // THE SOFTWARE.
20 
21 
22 #ifndef linear_scalar_shapeset_barycentric_hpp
23 #define linear_scalar_shapeset_barycentric_hpp
24 
25 #include "basis.hpp"
26 #include "basis_data.hpp"
27 #include "piecewise_linear_continuous_scalar_basis.hpp"
28 
29 namespace Fiber
30 {
31 
32 template <typename ValueType>
33 class LinearScalarShapesetBarycentric : public Basis<ValueType>
34 {
35 public:
36  typedef typename Basis<ValueType>::CoordinateType CoordinateType;
37  enum BasisType {TYPE1,TYPE2};
38 
39 public:
40 
41  LinearScalarShapesetBarycentric(BasisType type) :
42  m_type(type)
43  {
44  }
45 
46  virtual int size() const {
47  return 3;
48  }
49 
50  virtual int order() const {
51  return 1;
52  }
53 
54  virtual void evaluate(size_t what,
55  const arma::Mat<CoordinateType>& points,
56  LocalDofIndex localDofIndex,
57  BasisData<ValueType>& data) const {
58 
60 
61  linearBasis.evaluate(what,points,ALL_DOFS,temp);
62 
63  ValueType coeffs[2][3][3] =
64  {
65  {{1,1./3,1./2},{0,1./3,0},{0,1./3,1./2}},
66  {{1,1./2,1./3},{0,1./2,1./3},{0,0,1./3}}
67  };
68 
69  int type_index;
70  if (m_type == TYPE1)
71  type_index = 0;
72  else
73  type_index = 1;
74 
75  if (localDofIndex != ALL_DOFS){
76 
77  if (what & VALUES){
78  data.values.set_size(1,1,temp.values.extent(2));
79  for (int i=0;i<data.values.extent(2);++i){
80  data.values(0,0,i) = 0;
81  for (int j=0;j<3;++j)
82  data.values(0,0,i) +=coeffs[type_index][localDofIndex][j]*
83  temp.values(0,j,i);
84  }
85  }
86  if (what & DERIVATIVES){
87  data.derivatives.set_size(1,temp.derivatives.extent(1),1,temp.derivatives.extent(3));
88  for (int i=0;i<data.derivatives.extent(1);++i)
89  for (int j=0;j<data.derivatives.extent(3);++j){
90  data.derivatives(0,i,0,j) = 0;
91  for (int k=0;k<3;++k)
92  data.derivatives(0,i,0,j)+=
93  coeffs[type_index][localDofIndex][k]*
94  temp.derivatives(0,i,k,j);
95  }
96  }
97 
98  }
99  else {
100  if (what & VALUES){
101  data.values.set_size(1,3,temp.values.extent(2));
102  for (int dofIndex=0;dofIndex<3;++dofIndex){
103  for (int i=0;i<data.values.extent(2);++i){
104  data.values(0,dofIndex,i) = 0;
105  for (int j=0;j<3;++j)
106  data.values(0,dofIndex,i) +=coeffs[type_index][dofIndex][j]*
107  temp.values(0,j,i);
108 
109  }
110  }
111  }
112 
113  if (what & DERIVATIVES){
114  data.derivatives.set_size(1,temp.derivatives.extent(1),3,temp.derivatives.extent(3));
115  for (int dofIndex=0;dofIndex<3;++dofIndex)
116  for (int i=0;i<data.derivatives.extent(1);++i)
117  for (int j=0;j<data.derivatives.extent(3);++j){
118  data.derivatives(0,i,dofIndex,j) = 0;
119  for (int k=0;k<3;++k)
120  data.derivatives(0,i,dofIndex,j)+=
121  coeffs[type_index][dofIndex][k]*
122  temp.derivatives(0,i,k,j);
123  }
124 
125 
126  }
127  }
128 
129  }
130 
131 
132 
133  virtual std::pair<const char*,int> clCodeString (bool isTestBasis) const {
134  throw std::runtime_error("PiecewiseLinearContinuousScalarBasisBarycentric::clCodeString():"
135  "OpenCL not supported for this basis type.");
136 
137  }
138 
139 private:
141  BasisType m_type;
142 
143 
144 };
145 
146 } // namespace Fiber
147 
148 
149 
150 #endif
151 
virtual void evaluate(size_t what, const arma::Mat< CoordinateType > &points, LocalDofIndex localDofIndex, BasisData< ValueType > &data) const
Evaluate the shape functions making up this shapeset and/or their derivatives at specified points...
Definition: linear_scalar_shapeset_barycentric.hpp:54
Definition: linear_scalar_shapeset_barycentric.hpp:33
_4dArray< ValueType > derivatives
Derivatives of shape functions.
Definition: basis_data.hpp:59
virtual int size() const
Return the number of shape functions.
Definition: linear_scalar_shapeset_barycentric.hpp:46
virtual int order() const
Return the maximum polynomial order of shape functions.
Definition: linear_scalar_shapeset_barycentric.hpp:50
_3dArray< ValueType > values
Values of shape functions.
Definition: basis_data.hpp:51
virtual std::pair< const char *, int > clCodeString(bool isTestBasis) const
Returns an OpenCL code snippet for shape function evaluation.
Definition: linear_scalar_shapeset_barycentric.hpp:133
Collection of shape functions defined on a reference element.
Definition: basis.hpp:36
Storage of values and/or derivatives of shape functions.
Definition: basis_data.hpp:44
virtual void evaluate(size_t what, const arma::Mat< CoordinateType > &points, LocalDofIndex localDofIndex, BasisData< ValueType > &data) const
Evaluate the shape functions making up this shapeset and/or their derivatives at specified points...
Definition: linear_scalar_shapeset.hpp:88