BEM++  2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
default_test_kernel_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_kernel_trial_integral_imp_hpp
22 #define fiber_default_test_kernel_trial_integral_imp_hpp
23 
24 #include "default_test_kernel_trial_integral.hpp"
25 
26 #include "geometrical_data.hpp"
27 
28 #include <cassert>
29 
30 namespace Fiber
31 {
32 
33 template <typename IntegrandFunctor>
34 void DefaultTestKernelTrialIntegral<IntegrandFunctor>::
35 addGeometricalDependencies(size_t& testGeomDeps, size_t& trialGeomDeps) const
36 {
37  testGeomDeps |= INTEGRATION_ELEMENTS;
38  trialGeomDeps |= INTEGRATION_ELEMENTS;
39 
40  m_functor.addGeometricalDependencies(testGeomDeps, trialGeomDeps);
41 }
42 
43 template <typename IntegrandFunctor>
44 void DefaultTestKernelTrialIntegral<IntegrandFunctor>::
45 evaluateWithTensorQuadratureRule(
46  const GeometricalData<CoordinateType>& testGeomData,
47  const GeometricalData<CoordinateType>& trialGeomData,
48  const CollectionOf3dArrays<BasisFunctionType>& testValues,
49  const CollectionOf3dArrays<BasisFunctionType>& trialValues,
50  const CollectionOf4dArrays<KernelType>& kernelValues,
51  const std::vector<CoordinateType>& testQuadWeights,
52  const std::vector<CoordinateType>& trialQuadWeights,
53  arma::Mat<ResultType>& result) const
54 {
55  // Evaluate constants
56 
57  const size_t testDofCount = testValues[0].extent(1);
58  const size_t trialDofCount = trialValues[0].extent(1);
59 
60  const size_t testPointCount = testQuadWeights.size();
61  const size_t trialPointCount = trialQuadWeights.size();
62 
63  // Assert that array dimensions are correct
64 
65  for (size_t i = 0; i < kernelValues.size(); ++i) {
66  assert(kernelValues[i].extent(2) == testPointCount);
67  assert(kernelValues[i].extent(3) == trialPointCount);
68  }
69  for (size_t i = 0; i < testValues.size(); ++i)
70  assert(testValues[i].extent(2) == testPointCount);
71  for (size_t i = 0; i < trialValues.size(); ++i)
72  assert(trialValues[i].extent(2) == trialPointCount);
73 
74  assert(result.n_rows == testDofCount);
75  assert(result.n_cols == trialDofCount);
76 
77  // Integrate
78 
79  for (size_t trialDof = 0; trialDof < trialDofCount; ++trialDof)
80  for (size_t testDof = 0; testDof < testDofCount; ++testDof)
81  {
82  ResultType sum = 0.;
83  for (size_t trialPoint = 0; trialPoint < trialPointCount; ++trialPoint) {
84  const CoordinateType trialWeight =
85  trialGeomData.integrationElements(trialPoint) *
86  trialQuadWeights[trialPoint];
87  ResultType partialSum = 0.;
88  for (size_t testPoint = 0; testPoint < testPointCount; ++testPoint) {
89  const CoordinateType testWeight =
90  testGeomData.integrationElements(testPoint) *
91  testQuadWeights[testPoint];
92  partialSum += m_functor.evaluate(
93  testGeomData.const_slice(testPoint),
94  trialGeomData.const_slice(trialPoint),
95  testValues.const_slice(testDof, testPoint),
96  trialValues.const_slice(trialDof, trialPoint),
97  kernelValues.const_slice(testPoint, trialPoint)) *
98  testWeight;
99  }
100  sum += partialSum * trialWeight;
101  }
102  result(testDof, trialDof) = sum;
103  }
104 }
105 
106 template <typename IntegrandFunctor>
107 void DefaultTestKernelTrialIntegral<IntegrandFunctor>::
108 evaluateWithNontensorQuadratureRule(
109  const GeometricalData<CoordinateType>& testGeomData,
110  const GeometricalData<CoordinateType>& trialGeomData,
111  const CollectionOf3dArrays<BasisFunctionType>& testValues,
112  const CollectionOf3dArrays<BasisFunctionType>& trialValues,
113  const CollectionOf3dArrays<KernelType>& kernelValues,
114  const std::vector<CoordinateType>& quadWeights,
115  arma::Mat<ResultType>& result) const
116 {
117  // Evaluate constants
118 
119  const size_t testDofCount = testValues[0].extent(1);
120  const size_t trialDofCount = trialValues[0].extent(1);
121 
122  const size_t pointCount = quadWeights.size();
123 
124  // Assert that array dimensions are correct
125 
126  for (size_t i = 0; i < kernelValues.size(); ++i)
127  assert(kernelValues[i].extent(2) == pointCount);
128  for (size_t i = 0; i < testValues.size(); ++i)
129  assert(testValues[i].extent(2) == pointCount);
130  for (size_t i = 0; i < trialValues.size(); ++i)
131  assert(trialValues[i].extent(2) == pointCount);
132 
133  // Integrate
134 
135  for (size_t trialDof = 0; trialDof < trialDofCount; ++trialDof)
136  for (size_t testDof = 0; testDof < testDofCount; ++testDof)
137  {
138  ResultType sum = 0.;
139  for (size_t point = 0; point < pointCount; ++point)
140  sum += m_functor.evaluate(
141  testGeomData.const_slice(point),
142  trialGeomData.const_slice(point),
143  testValues.const_slice(testDof, point),
144  trialValues.const_slice(trialDof, point),
145  kernelValues.const_slice(point)) *
146  (testGeomData.integrationElements(point) *
147  trialGeomData.integrationElements(point) *
148  quadWeights[point]);
149  result(testDof, trialDof) = sum;
150  }
151 }
152 
153 } // namespace Fiber
154 
155 #endif