WaveBlocksND
Public Types | Public Member Functions | Private Attributes | List of all members
waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N > Class Template Reference

Evaluates a wavepacket slice by slice. More...

#include <hawp_evaluator.hpp>

Public Types

typedef Eigen::Matrix< complex_t, D, D > CMatrixDD
 
typedef Eigen::Array< complex_t, 1, N > CArray1N
 
typedef Eigen::Matrix< complex_t, D, N > CMatrixDN
 
typedef Eigen::Matrix< real_t, D, N > RMatrixDN
 

Public Member Functions

 HaWpEvaluator (real_t eps, const HaWpParamSet< D > *parameters, const ShapeEnum< D, MultiIndex > *enumeration, const CMatrixDN &x)
 
CArray1N seed () const
 Evaluates basis function on node \( \underline{0} \) (ground-state). More...
 
HaWpBasisVector< N > step (std::size_t islice, const HaWpBasisVector< N > &prev_basis, const HaWpBasisVector< N > &curr_basis) const
 Having basis function values on previous and current slice, this member function computes basis function values on the next slice (using recursive evaluation formula). More...
 
HaWpBasisVector< N > all () const
 Evaluates all basis functions. More...
 
CArray< 1, N > reduce (const Coefficients &coefficients) const
 Evaluates wavepacket in a memory efficient manner. More...
 
CArray< Eigen::Dynamic, N > vector_reduce (ShapeEnum< D, MultiIndex > **subset_enums, complex_t const **subset_coeffs, std::size_t n_components) const
 Efficiently evaluates a vectorial wavepacket with shared Hagedorn parameter set, but different basis shapes. More...
 
CArray< Eigen::Dynamic, N > vector_reduce (complex_t const **coefficients, std::size_t n_components) const
 Efficiently evaluates a vectorial wavepacket with shared Hagedorn parameter set and shared basis shapes. More...
 

Private Attributes

real_t eps_
 
const HaWpParamSet< D > * parameters_
 
const ShapeEnum< D, MultiIndex > * enumeration_
 
int npts_
 
CMatrixDN dx_
 
CMatrixDD Qinv_
 
CMatrixDD Qh_Qinvt_
 
CMatrixDN Qinv_dx_
 
std::vector< real_tsqrt_
 

Detailed Description

template<dim_t D, class MultiIndex, int N>
class waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >

Evaluates a wavepacket slice by slice.

This class is low-level. You should not use it directly. You should use the high-level member functions AbstractScalarHaWp::evaluate() and AbstractScalarHaWpBasis::evaluate_basis().

The only reason you may want to use HaWpEvaluator directly is when you gain an advantage by evaluating a wavepacket slice-by-slice. The slice-by-slice evaluation reduces memory since you don't have to store all basis function values, but only the 'active' ones. Take a look at the implementation of the member functions all() or reduce() to learn, how to evaluate a wave packet slice-by-slice.

Template Parameters
Ddimensionality of wavepacket
MulitIndexThe type used to represent multi-indices.
NNumber of quadrature points. Don't use Eigen::Dynamic. It works, but performance is bad.

Member Typedef Documentation

template<dim_t D, class MultiIndex , int N>
typedef Eigen::Array<complex_t,1,N> waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::CArray1N
template<dim_t D, class MultiIndex , int N>
typedef Eigen::Matrix<complex_t,D,D> waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::CMatrixDD
template<dim_t D, class MultiIndex , int N>
typedef Eigen::Matrix<complex_t,D,N> waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::CMatrixDN
template<dim_t D, class MultiIndex , int N>
typedef Eigen::Matrix<real_t,D,N> waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::RMatrixDN

Constructor & Destructor Documentation

template<dim_t D, class MultiIndex , int N>
waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::HaWpEvaluator ( real_t  eps,
const HaWpParamSet< D > *  parameters,
const ShapeEnum< D, MultiIndex > *  enumeration,
const CMatrixDN x 
)
inline
Parameters
[in]epsThe semi-classical scaling parameter \( \varepsilon \) of the wavepacket.
[in]parametersThe Hagedorn parameter set \( \Pi \) of the wavepacket.
[in]enumerationThe basis shape \( \mathfrak{K} \) of the wavepacket. Or the union of all component's basis shapes, if you want to evaluate a vectorial wavepacket.
[in]xQuadrature points: Complex matrix of shape \( (D \times N) \), where \( N \) number of quadrature points).

Member Function Documentation

template<dim_t D, class MultiIndex , int N>
HaWpBasisVector<N> waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::all ( ) const
inline

Evaluates all basis functions.

Returns
Complex 2D-Array of shape \( (|\mathfrak{K}| \times N) \), where \( N \) is the number of quadrature points.
template<dim_t D, class MultiIndex , int N>
CArray<1,N> waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::reduce ( const Coefficients coefficients) const
inline

Evaluates wavepacket in a memory efficient manner.

This function computes the dot product (thus the name 'reduce') of the wavepacket coefficients with the wavepacket basis. This is done by evaluating the wavepacket basis slice by slice and multiplying the basis values with the coefficients on the fly. Computed basis values are discarded, once they are not needed any more.

Parameters
[in]coefficientsVector of wavepacket coefficients, length is \( |\mathfrak{K}| \).
Returns
Complex 2D-Array of shape \( (1 \times N) \), where \( N \) is the number of quadrature points.
template<dim_t D, class MultiIndex , int N>
CArray1N waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::seed ( ) const
inline

Evaluates basis function on node \( \underline{0} \) (ground-state).

Returns
Complex 2D-Array of shape \( (D \times N) \), where \( N \) is the number of quadrature points.
template<dim_t D, class MultiIndex , int N>
HaWpBasisVector<N> waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::step ( std::size_t  islice,
const HaWpBasisVector< N > &  prev_basis,
const HaWpBasisVector< N > &  curr_basis 
) const
inline

Having basis function values on previous and current slice, this member function computes basis function values on the next slice (using recursive evaluation formula).

Hint: Use function seed() to bootstrap recursion.

Parameters
[in]isliceOrdinal of current slice.
[in]prev_basisBasis values on previous slice. Type: Complex 2D-array of shape \( (S \times N) \), where \( S \) is the number of nodes in the current slice and \( N \) is the number of quadrature points.
[in]curr_basisB Basis values on current slice. Type: Complex 2D-array of shape \( (S \times N) \)
Returns
Computed basis values on next slice. Type: Complex 2D-Array of shape \( (S^+ \times N) \), where \( S^+ \) is the number of nodes in the next slice.
template<dim_t D, class MultiIndex , int N>
CArray<Eigen::Dynamic,N> waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::vector_reduce ( ShapeEnum< D, MultiIndex > **  subset_enums,
complex_t const **  subset_coeffs,
std::size_t  n_components 
) const
inline

Efficiently evaluates a vectorial wavepacket with shared Hagedorn parameter set, but different basis shapes.

The union of all component's basis shapes is passed to the constructor.

This function is used to evaluate a homogeneous wavepacket (see HomogeneousHaWp::evaluate())

Parameters
[in]subset_enumsThe shape enumerations of each wavepacket component.
[in]subset_coeffsThe coefficients of each wavepacket component.
[in]n_componentsThe number of wavepacket components \( C \).
Returns
Complex 2D-Array of shape \( (C \times N) \), where \( N \) is the number of quadrature points.
template<dim_t D, class MultiIndex , int N>
CArray<Eigen::Dynamic,N> waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::vector_reduce ( complex_t const **  coefficients,
std::size_t  n_components 
) const
inline

Efficiently evaluates a vectorial wavepacket with shared Hagedorn parameter set and shared basis shapes.

This function is used to evaluate a wavepacket gradient (see HaWpGradient::evaluate()).

Parameters
[in]coefficientsThe coefficients of each wavepacket component.
[in]n_componentsThe number of wavepacket components.
[in]n_componentsThe number of wavepacket components \( C \).
Returns
Complex 2D-Array of shape \( (C \times N) \), where \( N \) is the number of quadrature points.

Member Data Documentation

template<dim_t D, class MultiIndex , int N>
CMatrixDN waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::dx_
private

precomputed expression: x - q

template<dim_t D, class MultiIndex , int N>
const ShapeEnum<D,MultiIndex>* waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::enumeration_
private
template<dim_t D, class MultiIndex , int N>
real_t waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::eps_
private
template<dim_t D, class MultiIndex , int N>
int waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::npts_
private

number of quadrature points

template<dim_t D, class MultiIndex , int N>
const HaWpParamSet<D>* waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::parameters_
private
template<dim_t D, class MultiIndex , int N>
CMatrixDD waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::Qh_Qinvt_
private

precomputed expression: Q^H * Q^{-T}

template<dim_t D, class MultiIndex , int N>
CMatrixDD waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::Qinv_
private

precomputed expression: Q^{-1}

template<dim_t D, class MultiIndex , int N>
CMatrixDN waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::Qinv_dx_
private

precomputed expression: Q^{-1} * (x - q)

template<dim_t D, class MultiIndex , int N>
std::vector<real_t> waveblocks::wavepackets::HaWpEvaluator< D, MultiIndex, N >::sqrt_
private

lookup-table for sqrt


The documentation for this class was generated from the following file: