HagedornWavepacketInhomogeneous¶
About the HagedornWavepacketInhomogeneous
class¶
The WaveBlocks Project
@author: R. Bourquin @copyright: Copyright (C) 2010, 2011, 2012, 2013, 2014, 2015, 2016 R. Bourquin @license: Modified BSD License
Inheritance diagram¶
Class documentation¶
-
class
WaveBlocksND.
HagedornWavepacketInhomogeneous
(dimension, ncomponents, eps)[source]¶ This class represents inhomogeneous vector valued Hagedorn wavepackets with components in space dimensions.
-
__init__
(dimension, ncomponents, eps)[source]¶ Initialize a new in homogeneous Hagedorn wavepacket.
Parameters: - dimension – The space dimension the packet has.
- ncomponents – The number of components the packet has.
- eps – The semi-classical scaling parameter of the basis functions.
Returns: An instance of
HagedornWavepacketInhomogeneous
.
-
evaluate_at
(grid, *, component=None, prefactor=False)¶ Evaluate the Hagedorn wavepacket at the given nodes .
Parameters: - grid (A class having a
get_nodes(...)()
method.) – The grid containing the nodes . - component – The index of a single component to evaluate.
(Defaults to
None
for evaluating all components.) - prefactor (Boolean, default is
False
.) – Whether to include a factor of .
Returns: A list of arrays or a single array containing the values of the at the nodes .
- grid (A class having a
-
evaluate_basis_at
(grid, component, *, prefactor=False)¶ Evaluate the basis functions recursively at the given nodes .
Parameters: - grid (A class having a
get_nodes(...)()
method.) – The grid containing the nodes . - component – The index of a single component to evaluate.
- prefactor (Boolean, default is
False
.) – Whether to include a factor of .
Returns: A two-dimensional ndarray of shape where the entry is the value of .
- grid (A class having a
-
gen_id
()¶ Generate an (unique) ID per wavepacket instance.
Note
The
packet id
is a string of length 32 because this is exactly the length of an ‘md5’ digest in hex representation.
-
get_basis_shapes
(*, component=None)¶ Retrieve the basis shapes for each component .
Parameters: component (int) – The component whose basis shape we request. (Default is None
which means to return the basis shapes for all components.Returns: The basis shape for an individual component or a list with all shapes.
-
get_coefficient
(component, index)¶ Retrieve a single coefficient of the specified component of .
Parameters: - component – The index of the component we want to update.
- index (A tuple of integers.) – The multi-index of the coefficient we want to update.
Returns: A single complex number.
Raise: ValueError
For invalid indices or .
-
get_coefficient_vector
(*, component=None)¶ Retrieve the coefficients for all components simultaneously.
Warning
This function does not copy the input data! This is for efficiency as this routine is used in the innermost loops.
Parameters: component (int) – The component whose coefficients we request. (Default is None
which means to return the coefficients for all components.Returns: The coefficients of all components stacked into a single long column vector.
-
get_coefficients
(*, component=None)¶ Returns the coefficients for some component of or all the coefficients of all components.
Note: this method copies the data arrays.
Parameters: component (int (Default is None
meaning all)) – The index of the component we want to retrieve.Returns: A single ndarray with the coefficients of the given component or a list containing the ndarrays for each component. Each ndarray is two-dimensional with a shape of . Raise: ValueError
For invalid component indices .
-
get_description
()[source]¶ Return a description of this wavepacket object. A description is a
dict
containing all key-value pairs necessary to reconstruct the current instance. A description never contains any data.
-
get_dimension
()¶ Returns: The space dimension of the wavepacket .
-
get_eps
()¶ Retrieve the semi-classical scaling parameter of the wavepacket.
Returns: The value of .
-
get_gradient_operator
()¶ Return the
Gradient
subclass suitable for computing gradients of this wavepacket.Returns: A GradientHAWP
instance.
-
get_id
()¶ Return the packet ID of this wavepacket instance. The ID may be used for storing packets in associative lists.
Returns: The ID of the current instance.
-
get_innerproduct
()¶ Return the
InnerProduct
subclass instance used computing inner products and evaluating brakets.Returns: The current InnerProduct
subclass instance.
-
get_number_components
()¶ Returns: The number of components the wavepacket has.
-
get_parameters
(component=None, aslist=False, key=('q', 'p', 'Q', 'P', 'S'))[source]¶ Get the Hagedorn parameter set of each component :math`Phi_i` of the wavepacket .
Parameters: - component – The index of the component whose parameters we want to get.
- aslist – Dummy parameter for API compatibility with the homogeneous packets.
Returns: A list with all parameter sets or a single parameter set. The parameters are always in this order.
-
norm
(component=None, summed=False)¶ Calculate the norm of the wavepacket .
Parameters: - component (int or
None
.) – The index of the component whose norm is calculated. The default value isNone
which means to compute the norms of all components. - summed (Boolean, default is
False
.) – Whether to sum up the norms of the individual components .
Returns: The norm of or the norm of or a list with the norms of all components. Depending on the values of
component
andsummed
.- component (int or
-
set_basis_shapes
(basis_shape, *, component=None)¶ Set the basis shape of a given component or for all components.
Parameters: - basis_shape (A subclass of
BasisShape
.) – The basis shape for an individual component or a list with all shapes. - component (int) – The component whose basis shape we want to set. (Default is
None
which means to set the basis shapes for all components.
- basis_shape (A subclass of
-
set_coefficient
(component, index, value)¶ Set a single coefficient of the specified component of .
Parameters: - component – The index of the component we want to update.
- index (A tuple of integers.) – The multi-index of the coefficient we want to update.
- value – The new value of the coefficient .
Raise: ValueError
For invalid indices or .
-
set_coefficient_vector
(vector)¶ Set the coefficients for all components simultaneously.
Warning
This function does not copy the input data! This is for efficiency as this routine is used in the innermost loops.
Parameters: vector (A two-dimensional ndarray of appropriate shape.) – The coefficients of all components as a single long column vector.
-
set_coefficients
(values, *, component=None)¶ Update all the coefficients of or update the coefficients of the components only.
Note: this method copies the data arrays.
Parameters: - values (An ndarray of suitable shape or a list of ndarrays.) – The new values of the coefficients of .
- component (int (Default is
None
meaning all)) – The index of the component we want to update with new coefficients.
Raise: ValueError
For invalid component indices .
-
set_id
(anid)¶ Manually set a new ID for the current wavepacket instance.
Parameters: anid (str) – The new ID.
-
set_innerproduct
(innerproduct)¶ Set the
InnerProduct
subclass instance used for computing inner products and evaluating brakets.Parameters: innerproduct – The new InnerProduct
subclass instance.
-
set_parameters
(Pi, component=None, key=('q', 'p', 'Q', 'P', 'S'))[source]¶ Set the Hagedorn parameter set of each component :math`Phi_i` of the wavepacket .
Parameters: - Pi (A single tuple or a list of tuples) – The parameter sets with its values in this order.
- component – The index of the component whose parameters we want to update.
-
slim_recursion
(grid, component, *, prefactor=False)¶ Evaluate the Hagedorn wavepacket at the given nodes . This routine is a slim version compared to the full basis evaluation. At every moment we store only the data we really need to compute the next step until we hit the highest order basis functions.
Parameters: - grid (A class having a
get_nodes(...)()
method.) – The grid containing the nodes . - component – The index of a single component to evaluate.
- prefactor (Boolean, default is
False
.) – Whether to include a factor of .
Returns: A list of arrays or a single array containing the values of the at the nodes .
Note that this function does not include the global phase .
- grid (A class having a
-