BEM++  2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
default_local_assembler_for_integral_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_integral_operators_on_surfaces_hpp
22 #define fiber_default_local_assembler_for_integral_operators_on_surfaces_hpp
23 
24 #include "../common/common.hpp"
25 
26 #include "local_assembler_for_integral_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 "test_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 <tbb/mutex.h>
42 #include <cstring>
43 #include <climits>
44 #include <set>
45 #include <utility>
46 #include <vector>
47 
48 namespace Fiber
49 {
50 
52 class OpenClHandler;
53 template <typename CoordinateType> class CollectionOfShapesetTransformations;
54 template <typename ValueType> class CollectionOfKernels;
55 template <typename BasisFunctionType, typename KernelType, typename ResultType>
56 class TestKernelTrialIntegral;
57 
58 template <typename CoordinateType> class RawGridGeometry;
59 
60 template <typename CoordinateType>
61 class QuadratureDescriptorSelectorForIntegralOperators;
62 template <typename CoordinateType> class DoubleQuadratureRuleFamily;
65 template <typename BasisFunctionType, typename KernelType,
66  typename ResultType, typename GeometryFactory>
68  public LocalAssemblerForIntegralOperators<ResultType>
69 {
70 public:
71  typedef typename ScalarTraits<ResultType>::RealType CoordinateType;
72 
74  const shared_ptr<const GeometryFactory>& testGeometryFactory,
75  const shared_ptr<const GeometryFactory>& trialGeometryFactory,
76  const shared_ptr<const RawGridGeometry<CoordinateType> >& testRawGeometry,
77  const shared_ptr<const RawGridGeometry<CoordinateType> >& trialRawGeometry,
78  const shared_ptr<const std::vector<const Shapeset<BasisFunctionType>*> >& testShapesets,
79  const shared_ptr<const std::vector<const Shapeset<BasisFunctionType>*> >& trialShapesets,
80  const shared_ptr<const CollectionOfShapesetTransformations<CoordinateType> >& testTransformations,
81  const shared_ptr<const CollectionOfKernels<KernelType> >& kernel,
82  const shared_ptr<const CollectionOfShapesetTransformations<CoordinateType> >& trialTransformations,
84  const shared_ptr<const OpenClHandler>& openClHandler,
85  const ParallelizationOptions& parallelizationOptions,
86  VerbosityLevel::Level verbosityLevel,
87  bool cacheSingularIntegrals,
88  const shared_ptr<const QuadratureDescriptorSelectorForIntegralOperators<CoordinateType> >& quadDescSelector,
89  const shared_ptr<const DoubleQuadratureRuleFamily<CoordinateType> >& quadRuleFamily);
91 
92 public:
93  virtual void evaluateLocalWeakForms(
94  CallVariant callVariant,
95  const std::vector<int>& elementIndicesA,
96  int elementIndexB,
97  LocalDofIndex localDofIndexB,
98  std::vector<arma::Mat<ResultType> >& result,
99  CoordinateType nominalDistance = -1.);
100 
101  virtual void evaluateLocalWeakForms(
102  const std::vector<int>& testElementIndices,
103  const std::vector<int>& trialElementIndices,
104  Fiber::_2dArray<arma::Mat<ResultType> >& result,
105  CoordinateType nominalDistance = -1.);
106 
107  virtual CoordinateType estimateRelativeScale(CoordinateType minDist) const;
108 
109 private:
112  typedef typename Integrator::ElementIndexPair ElementIndexPair;
114  BasisFunctionType> Utilities;
115 
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);
125  }
126  };
127 
133  typedef alternative_less<
134  typename ElementIndexPair::first_type,
135  typename ElementIndexPair::second_type> ElementIndexPairCompare;
136 
145  typedef std::set<ElementIndexPair, ElementIndexPairCompare> ElementIndexPairSet;
146 
147  bool testAndTrialGridsAreIdentical() const;
148 
149  void cacheSingularLocalWeakForms();
150  void findPairsOfAdjacentElements(ElementIndexPairSet& pairs) const;
151  void cacheLocalWeakForms(const ElementIndexPairSet& elementIndexPairs);
152 
153  const Integrator& selectIntegrator(
154  int testElementIndex, int trialElementIndex,
155  CoordinateType nominalDistance = -1.);
156 
157  const Integrator& getIntegrator(const DoubleQuadratureDescriptor& index);
158 
159 private:
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;
171  ParallelizationOptions m_parallelizationOptions;
172  VerbosityLevel::Level m_verbosityLevel;
173  shared_ptr<const QuadratureDescriptorSelectorForIntegralOperators<CoordinateType> > m_quadDescSelector;
174  shared_ptr<const DoubleQuadratureRuleFamily<CoordinateType> > m_quadRuleFamily;
175 
176  typedef tbb::concurrent_unordered_map<DoubleQuadratureDescriptor,
177  Integrator*> IntegratorMap;
178  IntegratorMap m_testKernelTrialIntegrators;
179  mutable tbb::mutex m_integratorCreationMutex;
180 
181  enum { INVALID_INDEX = INT_MAX };
193  Cache m_cache;
195 };
196 
197 } // namespace Fiber
198 
199 #include "default_local_assembler_for_integral_operators_on_surfaces_imp.hpp"
200 
201 #endif
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