OpenMCTDHB v2.3
Functions/Subroutines

modulegrid Module Reference

This module contains routines that return position and momentum space related arrays. More...

Functions/Subroutines

subroutine get_positionGrid (pos_i, pos_f, gdim, grid)
 This routine returns position grid in ascending order.
subroutine get_momentumGrid (xl, xu, gdim, momentumGrid)
 This routine returns momentum grid in ascending order, i.e. in the order pmin,...,pmax.
subroutine get_FFTmomentumGrid (pos_i, pos_f, gdim, FFTmomentumGrid)
 This routine returns momentum grid in fft order, i.e. in the order 0,dp,2dp,...,pmax,pmin,...,-dp.
subroutine get_p2over2m (NDX, NDY, NDZ, xi, xf, yi, yf, zi, zf, Dims, p2over2m)
 This routine computes a vector for all gridpoints containing the values of p^2/2m in FFT order This is useful for the computation of the kinetic energy.
subroutine get_ijk_from_m (m, NDX, NDY, i, j, k)
 This routine computes the indices i,j,k of the grid point r=(x_i,y_j,z_k) from the one dimenional array index m = i + NDX*(j-1) + (NDX*NDY)(k-1) that runs over all grid points. This is the order in which all grid points in 1D, 2D and 3D are stored. i,j,k all start at 1.
subroutine check_grid (Dims, NDX, NDY, NDZ)
 This routine checks for consistency of the grid with dimensions.
subroutine get_defaultOrbitals (NDX, NDY, NDZ, Dims, gridx, gridy, gridz, MOrb, Psi)
 This routine computes (unnormalized) gaussians as initial orbitals. It is useful in relaxation runs. The sigma of the gaussians is one.

Detailed Description

This module contains routines that return position and momentum space related arrays.

Authors:
Kaspar Sakmann (2007)
Since:
2007

Function/Subroutine Documentation

subroutine modulegrid::check_grid ( INTEGER,intent(in)  Dims,
INTEGER,intent(in)  NDX,
INTEGER,intent(in)  NDY,
INTEGER,intent(in)  NDZ 
)

This routine checks for consistency of the grid with dimensions.

Authors:
Kaspar Sakmann (2007)
Since:
2007
Todo:
there should be a check to see whether the DFT is an accurate enough approximation to the Fourier transform.
Parameters:
Dimsnumber of dimensions used in the current problem
NDXnumber of dimensions used in direction x
NDYnumber of dimensions used in direction y
NDZnumber of dimensions used in direction z

Definition at line 245 of file modulegrid.F90.

Here is the caller graph for this function:

subroutine modulegrid::get_defaultOrbitals ( INTEGER,intent(in)  NDX,
INTEGER,intent(in)  NDY,
INTEGER,intent(in)  NDZ,
INTEGER,intent(in)  Dims,
REAL*8,dimension(ndx),intent(in)  gridx,
REAL*8,dimension(ndy),intent(in)  gridy,
REAL*8,dimension(ndz),intent(in)  gridz,
INTEGER,intent(in)  MOrb,
COMPLEX*16,dimension(ndx*ndy*ndz,morb),intent(out)  Psi 
)

This routine computes (unnormalized) gaussians as initial orbitals. It is useful in relaxation runs. The sigma of the gaussians is one.

Authors:
Kaspar Sakmann (2007)
Since:
2007
Parameters:
NDXnumber of dimensions used in direction x
NDYnumber of dimensions used in direction y
NDZnumber of dimensions used in direction z
Dimsnumber of dimensions used in the current problem
gridxx-position grid
gridyy-position grid
gridzz-position grid
MOrbnumber of orbitals used in the current problem
Psiarray of default orbitals (gaussians)

Definition at line 299 of file modulegrid.F90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine modulegrid::get_FFTmomentumGrid ( REAL*8,intent(in)  pos_i,
REAL*8,intent(in)  pos_f,
INTEGER,intent(in)  gdim,
REAL*8,dimension(gdim),intent(out)  FFTmomentumGrid 
)

This routine returns momentum grid in fft order, i.e. in the order 0,dp,2dp,...,pmax,pmin,...,-dp.

For example: gdim=5: (0,dp,2dp,-2dp,-dp)

gdim=6: (0,dp,2dp,3dp,-2dp,-dp)

Authors:
Kaspar Sakmann (2007)
Since:
2007
Parameters:
pos_ileft border of the box
pos_fright border of the box ( =last grid point + dx)
gdimthe number of grid points in this direction
FFTmomentumGridarray of momentum space grid points in FFT order

Definition at line 102 of file modulegrid.F90.

Here is the caller graph for this function:

subroutine modulegrid::get_ijk_from_m ( INTEGER,intent(in)  m,
INTEGER,intent(in)  NDX,
INTEGER,intent(in)  NDY,
INTEGER,intent(out)  i,
INTEGER,intent(out)  j,
INTEGER,intent(out)  k 
)

This routine computes the indices i,j,k of the grid point r=(x_i,y_j,z_k) from the one dimenional array index m = i + NDX*(j-1) + (NDX*NDY)(k-1) that runs over all grid points. This is the order in which all grid points in 1D, 2D and 3D are stored. i,j,k all start at 1.

Authors:
Kaspar Sakmann (2007)
Since:
2007
Parameters:
mthe index that runs over all grid points in 1d,2d and 3d
NDXnumber of grid points in x-direction
NDYnumber of grid points in y-direction
iindex of grid point x_i = positionGridx(i)
jindex of grid point y_j = positionGridx(j)
kindex of grid point z_k = positionGridx(k)

Definition at line 215 of file modulegrid.F90.

Here is the caller graph for this function:

subroutine modulegrid::get_momentumGrid ( REAL*8,intent(in)  xl,
REAL*8,intent(in)  xu,
INTEGER,intent(in)  gdim,
REAL*8,dimension(gdim),intent(out)  momentumGrid 
)

This routine returns momentum grid in ascending order, i.e. in the order pmin,...,pmax.

For example: gdim=5: (-2dp,-dp,0,dp,2dp)

gdim=6: (-2dp,-dp,0,dp,2dp,3dp)

Authors:
Kaspar Sakmann (2007)
Since:
2007
Parameters:
xlleft border of box
xuright border of box( =last grid point + dx)
gdimnumber of grid points in this direction
momentumGridarray of momentum space grid points in ascending order

Definition at line 47 of file modulegrid.F90.

Here is the caller graph for this function:

subroutine modulegrid::get_p2over2m ( INTEGER,intent(in)  NDX,
INTEGER,intent(in)  NDY,
INTEGER,intent(in)  NDZ,
REAL*8,intent(in)  xi,
REAL*8,intent(in)  xf,
REAL*8,intent(in)  yi,
REAL*8,intent(in)  yf,
REAL*8,intent(in)  zi,
REAL*8,intent(in)  zf,
INTEGER,intent(in)  Dims,
REAL*8,dimension(ndx*ndy*ndz),intent(out)  p2over2m 
)

This routine computes a vector for all gridpoints containing the values of p^2/2m in FFT order This is useful for the computation of the kinetic energy.

Authors:
Kaspar Sakmann (2007)
Since:
2007
Parameters:
NDXnumber of grid points in x
NDYnumber of grid points in y
NDZnumber of grid points in z
xileft border of box
xfright border of box( =last grid point + dx)
yileft border of box
yfright border of box( =last grid point + dy)
zileft border of box
zfright border of box( =last grid point + dz)
Dimsnumber of dimensions
p2over2mthe array p^2/2m in FFT order

Definition at line 134 of file modulegrid.F90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine modulegrid::get_positionGrid ( REAL*8,intent(in)  pos_i,
REAL*8,intent(in)  pos_f,
INTEGER,intent(in)  gdim,
REAL*8,dimension(gdim),intent(out)  grid 
)

This routine returns position grid in ascending order.

Authors:
Kaspar Sakmann (2007)
Since:
2007
Parameters:
pos_ileft border of box
pos_fright border of box( =last grid point + dx)
gdimnumber of grid points in this direction
gridthe position grid

Definition at line 18 of file modulegrid.F90.

Here is the caller graph for this function:

 All Namespaces Files Functions Variables