molecular_grid Module

MOLECULAR_GRID module to define a multicenter molecular grid system written by Son, Sang-Kil in Mar. 2014 author: Sang-Kil Son USAGE: type(Grid3D) :: G call construct_Grid3D( P, A, G )



Derived Types

type, public ::  grid3d

molecular grid information

Components

Type Visibility Attributes Name Initial
integer, public :: n_grid
integer, public :: n_nuc
real(kind=long), public, allocatable :: x(:)
real(kind=long), public, allocatable :: y(:)
real(kind=long), public, allocatable :: z(:)
real(kind=long), public, allocatable :: weight(:)
integer, public, allocatable :: i_nuc(:)
integer, public, allocatable :: i_r(:)
integer, public, allocatable :: i_ang(:)
integer, public, allocatable :: i_start(:)
integer, public, allocatable :: i_end(:)
real(kind=long), public, allocatable :: weight2c(:,:)
real(kind=long), public, allocatable :: weight_nuc(:,:)
real(kind=long), public, allocatable :: w_decomp(:)
real(kind=long), public, allocatable :: r_a(:)
real(kind=long), public, allocatable :: y_lm_prefactor(:)
logical, public :: yes_periodic
logical, public :: yes_periodic_fuzzy
real(kind=long), public :: unitcell_x
real(kind=long), public :: unitcell_y
real(kind=long), public :: unitcell_z
logical, public :: yes_cutoff
real(kind=long), public :: r_cutoff
real(kind=long), public, allocatable :: z_a(:)
real(kind=long), public, allocatable :: xx(:)
real(kind=long), public, allocatable :: yy(:)
real(kind=long), public, allocatable :: zz(:)
real(kind=long), public, allocatable :: rr(:)
real(kind=long), public, allocatable :: rp(:,:)
integer, public :: n_tr
real(kind=long), public, allocatable :: rr_tr(:,:,:)
integer, public, allocatable :: ix_tr(:,:,:)
integer, public, allocatable :: iy_tr(:,:,:)
integer, public, allocatable :: iz_tr(:,:,:)

Functions

public function integral(g, f) result(value)

integral on the grid

Arguments

Type IntentOptional Attributes Name
type(grid3d), intent(in) :: g
real(kind=long), intent(in) :: f(g%n_grid)

Return Value real(kind=long)

public function integral_single_center(g, f, i_nuc) result(value)

single center grid integral

Arguments

Type IntentOptional Attributes Name
type(grid3d), intent(in) :: g
real(kind=long), intent(in) :: f(g%n_grid)
integer, intent(in) :: i_nuc

Return Value real(kind=long)

public function integral_double_center(g, f, i_nuc, j_nuc) result(value)

double center grid integral

Arguments

Type IntentOptional Attributes Name
type(grid3d), intent(in) :: g
real(kind=long), intent(in) :: f(g%n_grid)
integer, intent(in) :: i_nuc
integer, intent(in) :: j_nuc

Return Value real(kind=long)

public function obtain_grid_position(g, i, j_nuc, k_tr) result(r)

Arguments

Type IntentOptional Attributes Name
type(grid3d), intent(in) :: g
integer, intent(in) :: i

i-th grid point

integer, intent(in) :: j_nuc

j_nuc-th atom

integer, intent(in), optional :: k_tr

k-th closest point if periodic

Return Value real(kind=long), (3)

(x, y, z) of i-th grid point

public function obtain_nuclear_position(g, i_nuc, j_nuc, k_tr) result(r)

Arguments

Type IntentOptional Attributes Name
type(grid3d), intent(in) :: g
integer, intent(in) :: i_nuc
integer, intent(in) :: j_nuc
integer, intent(in), optional :: k_tr

Return Value real(kind=long), (3)

public function bond_distance(g, i_nuc, j_nuc, k_tr) result(r_ab)

calculated distances between atoms (employing periodic boundary conditions)

Arguments

Type IntentOptional Attributes Name
type(grid3d), intent(in) :: g
integer, intent(in) :: i_nuc
integer, intent(in) :: j_nuc
integer, intent(in), optional :: k_tr

Return Value real(kind=long)

public function translation_in_crystal_fuzzy(g, x, y, z, n_fuzzy) result(ret4x8)

Arguments

Type IntentOptional Attributes Name
type(grid3d), intent(in) :: g
real(kind=long), intent(in) :: x
real(kind=long), intent(in) :: y
real(kind=long), intent(in) :: z
integer, intent(in), optional :: n_fuzzy

Return Value real(kind=long), (4,8)


Subroutines

public subroutine construct_grid3d(p, a, g)

creates the molecular Grid G

Arguments

Type IntentOptional Attributes Name
type(param), intent(in) :: p
type(atom_set), intent(in) :: a(:)
type(grid3d), intent(out) :: g

public subroutine reconstruct_grid3d(a, g)

Read more…

Arguments

Type IntentOptional Attributes Name
type(atom_set), intent(in) :: a(:)
type(grid3d), intent(inout) :: g

public subroutine update_grid3d_with_displacement(a, g, displacement)

update Grid3D with displacement use this subroutine if only molecular geometry in A has been changed. note that A should be already updated with update_Atom_set_with_displacement() in atom_in_mol.f90.

Read more…

Arguments

Type IntentOptional Attributes Name
type(atom_set), intent(in) :: a(:)
type(grid3d), intent(inout) :: g
real(kind=long), intent(in) :: displacement(:)

public subroutine place_ag_on_grid3d(a, g, i_start)

Read more…

Arguments

Type IntentOptional Attributes Name
type(atom_set), intent(in) :: a
type(grid3d), intent(inout) :: g
integer, intent(in) :: i_start

public subroutine calculate_weight_grid3d(a, g)

calculate all integration weights on Grid3D and

Arguments

Type IntentOptional Attributes Name
type(atom_set), intent(in) :: a(:)
type(grid3d), intent(inout) :: g

public subroutine purge_grid3d(g)

Arguments

Type IntentOptional Attributes Name
type(grid3d), intent(inout) :: g

public subroutine translation_in_crystal(g, x, y, z)

Arguments

Type IntentOptional Attributes Name
type(grid3d), intent(in) :: g
real(kind=long), intent(inout) :: x
real(kind=long), intent(inout) :: y
real(kind=long), intent(inout) :: z

public subroutine find_neighboring_atoms(g, a, output)

count how many neighboring atoms for each atom and find the bond distances

Arguments

Type IntentOptional Attributes Name
type(grid3d), intent(in) :: g
type(atom_set), intent(in) :: a(:)
integer :: output