MCTDH-X  v2.22
 All Classes Files Functions Variables Pages
dav_integrators Module Reference

Module for the davidson integrators. More...

Public Member Functions

subroutine davstep (Psi, DtPsi, PsiDim, IntOrder, TolError,
 SUBROUTINE DAVSTEP *

  • Diagonalises the Hamiltonian represented in the SPFs. * Improved relaxation for Hermitian (complex) Hamiltonians. *
  • Input parameters: * Psi: The (complex) initial-value vector. * DtPsi: Action of Hamiltonian on the initial-value vector, * i.e. H|Psi> . * PsiDim Length of Psi and DtPsi vectors (i.e. adim). * IntOrder: Maximum interation order. * TolError: Maximum error that is tolerated. * Steps: Number of steps made so far (Steps is passed to * "WriteStepDav"). * Krylov: Matrix with minimum size PsiDim*IntOrder with * columns containing the Davison vectors. * HKrylov: Matrix containing H*Krylov * rlxemin: Lower energy to restrict the search for max. ovlp * rlxemax: Upper energy to restrict the search for max. ovlp *
  • Output parameters: * Psi: Propagated Psi. * Krylov: Matrix with minimum size PsiDim*(IntOrder-1) with * columns containing the Krylov vectors H^2|psi> to * H^TrueOrder|psi>. * TrueOrder: The order that has actually been used (may be less * than IntOrder). * Steps: Same as on entry. * ErrorCode: Error code having the following meaning: * 0: everything was o. k., * 1: illegal integration order, * 3: diagonalization failed, *
  • External routines: * Func: Computes the action of the Hamiltonian on Psi. * Called as * Func(Abs(Time),Psi,DtPsi). * WriteStep: Writes the stepsize and error to a file. * Called as WriteStepDav(Steps,Order,E,beta,Error,Q,t). * wrtrlxdav: Writes info to rlx_info file. *
  • Some new Variables: * rlaxnum : This is the argument of the keyword relaxation. * If relaxation has no argument, rlaxnum=-1 . * rlaxnum=0 : relaxation to ground-state. * rlaxnum=n : relaxation to n-th state (n=0...900). * rlaxnum=920 : relaxation=follow. * rlaxnum=930 : relaxation=lock. *
  • V8.3 06/2003 HDM *.
More...
 
subroutine bdavstep (Psi, DtPsi, PsiDim, IntOrder, TolError,
 
More...
 
subroutine orthonc (Krylov, v, ovl, n, PsiDim)
 

Detailed Description

Module for the davidson integrators.

Member Function/Subroutine Documentation

subroutine dav_integrators::bdavstep ( complex*16, dimension(adim)  Psi,
complex*16, dimension(adim)  DtPsi,
integer  PsiDim,
integer  IntOrder,
real*8  TolError 
)


  • *
  • SUBROUTINE BDAVSTEP *
  • *
  • Block-Davidson diagonalisation. (single-set multi-packet run). *
  • Diagonalises the Hamiltonian represented in the SPFs. *
  • Assuming hermitian Hamiltonian and complex WF. *
  • *
  • *
  • *
  • Input parameters: *
  • Psi: The block of (complex) initial-value vectors. *
  • DtPsi: Action of Hamiltonian on the initial-value vector, *
  • i.e. H|Psi> . *
  • PsiDim Length of Psi and DtPsi vectors (i.e. adim/npacket) *
  • BDim Block size. BDim=npacket (=gdim(feb) *
  • IntOrder: Maximum interation order. *
  • TolError: Maximum error that is tolerated. *
  • Krylov: Matrix with minimum size PsiDim*IntOrder with *
  • columns containing the Davison vectors. (complex) *
  • HKrylov: Matrix containing H*Krylov (complex) *
  • *
  • Output parameters: *
  • Psi: Propagated Psi. *
  • Krylov: Matrix with minimum size PsiDim*(IntOrder-1) with *
  • columns containing the Krylov vectors H^2|psi> to *
  • H^TrueOrder|psi>. *
  • TrueOrder: The order that has actually been used (may be less *
  • than IntOrder). *
  • ErrorCode: Error code having the following meaning: *
  • 0: everything was o. k., *
  • 1: illegal integration order, *
  • 3: diagonalization failed, *
  • *
  • Other parameters: *
  • CData: Complex*16 data needed in the Func routine. *
  • RData: Real*8 data needed in the Func routine. *
  • IData: Integer data needed in the Func routine. *
  • LData: Logical data needed in the Func routine. *
  • *
  • External routines: *
  • Func: Computes the action of the Hamiltonian on Psi. *
  • Called as *
  • Func(Abs(Time),Psi,DtPsi,CData,RData,IData,LData). *
  • WriteStepBDav: Writes the energy,beta, and error to a file. *

Called as WriteStepBDav(Steps,TrueOrder,E,beta,Error,Q,QQ,Time) *

  • *
  • Some new Variables: *
  • rlaxnum : This is the argument of the keyword relaxation. *
  • If relaxation has no argument, rlaxnum=-1 . *
  • rlaxnum=0 : relaxation to ground-state. *
  • rlaxnum=n : relaxation to n-th state (n=0...900). *
  • rlaxnum=920 : relaxation=follow. *
  • rlaxnum=930 : relaxation=lock. *
  • *
  • V8.3 03/2010 HDM *
    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! If you do some changes in bdavstep, please do them in mpibdavstep too! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

References function_library::beta(), orbital_equationofmotion::func(), hamiltonianaction_coefficients::get_hamiltoniandiagonal(), orthonc(), and function_library::zerovxd().

Referenced by integration::integrator_block_ci().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine dav_integrators::davstep ( complex*16, dimension(psidim)  Psi,
complex*16, dimension(psidim)  DtPsi,
integer  PsiDim,
integer  IntOrder,
real*8  TolError 
)

SUBROUTINE DAVSTEP *

  • Diagonalises the Hamiltonian represented in the SPFs. * Improved relaxation for Hermitian (complex) Hamiltonians. *
  • Input parameters: * Psi: The (complex) initial-value vector. * DtPsi: Action of Hamiltonian on the initial-value vector, * i.e. H|Psi> . * PsiDim Length of Psi and DtPsi vectors (i.e. adim). * IntOrder: Maximum interation order. * TolError: Maximum error that is tolerated. * Steps: Number of steps made so far (Steps is passed to * "WriteStepDav"). * Krylov: Matrix with minimum size PsiDim*IntOrder with * columns containing the Davison vectors. * HKrylov: Matrix containing H*Krylov * rlxemin: Lower energy to restrict the search for max. ovlp * rlxemax: Upper energy to restrict the search for max. ovlp *
  • Output parameters: * Psi: Propagated Psi. * Krylov: Matrix with minimum size PsiDim*(IntOrder-1) with * columns containing the Krylov vectors H^2|psi> to * H^TrueOrder|psi>. * TrueOrder: The order that has actually been used (may be less * than IntOrder). * Steps: Same as on entry. * ErrorCode: Error code having the following meaning: * 0: everything was o. k., * 1: illegal integration order, * 3: diagonalization failed, *
  • External routines: * Func: Computes the action of the Hamiltonian on Psi. * Called as * Func(Abs(Time),Psi,DtPsi). * WriteStep: Writes the stepsize and error to a file. * Called as WriteStepDav(Steps,Order,E,beta,Error,Q,t). * wrtrlxdav: Writes info to rlx_info file. *
  • Some new Variables: * rlaxnum : This is the argument of the keyword relaxation. * If relaxation has no argument, rlaxnum=-1 . * rlaxnum=0 : relaxation to ground-state. * rlaxnum=n : relaxation to n-th state (n=0...900). * rlaxnum=920 : relaxation=follow. * rlaxnum=930 : relaxation=lock. *
  • V8.3 06/2003 HDM *.

References function_library::beta(), orbital_equationofmotion::func(), and hamiltonianaction_coefficients::get_hamiltoniandiagonal().

Referenced by integration::integrator_ci().

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine dav_integrators::orthonc ( complex*16, dimension(psidim,0:n)  Krylov,
complex*16, dimension(psidim)  v,
real*8  ovl,
integer  n,
integer  PsiDim 
)

Referenced by bdavstep().

Here is the caller graph for this function:


The documentation for this module was generated from the following file: