BEM++  2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
default_local_assembler_for_local_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_local_operators_on_surfaces_hpp
22 #define fiber_default_local_assembler_for_local_operators_on_surfaces_hpp
23 
24 #include "../common/common.hpp"
25 
26 #include "local_assembler_for_local_operators.hpp"
27 
28 #include "numerical_quadrature.hpp"
29 #include "numerical_test_trial_integrator.hpp"
30 #include "shared_ptr.hpp"
31 
32 #include "../common/armadillo_fwd.hpp"
33 #include "../common/boost_ptr_map_fwd.hpp"
34 #include <cstring>
35 #include <iostream>
36 #include <map>
37 #include <set>
38 #include <utility>
39 #include <vector>
40 
41 namespace Fiber
42 {
43 
45 class OpenClHandler;
46 
47 template <typename CoordinateType>
48 class QuadratureDescriptorSelectorForLocalOperators;
49 template <typename CoordinateType> class SingleQuadratureRuleFamily;
52 template <typename BasisFunctionType, typename ResultType, typename GeometryFactory>
54  public LocalAssemblerForLocalOperators<ResultType>
55 {
56 public:
57  typedef typename ScalarTraits<ResultType>::RealType CoordinateType;
58 
60  const shared_ptr<const GeometryFactory>& geometryFactory,
61  const shared_ptr<const RawGridGeometry<CoordinateType> >& rawGeometry,
62  const shared_ptr<const std::vector<const Shapeset<BasisFunctionType>*> >& testShapesets,
63  const shared_ptr<const std::vector<const Shapeset<BasisFunctionType>*> >& trialShapesets,
64  const shared_ptr<const CollectionOfShapesetTransformations<CoordinateType> >& testTransformations,
65  const shared_ptr<const CollectionOfShapesetTransformations<CoordinateType> >& trialTransformations,
66  const shared_ptr<const TestTrialIntegral<BasisFunctionType, ResultType> >& integral,
67  const shared_ptr<const OpenClHandler>& openClHandler,
68  const shared_ptr<const QuadratureDescriptorSelectorForLocalOperators<CoordinateType> >& quadDescSelector,
69  const shared_ptr<const SingleQuadratureRuleFamily<CoordinateType> >& quadRuleFamily);
70 
71  virtual void evaluateLocalWeakForms(
72  const std::vector<int>& elementIndices,
73  std::vector<arma::Mat<ResultType> >& result);
74 
75 private:
78  selectIntegrator(int elementIndex);
79 
81  const SingleQuadratureDescriptor& desc);
82 
83 private:
84  typedef boost::ptr_map<SingleQuadratureDescriptor,
87  BasisFunctionType> Utilities;
88 
89 private:
90  shared_ptr<const GeometryFactory> m_geometryFactory;
91  shared_ptr<const RawGridGeometry<CoordinateType> > m_rawGeometry;
92  shared_ptr<const std::vector<const Shapeset<BasisFunctionType>*> > m_testShapesets;
93  shared_ptr<const std::vector<const Shapeset<BasisFunctionType>*> > m_trialShapesets;
94  shared_ptr<const CollectionOfShapesetTransformations<CoordinateType> > m_testTransformations;
95  shared_ptr<const CollectionOfShapesetTransformations<CoordinateType> > m_trialTransformations;
96  shared_ptr<const TestTrialIntegral<BasisFunctionType, ResultType> > m_integral;
97  shared_ptr<const OpenClHandler> m_openClHandler;
98  shared_ptr<const QuadratureDescriptorSelectorForLocalOperators<CoordinateType> > m_quadDescSelector;
99  shared_ptr<const SingleQuadratureRuleFamily<CoordinateType> > m_quadRuleFamily;
100 
101  IntegratorMap m_testTrialIntegrators;
103 };
104 
105 } // namespace Fiber
106 
107 #include "default_local_assembler_for_local_operators_on_surfaces_imp.hpp"
108 
109 #endif
Traits of scalar types.
Definition: scalar_traits.hpp:40
Abstract interface of a local assembler for local operators.
Definition: local_assembler_for_local_operators.hpp:48
Integration of products of test and trial functions over elements.
Definition: test_trial_integrator.hpp:40
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
Quadrature descriptor selector used during the discretization of local boundary operators.
Definition: quadrature_descriptor_selector_factory.hpp:38
Definition: default_local_assembler_for_local_operators_on_surfaces.hpp:53
Definition: default_local_assembler_for_operators_on_surfaces_utilities.hpp:35
Family of quadrature rules over single elements.
Definition: numerical_quadrature_strategy.hpp:35
virtual void evaluateLocalWeakForms(const std::vector< int > &elementIndices, std::vector< arma::Mat< ResultType > > &result)
Assemble local weak forms.
Definition: default_local_assembler_for_local_operators_on_surfaces_imp.hpp:65
Parameters of a quadrature rule used in the evaluation of integrals over single elements.
Definition: single_quadrature_descriptor.hpp:34