HagedornWavepacket

About the HagedornWavepacket class

The WaveBlocks Project

This file contains the class which represents a homogeneous Hagedorn wavepacket.

@author: R. Bourquin @copyright: Copyright (C) 2010, 2011 R. Bourquin @license: Modified BSD License

Inheritance diagram

Inheritance diagram of HagedornWavepacket

Class documentation

class WaveBlocks.HagedornWavepacket.HagedornWavepacket(parameters)[source]

This class represents homogeneous vector valued wavepackets |\Psi\rangle.

coefficients = None

The coefficients c^i of the linear combination for each component \Phi_k.

evaluate_at(nodes, component=None, prefactor=False)[source]

Evaluete the Hagedorn wavepacket \Psi at the given nodes \gamma.

Parameters:
  • nodes – The nodes \gamma at which the Hagedorn wavepacket gets evaluated.
  • component – The index i of a single component \Phi_i to evaluate. (Defaults to ‘None’ for evaluating all components.)
  • prefactor – Whether to include a factor of \left(\det\left(Q\right)\right)^{-\frac{1}{2}}.
Returns:

A list of arrays or a single array containing the values of the \Phi_i at the nodes \gamma.

evaluate_basis_at(nodes, component=None, prefactor=False)[source]

Evaluate the Hagedorn functions \phi_k recursively at the given nodes \gamma.

Parameters:
  • nodes – The nodes \gamma at which the Hagedorn functions are evaluated.
  • component – Takes the basis size K_i of this component i as upper bound for K.
  • prefactor – Whether to include a factor of \left(\det\left(Q\right)\right)^{-\frac{1}{2}}.
Returns:

Returns a twodimensional K times #nodes array H where the entry H[k,i] is the value of the k-th Hagedorn function evaluated at the node i.

gen_id()

Generate an (unique) ID per wavepacket instance.

get_basis_size(component=None)
Returns:The size of the basis, i.e. the number K of {\phi_k}_{k=1}^K.
get_coefficient_vector()
Returns:The coefficients c^i of all components \Phi_i as a single long column vector.
get_coefficients(component=None)

Returns the coefficients c^i for some components \Phi_i of |\Psi\rangle.

Parameters:component – The index i of the coefficients c^i we want to get.
Returns:The coefficients c^i either for all components \Phi_i or for a specified one.
get_id()

Return the packet ID of this wavepacket instance. The ID may be used for storing packets in associative lists.

get_norm(component=None, summed=False)[source]

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

Parameters:
  • component – The component \Phi_i of which the norm is calculated.
  • summed – Whether to sum up the norms of the individual components \Phi_i.
Returns:

A list containing the norms of all components \Phi_i or the overall norm of \Psi.

get_number_components()
Returns:The number N of components the wavepacket \Psi has.
get_parameters(component=None, aslist=False)[source]

Get the Hagedorn parameters \Pi of the wavepacket \Psi.

Parameters:
  • component – Dummy parameter for API compatibility with the inhomogeneous packets.
  • aslist – Return a list of N parameter tuples. This is for API compatibility with inhomogeneous packets.
Returns:

The Hagedorn parameters P, Q, S, p, q of \Psi in this order.

get_quadrature()[source]

Return the HomogeneousQuadrature instance used for evaluating brakets.

Returns:The current instance HomogeneousQuadrature.
grady(component)[source]

Compute the effect of the operator -i \varepsilon^2 \frac{\partial}{\partial x} on the basis functions of a component \Phi_i of the Hagedorn wavepacket \Psi.

Parameters:component – The index i of the component \Phi_i on which we apply the above operator.
Returns:The modified coefficients.
kinetic_energy(summed=False)[source]

Calculate the kinetic energy \langle\Psi|T|\Psi\rangle of the wavepacket componentwise.

Parameters:summed – Wheter to sum up the individual integrals \langle\Phi_i|T_{i,j}|\Phi_j\rangle.
Returns:The kinetic energy of the wavepacket’s components \Phi_i or the overall kinetic energy of \Psi.
number_components = None

Number of components \Phi_i the wavepacket |\Psi\rangle has got.

potential_energy(potential, summed=False)[source]

Calculate the potential energy \langle\Psi|V|\Psi\rangle of the wavepacket componentwise.

Parameters:
  • potential – The potential energy operator V as function.
  • summed – Wheter to sum up the individual integrals \langle\Phi_i|V_{i,j}|\Phi_j\rangle.
Returns:

The potential energy of the wavepacket’s components \Phi_i or the overall potential energy of \Psi.

project_to_canonical(potential, assign=True)[source]

Project the Hagedorn wavepacket into the canonical basis.

Parameters:
  • potential – The potential V whose eigenvectors nu_l are used for the transformation.
  • assign – Whether to assign the new coefficient values to the wavepacket. Default true.

Note

This function is expensive and destructive! It modifies the coefficients of the self instance if the assign parameter is True (default).

project_to_eigen(potential, assign=True)[source]

Project the Hagedorn wavepacket into the eigenbasis of a given potential V.

Parameters:
  • potential – The potential V whose eigenvectors nu_l are used for the transformation.
  • assign – Whether to assign the new coefficient values to the wavepacket. Default true.

Note

This function is expensive and destructive! It modifies the coefficients of the self instance if the assign parameter is True (default).

quadrature = None

An object that can compute brakets via quadrature.

set_basis_size(basis_size, component=None)

Set the size of the basis of a given component or all components.

Parameters:
  • basis_size – An single positive integer or a list of N positive integers.
  • component – The component for which we want to set the basis size. Default is None which means ‘all’.
set_coefficient(component, index, value)

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

Parameters:
  • component – The index i of the component \Phi_i we want to update.
  • index – The index k of the coefficient c^i_k we want to update.
  • value – The new value of the coefficient c^i_k.
Raises ValueError:
 

For invalid indices i or k.

set_coefficient_vector(vector)

Set the coefficients for all components \Phi_i simultaneously.

Parameters:vector – The coefficients of all components as a single long column vector.

Note

This function does NOT copy the input data! This is for efficiency as this routine is used in the innermost loops.

set_coefficients(values, component=None)

Update the coefficients c of \Psi.

Parameters:
  • values – The new values of the coefficients c^i of \Phi_i.
  • component – The index i of the component we want to update with new coefficients.
Raises ValueError:
 

For invalid indices i.

Note

This function can either set new coefficients for a single component \Phi_i only if the component attribute is set or for all components simultaneously if values is a list of arrays.

set_id(anid)

Manually set an ID for the current wavepacket instance.

set_parameters(parameters, component=None)[source]

Set the Hagedorn parameters \Pi of the wavepacket \Psi.

Parameters:
  • parameters – The Hagedorn parameters P, Q, S, p, q of \Psi in this order.
  • component – Dummy parameter for API compatibility with the inhomogeneous packets.
set_quadrature(quadrature)[source]

Set the HomogeneousQuadrature instance used for evaluating brakets.

Parameters:quadrature – The new HomogeneousQuadrature instance. May be None to use a dafault one with a quadrature rule of order K+4.
to_fourier_space(assign=True)[source]

Transform the wavepacket to Fourier space.

Parameters:assign – Whether to assign the transformation to this packet or return a cloned packet.

Note

This is the inverse of the method to_real_space().

to_real_space(assign=True)[source]

Transform the wavepacket to real space.

Parameters:assign – Whether to assign the transformation to this packet or return a cloned packet.

Note

This is the inverse of the method to_fourier_space().

Table Of Contents

Previous topic

Wavepacket

Next topic

HagedornWavepacketInhomogeneous

This Page