BEM++  2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Groups Pages
default_local_assembler_for_grid_functions_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_grid_functions_on_surfaces_hpp
22 #define fiber_default_local_assembler_for_grid_functions_on_surfaces_hpp
23 
24 #include "../common/common.hpp"
25 
26 #include "local_assembler_for_grid_functions.hpp"
27 #include "numerical_quadrature.hpp"
28 #include "quadrature_options.hpp"
29 #include "test_function_integrator.hpp"
30 #include "scalar_traits.hpp"
31 
32 #include <tbb/concurrent_unordered_map.h>
33 #include <map>
34 #include <vector>
35 
36 namespace Fiber
37 {
38 
40 class OpenClHandler;
41 template <typename CoordinateType> class CollectionOfShapesetTransformations;
42 template <typename ValueType> class Function;
43 template <typename CoordinateType> class RawGridGeometry;
44 
45 template <typename CoordinateType> class QuadratureDescriptorSelectorForGridFunctions;
46 template <typename CoordinateType> class SingleQuadratureRuleFamily;
49 template <typename BasisFunctionType, typename UserFunctionType,
50  typename ResultType, typename GeometryFactory>
52  public LocalAssemblerForGridFunctions<ResultType>
53 {
54 public:
55  typedef typename ScalarTraits<ResultType>::RealType CoordinateType;
56 
58  const shared_ptr<const GeometryFactory>& geometryFactory,
59  const shared_ptr<const RawGridGeometry<CoordinateType> >& rawGeometry,
60  const shared_ptr<const std::vector<
61  const Shapeset<BasisFunctionType>*> >& testShapesets,
62  const shared_ptr<const CollectionOfShapesetTransformations<
63  CoordinateType> >& testTransformations,
64  const shared_ptr<const Function<UserFunctionType> >& function,
65  const shared_ptr<const OpenClHandler>& openClHandler,
66  const shared_ptr<const QuadratureDescriptorSelectorForGridFunctions<
67  CoordinateType> >& quadDescSelector,
68  const shared_ptr<const SingleQuadratureRuleFamily<
69  CoordinateType> >& quadRuleFamily);
71 
72 public:
73  virtual void evaluateLocalWeakForms(
74  const std::vector<int>& elementIndices,
75  std::vector<arma::Col<ResultType> >& result);
76 
77 private:
79 
80  const Integrator& selectIntegrator(int elementIndex);
81 
82  const Integrator& getIntegrator(const SingleQuadratureDescriptor& index);
83 
84 private:
85  typedef tbb::concurrent_unordered_map<SingleQuadratureDescriptor,
86  Integrator*> IntegratorMap;
87 
88 private:
89  shared_ptr<const GeometryFactory> m_geometryFactory;
90  shared_ptr<const RawGridGeometry<CoordinateType> > m_rawGeometry;
91  shared_ptr<const std::vector<const Shapeset<BasisFunctionType>*> > m_testShapesets;
92  shared_ptr<const CollectionOfShapesetTransformations<CoordinateType> > m_testTransformations;
93  shared_ptr<const Function<UserFunctionType> > m_function;
94  shared_ptr<const OpenClHandler> m_openClHandler;
95  shared_ptr<const QuadratureDescriptorSelectorForGridFunctions<CoordinateType> > m_quadDescSelector;
96  shared_ptr<const SingleQuadratureRuleFamily<CoordinateType> > m_quadRuleFamily;
97 
98  IntegratorMap m_testFunctionIntegrators;
99 };
100 
101 } // namespace Fiber
102 
103 #include "default_local_assembler_for_grid_functions_on_surfaces_imp.hpp"
104 
105 #endif
Traits of scalar types.
Definition: scalar_traits.hpp:40
Integration of products of test functions and arbitrary functions over elements.
Definition: test_function_integrator.hpp:41
Definition: default_local_assembler_for_grid_functions_on_surfaces.hpp:51
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 functions.
Definition: quadrature_descriptor_selector_factory.hpp:36
Definition: default_local_assembler_for_operators_on_surfaces_utilities.hpp:35
Family of quadrature rules over single elements.
Definition: numerical_quadrature_strategy.hpp:35
Collection of shape function transformations.
Definition: collection_of_shapeset_transformations.hpp:64
Local assembler for grid functions.
Definition: local_assembler_for_grid_functions.hpp:41
virtual void evaluateLocalWeakForms(const std::vector< int > &elementIndices, std::vector< arma::Col< ResultType > > &result)
Assemble local weak forms of a source term on specified elements.
Definition: default_local_assembler_for_grid_functions_on_surfaces_imp.hpp:80
Parameters of a quadrature rule used in the evaluation of integrals over single elements.
Definition: single_quadrature_descriptor.hpp:34