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.

std::vector<std::vector<double>> const &GetDiscreteToMomentOperator() const

Return a reference to the precomputed discrete-to-moment operator.

The operator is accessed as [d][m], where m is the moment index and d is the direction index.

std::vector<std::vector<double>> const &GetMomentToDiscreteOperator() const

Return a reference to the precomputed moment-to-discrete operator.

The operator is accessed as [d][m], where m is the moment index and d is the direction 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

Public Members

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.

Protected Functions

inline explicit AngularQuadrature(AngularQuadratureType type, unsigned int dimension, unsigned int scattering_order)
void MakeHarmonicIndices()

Populate the map of moment index to spherical harmonic indices.

Protected Attributes

std::vector<std::vector<double>> d2m_op_

Discrete-to-moment operator matrix.

std::vector<std::vector<double>> m2d_op_

Moment-to-discrete operator matrix.

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.

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.