BEM++  2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
test_kernel_trial_integral.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_test_kernel_trial_integral_hpp
22 #define fiber_test_kernel_trial_integral_hpp
23 
24 #include "../common/common.hpp"
25 
26 #include "scalar_traits.hpp"
27 
28 #include "../common/armadillo_fwd.hpp"
29 #include <vector>
30 
31 namespace Fiber
32 {
33 
35 template <typename T> class CollectionOf3dArrays;
36 template <typename T> class CollectionOf4dArrays;
37 template <typename CoordinateType> class GeometricalData;
68 template <typename BasisFunctionType_, typename KernelType_,
69  typename ResultType_>
71 {
72 public:
73  typedef BasisFunctionType_ BasisFunctionType;
74  typedef KernelType_ KernelType;
75  typedef ResultType_ ResultType;
76  typedef typename ScalarTraits<ResultType>::RealType CoordinateType;
77 
80  {}
81 
110  virtual void addGeometricalDependencies(
111  size_t& testGeomDeps, size_t& trialGeomDeps) const = 0;
112 
174  const GeometricalData<CoordinateType>& testGeomData,
175  const GeometricalData<CoordinateType>& trialGeomData,
176  const CollectionOf3dArrays<BasisFunctionType>& testTransformations,
177  const CollectionOf3dArrays<BasisFunctionType>& trialTransformations,
178  const CollectionOf4dArrays<KernelType>& kernels,
179  const std::vector<CoordinateType>& testQuadWeights,
180  const std::vector<CoordinateType>& trialQuadWeights,
181  arma::Mat<ResultType>& result) const = 0;
182 
241  const GeometricalData<CoordinateType>& testGeomData,
242  const GeometricalData<CoordinateType>& trialGeomData,
243  const CollectionOf3dArrays<BasisFunctionType>& testTransformations,
244  const CollectionOf3dArrays<BasisFunctionType>& trialTransformations,
245  const CollectionOf3dArrays<KernelType>& kernels,
246  const std::vector<CoordinateType>& quadWeights,
247  arma::Mat<ResultType>& result) const = 0;
248 };
249 
250 } // namespace Fiber
251 
252 #endif
Traits of scalar types.
Definition: scalar_traits.hpp:40
An integral representing the weak form of an integral operator.
Definition: test_kernel_trial_integral.hpp:70
virtual void evaluateWithNontensorQuadratureRule(const GeometricalData< CoordinateType > &testGeomData, const GeometricalData< CoordinateType > &trialGeomData, const CollectionOf3dArrays< BasisFunctionType > &testTransformations, const CollectionOf3dArrays< BasisFunctionType > &trialTransformations, const CollectionOf3dArrays< KernelType > &kernels, const std::vector< CoordinateType > &quadWeights, arma::Mat< ResultType > &result) const =0
Evaluate the integral using a non-tensor-product quadrature rule.
Storage of geometrical data.
Definition: geometrical_data.hpp:54
Definition: collection_of_3d_arrays.hpp:39
virtual void addGeometricalDependencies(size_t &testGeomDeps, size_t &trialGeomDeps) const =0
Retrieve types of geometrical data on which the integrand of this integral depends explicitly...
virtual void evaluateWithTensorQuadratureRule(const GeometricalData< CoordinateType > &testGeomData, const GeometricalData< CoordinateType > &trialGeomData, const CollectionOf3dArrays< BasisFunctionType > &testTransformations, const CollectionOf3dArrays< BasisFunctionType > &trialTransformations, const CollectionOf4dArrays< KernelType > &kernels, const std::vector< CoordinateType > &testQuadWeights, const std::vector< CoordinateType > &trialQuadWeights, arma::Mat< ResultType > &result) const =0
Evaluate the integral using a tensor-product quadrature rule.
Definition: collection_of_4d_arrays.hpp:41
virtual ~TestKernelTrialIntegral()
Destructor.
Definition: test_kernel_trial_integral.hpp:79