Class AngularQuadrature

Nested Relationships

Nested Types

Inheritance Relationships

Derived Types

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 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 size_t GetNumAngles() const
inline bool Empty() const
inline AngularQuadratureType GetType() const
virtual std::string GetName() const = 0

Return a user-facing quadrature name.

virtual double GetWeightSum() const

Return the sum of all quadrature weights.

inline const std::vector<QuadraturePointPhiTheta> &GetAbscissae() const
inline const QuadraturePointPhiTheta &GetAbscissa(size_t angle_index) const
inline const std::vector<double> &GetWeights() const
inline double GetWeight(size_t angle_index) const
inline const std::vector<Vector3> &GetOmegas() const
inline const Vector3 &GetOmega(size_t angle_index) const

Protected Functions

inline explicit AngularQuadrature(AngularQuadratureType type, unsigned int dimension, unsigned int scattering_order, OperatorConstructionMethod method = OperatorConstructionMethod::STANDARD)
inline void SetOperatorConstructionMethod(OperatorConstructionMethod method)

Sets the operator construction method.

void MakeHarmonicIndices()

Populate the map of moment index to spherical harmonic indices.

inline virtual unsigned int GetQuadratureOrder() const

Quadrature-specific order parameter used by harmonic-selection rules.

inline virtual unsigned int GetNumPolarAngles() const

Number of polar angles used by harmonic-selection rules.

inline virtual unsigned int GetNumAzimuthalAngles() const

Number of azimuthal angles used by harmonic-selection rules.

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_
std::vector<QuadraturePointPhiTheta> abscissae_ = {}

Quadrature point abscissae in spherical coordinates.

std::vector<double> weights_ = {}

Quadrature weights.

std::vector<Vector3> omegas_ = {}

Quadrature point direction vectors in Cartesian coordinates.

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

Public Members

unsigned int ell = 0

Degree of the spherical harmonic.

int m = 0

Order of the spherical harmonic.