MPI-SCATCI 2.0
An MPI version of SCATCI
Loading...
Searching...
No Matches
Hamiltonian_module::BaseHamiltonian Type Referenceabstract

This is an abstract class that contains the majority of functionality required to construct hamiltonians. More...

Public Member Functions

procedure, public construct (this, option, csfs, orbitals, integral)
 Constructs the hamiltonain object and provides all the various componenets required for building.
procedure(generic_build), deferred build_hamiltonian (this, matrix_elements)
procedure, public slater_rules (this, csf_a, csf_b, sym_elements, flag, job_)
 Performs Slater rules and returns symbolic elements.
procedure, public evaluate_integrals (this, symbolic_elements, dummy_orb)
 Evaluates all integrals from a list of symbolic elements.
procedure, public evaluate_integrals_singular (this, label, coeff, dummy_orb)
 Evaluates a single integral from labels and also checks for dummy orbitals.
procedure, public my_job (this)
 Used by the easy MPI parallelization of slater loops.

Public Attributes

class(orbitaltable), pointer orbitals
 Our orbitals required to generate symblic elements.
class(options), pointer options
 Scatci program settings.
class(baseintegral), pointer integral
 The integrals we are using.
class(csfobject), dimension(:), pointer csfs
 Our configuration state functions.
integer nflg = 0
integer diagonal_flag
integer positron_flag
 Positron aware flag.
integer phase_flag
 whether to evaluate integrals whilst dealing with phase
logical constructed = .false.
 Has the hamiltonain been constructed.
logical initialized = .false.
 Has the hamiltonian been initialized.
integer job_id = 0
 Whose job it is to (soon to be deprecated)
integer number_of_integrals = 0
 How many integrals have been evaluated?
real(wp) element_one = 0.0
 First element for idiag = 0.
real(wp) enh_factor
 Enhancement factor set in options.
integer enh_factor_type
 Enhancement factor type set in options.
type(symbolicelementvector) reference_symbol
 Symbols for idiag = 0.

Private Member Functions

procedure, private handle_two_pair (this, dtrs, coeff, sym_element, flag)
 Handles if there are two pairs of spin orbitals that differ.
procedure, private handle_one_pair (this, dtrs, coeff, sym_element, csf, flag)
 Handles if there are one pair of spin orbitals that differ.
procedure, private handle_same (this, dtrs, coeff, sym_element, csf, flag)
 Handles if there are no differences between spin orbitals.

Detailed Description

This is an abstract class that contains the majority of functionality required to construct hamiltonians.

Authors
A Al-Refaie
Date
2017

This class provides an abstraction of a lot of the components required to build the hamiltonian. For example, the user does not need to worry about the specific implementation of Slaters rules or evaluating integrals this class wraps these features for you and should allow one to implement hamiltonians closer to those described in papers. This is not a Matrix class. BaseHamiltonain deals with building the matrix whilst BaseMatrix deals with storing the matrix.

Definition at line 58 of file Hamiltonian_module.f90.

Member Function/Subroutine Documentation

◆ build_hamiltonian()

procedure(generic_build), deferred Hamiltonian_module::BaseHamiltonian::build_hamiltonian ( class(basehamiltonian) this,
class(basematrix), intent(inout) matrix_elements )
pure virtual

Definition at line 80 of file Hamiltonian_module.f90.

◆ construct()

procedure, public Hamiltonian_module::BaseHamiltonian::construct ( 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 79 of file Hamiltonian_module.f90.

◆ evaluate_integrals()

procedure, public Hamiltonian_module::BaseHamiltonian::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 85 of file Hamiltonian_module.f90.

◆ evaluate_integrals_singular()

procedure, public Hamiltonian_module::BaseHamiltonian::evaluate_integrals_singular ( class(basehamiltonian), intent(inout) this,
integer(longint), dimension(nidx), 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 86 of file Hamiltonian_module.f90.

◆ handle_one_pair()

procedure, private Hamiltonian_module::BaseHamiltonian::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 )
private

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 91 of file Hamiltonian_module.f90.

◆ handle_same()

procedure, private Hamiltonian_module::BaseHamiltonian::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 )
private

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 92 of file Hamiltonian_module.f90.

◆ handle_two_pair()

procedure, private Hamiltonian_module::BaseHamiltonian::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 )
private

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 90 of file Hamiltonian_module.f90.

◆ my_job()

procedure, public Hamiltonian_module::BaseHamiltonian::my_job ( class(basehamiltonian) this)

Used by the easy MPI parallelization of slater loops.

Definition at line 93 of file Hamiltonian_module.f90.

◆ slater_rules()

procedure, public Hamiltonian_module::BaseHamiltonian::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 83 of file Hamiltonian_module.f90.

Member Data Documentation

◆ constructed

logical Hamiltonian_module::BaseHamiltonian::constructed = .false.

Has the hamiltonain been constructed.

Definition at line 68 of file Hamiltonian_module.f90.

◆ csfs

class(csfobject), dimension(:), pointer Hamiltonian_module::BaseHamiltonian::csfs

Our configuration state functions.

Definition at line 62 of file Hamiltonian_module.f90.

◆ diagonal_flag

integer Hamiltonian_module::BaseHamiltonian::diagonal_flag

Definition at line 65 of file Hamiltonian_module.f90.

◆ element_one

real(wp) Hamiltonian_module::BaseHamiltonian::element_one = 0.0

First element for idiag = 0.

Definition at line 72 of file Hamiltonian_module.f90.

◆ enh_factor

real(wp) Hamiltonian_module::BaseHamiltonian::enh_factor

Enhancement factor set in options.

Definition at line 73 of file Hamiltonian_module.f90.

◆ enh_factor_type

integer Hamiltonian_module::BaseHamiltonian::enh_factor_type

Enhancement factor type set in options.

Definition at line 74 of file Hamiltonian_module.f90.

◆ initialized

logical Hamiltonian_module::BaseHamiltonian::initialized = .false.

Has the hamiltonian been initialized.

Definition at line 69 of file Hamiltonian_module.f90.

◆ integral

class(baseintegral), pointer Hamiltonian_module::BaseHamiltonian::integral

The integrals we are using.

Definition at line 61 of file Hamiltonian_module.f90.

◆ job_id

integer Hamiltonian_module::BaseHamiltonian::job_id = 0

Whose job it is to (soon to be deprecated)

Definition at line 70 of file Hamiltonian_module.f90.

◆ nflg

integer Hamiltonian_module::BaseHamiltonian::nflg = 0

Definition at line 64 of file Hamiltonian_module.f90.

◆ number_of_integrals

integer Hamiltonian_module::BaseHamiltonian::number_of_integrals = 0

How many integrals have been evaluated?

Definition at line 71 of file Hamiltonian_module.f90.

◆ options

class(options), pointer Hamiltonian_module::BaseHamiltonian::options

Scatci program settings.

Definition at line 60 of file Hamiltonian_module.f90.

◆ orbitals

class(orbitaltable), pointer Hamiltonian_module::BaseHamiltonian::orbitals

Our orbitals required to generate symblic elements.

Definition at line 59 of file Hamiltonian_module.f90.

◆ phase_flag

integer Hamiltonian_module::BaseHamiltonian::phase_flag

whether to evaluate integrals whilst dealing with phase

Definition at line 67 of file Hamiltonian_module.f90.

◆ positron_flag

integer Hamiltonian_module::BaseHamiltonian::positron_flag

Positron aware flag.

Definition at line 66 of file Hamiltonian_module.f90.

◆ reference_symbol

type(symbolicelementvector) Hamiltonian_module::BaseHamiltonian::reference_symbol

Symbols for idiag = 0.

Definition at line 76 of file Hamiltonian_module.f90.


The documentation for this type was generated from the following file: