21 #ifndef fiber_default_local_assembler_for_integral_operators_on_surfaces_hpp
22 #define fiber_default_local_assembler_for_integral_operators_on_surfaces_hpp
24 #include "../common/common.hpp"
26 #include "local_assembler_for_integral_operators.hpp"
28 #include "_2d_array.hpp"
29 #include "accuracy_options.hpp"
30 #include "default_local_assembler_for_operators_on_surfaces_utilities.hpp"
31 #include "element_pair_topology.hpp"
33 #include "parallelization_options.hpp"
34 #include "shared_ptr.hpp"
35 #include "test_kernel_trial_integrator.hpp"
36 #include "verbosity_level.hpp"
38 #include <boost/static_assert.hpp>
39 #include <boost/tuple/tuple_comparison.hpp>
40 #include <tbb/concurrent_unordered_map.h>
41 #include <tbb/mutex.h>
53 template <
typename CoordinateType>
class CollectionOfShapesetTransformations;
54 template <
typename ValueType>
class CollectionOfKernels;
55 template <
typename BasisFunctionType,
typename KernelType,
typename ResultType>
56 class TestKernelTrialIntegral;
58 template <
typename CoordinateType>
class RawGridGeometry;
60 template <
typename CoordinateType>
61 class QuadratureDescriptorSelectorForIntegralOperators;
62 template <
typename CoordinateType>
class DoubleQuadratureRuleFamily;
65 template <
typename BasisFunctionType,
typename KernelType,
66 typename ResultType,
typename GeometryFactory>
74 const shared_ptr<const GeometryFactory>& testGeometryFactory,
75 const shared_ptr<const GeometryFactory>& trialGeometryFactory,
84 const shared_ptr<const OpenClHandler>& openClHandler,
86 VerbosityLevel::Level verbosityLevel,
87 bool cacheSingularIntegrals,
94 CallVariant callVariant,
95 const std::vector<int>& elementIndicesA,
97 LocalDofIndex localDofIndexB,
98 std::vector<arma::Mat<ResultType> >& result,
99 CoordinateType nominalDistance = -1.);
102 const std::vector<int>& testElementIndices,
103 const std::vector<int>& trialElementIndices,
105 CoordinateType nominalDistance = -1.);
112 typedef typename Integrator::ElementIndexPair ElementIndexPair;
114 BasisFunctionType> Utilities;
120 template <
typename T1,
typename T2>
121 struct alternative_less {
122 bool operator() (
const std::pair<T1, T2>& a,
const std::pair<T1, T2>& b)
const {
123 return a.second < b.second
124 || (!(b.second < a.second) && a.first < b.first);
133 typedef alternative_less<
134 typename ElementIndexPair::first_type,
135 typename ElementIndexPair::second_type> ElementIndexPairCompare;
145 typedef std::set<ElementIndexPair, ElementIndexPairCompare> ElementIndexPairSet;
147 bool testAndTrialGridsAreIdentical()
const;
149 void cacheSingularLocalWeakForms();
150 void findPairsOfAdjacentElements(ElementIndexPairSet& pairs)
const;
151 void cacheLocalWeakForms(
const ElementIndexPairSet& elementIndexPairs);
153 const Integrator& selectIntegrator(
154 int testElementIndex,
int trialElementIndex,
155 CoordinateType nominalDistance = -1.);
160 shared_ptr<const GeometryFactory> m_testGeometryFactory;
161 shared_ptr<const GeometryFactory> m_trialGeometryFactory;
162 shared_ptr<const RawGridGeometry<CoordinateType> > m_testRawGeometry;
163 shared_ptr<const RawGridGeometry<CoordinateType> > m_trialRawGeometry;
164 shared_ptr<const std::vector<const Shapeset<BasisFunctionType>*> > m_testShapesets;
165 shared_ptr<const std::vector<const Shapeset<BasisFunctionType>*> > m_trialShapesets;
166 shared_ptr<const CollectionOfShapesetTransformations<CoordinateType> > m_testTransformations;
167 shared_ptr<const CollectionOfKernels<KernelType> > m_kernels;
168 shared_ptr<const CollectionOfShapesetTransformations<CoordinateType> > m_trialTransformations;
169 shared_ptr<const TestKernelTrialIntegral<BasisFunctionType, KernelType, ResultType> > m_integral;
170 shared_ptr<const OpenClHandler> m_openClHandler;
172 VerbosityLevel::Level m_verbosityLevel;
173 shared_ptr<const QuadratureDescriptorSelectorForIntegralOperators<CoordinateType> > m_quadDescSelector;
174 shared_ptr<const DoubleQuadratureRuleFamily<CoordinateType> > m_quadRuleFamily;
177 Integrator*> IntegratorMap;
178 IntegratorMap m_testKernelTrialIntegrators;
179 mutable tbb::mutex m_integratorCreationMutex;
181 enum { INVALID_INDEX = INT_MAX };
199 #include "default_local_assembler_for_integral_operators_on_surfaces_imp.hpp"
Traits of scalar types.
Definition: scalar_traits.hpp:40
Parallel operation settings.
Definition: parallelization_options.hpp:32
Integration over pairs of elements.
Definition: test_kernel_trial_integrator.hpp:42
Definition: default_local_assembler_for_operators_on_surfaces_utilities.hpp:42
virtual CoordinateType estimateRelativeScale(CoordinateType minDist) const
Estimate how fast the entries in the matrix of this operator decay with interelement distance...
Definition: default_local_assembler_for_integral_operators_on_surfaces_imp.hpp:382
virtual void evaluateLocalWeakForms(CallVariant callVariant, const std::vector< int > &elementIndicesA, int elementIndexB, LocalDofIndex localDofIndexB, std::vector< arma::Mat< ResultType > > &result, CoordinateType nominalDistance=-1.)
Assemble local weak forms.
Definition: default_local_assembler_for_integral_operators_on_surfaces_imp.hpp:160
Parameters of a quadrature rule used in the evaluation of integrals over pairs of elements...
Definition: double_quadrature_descriptor.hpp:36
Collection of shape functions defined on a reference element.
Definition: default_local_assembler_for_operators_on_surfaces_utilities.hpp:34
Abstract interface of a local assembler for integral operators.
Definition: local_assembler_for_integral_operators.hpp:48
Definition: default_local_assembler_for_integral_operators_on_surfaces.hpp:67
Definition: default_local_assembler_for_operators_on_surfaces_utilities.hpp:35
Family of quadrature rules over pairs of elements.
Definition: double_quadrature_rule_family.hpp:52
Quadrature descriptor selector used during the discretization of boundary integral operators...
Definition: quadrature_descriptor_selector_factory.hpp:37
Simple implementation of a 2D Fortran-ordered array.
Definition: _2d_array.hpp:41