BEM++  2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
constant_scalar_shapeset.hpp
1 // Copyright (C) 2011-2012 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 #ifndef fiber_constant_scalar_shapeset_hpp
22 #define fiber_constant_scalar_shapeset_hpp
23 
24 #include "../common/common.hpp"
25 
26 #include "basis.hpp"
27 #include "basis_data.hpp"
28 
29 #include <algorithm>
30 
31 namespace Fiber
32 {
33 
35 template <typename ValueType>
36 class ConstantScalarShapeset : public Basis<ValueType>
37 {
38 public:
39  typedef typename Basis<ValueType>::CoordinateType CoordinateType;
40 
41  virtual int size() const {
42  return 1;
43  }
44 
45  virtual int order() const {
46  return 0;
47  }
48 
49  virtual void evaluate(size_t what,
50  const arma::Mat<CoordinateType>& points,
51  LocalDofIndex localDofIndex,
52  BasisData<ValueType>& data) const {
53  if (localDofIndex != ALL_DOFS && localDofIndex != 0)
54  throw std::invalid_argument("ConstantScalarShapeset::evaluate(): "
55  "Invalid localDofIndex");
56  // Since there is only one shape function, there is no difference
57  // between calculating all shape functions and just one.
58 
59  const int componentCount = 1;
60  const int functionCount = 1;
61  const int pointCount = points.n_cols;
62  if (what & VALUES)
63  {
64  data.values.set_size(componentCount, functionCount, pointCount);
65  std::fill(data.values.begin(), data.values.end(), 1.);
66  }
67  if (what & DERIVATIVES)
68  {
69  const int coordCount = points.n_rows;
70  data.derivatives.set_size(componentCount, coordCount,
71  functionCount, pointCount);
72  std::fill(data.derivatives.begin(), data.derivatives.end(), 0.);
73  }
74  }
75 };
76 
77 } // namespace Fiber
78 
79 #endif
_4dArray< ValueType > derivatives
Derivatives of shape functions.
Definition: basis_data.hpp:59
virtual int order() const
Return the maximum polynomial order of shape functions.
Definition: constant_scalar_shapeset.hpp:45
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: constant_scalar_shapeset.hpp:49
virtual int size() const
Return the number of shape functions.
Definition: constant_scalar_shapeset.hpp:41
_3dArray< ValueType > values
Values of shape functions.
Definition: basis_data.hpp:51
A shapeset containing only one function: the constant function.
Definition: constant_scalar_shapeset.hpp:36
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