BEM++  2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
default_local_assembler_for_potential_operators_on_surfaces.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_default_local_assembler_for_potential_operators_on_surfaces_hpp
22 #define fiber_default_local_assembler_for_potential_operators_on_surfaces_hpp
23 
24 #include "../common/common.hpp"
25 
26 #include "local_assembler_for_potential_operators.hpp"
27 
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"
32 #include "numerical_quadrature.hpp"
33 #include "parallelization_options.hpp"
34 #include "shared_ptr.hpp"
35 #include "kernel_trial_integrator.hpp"
36 #include "verbosity_level.hpp"
37 
38 #include <boost/static_assert.hpp>
39 #include <boost/tuple/tuple_comparison.hpp>
40 #include <tbb/concurrent_unordered_map.h>
41 #include <cstring>
42 #include <climits>
43 #include <set>
44 #include <utility>
45 #include <vector>
46 
47 namespace Fiber
48 {
49 
51 class OpenClHandler;
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>
65  public LocalAssemblerForPotentialOperators<ResultType>
66 {
67 public:
68  typedef typename ScalarTraits<ResultType>::RealType CoordinateType;
69 
71  const arma::Mat<CoordinateType>& points,
72  const shared_ptr<const GeometryFactory>& geometryFactory,
73  const shared_ptr<const RawGridGeometry<CoordinateType> >& rawGeometry,
74  const shared_ptr<const std::vector<const Shapeset<BasisFunctionType>*> >& trialShapesets,
75  const shared_ptr<const CollectionOfKernels<KernelType> >& kernels,
76  const shared_ptr<const CollectionOfShapesetTransformations<CoordinateType> >& trialTransformations,
77  const shared_ptr<const KernelTrialIntegral<BasisFunctionType, KernelType, ResultType> >& integral,
78  const ParallelizationOptions& parallelizationOptions,
79  VerbosityLevel::Level verbosityLevel,
81  BasisFunctionType> >& quadDescSelector,
82  const shared_ptr<const SingleQuadratureRuleFamily<CoordinateType> >& quadRuleFamily);
84 
85  virtual void evaluateLocalContributions(
86  const std::vector<int>& pointIndices,
87  int trialElementIndex,
88  LocalDofIndex localTrialDofIndex,
89  std::vector<arma::Mat<ResultType> >& result,
90  CoordinateType nominalDistance = -1.);
91 
92  virtual void evaluateLocalContributions(
93  int pointIndex,
94  int componentIndex,
95  const std::vector<int>& trialElementIndices,
96  std::vector<arma::Mat<ResultType> >& result,
97  CoordinateType nominalDistance = -1.);
98 
99  virtual void evaluateLocalContributions(
100  const std::vector<int>& pointIndices,
101  const std::vector<int>& trialElementIndices,
102  Fiber::_2dArray<arma::Mat<ResultType> >& result,
103  CoordinateType nominalDistance = -1.);
104 
105  virtual int resultDimension() const;
106 
107  virtual CoordinateType estimateRelativeScale(CoordinateType minDist) const;
108 
109 private:
112  Integrator;
114  BasisFunctionType> Utilities;
115 
116  const Integrator& selectIntegrator(
117  int testElementIndex, int trialElementIndex,
118  CoordinateType nominalDistance = -1.);
119 
120  int order(int pointIndex, int trialElementIndex,
121  CoordinateType nominalDistance) const;
122 
123  const Integrator& getIntegrator(const SingleQuadratureDescriptor& index);
124 
125  CoordinateType pointElementDistanceSquared(
126  int pointIndex, int trialElementIndex) const;
127 
128  void precalculateElementSizesAndCenters();
129 
130 private:
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;
138  ParallelizationOptions m_parallelizationOptions;
139  VerbosityLevel::Level m_verbosityLevel;
140  shared_ptr<const QuadratureDescriptorSelectorForPotentialOperators<BasisFunctionType> >
141  m_quadDescSelector;
142  shared_ptr<const SingleQuadratureRuleFamily<CoordinateType> > m_quadRuleFamily;
143 
144 
145  typedef tbb::concurrent_unordered_map<SingleQuadratureDescriptor,
146  Integrator*> IntegratorMap;
147  IntegratorMap m_kernelTrialIntegrators;
148 
149  enum { INVALID_INDEX = INT_MAX };
151 };
152 
153 } // namespace Fiber
154 
155 #include "default_local_assembler_for_potential_operators_on_surfaces_imp.hpp"
156 
157 #endif
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