Chi-Tech
angular_quadrature_base.h
Go to the documentation of this file.
1#ifndef CHI_MATH_ANGULAR_QUADRATURE_H
2#define CHI_MATH_ANGULAR_QUADRATURE_H
3
4#include <vector>
5
6#include "mesh/chi_mesh.h"
7
8namespace chi_math
9{
10struct QuadraturePointPhiTheta;
11
13{
14 Arbitrary = 1,
16 SLDFESQ = 3
17};
18class AngularQuadrature;
19class AngularQuadratureCustom;
20} // namespace chi_math
21
22/**Simple structure to add names to the angle components.*/
24{
25 double phi = 0.0;
26 double theta = 0.0;
27 QuadraturePointPhiTheta(const double phi, const double theta)
28 : phi(phi), theta(theta)
29 {
30 }
31};
32
33// ################################################################### Class def
34/**Base class for angular quadratures.*/
36{
37public:
39
40public:
41 std::vector<chi_math::QuadraturePointPhiTheta> abscissae_;
42 std::vector<double> weights_;
43 std::vector<chi_mesh::Vector3> omegas_;
44
46 {
47 unsigned int ell = 0;
48 int m = 0;
49
50 HarmonicIndices() = default;
51 HarmonicIndices(unsigned int in_ell, int in_m) : ell(in_ell), m(in_m) {}
52
53 bool operator==(const HarmonicIndices& other) const
54 {
55 return (ell == other.ell and m == other.m);
56 }
57 };
58
59protected:
60 std::vector<std::vector<double>> d2m_op_;
61 std::vector<std::vector<double>> m2d_op_;
62 std::vector<HarmonicIndices> m_to_ell_em_map_;
63 bool d2m_op_built_ = false;
64 bool m2d_op_built_ = false;
65
66public:
68
70 : type_(in_type)
71 {
72 }
73
74 virtual ~AngularQuadrature() = default;
75
76 virtual void OptimizeForPolarSymmetry(double normalization);
77
78protected:
79 virtual void MakeHarmonicIndices(unsigned int scattering_order,
80 int dimension);
81
82public:
83 virtual void BuildDiscreteToMomentOperator(unsigned int scattering_order,
84 int dimension);
85 virtual void BuildMomentToDiscreteOperator(unsigned int scattering_order,
86 int dimension);
87
88 std::vector<std::vector<double>> const& GetDiscreteToMomentOperator() const;
89
90 std::vector<std::vector<double>> const& GetMomentToDiscreteOperator() const;
91
92 const std::vector<HarmonicIndices>& GetMomentToHarmonicsIndexMap() const;
93};
94
96{
97public:
98 AngularQuadratureCustom(std::vector<double>& azimuthal,
99 std::vector<double>& polar,
100 std::vector<double>& in_weights,
101 bool verbose);
102};
103
104#endif // CHI_MATH_ANGULAR_QUADRATURE_H
AngularQuadratureCustom(std::vector< double > &azimuthal, std::vector< double > &polar, std::vector< double > &in_weights, bool verbose)
std::vector< std::vector< double > > const & GetMomentToDiscreteOperator() const
const std::vector< HarmonicIndices > & GetMomentToHarmonicsIndexMap() const
virtual ~AngularQuadrature()=default
const chi_math::AngularQuadratureType type_
virtual void MakeHarmonicIndices(unsigned int scattering_order, int dimension)
std::vector< std::vector< double > > m2d_op_
std::vector< chi_math::QuadraturePointPhiTheta > abscissae_
virtual void OptimizeForPolarSymmetry(double normalization)
std::vector< std::vector< double > > const & GetDiscreteToMomentOperator() const
std::vector< std::vector< double > > d2m_op_
virtual void BuildDiscreteToMomentOperator(unsigned int scattering_order, int dimension)
std::vector< HarmonicIndices > m_to_ell_em_map_
virtual void BuildMomentToDiscreteOperator(unsigned int scattering_order, int dimension)
AngularQuadrature(chi_math::AngularQuadratureType in_type)
std::vector< chi_mesh::Vector3 > omegas_
bool operator==(const HarmonicIndices &other) const
QuadraturePointPhiTheta(const double phi, const double theta)