HagedornWavepacketBase

About the HagedornWavepacketBase 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

Inheritance diagram of HagedornWavepacketBase

Class documentation

class WaveBlocksND.HagedornWavepacketBase(parameters)[source]

This class implements the abstract Wavepacket interface and contains code common to all types of Hagedorn wavepackets.

clone()

Clone the wavepacket. Return a new copy of the wavepacket and make sure that all references between the two wavepackets get broken.

Raise:NotImplementedError Abstract interface.
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)[source]

Retrieve the basis shapes \mathfrak{K}_i for each component i.

Parameters:component (int) – The component i 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)[source]

Retrieve a single coefficient c^i_k of the specified component \Phi_i of \Psi.

Parameters:
  • component – The index i of the component \Phi_i we want to update.
  • index (A tuple of D integers.) – The multi-index k of the coefficient c^i_k we want to update.
Returns:

A single complex number.

Raise:

ValueError For invalid indices i or k.

get_coefficient_vector(*, component=None)[source]

Retrieve the coefficients for all components \Phi_i 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 i whose coefficients we request. (Default is None which means to return the coefficients for all components.
Returns:The coefficients c^i of all components \Phi_i stacked into a single long column vector.
get_coefficients(*, component=None)[source]

Returns the coefficients c^i for some component \Phi_i of \Psi or all the coefficients c of all components.

Note: this method copies the data arrays.

Parameters:component (int (Default is None meaning all)) – The index i 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 (|\mathfrak{K}_i|, 1).
Raise:ValueError For invalid component indices i.
get_dimension()
Returns:The space dimension D of the wavepacket \Psi.
get_eps()[source]

Retrieve the semi-classical scaling parameter \varepsilon of the wavepacket.

Returns:The value of \varepsilon.
get_gradient_operator()[source]

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()[source]

Return the InnerProduct subclass instance used computing inner products and evaluating brakets.

Returns:The current InnerProduct subclass instance.
get_number_components()
Returns:The number N of components the wavepacket \Psi has.
norm(component=None, summed=False)[source]

Calculate the L^2 norm \langle\Psi|\Psi\rangle of the wavepacket \Psi.

Parameters:
  • component (int or None.) – The index i of the component \Phi_i whose norm is calculated. The default value is None which means to compute the norms of all N components.
  • summed (Boolean, default is False.) – Whether to sum up the norms \langle\Phi_i|\Phi_i\rangle of the individual components \Phi_i.
Returns:

The norm of \Psi or the norm of \Phi_i or a list with the N norms of all components. Depending on the values of component and summed.

set_basis_shapes(basis_shape, *, component=None)[source]

Set the basis shape \mathfrak{K} 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 N shapes.
  • component (int) – The component i whose basis shape we want to set. (Default is None which means to set the basis shapes for all components.
set_coefficient(component, index, value)[source]

Set a single coefficient c^i_k of the specified component \Phi_i of \Psi.

Parameters:
  • component – The index i of the component \Phi_i we want to update.
  • index (A tuple of D integers.) – The multi-index k of the coefficient c^i_k we want to update.
  • value – The new value of the coefficient c^i_k.
Raise:

ValueError For invalid indices i or k.

set_coefficient_vector(vector)[source]

Set the coefficients for all components \Phi_i 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)[source]

Update all the coefficients c of \Psi or update the coefficients c^i of the components \Phi_i 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 c^i of \Phi_i.
  • component (int (Default is None meaning all)) – The index i of the component we want to update with new coefficients.
Raise:

ValueError For invalid component indices i.

set_id(anid)

Manually set a new ID for the current wavepacket instance.

Parameters:anid (str) – The new ID.
set_innerproduct(innerproduct)[source]

Set the InnerProduct subclass instance used for computing inner products and evaluating brakets.

Parameters:innerproduct – The new InnerProduct subclass instance.