BEM++  2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
default_test_trial_integral_imp.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_default_test_trial_integral_imp_hpp
22 #define fiber_default_test_trial_integral_imp_hpp
23 
24 #include "default_test_trial_integral.hpp"
25 
26 #include "geometrical_data.hpp"
27 
28 #include <algorithm>
29 
30 namespace Fiber
31 {
32 
33 /*
34 template <typename BasisFunctionType_, typename KernelType_,
35  typename ResultType_>
36 class IntegrandFunctor
37 {
38 public:
39  typedef BasisFunctionType_ BasisFunctionType;
40  typedef KernelType_ KernelType;
41  typedef ResultType_ ResultType;
42  typedef typename ScalarTraits<ResultType>::RealType CoordinateType;
43 
44  void addGeometricalDependencies(int& geomDeps) const;
45 
46  ResultType evaluate(
47  const ConstGeometricalDataSlice<CoordinateType>& geomData,
48  const CollectionOf1dSlicesOfConst3dArrays<BasisFunctionType>&
49  testTransformations,
50  const CollectionOf1dSlicesOfConst3dArrays<BasisFunctionType>&
51  trialTransformations) const;
52 };
53 */
54 
55 template <typename IntegrandFunctor>
56 void DefaultTestTrialIntegral<IntegrandFunctor>::
57 addGeometricalDependencies(size_t& geomDeps) const
58 {
59  geomDeps |= INTEGRATION_ELEMENTS;
60  m_functor.addGeometricalDependencies(geomDeps);
61 }
62 
63 template <typename IntegrandFunctor>
64 void DefaultTestTrialIntegral<IntegrandFunctor>::evaluate(
65  const GeometricalData<CoordinateType>& geomData,
66  const CollectionOf3dArrays<BasisFunctionType>& testValues,
67  const CollectionOf3dArrays<BasisFunctionType>& trialValues,
68  const std::vector<CoordinateType>& weights,
69  arma::Mat<ResultType>& result) const
70 {
71  const size_t pointCount = weights.size();
72  assert(testValues.size() == 1);
73  assert(trialValues.size() == 1);
74 
75  // We don't care about the number of rows of testValues[0] and trialValues[0]
76  // -- it's up to the integrand functor
77  const size_t testDofCount = testValues[0].extent(1);
78  const size_t trialDofCount = trialValues[0].extent(1);
79  assert(testValues[0].extent(2) == pointCount);
80  assert(trialValues[0].extent(2) == pointCount);
81  for (size_t trialDof = 0; trialDof < trialDofCount; ++trialDof)
82  for (size_t testDof = 0; testDof < testDofCount; ++testDof) {
83  ResultType sum = 0.;
84  for (size_t point = 0; point < pointCount; ++point)
85  sum += m_functor.evaluate(
86  geomData.const_slice(point),
87  testValues.const_slice(testDof, point),
88  trialValues.const_slice(trialDof, point)) *
89  geomData.integrationElements(point) *
90  weights[point];
91  result(testDof, trialDof) = sum;
92  }
93 }
94 
95 } // namespace Fiber
96 
97 #endif