Celeritas 0.7.0-dev.246+develop.67e9bdac
Loading...
Searching...
No Matches
Functions
ArrayUtils.hh File Reference

Math functions using celeritas::Array. More...

#include <cmath>
#include "corecel/Assert.hh"
#include "corecel/Types.hh"
#include "corecel/cont/Array.hh"
#include "Algorithms.hh"
#include "ArraySoftUnit.hh"
#include "SoftEqual.hh"
#include "detail/ArrayUtilsImpl.hh"

Functions

template<class T , size_type N>
void celeritas::axpy (T a, Array< T, N > const &x, Array< T, N > *y)
 Increment a vector by another vector multiplied by a scalar.
 
template<class T , size_type N>
T celeritas::dot_product (Array< T, N > const &x, Array< T, N > const &y)
 Dot product of two vectors.
 
template<class T >
Array< T, 3 > celeritas::cross_product (Array< T, 3 > const &x, Array< T, 3 > const &y)
 Cross product of two space vectors.
 
template<class T , size_type N>
T celeritas::norm (Array< T, N > const &v)
 Calculate the Euclidean (2) norm of a vector.
 
template<class T , size_type N>
Array< T, N > celeritas::make_unit_vector (Array< T, N > const &v)
 Construct a unit vector.
 
template<class T , size_type N>
Array< T, N > celeritas::make_orthogonal (Array< T, N > const &x, Array< T, N > const &y)
 Return the component of x that is orthogonal to the unit vector y.
 
template<class T , size_type N>
bool celeritas::is_soft_orthogonal (Array< T, N > const &x, Array< T, N > const &y)
 Check whether two unit vectors are approximately orthogonal.
 
template<class T , size_type N>
bool celeritas::is_soft_collinear (Array< T, N > const &x, Array< T, N > const &y)
 Check whether two vectors are approximately collinear.
 
template<class T , size_type N>
T celeritas::distance (Array< T, N > const &x, Array< T, N > const &y)
 Calculate the Euclidean (2) distance between two points.
 
template<class T >
Array< T, 3 > celeritas::from_spherical (T costheta, T phi)
 Calculate a Cartesian vector from spherical coordinates.
 
template<class T >
Array< T, 3 > celeritas::rotate (Array< T, 3 > const &dir, Array< T, 3 > const &rot)
 Rotate a direction about the given scattering direction.
 
template<class T1 , class T2 , size_type N>
Array< T1, N > celeritas::array_cast (Array< T2, N > const &x)
 Convert an array from type T2 to T1.
 

Detailed Description

Math functions using celeritas::Array.

Function Documentation

◆ axpy()

template<class T , size_type N>
void celeritas::axpy ( T  a,
Array< T, N > const x,
Array< T, N > *  y 
)
inline

Increment a vector by another vector multiplied by a scalar.

\[ y \gets \alpha x + y \]

Note that this uses celeritas::fma which supports types other than floating point.

◆ dot_product()

template<class T , size_type N>
T celeritas::dot_product ( Array< T, N > const x,
Array< T, N > const y 
)
inline

Dot product of two vectors.

Note that this uses celeritas::fma which supports types other than floating point.

◆ from_spherical()

template<class T >
Array< T, 3 > celeritas::from_spherical ( T  costheta,
T  phi 
)
inline

Calculate a Cartesian vector from spherical coordinates.

Theta is the angle between the z axis and the outgoing vector, and phi is the angle between the x axis and the projection of the vector onto the x-y plane.

◆ is_soft_collinear()

template<class T , size_type N>
bool celeritas::is_soft_collinear ( Array< T, N > const x,
Array< T, N > const y 
)
inline

Check whether two vectors are approximately collinear.

Precondition
Vectors must be normalized, i.e., have unit magnitude.

◆ is_soft_orthogonal()

template<class T , size_type N>
bool celeritas::is_soft_orthogonal ( Array< T, N > const x,
Array< T, N > const y 
)
inline

Check whether two unit vectors are approximately orthogonal.

Note that the test for orthogonality should use relative tolerance, not absolute.

◆ make_orthogonal()

template<class T , size_type N>
Array< T, N > celeritas::make_orthogonal ( Array< T, N > const x,
Array< T, N > const y 
)
inline

Return the component of x that is orthogonal to the unit vector y.

In this implementation, y must be normalized, and the result is not normalized.

\[ \mathbf{x}' \gets \mathbf{x} - (\mathbf{x} \cdot \mathbf{y}) \mathbf{y} \, , \quad \|\mathbf{y}\| = 1 \]

◆ make_unit_vector()

template<class T , size_type N>
Array< T, N > celeritas::make_unit_vector ( Array< T, N > const v)
inline

Construct a unit vector.

Unit vectors have an Euclidean norm magnitude of 1.

◆ rotate()

template<class T >
Array< T, 3 > celeritas::rotate ( Array< T, 3 > const dir,
Array< T, 3 > const rot 
)
inline

Rotate a direction about the given scattering direction.

The equivalent to calling the Shift transport code's

void cartesian_vector_transform(
double costheta,
double phi,
Vector_View vector);

is the call

vector = rotate(from_spherical(costheta, phi), vector);

This code effectively decomposes the given rotation vector rot into two sequential transform matrices, one with an angle theta about the y axis and one about phi rotating around the z axis. These two angles are the spherical coordinate transform of the given rot cartesian direction vector.

There is some extra code in here to deal with loss of precision when the incident direction is along the z axis. As rot approaches z, the azimuthal angle \( \phi \) must be calculated carefully from both the x and y components of the vector, not independently. If rot actually equals z then the azimuthal angle is completely indeterminate so we arbitrarily choose \( \phi = 0 \).

This function is often used for calculating exiting scattering angles. In that case, dir is the exiting angle from the scattering calculation, and rot is the original direction of the particle. The direction vectors are defined as

\[ \vec{\Omega} = \sin\theta\cos\phi\vec{i} + \sin\theta\sin\phi\vec{j} + \cos\theta\vec{k} \,. \]