21 #ifndef fiber_test_kernel_trial_integral_hpp
22 #define fiber_test_kernel_trial_integral_hpp
24 #include "../common/common.hpp"
26 #include "scalar_traits.hpp"
28 #include "../common/armadillo_fwd.hpp"
35 template <
typename T>
class CollectionOf3dArrays;
36 template <
typename T>
class CollectionOf4dArrays;
37 template <
typename CoordinateType>
class GeometricalData;
68 template <
typename BasisFunctionType_,
typename KernelType_,
73 typedef BasisFunctionType_ BasisFunctionType;
74 typedef KernelType_ KernelType;
75 typedef ResultType_ ResultType;
111 size_t& testGeomDeps,
size_t& trialGeomDeps)
const = 0;
179 const std::vector<CoordinateType>& testQuadWeights,
180 const std::vector<CoordinateType>& trialQuadWeights,
181 arma::Mat<ResultType>& result)
const = 0;
246 const std::vector<CoordinateType>& quadWeights,
247 arma::Mat<ResultType>& result)
const = 0;
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