Chi-Tech
angular_product_quadrature.h
Go to the documentation of this file.
1#ifndef _product_quadrature_h
2#define _product_quadrature_h
3
4#include <map>
5#include <vector>
6
8
9namespace chi_math
10{
12 {
13 UNKNOWN = 0,
19 };
20
21 class ProductQuadrature;
22 class AngularQuadratureProdGL;
23 class AngularQuadratureProdGLL;
24 class AngularQuadratureProdGLC;
25 class AngularQuadratureProdCustom;
26}
27
28//######################################################### Class def
29/** Class for product quadratures*/
31{
32public:
33 std::vector<double> polar_ang_;
34 std::vector<double> azimu_ang_;
35protected:
36 /** Linear indices of ordered directions mapped to polar level. */
37 std::map<unsigned int, std::vector<unsigned int>> map_directions_;
38
39protected:
42 {}
43
44public:
45 ~ProductQuadrature() override = default;
46
47 void AssembleCosines(const std::vector<double>& azimuthal,
48 const std::vector<double>& polar,
49 const std::vector<double>& in_weights,
50 bool verbose);
51
52 void OptimizeForPolarSymmetry(double normalization) override;
53 /**Obtains the abscissae index given the indices of the
54 * polar angle index and the azimuthal angle index.*/
55 unsigned int GetAngleNum(const unsigned int polar_angle_index,
56 const unsigned int azimu_angle_index) const
57 { return map_directions_.at(polar_angle_index)[azimu_angle_index]; }
58 /** Return constant reference to map_directions. */
59 const std::map<unsigned int,
60 std::vector<unsigned int>>& GetDirectionMap() const
61 { return map_directions_; }
62};
63
64//######################################################### Class def
66{
67public:
68 explicit AngularQuadratureProdGL(int Np, bool verbose=false);
69};
70
71//######################################################### Class def
73{
74public:
75 explicit
76 AngularQuadratureProdGLL(int Na, int Np, bool verbose=false);
77};
78
79//######################################################### Class def
81{
82public:
83 explicit
84 AngularQuadratureProdGLC(int Na, int Np, bool verbose=false);
85};
86
87//######################################################### Class def
89{
90public:
91 AngularQuadratureProdCustom(const std::vector<double>& azimuthal,
92 const std::vector<double>& polar,
93 const std::vector<double>& in_weights, bool verbose);
94};
95
96#endif
AngularQuadratureProdCustom(const std::vector< double > &azimuthal, const std::vector< double > &polar, const std::vector< double > &in_weights, bool verbose)
AngularQuadratureProdGLC(int Na, int Np, bool verbose=false)
AngularQuadratureProdGL(int Np, bool verbose=false)
AngularQuadratureProdGLL(int Na, int Np, bool verbose=false)
void AssembleCosines(const std::vector< double > &azimuthal, const std::vector< double > &polar, const std::vector< double > &in_weights, bool verbose)
~ProductQuadrature() override=default
const std::map< unsigned int, std::vector< unsigned int > > & GetDirectionMap() const
unsigned int GetAngleNum(const unsigned int polar_angle_index, const unsigned int azimu_angle_index) const
void OptimizeForPolarSymmetry(double normalization) override
std::map< unsigned int, std::vector< unsigned int > > map_directions_