21 #ifndef fiber_default_local_assembler_for_potential_operators_on_surfaces_hpp
22 #define fiber_default_local_assembler_for_potential_operators_on_surfaces_hpp
24 #include "../common/common.hpp"
26 #include "local_assembler_for_potential_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 "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>
52 template <
typename CoordinateType>
class CollectionOfShapesetTransformations;
53 template <
typename ValueType>
class CollectionOfKernels;
54 template <
typename BasisFunctionType,
typename KernelType,
typename ResultType>
55 class KernelTrialIntegral;
56 template <
typename CoordinateType>
class RawGridGeometry;
57 template <
typename BasisFunctionType>
58 class QuadratureDescriptorSelectorForPotentialOperators;
59 template <
typename CoordinateType>
class SingleQuadratureRuleFamily;
62 template <
typename BasisFunctionType,
typename KernelType,
63 typename ResultType,
typename GeometryFactory>
71 const arma::Mat<CoordinateType>& points,
72 const shared_ptr<const GeometryFactory>& geometryFactory,
79 VerbosityLevel::Level verbosityLevel,
81 BasisFunctionType> >& quadDescSelector,
86 const std::vector<int>& pointIndices,
87 int trialElementIndex,
88 LocalDofIndex localTrialDofIndex,
89 std::vector<arma::Mat<ResultType> >& result,
90 CoordinateType nominalDistance = -1.);
95 const std::vector<int>& trialElementIndices,
96 std::vector<arma::Mat<ResultType> >& result,
97 CoordinateType nominalDistance = -1.);
100 const std::vector<int>& pointIndices,
101 const std::vector<int>& trialElementIndices,
103 CoordinateType nominalDistance = -1.);
105 virtual int resultDimension()
const;
107 virtual CoordinateType estimateRelativeScale(CoordinateType minDist)
const;
114 BasisFunctionType> Utilities;
116 const Integrator& selectIntegrator(
117 int testElementIndex,
int trialElementIndex,
118 CoordinateType nominalDistance = -1.);
120 int order(
int pointIndex,
int trialElementIndex,
121 CoordinateType nominalDistance)
const;
125 CoordinateType pointElementDistanceSquared(
126 int pointIndex,
int trialElementIndex)
const;
128 void precalculateElementSizesAndCenters();
131 arma::Mat<CoordinateType> m_points;
132 shared_ptr<const GeometryFactory> m_geometryFactory;
133 shared_ptr<const RawGridGeometry<CoordinateType> > m_rawGeometry;
134 shared_ptr<const std::vector<const Shapeset<BasisFunctionType>*> > m_trialShapesets;
135 shared_ptr<const CollectionOfKernels<KernelType> > m_kernels;
136 shared_ptr<const CollectionOfShapesetTransformations<CoordinateType> > m_trialTransformations;
137 shared_ptr<const KernelTrialIntegral<BasisFunctionType, KernelType, ResultType> > m_integral;
139 VerbosityLevel::Level m_verbosityLevel;
140 shared_ptr<const QuadratureDescriptorSelectorForPotentialOperators<BasisFunctionType> >
142 shared_ptr<const SingleQuadratureRuleFamily<CoordinateType> > m_quadRuleFamily;
146 Integrator*> IntegratorMap;
147 IntegratorMap m_kernelTrialIntegrators;
149 enum { INVALID_INDEX = INT_MAX };
155 #include "default_local_assembler_for_potential_operators_on_surfaces_imp.hpp"
Traits of scalar types.
Definition: scalar_traits.hpp:40
virtual void evaluateLocalContributions(const std::vector< int > &pointIndices, int trialElementIndex, LocalDofIndex localTrialDofIndex, std::vector< arma::Mat< ResultType > > &result, CoordinateType nominalDistance=-1.)
Assemble local contributions.
Definition: default_local_assembler_for_potential_operators_on_surfaces_imp.hpp:94
Integration over pairs of elements.
Definition: kernel_trial_integrator.hpp:42
Parallel operation settings.
Definition: parallelization_options.hpp:32
Definition: default_local_assembler_for_potential_operators_on_surfaces.hpp:64
Abstract interface of a local assembler for potential operators.
Definition: local_assembler_for_potential_operators.hpp:54
Definition: default_local_assembler_for_operators_on_surfaces_utilities.hpp:42
Collection of shape functions defined on a reference element.
Definition: default_local_assembler_for_operators_on_surfaces_utilities.hpp:34
Definition: default_local_assembler_for_operators_on_surfaces_utilities.hpp:35
Family of quadrature rules over single elements.
Definition: numerical_quadrature_strategy.hpp:35
Definition: kernel_trial_integral.hpp:42
Simple implementation of a 2D Fortran-ordered array.
Definition: _2d_array.hpp:41
Parameters of a quadrature rule used in the evaluation of integrals over single elements.
Definition: single_quadrature_descriptor.hpp:34
Quadrature descriptor selector used during the evaluation on potentials.
Definition: quadrature_descriptor_selector_factory.hpp:39