https://mooseframework.inl.gov
RayTracingAngularQuadrature.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 // MOOSE includes
13 #include "MooseTypes.h"
14 
15 // libMesh includes
16 #include "libmesh/dense_matrix.h"
17 #include "libmesh/point.h"
18 
20 {
21 public:
31  RayTracingAngularQuadrature(const unsigned int dim,
32  const unsigned int polar_order,
33  const unsigned int azimuthal_order,
34  const Real mu_min,
35  const Real mu_max);
36 
45  static void gaussLegendre(const unsigned int order, std::vector<Real> & x, std::vector<Real> & w);
46 
54  static void chebyshev(const unsigned int order, std::vector<Real> & x, std::vector<Real> & w);
55 
59  static void rotationMatrix(const libMesh::Point & direction, libMesh::DenseMatrix<Real> & matrix);
64 
68  void rotate(const libMesh::Point & rotation_direction);
69 
76  const libMesh::Point & getDirection(const unsigned int l) const;
87  const std::vector<Real> & getWeights(const unsigned int l) const;
99  Real getTotalWeight(const unsigned int l) const;
108  const std::vector<Real> & getPolarSins(const unsigned int l) const;
109 
117  std::size_t numDirections() const { return _current_directions.size(); }
121  bool hasDirection(const unsigned int l) const { return l < _current_directions.size(); }
122 
126  void checkDirection(const unsigned int l) const;
127 
137  std::size_t numPolar(const unsigned int l) const;
138 
142  unsigned int polarOrder() const { return _polar_order; }
146  unsigned int azimuthalOrder() const { return _azimuthal_order; }
150  Real muMin() const { return _mu_min; }
154  Real muMax() const { return _mu_max; }
155 
160 
164  unsigned int dim() const { return _dim; }
165 
166 private:
170  void build();
171 
173  const unsigned int _dim;
175  const unsigned int _polar_order;
177  const unsigned int _azimuthal_order;
179  const Real _mu_min;
181  const Real _mu_max;
182 
184  std::vector<Real> _phi;
186  std::vector<Real> _mu;
188  std::vector<Real> _polar_sin;
190  std::vector<Real> _w;
191 
201  std::vector<libMesh::Point> _current_directions;
203  std::vector<std::vector<Real>> _current_weights;
204  std::vector<std::vector<Real>> _current_polar_sins;
206 
209 };
static void gaussLegendre(const unsigned int order, std::vector< Real > &x, std::vector< Real > &w)
Builds Gauss-Legendre quadrature on 0, 1, with weights that sum to 1.
Real muMax() const
Get the maximum mu.
std::vector< std::vector< Real > > _current_polar_sins
const unsigned int _polar_order
The polar order.
static void chebyshev(const unsigned int order, std::vector< Real > &x, std::vector< Real > &w)
Builds Chebyshev quadrature on [0, 2] with weights that sum to 2.
std::vector< Real > _phi
Quadrature phi.
std::vector< Real > _mu
Quadrature mu.
const libMesh::Point & getDirection(const unsigned int l) const
Get the direction associated with direction l.
std::size_t numDirections() const
Get the number of directions in the rotated and projected quadrature.
const Real _mu_min
The minimum mu.
std::vector< Real > _w
Quadrature combined weights.
void build()
Build the quadrature.
const std::vector< double > x
std::vector< libMesh::Point > _current_directions
The current quadrature information.
libMesh::Point _current_rotation_direction
The current rotation direction.
std::vector< std::vector< Real > > _current_weights
std::size_t numPolar(const unsigned int l) const
The number of polar directions associated with the given direction.
const unsigned int _azimuthal_order
The azimuthal order.
unsigned int polarOrder() const
Get the polar order.
static void rotationMatrix(const libMesh::Point &direction, libMesh::DenseMatrix< Real > &matrix)
Builds the rotation matrix for direction direction into matrix.
const Real _mu_max
The maximum mu.
unsigned int dim() const
Get the quadrature dimension.
bool hasDirection(const unsigned int l) const
Whether or not the angular quadrature has direction l.
Real getTotalWeight(const unsigned int l) const
Gets the total of the weights associated with the direction l.
void rotate(const libMesh::Point &rotation_direction)
Rotates the quadrature to a given direction.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static const std::string v
Definition: NS.h:84
const std::vector< Real > & getPolarSins(const unsigned int l) const
Gets the polar sins for the direction l.
static libMesh::Point orthonormalVector(const libMesh::Point &v)
Gets the vector that is orthonormal to v.
void checkDirection(const unsigned int l) const
Throws a MooseError if the angular quadrature does not have direction l.
const std::vector< Real > & getWeights(const unsigned int l) const
Get the weights associated with the direction l.
unsigned int azimuthalOrder() const
Get the azimuthal order.
RayTracingAngularQuadrature(const unsigned int dim, const unsigned int polar_order, const unsigned int azimuthal_order, const Real mu_min, const Real mu_max)
Constructor.
const unsigned int _dim
The dimension.
const libMesh::Point & currentRotationDirection() const
Get the current rotation direction.
std::vector< Real > _polar_sin
Quadrature polar sin.
Real muMin() const
Get the minimum mu.