BEM++  2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
separable_numerical_test_kernel_trial_integrator.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_separable_numerical_test_kernel_trial_integrator_hpp
22 #define fiber_separable_numerical_test_kernel_trial_integrator_hpp
23 
24 #include "../common/common.hpp"
25 
26 #include "bempp/common/config_opencl.hpp"
27 
28 #include "test_kernel_trial_integrator.hpp"
29 
30 #include <tbb/enumerable_thread_specific.h>
31 
32 namespace Fiber
33 {
34 
36 class OpenClHandler;
37 template <typename CoordinateType> class CollectionOfShapesetTransformations;
38 template <typename ValueType> class CollectionOfKernels;
39 template <typename CoordinateType> class RawGridGeometry;
40 template <typename BasisFunctionType, typename KernelType, typename ResultType>
41 class TestKernelTrialIntegral;
45 template <typename BasisFunctionType, typename KernelType,
46  typename ResultType, typename GeometryFactory>
48  public TestKernelTrialIntegrator<BasisFunctionType, KernelType, ResultType>
49 {
50 public:
52  typedef typename Base::CoordinateType CoordinateType;
53  typedef typename Base::ElementIndexPair ElementIndexPair;
54 
56  const arma::Mat<CoordinateType>& localTestQuadPoints,
57  const arma::Mat<CoordinateType>& localTrialQuadPoints,
58  const std::vector<CoordinateType>& testQuadWeights,
59  const std::vector<CoordinateType>& trialQuadWeights,
60  const GeometryFactory& testGeometryFactory,
61  const GeometryFactory& trialGeometryFactory,
62  const RawGridGeometry<CoordinateType>& testRawGeometry,
63  const RawGridGeometry<CoordinateType>& trialRawGeometry,
64  const CollectionOfShapesetTransformations<CoordinateType>& testTransformations,
65  const CollectionOfKernels<KernelType>& kernels,
66  const CollectionOfShapesetTransformations<CoordinateType>& trialTransformations,
68  const OpenClHandler& openClHandler,
69  bool cacheGeometricalData = true);
70 
72 
73  virtual void integrate(
74  CallVariant callVariant,
75  const std::vector<int>& elementIndicesA,
76  int elementIndexB,
77  const Shapeset<BasisFunctionType>& basisA,
78  const Shapeset<BasisFunctionType>& basisB,
79  LocalDofIndex localDofIndexB,
80  const std::vector<arma::Mat<ResultType>*>& result) const;
81 
82  virtual void integrate(
83  const std::vector<ElementIndexPair>& elementIndexPairs,
84  const Shapeset<BasisFunctionType>& testShapeset,
85  const Shapeset<BasisFunctionType>& trialShapeset,
86  const std::vector<arma::Mat<ResultType>*>& result) const;
87 
88 private:
89  void integrateCpu(
90  CallVariant callVariant,
91  const std::vector<int>& elementIndicesA,
92  int elementIndexB,
93  const Shapeset<BasisFunctionType>& basisA,
94  const Shapeset<BasisFunctionType>& basisB,
95  LocalDofIndex localDofIndexB,
96  const std::vector<arma::Mat<ResultType>*>& result) const;
97 
98  void integrateCl(
99  CallVariant callVariant,
100  const std::vector<int>& elementIndicesA,
101  int elementIndexB,
102  const Shapeset<BasisFunctionType>& basisA,
103  const Shapeset<BasisFunctionType>& basisB,
104  LocalDofIndex localDofIndexB,
105  const std::vector<arma::Mat<ResultType>*>& result) const;
106 
107  void integrateCpu(
108  const std::vector<ElementIndexPair>& elementIndexPairs,
109  const Shapeset<BasisFunctionType>& testShapeset,
110  const Shapeset<BasisFunctionType>& trialShapeset,
111  const std::vector<arma::Mat<ResultType>*>& result) const;
112 
113  void integrateCl(
114  const std::vector<ElementIndexPair>& elementIndexPairs,
115  const Shapeset<BasisFunctionType>& testShapeset,
116  const Shapeset<BasisFunctionType>& trialShapeset,
117  const std::vector<arma::Mat<ResultType>*>& result) const;
118 
119  void precalculateGeometricalData();
120  void precalculateGeometricalDataOnSingleGrid(
121  const arma::Mat<CoordinateType>& localQuadPoints,
122  const GeometryFactory& geometryFactory,
123  const RawGridGeometry<CoordinateType>& rawGeometry,
124  size_t geomDeps,
125  std::vector<GeometricalData<CoordinateType> >& geomData);
126 
131  const std::pair<const char*, int> clStrIntegrateRowOrCol () const;
132 
133  arma::Mat<CoordinateType> m_localTestQuadPoints;
134  arma::Mat<CoordinateType> m_localTrialQuadPoints;
135  std::vector<CoordinateType> m_testQuadWeights;
136  std::vector<CoordinateType> m_trialQuadWeights;
137 
138  const GeometryFactory& m_testGeometryFactory;
139  const GeometryFactory& m_trialGeometryFactory;
140  const RawGridGeometry<CoordinateType>& m_testRawGeometry;
141  const RawGridGeometry<CoordinateType>& m_trialRawGeometry;
142 
143  const CollectionOfShapesetTransformations<CoordinateType>& m_testTransformations;
144  const CollectionOfKernels<KernelType>& m_kernels;
145  const CollectionOfShapesetTransformations<CoordinateType>& m_trialTransformations;
147 
148  const OpenClHandler& m_openClHandler;
149  bool m_cacheGeometricalData;
150 
151  std::vector<GeometricalData<CoordinateType> > m_cachedTestGeomData;
152  std::vector<GeometricalData<CoordinateType> > m_cachedTrialGeomData;
153  mutable tbb::enumerable_thread_specific<GeometricalData<CoordinateType> >
154  m_testGeomData, m_trialGeomData;
155 
156 
157 #ifdef WITH_OPENCL
158  cl::Buffer *clTestQuadPoints;
159  cl::Buffer *clTrialQuadPoints;
160  cl::Buffer *clTestQuadWeights;
161  cl::Buffer *clTrialQuadWeights;
162 #endif
163 };
164 
165 } // namespace Fiber
166 
167 #include "separable_numerical_test_kernel_trial_integrator_imp.hpp"
168 
169 #endif
Integration over pairs of elements.
Definition: test_kernel_trial_integrator.hpp:42
Storage of geometrical data.
Definition: geometrical_data.hpp:54
Collection of shape functions defined on a reference element.
Definition: default_local_assembler_for_operators_on_surfaces_utilities.hpp:34
Integration over pairs of elements on tensor-product point grids.
Definition: separable_numerical_test_kernel_trial_integrator.hpp:47
Definition: default_local_assembler_for_operators_on_surfaces_utilities.hpp:35
const std::pair< const char *, int > clStrIntegrateRowOrCol() const
Returns an OpenCL code snippet containing the clIntegrate kernel function for integrating a single ro...
Definition: separable_numerical_test_kernel_trial_integrator_imp.hpp:963
Definition: opencl_handler.hpp:208