(git:6a2e663)
qs_collocate_density Module Reference

Calculate the plane wave density by collocating the primitive Gaussian functions (pgf). More...

Functions/Subroutines

subroutine, public calculate_rho_nlcc (rho_nlcc, qs_env)
 computes the density of the non-linear core correction on the grid More...
 
subroutine, public calculate_ppl_grid (vppl, qs_env)
 computes the local pseudopotential (without erf term) on the grid More...
 
subroutine, public calculate_lri_rho_elec (lri_rho_g, lri_rho_r, qs_env, lri_coef, total_rho, basis_type, exact_1c_terms, pmat, atomlist)
 Collocates the fitted lri density on a grid. More...
 
subroutine, public calculate_drho_core (drho_core, qs_env, beta, lambda)
 Computes the derivative of the density of the core charges with respect to the nuclear coordinates on the grid. More...
 
subroutine, public calculate_rho_single_gaussian (rho_gb, qs_env, iatom_in)
 collocate a single Gaussian on the grid More...
 
subroutine, public calculate_rho_metal (rho_metal, coeff, total_rho_metal, qs_env)
 computes the image charge density on the grid (including coeffcients) More...
 
subroutine, public calculate_rho_resp_single (rho_gb, qs_env, eta, iatom_in)
 collocate a single Gaussian on the grid for periodic RESP fitting More...
 
subroutine, public calculate_rho_elec (matrix_p, matrix_p_kp, rho, rho_gspace, total_rho, ks_env, soft_valid, compute_tau, compute_grad, basis_type, der_type, idir, task_list_external, pw_env_external)
 computes the density corresponding to a given density matrix on the grid More...
 
subroutine, public calculate_drho_elec (matrix_p, matrix_p_kp, drho, drho_gspace, qs_env, soft_valid, basis_type)
 computes the gradient of the density corresponding to a given density matrix on the grid More...
 
subroutine, public calculate_drho_elec_dr (matrix_p, matrix_p_kp, drho, drho_gspace, qs_env, soft_valid, basis_type, beta, lambda)
 Computes the gradient wrt. nuclear coordinates of a density on the grid The density is given in terms of the density matrix_p. More...
 
subroutine, public collocate_single_gaussian (rho, rho_gspace, atomic_kind_set, qs_kind_set, cell, dft_control, particle_set, pw_env, required_function, basis_type)
 maps a single gaussian on the grid More...
 
subroutine, public calculate_wavefunction (mo_vectors, ivector, rho, rho_gspace, atomic_kind_set, qs_kind_set, cell, dft_control, particle_set, pw_env, basis_type, external_vector)
 maps a given wavefunction on the grid More...
 

Detailed Description

Calculate the plane wave density by collocating the primitive Gaussian functions (pgf).

History
  • rewrote collocate for increased accuracy and speed
  • introduced the PGI hack for increased speed with that compiler (22.02.02)
  • Added Multiple Grid feature
  • new way to get over the grid (01.03.02)
  • removed timing calls since they were getting expensive
  • Updated with the new QS data structures (09.04.02,MK)
  • introduction of the real space grid type ( prelim. version JVdV 05.02)
  • parallel FFT (JGH 22.05.02)
  • multigrid arrays independent from density (JGH 30.08.02)
  • old density stored in g space (JGH 30.08.02)
  • distributed real space code (JGH 17.07.03)
  • refactoring and new loop ordering (JGH 23.11.03)
  • OpenMP parallelization (JGH 03.12.03)
  • Modified to compute tau (Joost 12.03)
  • removed the incremental density rebuild (Joost 01.04)
  • introduced realspace multigridding (Joost 02.04)
  • introduced map_consistent (Joost 02.04)
  • Addition of the subroutine calculate_atomic_charge_density (TdK, 08.05)
  • rewrite of the collocate/integrate kernels (Joost VandeVondele, 03.07)
  • Extended by the derivatives for DFPT [Sandra Luber, Edward Ditler, 2021]
Author
Matthias Krack (03.04.2001) 1) Joost VandeVondele (01.2002) Thomas D. Kuehne (04.08.2005) Ole Schuett (2020)

Function/Subroutine Documentation

◆ calculate_rho_nlcc()

subroutine, public qs_collocate_density::calculate_rho_nlcc ( type(pw_r3d_rs_type), intent(inout)  rho_nlcc,
type(qs_environment_type), pointer  qs_env 
)

computes the density of the non-linear core correction on the grid

Parameters
rho_nlcc...
qs_env...

Definition at line 153 of file qs_collocate_density.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ calculate_ppl_grid()

subroutine, public qs_collocate_density::calculate_ppl_grid ( type(pw_r3d_rs_type), intent(inout)  vppl,
type(qs_environment_type), pointer  qs_env 
)

computes the local pseudopotential (without erf term) on the grid

Parameters
vppl...
qs_env...

Definition at line 332 of file qs_collocate_density.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ calculate_lri_rho_elec()

subroutine, public qs_collocate_density::calculate_lri_rho_elec ( type(pw_c1d_gs_type), intent(inout)  lri_rho_g,
type(pw_r3d_rs_type), intent(inout)  lri_rho_r,
type(qs_environment_type), pointer  qs_env,
type(lri_kind_type), dimension(:), pointer  lri_coef,
real(kind=dp), intent(out)  total_rho,
character(len=*), intent(in)  basis_type,
logical, intent(in)  exact_1c_terms,
type(dbcsr_type), optional  pmat,
integer, dimension(:), optional  atomlist 
)

Collocates the fitted lri density on a grid.

Parameters
lri_rho_g...
lri_rho_r...
qs_env...
lri_coef...
total_rho...
basis_type...
exact_1c_terms...
pmatreplicated block diagonal density matrix (optional)
atomlistlist of atoms to be included (optional)
History
04.2013
Author
Dorothea Golze

Definition at line 515 of file qs_collocate_density.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ calculate_drho_core()

subroutine, public qs_collocate_density::calculate_drho_core ( type(pw_c1d_gs_type), intent(inout)  drho_core,
type(qs_environment_type), pointer  qs_env,
integer, intent(in)  beta,
integer, intent(in)  lambda 
)

Computes the derivative of the density of the core charges with respect to the nuclear coordinates on the grid.

Parameters
drho_coreThe resulting density derivative
qs_env...
betaDerivative direction
lambdaAtom index
Note
SL November 2014, ED 2021

Definition at line 1069 of file qs_collocate_density.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ calculate_rho_single_gaussian()

subroutine, public qs_collocate_density::calculate_rho_single_gaussian ( type(pw_c1d_gs_type), intent(inout)  rho_gb,
type(qs_environment_type), pointer  qs_env,
integer, intent(in)  iatom_in 
)

collocate a single Gaussian on the grid

Parameters
rho_gbcharge density generated by a single gaussian
qs_envqs environment
iatom_inatom index
History
12.2011 created
Author
Dorothea Golze

Definition at line 1201 of file qs_collocate_density.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ calculate_rho_metal()

subroutine, public qs_collocate_density::calculate_rho_metal ( type(pw_c1d_gs_type), intent(inout)  rho_metal,
real(kind=dp), dimension(:), pointer  coeff,
real(kind=dp), intent(out), optional  total_rho_metal,
type(qs_environment_type), pointer  qs_env 
)

computes the image charge density on the grid (including coeffcients)

Parameters
rho_metalimage charge density
coeffexpansion coefficients of the image charge density, i.e. rho_metal=sum_a c_a*g_a
total_rho_metaltotal induced image charge density
qs_envqs environment
History
01.2012 created
Author
Dorothea Golze

Definition at line 1291 of file qs_collocate_density.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ calculate_rho_resp_single()

subroutine, public qs_collocate_density::calculate_rho_resp_single ( type(pw_c1d_gs_type), intent(inout)  rho_gb,
type(qs_environment_type), pointer  qs_env,
real(kind=dp), intent(in)  eta,
integer, intent(in)  iatom_in 
)

collocate a single Gaussian on the grid for periodic RESP fitting

Parameters
rho_gbcharge density generated by a single gaussian
qs_envqs environment
etawidth of single Gaussian
iatom_inatom index
History
06.2012 created
Author
Dorothea Golze

Definition at line 1397 of file qs_collocate_density.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ calculate_rho_elec()

subroutine, public qs_collocate_density::calculate_rho_elec ( type(dbcsr_type), optional, target  matrix_p,
type(dbcsr_p_type), dimension(:), optional, pointer  matrix_p_kp,
type(pw_r3d_rs_type), intent(inout)  rho,
type(pw_c1d_gs_type), intent(inout)  rho_gspace,
real(kind=dp), intent(out), optional  total_rho,
type(qs_ks_env_type), pointer  ks_env,
logical, intent(in), optional  soft_valid,
logical, intent(in), optional  compute_tau,
logical, intent(in), optional  compute_grad,
character(len=*), intent(in), optional  basis_type,
integer, intent(in), optional  der_type,
integer, intent(in), optional  idir,
type(task_list_type), optional, pointer  task_list_external,
type(pw_env_type), optional, pointer  pw_env_external 
)

computes the density corresponding to a given density matrix on the grid

Parameters
matrix_p...
matrix_p_kp...
rho...
rho_gspace...
total_rho...
ks_env...
soft_valid...
compute_tau...
compute_grad...
basis_type...
der_type...
idir...
task_list_external...
pw_env_external...
History
IAB (15-Feb-2010): Added OpenMP parallelisation to task loop (c) The Numerical Algorithms Group (NAG) Ltd, 2010 on behalf of the HECToR project Anything that is not the default ORB basis_type requires an external_task_list 12.2019, (A.Bussy) Ole Schuett (2020): Migrated to C, see grid_api.F
Note
both rho and rho_gspace contain the new rho (in real and g-space respectively)

Definition at line 1707 of file qs_collocate_density.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ calculate_drho_elec()

subroutine, public qs_collocate_density::calculate_drho_elec ( type(dbcsr_type), optional, target  matrix_p,
type(dbcsr_p_type), dimension(:), optional, pointer  matrix_p_kp,
type(pw_r3d_rs_type), dimension(3), intent(inout)  drho,
type(pw_c1d_gs_type), dimension(3), intent(inout)  drho_gspace,
type(qs_environment_type), pointer  qs_env,
logical, intent(in), optional  soft_valid,
character(len=*), intent(in), optional  basis_type 
)

computes the gradient of the density corresponding to a given density matrix on the grid

Parameters
matrix_p...
matrix_p_kp...
drho...
drho_gspace...
qs_env...
soft_valid...
basis_type...
Note
this is an alternative to calculate the gradient through FFTs

Definition at line 1874 of file qs_collocate_density.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ calculate_drho_elec_dr()

subroutine, public qs_collocate_density::calculate_drho_elec_dr ( type(dbcsr_type), optional, target  matrix_p,
type(dbcsr_p_type), dimension(:), optional, pointer  matrix_p_kp,
type(pw_r3d_rs_type), intent(inout)  drho,
type(pw_c1d_gs_type), intent(inout)  drho_gspace,
type(qs_environment_type), pointer  qs_env,
logical, intent(in), optional  soft_valid,
character(len=*), intent(in), optional  basis_type,
integer, intent(in)  beta,
integer, intent(in)  lambda 
)

Computes the gradient wrt. nuclear coordinates of a density on the grid The density is given in terms of the density matrix_p.

Parameters
matrix_pDensity matrix
matrix_p_kp...
drhoDensity gradient on the grid
drho_gspaceDensity gradient on the reciprocal grid
qs_env...
soft_valid...
basis_type...
betaDerivative direction
lambdaAtom index
Note
SL, ED 2021 Adapted from calculate_drho_elec

Definition at line 2261 of file qs_collocate_density.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ collocate_single_gaussian()

subroutine, public qs_collocate_density::collocate_single_gaussian ( type(pw_r3d_rs_type), intent(inout)  rho,
type(pw_c1d_gs_type), intent(inout)  rho_gspace,
type(atomic_kind_type), dimension(:), pointer  atomic_kind_set,
type(qs_kind_type), dimension(:), pointer  qs_kind_set,
type(cell_type), pointer  cell,
type(dft_control_type), pointer  dft_control,
type(particle_type), dimension(:), pointer  particle_set,
type(pw_env_type), pointer  pw_env,
integer, intent(in)  required_function,
character(len=*), intent(in), optional  basis_type 
)

maps a single gaussian on the grid

Parameters
rho...
rho_gspace...
atomic_kind_set...
qs_kind_set...
cell...
dft_control...
particle_set...
pw_env...
required_function...
basis_type...
History
08.2022 created from calculate_wavefunction
Note
modified calculate_wave function assuming that the collocation of only a single Gaussian is required. chooses a basis function (in contrast to calculate_rho_core or calculate_rho_single_gaussian)

Definition at line 2667 of file qs_collocate_density.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ calculate_wavefunction()

subroutine, public qs_collocate_density::calculate_wavefunction ( type(cp_fm_type), intent(in)  mo_vectors,
integer, intent(in)  ivector,
type(pw_r3d_rs_type), intent(inout)  rho,
type(pw_c1d_gs_type), intent(inout)  rho_gspace,
type(atomic_kind_type), dimension(:), pointer  atomic_kind_set,
type(qs_kind_type), dimension(:), pointer  qs_kind_set,
type(cell_type), pointer  cell,
type(dft_control_type), pointer  dft_control,
type(particle_type), dimension(:), pointer  particle_set,
type(pw_env_type), pointer  pw_env,
character(len=*), intent(in), optional  basis_type,
real(kind=dp), dimension(:), optional  external_vector 
)

maps a given wavefunction on the grid

Parameters
mo_vectors...
ivector...
rho...
rho_gspace...
atomic_kind_set...
qs_kind_set...
cell...
dft_control...
particle_set...
pw_env...
basis_type...
external_vector...
History
08.2002 created [Joost VandeVondele] 03.2006 made independent of qs_env [Joost VandeVondele]
Note
modified calculate_rho_elec, should write the wavefunction represented by it's presumably dominated by the FFT and the rs->pw and back routines

BUGS ??? does it take correctly the periodic images of the basis functions into account i.e. is only correct if the basis functions are localised enough to be just in 1 cell ?

Definition at line 2857 of file qs_collocate_density.F.

Here is the call graph for this function:
Here is the caller graph for this function: