OpenMCTDHB v2.3
Functions/Subroutines

modulecorrelationfunctions Module Reference

This module contains routines to compute correlation function and writes them to files, e.g. the one- and two particle densities and momentum distributions. More...

Functions/Subroutines

subroutine write_correlations (Dims, printCorrelationFuncs, MOMSPACE2D, REALSPACE2D, x1const, x1slice, y1const, y1slice, x2const, x2slice, y2const, y2slice, AbsTime, Psi, rho_jk, rho_ijkl, Pnot, xstart, xstop)
 wrapper routine which manages output of correlations. compute the first and second order correlation functions in 1D calculations and the first and second order correlation functions in 2D calculations if desired.
subroutine get_integratedDensity (rho_jk, AbsTime, PSI, writemode)
 This routine integrates the density over one or more cartesian coordinates and writes the output to a file. Useful for visualizing 2D and 3D densities Which coordinates are integrated out is controlled by the integer writemode=0,...,7 the bit representation of writemode indicates which of x,y,z are integrated out, e.g. writemode=2 is interpreted as 010 and the y-coordinate is integrated out.
subroutine get_densityNonescape (xstart, xstop, rho_jk, time, PSI)
 This routine computes the nonescape probability or any desired integral of the diagonal of the reduced one-body density matrix.
subroutine get_maximum (Psi, rho_jk, maximum)
 This routine computes the maximum of either the one-particle density or the momentum distribution depending on whether th array of orbitals Psi or its Fourier transform FTPsi is given as the first argument.
subroutine get_oneBodyRDMDiagonal (time, Psi, rho_jk, mode, dilation)
 This routine computes the diagonal of the first order correlation function G1(x,x) if Psi is provided or G1(k,k) if FT(Psi) is provided and writes it to a file.
subroutine get_correlations (time, Psi, mode, dilation, rho_jk, rho_ijkl)
 This routine computes the correlation functions G1(r,r), G1(r,r'), G1(r',r') and G2(r,r',r',r), where r=(x,y,z), and writes them to a file. The first six columns contain x,y,z,x',y',z'. Currently only 1D has been tested.
subroutine get_dilation (dilation, Psi, rho_jk)
 This routine was written to deal with the following problem. If the x-domain has length L, then the k-domain has a grid spacing dp=2PI/L. It often turns out that this spacing is very wide. However, if the wave function in x-space is negligibly small at the box boundaries, it can be assumed to be zero outside the box boundaries. Thereby, the length L can be artificially increased and the momentum distribution can be evaluated using a smaller grid. The factor dilation is integer and describes how many times the k-space interval can be made smaller so that still 0.99% fit into it, using the same number of grid points. Using get_dilation allows to make nice plots of the momentum distributions. It should only be used to generate nicer output though. If this is unwanted please set dilation = 1 everywhere the variable appears. This is the default.
subroutine pedestrian_FT (Psi, FTPsi, dilation)
 This routine computes the Fourier integral FTPsi of the orbital array Psi and normalizes it to one. If dilation=1 the FFT grid is used, i.e, the spacing is dp=2pi/L where L is the length of the grid in x-space. Otherwise the integral is evaluated using te spacing 2pi/(dilation*L). Note that this is only OK, if it is clear for some physical reason (e.g. a trap potential) that the wave function is zero outside the intervall of length L.
subroutine get_corr_slice (time, PSI, mode, rho_jk, rho_ijkl, x1const, x1slice, y1const, y1slice, x2const, x2slice, y2const, y2slice)
 computes cuts of the correlation functions G1(r,r), G1(r,r'), G1(r',r') and G2(r,r',r',r), where r=(x,y), and writes them to a file.
subroutine FT_SORTED (PSI, FTPSI)
 -------------------------------------------------------------------------------------------- This subroutine uses the FFTW/ACML routines to Fourier-transform all the Orbitals and returns them in strictly ascending order for analysis purposes.
subroutine sort_FFT_to_ascending_2d (FFT, NDX, NDY)
 ------------------------------------------------------------------------------- sort an array representing a 2D object given in FFT order, compact storage to be in strictly ascending order (needed for gnuplots pm3d), compact storage
subroutine getsliceAsString (d, x, doubleString)
 ------------------------------------------------------------------------------- This subroutine is called to generate the names of the output files of the slices through the 2D correlation functions.

Detailed Description

This module contains routines to compute correlation function and writes them to files, e.g. the one- and two particle densities and momentum distributions.

Authors:
Kaspar Sakmann (2007), Axel Lode (2011)
Since:
2007

Function/Subroutine Documentation

subroutine modulecorrelationfunctions::FT_SORTED ( complex*16,dimension(ndx*ndy*ndz,morb),intent(in)  PSI,
complex*16,dimension(ndx*ndy*ndz,morb),intent(out)  FTPSI 
)

-------------------------------------------------------------------------------------------- This subroutine uses the FFTW/ACML routines to Fourier-transform all the Orbitals and returns them in strictly ascending order for analysis purposes.

Authors:
Axel Lode (2011)
Since:
2011
Todo:
Probably move this routine to the FFT module

Definition at line 1440 of file modulecorrelationfunctions.F90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine modulecorrelationfunctions::get_corr_slice ( real*8,intent(in)  time,
COMPLEX*16,dimension(ndx*ndy*ndz,morb)  PSI,
integer,intent(in)  mode,
complex*16,dimension(morb,morb),intent(in)  rho_jk,
complex*16,dimension(morb,morb,morb,morb),intent(in)  rho_ijkl,
logical,intent(in)  x1const,
real*8,intent(in)  x1slice,
logical,intent(in)  y1const,
real*8,intent(in)  y1slice,
logical,intent(in)  x2const,
real*8,intent(in)  x2slice,
logical,intent(in)  y2const,
real*8,intent(in)  y2slice 
)

computes cuts of the correlation functions G1(r,r), G1(r,r'), G1(r',r') and G2(r,r',r',r), where r=(x,y), and writes them to a file.

the user has to choose two of the coordinates to stay constant, to obtain a 3D plot. x1,y1,x2,y2.

time : current time PSI : orbital array mode : mode=1: x-space (input PSI=PSI) mode=2: k-space (input PSI=FT[PSI]) rho_jk : the 1-body density matrix where rho_jk(j,k) = rho_jk rho_ijkl : the 2-body density matrix where rho_ijkl(j,k) = rho_ijkl

Authors:
Axel Lode (2011), Kaspar Sakmann (2007)
Todo:
variable documentation, variable cleaning

Definition at line 942 of file modulecorrelationfunctions.F90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine modulecorrelationfunctions::get_correlations ( real*8,intent(in)  time,
complex*16,dimension(ndx*ndy*ndz,morb),intent(in)  Psi,
integer,intent(in)  mode,
integer,intent(in)  dilation,
complex*16,dimension(morb,morb),intent(in)  rho_jk,
complex*16,dimension(morb,morb,morb,morb),intent(in)  rho_ijkl 
)

This routine computes the correlation functions G1(r,r), G1(r,r'), G1(r',r') and G2(r,r',r',r), where r=(x,y,z), and writes them to a file. The first six columns contain x,y,z,x',y',z'. Currently only 1D has been tested.

time : current time Psi : orbital array mode : mode=1: x-space (input Psi=Psi) mode=2: k-space (input Psi=FT[Psi]) dilation : the integer dilation gives the factor by which the momentum intervall should be shrunk or equivalently the factor by which the sampling rate is increased. The momentum distribution at the outer ends is padded with zeros then. Normally dilation is computed by get_dilation. If unsure put dilation = 1, but the momentum grid will then usually have a bad resolution. rho_jk : the 1-body density matrix where rho_jk(j,k) = rho_jk rho_ijkl : the 2-body density matrix where rho_ijkl(j,k) = rho_ijkl

Authors:
Kaspar Sakmann (2007)
Since:
2007
Todo:
: reduce use associations, variable documentation

Definition at line 588 of file modulecorrelationfunctions.F90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine modulecorrelationfunctions::get_densityNonescape ( real*8,intent(in)  xstart,
real*8,intent(in)  xstop,
complex*16,dimension(morb,morb),intent(in)  rho_jk,
real*8  time,
COMPLEX*16,dimension(ndx*ndy*ndz,morb)  PSI 
)

This routine computes the nonescape probability or any desired integral of the diagonal of the reduced one-body density matrix.

Authors:
Axel Lode
Since:
(2011)
Todo:
generalize to more than 1D

Definition at line 314 of file modulecorrelationfunctions.F90.

Here is the caller graph for this function:

subroutine modulecorrelationfunctions::get_dilation ( integer,intent(inout)  dilation,
complex*16,dimension(ndx*ndy*ndz,morb),intent(in)  Psi,
complex*16,dimension(morb,morb),intent(in)  rho_jk 
)

This routine was written to deal with the following problem. If the x-domain has length L, then the k-domain has a grid spacing dp=2PI/L. It often turns out that this spacing is very wide. However, if the wave function in x-space is negligibly small at the box boundaries, it can be assumed to be zero outside the box boundaries. Thereby, the length L can be artificially increased and the momentum distribution can be evaluated using a smaller grid. The factor dilation is integer and describes how many times the k-space interval can be made smaller so that still 0.99% fit into it, using the same number of grid points. Using get_dilation allows to make nice plots of the momentum distributions. It should only be used to generate nicer output though. If this is unwanted please set dilation = 1 everywhere the variable appears. This is the default.

Authors:
Kaspar Sakmann (2007)
Since:
2007
Todo:
get rid of use associations, variable documentation

Definition at line 794 of file modulecorrelationfunctions.F90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine modulecorrelationfunctions::get_integratedDensity ( complex*16,dimension(morb,morb),intent(in)  rho_jk,
real*8,intent(in)  AbsTime,
COMPLEX*16,dimension(ndx*ndy*ndz,morb),intent(in)  PSI,
INTEGER,intent(in)  writemode 
)

This routine integrates the density over one or more cartesian coordinates and writes the output to a file. Useful for visualizing 2D and 3D densities Which coordinates are integrated out is controlled by the integer writemode=0,...,7 the bit representation of writemode indicates which of x,y,z are integrated out, e.g. writemode=2 is interpreted as 010 and the y-coordinate is integrated out.

Authors:
Kaspar Sakmann
Since:
(2012)

Definition at line 83 of file modulecorrelationfunctions.F90.

Here is the caller graph for this function:

subroutine modulecorrelationfunctions::get_maximum ( complex*16,dimension(ndx*ndy*ndz,morb),intent(in)  Psi,
complex*16,dimension(morb,morb),intent(in)  rho_jk,
complex*16,intent(out)  maximum 
)

This routine computes the maximum of either the one-particle density or the momentum distribution depending on whether th array of orbitals Psi or its Fourier transform FTPsi is given as the first argument.

Authors:
Kaspar Sakmann (2007)
Since:
2007
Todo:
make all use associated variables arguments of the routine, variable documentation

Definition at line 372 of file modulecorrelationfunctions.F90.

Here is the caller graph for this function:

subroutine modulecorrelationfunctions::get_oneBodyRDMDiagonal ( real*8  time,
complex*16,dimension(ndx*ndy*ndz,morb),intent(in)  Psi,
complex*16,dimension(morb,morb),intent(in)  rho_jk,
integer,intent(in)  mode,
integer,intent(in)  dilation 
)

This routine computes the diagonal of the first order correlation function G1(x,x) if Psi is provided or G1(k,k) if FT(Psi) is provided and writes it to a file.

time : current time Psi : orbital array rho_jk : the 1-body density matrix where rho_jk(j,k)=rho_jk mode : mode=1: G1(x,x) mode=2: G1(k,k) is computed dilation : the integer dilation gives the factor by which the momentum intervall should be shrunk or equivalently the factor by which the sampling rate is increased. The momentum distribution at the outer ends is padded with zeros then. Normally dilation is computed by get_dilation. If unsure put dilation = 1, but the momentum grid will then usually have a bad resolution.

Authors:
Kaspar Sakmann (2007)
Since:
2007
Todo:
reduce the use associations as much as possible, variable documentation

Definition at line 426 of file modulecorrelationfunctions.F90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine modulecorrelationfunctions::getsliceAsString ( integer  d,
real(kind=DBL)  x,
CHARACTER(len=7)  doubleString 
)

------------------------------------------------------------------------------- This subroutine is called to generate the names of the output files of the slices through the 2D correlation functions.

Authors:
Axel Lode (2011)
Since:
2011
Todo:
Probably move this routine to the Input/Output module

Definition at line 1598 of file modulecorrelationfunctions.F90.

subroutine modulecorrelationfunctions::pedestrian_FT ( complex*16,dimension(ndx*ndy*ndz,morb),intent(in)  Psi,
complex*16,dimension(ndx*ndy*ndz,morb),intent(out)  FTPsi,
integer,intent(in)  dilation 
)

This routine computes the Fourier integral FTPsi of the orbital array Psi and normalizes it to one. If dilation=1 the FFT grid is used, i.e, the spacing is dp=2pi/L where L is the length of the grid in x-space. Otherwise the integral is evaluated using te spacing 2pi/(dilation*L). Note that this is only OK, if it is clear for some physical reason (e.g. a trap potential) that the wave function is zero outside the intervall of length L.

Authors:
Kaspar Sakmann (2007)
Since:
2007
Todo:
variable documentation

Definition at line 872 of file modulecorrelationfunctions.F90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine modulecorrelationfunctions::sort_FFT_to_ascending_2d ( COMPLEX*16,dimension(ndx*ndy)  FFT,
integer  NDX,
integer  NDY 
)

------------------------------------------------------------------------------- sort an array representing a 2D object given in FFT order, compact storage to be in strictly ascending order (needed for gnuplots pm3d), compact storage

Authors:
Axel Lode (2011)
Since:
2011
Todo:
Probably move this routine to the FFT module

Definition at line 1533 of file modulecorrelationfunctions.F90.

Here is the caller graph for this function:

subroutine modulecorrelationfunctions::write_correlations ( INTEGER,intent(in)  Dims,
INTEGER,intent(in)  printCorrelationFuncs,
LOGICAL,intent(in)  MOMSPACE2D,
LOGICAL,intent(in)  REALSPACE2D,
LOGICAL,intent(in)  x1const,
REAL*8,intent(in)  x1slice,
LOGICAL,intent(in)  y1const,
REAL*8,intent(in)  y1slice,
LOGICAL,intent(in)  x2const,
REAL*8,intent(in)  x2slice,
LOGICAL,intent(in)  y2const,
REAL*8,intent(in)  y2slice,
REAL*8,intent(in)  AbsTime,
COMPLEX*16,dimension(ndx*ndy*ndz,morb),intent(in)  Psi,
COMPLEX*16,dimension(morb*morb),intent(in)  rho_jk,
COMPLEX*16,dimension(morb*morb*morb*morb),intent(in)  rho_ijkl,
LOGICAL,intent(in)  Pnot,
REAL*8,intent(in)  xstart,
REAL*8,intent(in)  xstop 
)

wrapper routine which manages output of correlations. compute the first and second order correlation functions in 1D calculations and the first and second order correlation functions in 2D calculations if desired.

Authors:
Axel Lode, Kaspar Sakmann
Since:
(2011)

Definition at line 21 of file modulecorrelationfunctions.F90.

Here is the call graph for this function:

Here is the caller graph for this function:

 All Namespaces Files Functions Variables