BEM++  2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
numerical_test_function_integrator_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 #include "numerical_test_function_integrator.hpp" // To keep IDEs happy
22 
23 #include "../common/common.hpp"
24 
25 #include "shapeset.hpp"
26 #include "basis_data.hpp"
27 #include "collection_of_shapeset_transformations.hpp"
28 #include "conjugate.hpp"
29 #include "function.hpp"
30 #include "geometrical_data.hpp"
31 #include "opencl_handler.hpp"
32 #include "raw_grid_geometry.hpp"
33 #include "types.hpp"
34 
35 #include <stdexcept>
36 #include <memory>
37 
38 namespace Fiber
39 {
40 
41 template <typename BasisFunctionType, typename UserFunctionType,
42  typename ResultType, typename GeometryFactory>
43 NumericalTestFunctionIntegrator<
44 BasisFunctionType, UserFunctionType, ResultType, GeometryFactory>::
45 NumericalTestFunctionIntegrator(
46  const arma::Mat<CoordinateType>& localQuadPoints,
47  const std::vector<CoordinateType> quadWeights,
48  const GeometryFactory& geometryFactory,
49  const RawGridGeometry<CoordinateType>& rawGeometry,
50  const CollectionOfShapesetTransformations<CoordinateType>& testTransformations,
51  const Function<UserFunctionType>& function,
52  const OpenClHandler& openClHandler) :
53  m_localQuadPoints(localQuadPoints),
54  m_quadWeights(quadWeights),
55  m_geometryFactory(geometryFactory),
56  m_rawGeometry(rawGeometry),
57  m_testTransformations(testTransformations),
58  m_function(function),
59  m_openClHandler(openClHandler)
60 {
61  if (localQuadPoints.n_cols != quadWeights.size())
62  throw std::invalid_argument("NumericalTestTrialIntegrator::"
63  "NumericalTestTrialIntegrator(): "
64  "numbers of points and weights do not match");
65 
66  if (testTransformations.transformationCount() != 1)
67  throw std::invalid_argument("NumericalTestTrialIntegrator::"
68  "NumericalTestTrialIntegrator(): "
69  "test transformation collection "
70  "must contain exactly one element");
71 
72 }
73 
74 template <typename BasisFunctionType, typename UserFunctionType,
75  typename ResultType, typename GeometryFactory>
76 void NumericalTestFunctionIntegrator<
77 BasisFunctionType, UserFunctionType, ResultType, GeometryFactory>::
78 integrate(
79  const std::vector<int>& elementIndices,
80  const Shapeset<BasisFunctionType>& testShapeset,
81  arma::Mat<ResultType>& result) const
82 {
83  const size_t pointCount = m_localQuadPoints.n_cols;
84  const size_t elementCount = elementIndices.size();
85 
86  if (pointCount == 0 || elementCount == 0)
87  return;
88  // TODO: in the (pathological) case that pointCount == 0 but
89  // elementCount != 0, set elements of result to 0.
90 
91  // Evaluate constants
92  const int componentCount = m_testTransformations.resultDimension(0);
93  const int testDofCount = testShapeset.size();
94 
95  if (m_function.codomainDimension() != componentCount)
96  throw std::runtime_error("NumericalTestFunctionIntegrator::integrate(): "
97  "test functions and the \"arbitrary\" function "
98  "must have the same number of components");
99 
100  BasisData<BasisFunctionType> testBasisData;
101  GeometricalData<CoordinateType> geomData;
102 
103  size_t testBasisDeps = 0;
104  size_t geomDeps = INTEGRATION_ELEMENTS;
105 
106  m_testTransformations.addDependencies(testBasisDeps, geomDeps);
107  m_function.addGeometricalDependencies(geomDeps);
108 
109  typedef typename GeometryFactory::Geometry Geometry;
110  std::auto_ptr<Geometry> geometry(m_geometryFactory.make());
111 
113  arma::Mat<UserFunctionType> functionValues;
114 
115  result.set_size(testDofCount, elementCount);
116 
117  testShapeset.evaluate(testBasisDeps, m_localQuadPoints, ALL_DOFS, testBasisData);
118 
119  // Iterate over the elements
120  for (size_t e = 0; e < elementCount; ++e)
121  {
122  m_rawGeometry.setupGeometry(elementIndices[e], *geometry);
123  geometry->getData(geomDeps, m_localQuadPoints, geomData);
124  m_testTransformations.evaluate(testBasisData, geomData, testValues);
125  m_function.evaluate(geomData, functionValues);
126 
127  for (int testDof = 0; testDof < testDofCount; ++testDof)
128  {
129  ResultType sum = 0.;
130  for (size_t point = 0; point < pointCount; ++point)
131  for (int dim = 0; dim < componentCount; ++dim)
132  sum += m_quadWeights[point] *
133  geomData.integrationElements(point) *
134  conjugate(testValues[0](dim, testDof, point)) *
135  functionValues(dim, point);
136  result(testDof, e) = sum;
137  }
138  }
139 }
140 
141 } // namespace Fiber
Definition: collection_of_3d_arrays.hpp:39