|
MPI-SCATCI 2.0
An MPI version of SCATCI
|
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. | |
This is an abstract class that contains the majority of functionality required to construct hamiltonians.
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.
|
pure virtual |
Definition at line 80 of file Hamiltonian_module.f90.
| 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.
| [in,out] | this | Hamiltonian object to update. |
| [in] | option | An Options class containing run constants |
| [in] | csfs | An array of configuration state functions |
| [in] | orbitals | An initialized Orbitals class |
| [in] | integral | A BaseIntegral class for evaluating symbols |
Definition at line 79 of file Hamiltonian_module.f90.
| 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.
| [in,out] | this | Hamiltonian object to update. |
| [in] | symbolic_elements | Symbolic elements to be evaluated. |
| [in] | dummy_orb | Check wheter to ignore the integral if it contains a dummy orbital of this value. |
Definition at line 85 of file Hamiltonian_module.f90.
| 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.
| [in,out] | this | Hamiltonian object to update. |
| [in] | label | A packed integral label |
| [in] | coeff | The computed coefficients from the Slater rules |
| [in] | dummy_orb | Check wheter to ignore the integral if it contains a dummy orbital of this value |
Definition at line 86 of file Hamiltonian_module.f90.
|
private |
Handles if there are one pair of spin orbitals that differ.
| [in,out] | this | Hamiltonian object to query. |
| [in] | dtrs | A 4 element array of slater determinants |
| [in] | coeff | The computed coefficients from the Slater rules |
| [in] | csf | The first CSFs determinant |
| [out] | sym_element | Resulting symbolic elements |
| [in] | flag | Whether we are doing the diagonal or not |
Definition at line 91 of file Hamiltonian_module.f90.
|
private |
Handles if there are no differences between spin orbitals.
| [in,out] | this | Hamiltonian object to query. |
| [in] | dtrs | A 4 element array of slater determinants |
| [in] | coeff | The computed coefficients from the Slater rules |
| [in] | csf | The first CSFs determinant |
| [out] | sym_element | Resulting symbolic elements |
| [in] | flag | Whether we are doing the diagonal or not |
Definition at line 92 of file Hamiltonian_module.f90.
|
private |
Handles if there are two pairs of spin orbitals that differ.
| [in,out] | this | Hamiltonian object to query. |
| [in] | dtrs | A 4 element array of slater determinants that differ. |
| [in] | coeff | The computed coefficients from the Slater rules. |
| [out] | sym_element | Resulting symbolic elements. |
| [in] | flag | Whether we are doing the diagonal or not. |
Definition at line 90 of file Hamiltonian_module.f90.
| 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.
| 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.
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
| [in,out] | this | Hamiltonian object to query. |
| [in] | csf_a | First configuration state function to slater |
| [in] | csf_b | Second configuration state function to slater |
| [out] | sym_elements | Resulting smbolic elements |
| [in] | flag | Whether 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.
| logical Hamiltonian_module::BaseHamiltonian::constructed = .false. |
Has the hamiltonain been constructed.
Definition at line 68 of file Hamiltonian_module.f90.
| class(csfobject), dimension(:), pointer Hamiltonian_module::BaseHamiltonian::csfs |
Our configuration state functions.
Definition at line 62 of file Hamiltonian_module.f90.
| integer Hamiltonian_module::BaseHamiltonian::diagonal_flag |
Definition at line 65 of file Hamiltonian_module.f90.
| real(wp) Hamiltonian_module::BaseHamiltonian::element_one = 0.0 |
First element for idiag = 0.
Definition at line 72 of file Hamiltonian_module.f90.
| real(wp) Hamiltonian_module::BaseHamiltonian::enh_factor |
Enhancement factor set in options.
Definition at line 73 of file Hamiltonian_module.f90.
| integer Hamiltonian_module::BaseHamiltonian::enh_factor_type |
Enhancement factor type set in options.
Definition at line 74 of file Hamiltonian_module.f90.
| logical Hamiltonian_module::BaseHamiltonian::initialized = .false. |
Has the hamiltonian been initialized.
Definition at line 69 of file Hamiltonian_module.f90.
| class(baseintegral), pointer Hamiltonian_module::BaseHamiltonian::integral |
The integrals we are using.
Definition at line 61 of file Hamiltonian_module.f90.
| 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.
| integer Hamiltonian_module::BaseHamiltonian::nflg = 0 |
Definition at line 64 of file Hamiltonian_module.f90.
| integer Hamiltonian_module::BaseHamiltonian::number_of_integrals = 0 |
How many integrals have been evaluated?
Definition at line 71 of file Hamiltonian_module.f90.
| class(options), pointer Hamiltonian_module::BaseHamiltonian::options |
Scatci program settings.
Definition at line 60 of file Hamiltonian_module.f90.
| class(orbitaltable), pointer Hamiltonian_module::BaseHamiltonian::orbitals |
Our orbitals required to generate symblic elements.
Definition at line 59 of file Hamiltonian_module.f90.
| integer Hamiltonian_module::BaseHamiltonian::phase_flag |
whether to evaluate integrals whilst dealing with phase
Definition at line 67 of file Hamiltonian_module.f90.
| integer Hamiltonian_module::BaseHamiltonian::positron_flag |
Positron aware flag.
Definition at line 66 of file Hamiltonian_module.f90.
| type(symbolicelementvector) Hamiltonian_module::BaseHamiltonian::reference_symbol |
Symbols for idiag = 0.
Definition at line 76 of file Hamiltonian_module.f90.