21 #ifndef fiber_typical_test_scalar_kernel_trial_integral_hpp
22 #define fiber_typical_test_scalar_kernel_trial_integral_hpp
24 #include "test_kernel_trial_integral.hpp"
26 #include "default_test_kernel_trial_integral.hpp"
27 #include "test_scalar_kernel_trial_integrand_functor.hpp"
29 #include <tbb/enumerable_thread_specific.h>
34 template <
typename BasisFunctionType_,
typename KernelType_,
typename ResultType_>
40 typedef typename Base::CoordinateType CoordinateType;
41 typedef typename Base::BasisFunctionType BasisFunctionType;
42 typedef typename Base::KernelType KernelType;
43 typedef typename Base::ResultType ResultType;
46 size_t& testGeomDeps,
size_t& trialGeomDeps)
const;
70 template <
typename BasisFunctionType_,
typename KernelType_,
typename ResultType_>
73 BasisFunctionType_, KernelType_, ResultType_>
80 template <
typename BasisFunctionType_,
typename ResultType_>
82 BasisFunctionType_, BasisFunctionType_, ResultType_> :
84 BasisFunctionType_, BasisFunctionType_, ResultType_>
87 BasisFunctionType_, BasisFunctionType_, ResultType_>
Base;
89 typedef typename Base::CoordinateType CoordinateType;
90 typedef typename Base::BasisFunctionType BasisFunctionType;
91 typedef typename Base::KernelType KernelType;
92 typedef typename Base::ResultType ResultType;
102 const std::vector<CoordinateType>& testQuadWeights,
103 const std::vector<CoordinateType>& trialQuadWeights,
104 arma::Mat<ResultType>& result)
const;
112 const std::vector<CoordinateType>& quadWeights,
113 arma::Mat<ResultType>& result)
const;
116 template <
typename CoordinateType_>
118 std::complex<CoordinateType_>, std::complex<CoordinateType_> > :
120 std::complex<CoordinateType_>, std::complex<CoordinateType_> >
123 std::complex<CoordinateType_>, std::complex<CoordinateType_> >
Base;
125 typedef typename Base::CoordinateType CoordinateType;
126 typedef typename Base::BasisFunctionType BasisFunctionType;
127 typedef typename Base::KernelType KernelType;
128 typedef typename Base::ResultType ResultType;
138 const std::vector<CoordinateType>& testQuadWeights,
139 const std::vector<CoordinateType>& trialQuadWeights,
140 arma::Mat<ResultType>& result)
const;
148 const std::vector<CoordinateType>& quadWeights,
149 arma::Mat<ResultType>& result)
const;
152 template <
typename CoordinateType_>
154 CoordinateType_, std::complex<CoordinateType_> > :
156 std::complex<CoordinateType_>, CoordinateType_, std::complex<CoordinateType_> >
159 std::complex<CoordinateType_>, CoordinateType_, std::complex<CoordinateType_> >
162 typedef typename Base::CoordinateType CoordinateType;
163 typedef typename Base::BasisFunctionType BasisFunctionType;
164 typedef typename Base::KernelType KernelType;
165 typedef typename Base::ResultType ResultType;
176 const std::vector<CoordinateType>& testQuadWeights,
177 const std::vector<CoordinateType>& trialQuadWeights,
178 arma::Mat<ResultType>& result)
const;
186 const std::vector<CoordinateType>& quadWeights,
187 arma::Mat<ResultType>& result)
const;
190 DefaultTestKernelTrialIntegral<
192 BasisFunctionType, KernelType, ResultType> > m_standardIntegral;
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
Retrieve types of geometrical data on which the integrand of this integral depends explicitly...
Definition: typical_test_scalar_kernel_trial_integral.cpp:376
Functor evaluating the integrand of most standard boundary integral operators.
Definition: test_scalar_kernel_trial_integrand_functor.hpp:43
Implementation of the TestKernelTrialIntegral interface for "typical" integrals, taking advantage of ...
Definition: typical_test_scalar_kernel_trial_integral.hpp:71
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: typical_test_scalar_kernel_trial_integral.hpp:35
Definition: collection_of_4d_arrays.hpp:41