|
| 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...
|
|
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
-
D | dimensionality of wavepacket |
MulitIndex | The type used to represent multi-indices. |
N | Number of quadrature points. Don't use Eigen::Dynamic. It works, but performance is bad. |
template<dim_t D, class MultiIndex , int N>
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] | coefficients | Vector 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>
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_enums | The shape enumerations of each wavepacket component. |
[in] | subset_coeffs | The coefficients of each wavepacket component. |
[in] | n_components | The number of wavepacket components \( C \). |
- Returns
- Complex 2D-Array of shape \( (C \times N) \), where \( N \) is the number of quadrature points.