Geophysical Inversion and Modelling Library v1.5.4
Loading...
Searching...
No Matches
integration.h
1/******************************************************************************
2 * Copyright (C) 2006-2024 by the GIMLi development team *
3 * Carsten Rücker carsten@resistivity.net *
4 * *
5 * Licensed under the Apache License, Version 2.0 (the "License"); *
6 * you may not use this file except in compliance with the License. *
7 * You may obtain a copy of the License at *
8 * *
9 * http://www.apache.org/licenses/LICENSE-2.0 *
10 * *
11 * Unless required by applicable law or agreed to in writing, software *
12 * distributed under the License is distributed on an "AS IS" BASIS, *
13 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. *
14 * See the License for the specific language governing permissions and *
15 * limitations under the License. *
16 * *
17 ******************************************************************************/
18
19#ifndef _GIMLI_INTEGRATION__H
20#define _GIMLI_INTEGRATION__H
21
22#include "gimli.h"
23
24namespace GIMLI{
25
27class DLLEXPORT IntegrationRules : public Singleton< IntegrationRules >{
28public:
29 friend class Singleton< IntegrationRules >;
30
32 const R3Vector & gauAbscissa(Index order) const;
33
35 const RVector & gauWeights(Index order) const;
36
43 const R3Vector & triGLAbscissa(Index order) const;
44
51 const RVector & triGLWeights(Index order) const;
52
53 const R3Vector & edgAbscissa(Index order) const;
54 const RVector & edgWeights(Index order) const;
55
56 const R3Vector & triAbscissa(Index order) const;
57 const RVector & triWeights(Index order) const;
58
59 const R3Vector & tetAbscissa(Index order) const;
60 const RVector & tetWeights(Index order) const;
61
62 const R3Vector & quaAbscissa(Index order) const;
63 const RVector & quaWeights(Index order) const;
64
65 const R3Vector & hexAbscissa(Index order) const;
66 const RVector & hexWeights(Index order) const;
67
68 const R3Vector & priAbscissa(Index order) const;
69 const RVector & priWeights(Index order) const;
73 const R3Vector & abscissa(const Shape & shape, uint order) const;
77 const RVector & weights(const Shape & shape, uint order) const;
78
81 inline void setTriUseGaussLegendre(bool use){ triUseGaussLegendre_ = use;}
82 inline bool triUseGaussLegendre() const { return triUseGaussLegendre_; }
83
84private:
86 IntegrationRules();
87
89 ~IntegrationRules();
90
91protected:
92
93 void initGau_();
94 void initTriGL_();
95
96 void initEdg_();
97 void initTri_();
98 void initTet_();
99 void initQua_();
100 void initHex_();
101 void initPri_();
102
103 bool triUseGaussLegendre_;
104
105 std::vector < R3Vector > gauAbscissa_;
106 std::vector < RVector > gauWeights_;
107
109 std::vector < R3Vector > triGLAbscissa_;
111 std::vector < RVector > triGLWeights_;
112
113 std::vector < R3Vector > edgAbscissa_;
114 std::vector < RVector > edgWeights_;
115
116 std::vector < R3Vector > triAbscissa_;
117 std::vector < RVector > triWeights_;
118
119 std::vector < R3Vector > tetAbscissa_;
120 std::vector < RVector > tetWeights_;
121
122 std::vector < R3Vector > quaAbscissa_;
123 std::vector < RVector > quaWeights_;
124
125 std::vector < R3Vector > hexAbscissa_;
126 std::vector < RVector > hexWeights_;
127
128 std::vector < R3Vector > priAbscissa_;
129 std::vector < RVector > priWeights_;
130
131};
132
133} // namespace GIMLI{
134
135#endif // _GIMLI_INTEGRATION__H
const RVector & weights(const Shape &shape, uint order) const
Definition integration.cpp:138
std::vector< R3Vector > triGLAbscissa_
Definition integration.h:109
const R3Vector & gauAbscissa(Index order) const
Definition integration.cpp:47
std::vector< RVector > triGLWeights_
Definition integration.h:111
const R3Vector & triGLAbscissa(Index order) const
Definition integration.cpp:57
const RVector & gauWeights(Index order) const
Definition integration.cpp:52
void setTriUseGaussLegendre(bool use)
Definition integration.h:81
const RVector & triGLWeights(Index order) const
Definition integration.cpp:62
const R3Vector & abscissa(const Shape &shape, uint order) const
Definition integration.cpp:116
A Shape defines a geometrical primitive.
Definition shape.h:181
Singleton()
Definition gimli.h:634
GIMLi main namespace for the Geophyiscal Inversion and Modelling Library.
Definition baseentity.h:24