Class AngularQuadrature
Defined in File angular_quadrature.h
Nested Relationships
Nested Types
Inheritance Relationships
Derived Types
public opensn::LebedevQuadrature2DXY(Class LebedevQuadrature2DXY)public opensn::LebedevQuadrature3DXYZ(Class LebedevQuadrature3DXYZ)public opensn::ProductQuadrature(Class ProductQuadrature)public opensn::SLDFEsqQuadrature(Class SLDFEsqQuadrature)public opensn::TriangularQuadrature(Class TriangularQuadrature)
Class Documentation
-
class AngularQuadrature
Base class for angular quadratures used in discrete ordinates transport calculations.
Subclassed by opensn::LebedevQuadrature2DXY, opensn::LebedevQuadrature3DXYZ, opensn::ProductQuadrature, opensn::SLDFEsqQuadrature, opensn::TriangularQuadrature
Public Functions
-
virtual ~AngularQuadrature() = default
-
void BuildDiscreteToMomentOperator()
Compute the discrete-to-moment operator.
-
void BuildMomentToDiscreteOperator()
Compute the moment-to-discrete operator.
-
inline void SetOperatorConstructionMethod(OperatorConstructionMethod method)
Sets the operator construction method.
-
inline OperatorConstructionMethod GetOperatorConstructionMethod() const
Gets the current operator construction method.
-
NDArray<double, 2> const &GetDiscreteToMomentOperator() const
Return a reference to the precomputed discrete-to-moment operator.
The operator is accessed as (d,m), where d is the direction index and m is the moment index.
-
NDArray<double, 2> const &GetMomentToDiscreteOperator() const
Return a reference to the precomputed moment-to-discrete operator.
The operator is accessed as (d,m), where d is the direction index and m is the moment index.
-
const std::vector<HarmonicIndices> &GetMomentToHarmonicsIndexMap() const
Return a reference to the precomputed harmonic index map.
-
inline unsigned int GetDimension() const
-
inline unsigned int GetScatteringOrder() const
-
inline unsigned int GetNumMoments() const
-
inline AngularQuadratureType GetType() const
-
inline void SetQuadratureOrder(unsigned int order)
Set the quadrature-specific order parameter.
For Lebedev: Lebedev order {3, 5, 7, …}; for SLDFESQ: uniform refinement level.
-
inline void SetNumberOfPolar(unsigned int num_polar)
Set the number of positive-mu polar points (product quadrature types).
-
inline void SetNumberOfAzimuthal(unsigned int num_azimu)
Set the total number of azimuthal points (product quadrature types).
Public Members
-
std::vector<QuadraturePointPhiTheta> abscissae = {}
Quadrature point abscissae in spherical coordinates.
-
std::vector<double> weights = {}
Quadrature weights.
Protected Functions
-
inline explicit AngularQuadrature(AngularQuadratureType type, unsigned int dimension, unsigned int scattering_order, OperatorConstructionMethod method = OperatorConstructionMethod::STANDARD)
-
void MakeHarmonicIndices()
Populate the map of moment index to spherical harmonic indices.
Protected Attributes
-
NDArray<double, 2> d2m_op_ = {}
Discrete-to-moment operator matrix.
OpenSn storage is indexed [angle][moment].
-
NDArray<double, 2> m2d_op_ = {}
Moment-to-discrete operator matrix.
OpenSn storage is indexed [angle][moment].
-
std::vector<HarmonicIndices> m_to_ell_em_map_ = {}
Mapping from moment index to spherical harmonic indices.
-
AngularQuadratureType type_
Quadrature type identifier.
-
unsigned int dimension_
Spatial dimension of the quadrature.
-
unsigned int scattering_order_
Maximum scattering order for moment calculations.
-
OperatorConstructionMethod construction_method_
-
unsigned int quadrature_order_ = 0
-
unsigned int n_polar_ = 0
-
unsigned int n_azimuthal_ = 0
-
struct HarmonicIndices
Spherical harmonic indices for moment-to-direction mappings.
Public Functions
-
HarmonicIndices() = default
-
inline HarmonicIndices(unsigned int ell, int m)
-
inline bool operator==(const HarmonicIndices &other) const
-
HarmonicIndices() = default
-
virtual ~AngularQuadrature() = default