Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Class List | Directories | File List | Namespace Members | Class Members | File Members

Methods for evolving states unitarily

A set of methods for evolving states according to Hamiltonian, using a Chebyshev polynomial expansion. See extensive documentation in chebyshev.cpp for further details. More...

Functions

void GenerateChebyshevMatrix (itype N, DenseMatrix< complex > &M)
 Computes an NxN matrix, where $ M_{m,n} $ is the coefficient of $ x^m $ in the $ n $ 'th Chebyshev polynomial.
void ResizeChebyshevMatrix (itype N, DenseMatrix< complex > &M)
 Efficiently resizes an existing Chebyshev polynomial matrix to NxN, using the coefficients already computed.
void GenerateChebyshevEvolutionVector (double time, itype N, DenseVector< complex > &V)
 Computes an N-element vector, where $ v_n $ is the coefficient of the $ n $ 'th Chebyshev polynomial in the expansion of $ e^{-ix} $ .
int ChebyshevPrep (double timestep, double prec, DenseVector< complex > &PowerCoeffs, DenseVector< double > &CoeffsBoundSeq)
 Precomputes the coefficients needed to evolve a state using ChebyshevStep(), for time timestep, with accuracy prec.
template<typename T>
int ChebyshevStep (const State &initial, State &final, DenseVector< complex > &PowerCoeffs, DenseVector< double > &CoeffsBoundSeq, const T &Hamiltonian, double E_bound, double prec)
 Transforms the quantum state initial by $ e^{-i\mathbf{H}} $ , placing the result in final. Requires precomputation by ChebyshevPrep().
template<typename T>
int ChebyshevEvolve (const State &initial, State &final, double time, const T &Hamiltonian, double E_bound, double prec=1e-10)
 Evolves the quantum state initial for a single step of duration time according to Hamiltonian, placing the result in final.
template<typename T>
int ChebyshevStepEvolve (const State &initial, State &final, double time, double maxtimestep, const T &Hamiltonian, double E_bound, double prec=1e-10)
 Evolves the quantum state initial ==> final, for duration time, in steps no longer than maxtimestep, according to Hamiltonian.

Detailed Description

A set of methods for evolving states according to Hamiltonian, using a Chebyshev polynomial expansion. See extensive documentation in chebyshev.cpp for further details.


Function Documentation

void GenerateChebyshevMatrix itype  N,
DenseMatrix< complex > &  M
 

Computes an NxN matrix, where $ M_{m,n} $ is the coefficient of $ x^m $ in the $ n $ 'th Chebyshev polynomial.

The matrix generated by GenerateChebyshevMatrix() represents all the Chebyshev polynomials of order N or less. $ T_n(x) $ is represented by the $ n $ th column of M. The entry in the $ m $ th row is the coefficient of $ x^m $ . M can therefore transform a column vector that represents $ f(x) $ as a weighted sum of Chebyshev polynomials (i.e., a Chebyshev approximation to some other function) into a column vector that represents the same function as a polynomial in $ x $ .

The coefficients are computed efficiently (in $ O(N^2) $ time) using the recursion relation $ T_{n+1}(x) = 2xT_n(x)-T_{n-1}(x) $ , with $ T_0(x) = 1 $ and $ T_1(x) = x $ .

void ResizeChebyshevMatrix itype  N,
DenseMatrix< complex > &  M
 

Efficiently resizes an existing Chebyshev polynomial matrix to NxN, using the coefficients already computed.

The structure of the Chebyshev polynomial recursion relation means that:

  • high order polynomials cannot be computed without first computing all the lower order ones,
  • given the lower order ones, the higher order polynomials can be computed very rapidly. Therefore, if we discover that we need the Chebyshev polynomials of up to Nth order, and we already have the ones of up to (N-k)th order, then ResizeChebyshevMatrix() can leverage the existing coefficients to compute the additional ones. This takes $ O(2Nk) $ time, as compared to $ O(N^2) $ to recompute from scratch.

void GenerateChebyshevEvolutionVector double  time,
itype  N,
DenseVector< complex > &  V
 

Computes an N-element vector, where $ v_n $ is the coefficient of the $ n $ 'th Chebyshev polynomial in the expansion of $ e^{-ix} $ .

Evolving a state vector using the Chebyshev polynomial approximation requires first expanding $ e^{-itx} $ in Chebyshev polynomials. The expansion is computed by GenerateChebyshevEvolutionVector(), and the coefficient of $ T_n(x) $ in the expansion is stored in $ V_n $ . The coefficients can be worked out analytically, and turn out to be Bessel $ J $ functions:

\[ V_n = (2i)^n J_n(t) \]

The Bessel functions are computed all at once using the Gnu Scientific Library's gsl_sf_bessel_Jn_array() function.

int ChebyshevPrep double  timestep,
double  prec,
DenseVector< complex > &  PowerCoeffs,
DenseVector< double > &  CoeffsBoundSeq
 

Precomputes the coefficients needed to evolve a state using ChebyshevStep(), for time timestep, with accuracy prec.

ChebyshevPrep() constructs the coefficients of a polynomial in a Hamiltonian $ \mathbf{\tilde{H}} $ that approximates $ e^{-it\mathbf{\tilde{H}}} $ (where $ t $ is given by timestep) with an error of no more than prec (provided that $ ||\mathbf{\tilde{H}}||_\infty \leq 1 $ . These coefficients are returned in the vector PowerCoeffs.

The order of approximation required (N) depends on both timestep and prec. We first compute an upper bound $ N_{max} \geq N $ , then generate the coefficients $ c_n $ of the $ N_{max} $ th order approximation, and determine the $ N $ such that $ c_n < \mathrm{prec}\ \forall\ n>N $ . We then terminate the expansion at $ N $ th order, generate the matrix of Chebyshev polynomial coefficients, and transform the series sum of Chebyshev polynomials into a sum of monomials $ \mathbf{\tilde{H}}^n $ .

Because the coefficient sequence ($ c_n $ ) is not necessarily decreasing in absolute value, we also compute a bounding sequence ($ b_n $ , or CoeffsBoundSeq). The bounding sequence is defined so that $ b_n \geq |c_m|\ \forall\ m\geq n $ .

The bounding sequence allows the routine that actually applies $ \mathbf{ \mathbf{\tilde{H}}} $ repeatedly to a state (ChebyshevStep()) to do dynamic error estimation. Since $ \mathbf{ \mathbf{\tilde{H}}} $ will typically be very contracting (all eigenvalues will be closer to zero than to -1 or 1), the norm of $ \mathbf{ \mathbf{\tilde{H}}}^n\left|\psi\right\rangle $ will decrease rapidly with $ n $ . When it decreases enough (relative to the bounding sequence of the coefficients), ChebyshevStep() can terminate the expansion.


Generated on Wed Jun 14 22:25:40 2006 for linalg by  doxygen 1.4.4