MPI-SCATCI  2.0
An MPI version of SCATCI
hamiltonian_module Module Reference

Base Hamiltonian module. More...

Data Types

type  basehamiltonian
 This is an abstract class that contains the majority of functionality required to construct hamiltonians. More...
 
interface  generic_build
 Main build routine of the hamiltonian. More...
 

Functions/Subroutines

subroutine construct_base_hamiltonian (this, option, csfs, orbitals, integral)
 Constructs the hamiltonain object and provides all the various componenets required for building. More...
 
subroutine slater_rules (this, csf_a, csf_b, sym_elements, flag, job_)
 Performs Slater rules and returns symbolic elements. More...
 
subroutine handle_two_pair (this, dtrs, coeff, sym_element, flag)
 Handles if there are two pairs of spin orbitals that differ. More...
 
subroutine handle_one_pair (this, dtrs, coeff, sym_element, csf, flag)
 Handles if there are one pair of spin orbitals that differ. More...
 
subroutine handle_same (this, dtrs, coeff, sym_element, csf, flag)
 Handles if there are no differences between spin orbitals. More...
 
real(wp) function evaluate_integrals_singular (this, label, coeff, dummy_orb)
 Evaluates a single integral from labels and also checks for dummy orbitals. More...
 
real(wp) function evaluate_integrals (this, symbolic_elements, dummy_orb)
 Evaluates all integrals from a list of symbolic elements. More...
 
logical function my_job (this)
 Used by the easy MPI parallelization of slater loops. More...
 
integer function fast_slat (rep1, rep2, fin1, fin2, num_rep1, num_rep2, spin_orbs, tot_num_spin_orbs)
 ? More...
 

Detailed Description

Base Hamiltonian module.

Authors
A Al-Refaie
Date
2017

Provides an abstract hamiltonian class that has built in functionalities for all inherited hamiltonians.

Todo:
Rewrite the idiag = 0 matrix evaluation into a more easily handled form and chenge variable naming to be less confusing.
Note
30/01/2017 - Ahmed Al-Refaie: Initial documentation version
16/01/2019 - Jakub Benda: Unifom coding style and expanded documentation.

Function/Subroutine Documentation

◆ construct_base_hamiltonian()

subroutine hamiltonian_module::construct_base_hamiltonian ( class(basehamiltonian), intent(inout)  this,
class(options), intent(in), target  option,
class(csfobject), dimension(:), intent(in), target  csfs,
class(orbitaltable), intent(in), target  orbitals,
class(baseintegral), intent(in), target  integral 
)

Constructs the hamiltonain object and provides all the various componenets required for building.

Authors
A Al-Refaie
Date
2017
Parameters
[in,out]thisHamiltonian object to update.
[in]optionAn Options class containing run constants
[in]csfsAn array of configuration state functions
[in]orbitalsAn initialized Orbitals class
[in]integralA BaseIntegral class for evaluating symbols

Definition at line 121 of file Hamiltonian_module.f90.

◆ evaluate_integrals()

real(wp) function hamiltonian_module::evaluate_integrals ( class(basehamiltonian), intent(inout)  this,
class(symbolicelementvector), intent(in)  symbolic_elements,
integer, intent(in)  dummy_orb 
)

Evaluates all integrals from a list of symbolic elements.

Authors
A Al-Refaie
Date
2017
Parameters
[in,out]thisHamiltonian object to update.
[in]symbolic_elementsSymbolic elements to be evaluated.
[in]dummy_orbCheck wheter to ignore the integral if it contains a dummy orbital of this value.
Returns
Total evaluated integrals from symbols.

Definition at line 457 of file Hamiltonian_module.f90.

◆ evaluate_integrals_singular()

real(wp) function hamiltonian_module::evaluate_integrals_singular ( class(basehamiltonian), intent(inout)  this,
integer(longint), dimension(2), intent(in)  label,
real(wp), intent(in)  coeff,
integer, intent(in)  dummy_orb 
)

Evaluates a single integral from labels and also checks for dummy orbitals.

Authors
A Al-Refaie
Date
2017
Parameters
[in,out]thisHamiltonian object to update.
[in]labelA packed integral label
[in]coeffThe computed coefficients from the Slater rules
[in]dummy_orbCheck wheter to ignore the integral if it contains a dummy orbital of this value
Returns
Evaluated integral

Definition at line 413 of file Hamiltonian_module.f90.

◆ fast_slat()

integer function hamiltonian_module::fast_slat ( integer, dimension(num_rep1), intent(in)  rep1,
integer, dimension(num_rep2), intent(in)  rep2,
integer, dimension(num_rep1), intent(in)  fin1,
integer, dimension(num_rep2), intent(in)  fin2,
integer, intent(in)  num_rep1,
integer, intent(in)  num_rep2,
integer, dimension(tot_num_spin_orbs), intent(inout)  spin_orbs,
integer, intent(in)  tot_num_spin_orbs 
)

?

Authors
A Al-Refaie
Date
2017
Deprecated:
No longer used.

Definition at line 504 of file Hamiltonian_module.f90.

◆ handle_one_pair()

subroutine hamiltonian_module::handle_one_pair ( class(basehamiltonian this,
integer, dimension(4), intent(inout)  dtrs,
real(wp), intent(in)  coeff,
class(symbolicelementvector sym_element,
class(csforbital), intent(in)  csf,
integer, intent(in)  flag 
)

Handles if there are one pair of spin orbitals that differ.

Authors
A Al-Refaie
Date
2017
Parameters
[in,out]thisHamiltonian object to query.
[in]dtrsA 4 element array of slater determinants
[in]coeffThe computed coefficients from the Slater rules
[in]csfThe first CSFs determinant
[out]sym_elementResulting symbolic elements
[in]flagWhether we are doing the diagonal or not

Definition at line 258 of file Hamiltonian_module.f90.

◆ handle_same()

subroutine hamiltonian_module::handle_same ( class(basehamiltonian this,
integer, dimension(4), intent(inout)  dtrs,
real(wp), intent(in)  coeff,
class(symbolicelementvector sym_element,
class(csforbital), intent(in)  csf,
integer, intent(in)  flag 
)

Handles if there are no differences between spin orbitals.

Authors
A Al-Refaie
Date
2017
Parameters
[in,out]thisHamiltonian object to query.
[in]dtrsA 4 element array of slater determinants
[in]coeffThe computed coefficients from the Slater rules
[in]csfThe first CSFs determinant
[out]sym_elementResulting symbolic elements
[in]flagWhether we are doing the diagonal or not

Definition at line 345 of file Hamiltonian_module.f90.

◆ handle_two_pair()

subroutine hamiltonian_module::handle_two_pair ( class(basehamiltonian this,
integer, dimension(4), intent(in)  dtrs,
real(wp), intent(in)  coeff,
class(symbolicelementvector sym_element,
integer, intent(in)  flag 
)

Handles if there are two pairs of spin orbitals that differ.

Authors
A Al-Refaie
Date
2017
Parameters
[in,out]thisHamiltonian object to query.
[in]dtrsA 4 element array of slater determinants that differ.
[in]coeffThe computed coefficients from the Slater rules.
[out]sym_elementResulting symbolic elements.
[in]flagWhether we are doing the diagonal or not.

Definition at line 232 of file Hamiltonian_module.f90.

◆ my_job()

logical function hamiltonian_module::my_job ( class(basehamiltonian this)

Used by the easy MPI parallelization of slater loops.

Definition at line 479 of file Hamiltonian_module.f90.

◆ slater_rules()

subroutine hamiltonian_module::slater_rules ( class(basehamiltonian), intent(in)  this,
class(csfobject), intent(in)  csf_a,
class(csfobject), intent(in)  csf_b,
class(symbolicelementvector), intent(inout)  sym_elements,
integer, intent(in)  flag,
logical, intent(in), optional  job_ 
)

Performs Slater rules and returns symbolic elements.

Authors
A Al-Refaie
Date
2017

Computes the Slater rules and adds them to the symbolic elements class that is provided. Additionaly has an automatic MPI mode by usage of the job_ = .true. keyword. Whilst not faster than manually unrolling a loop, it provides a quick and easy way to fairly effectively parallelize the loop. However this is not OpenMP compatible

Parameters
[in,out]thisHamiltonian object to query.
[in]csf_aFirst configuration state function to slater
[in]csf_bSecond configuration state function to slater
[out]sym_elementsResulting smbolic elements
[in]flagWhether we are doing the diagonal or not
[in]job_Whether to use the 'easy' MPI split method

Definition at line 159 of file Hamiltonian_module.f90.