MPI-SCATCI 2.0
An MPI version of SCATCI
Loading...
Searching...
No Matches
Uncontracted_Hamiltonian_module.f90
Go to the documentation of this file.
1! Copyright 2019
2!
3! For a comprehensive list of the developers that contributed to these codes
4! see the UK-AMOR website.
5!
6! This file is part of UKRmol-in (UKRmol+ suite).
7!
8! UKRmol-in is free software: you can redistribute it and/or modify
9! it under the terms of the GNU General Public License as published by
10! the Free Software Foundation, either version 3 of the License, or
11! (at your option) any later version.
12!
13! UKRmol-in is distributed in the hope that it will be useful,
14! but WITHOUT ANY WARRANTY; without even the implied warranty of
15! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16! GNU General Public License for more details.
17!
18! You should have received a copy of the GNU General Public License
19! along with UKRmol-in (in source/COPYING). Alternatively, you can also visit
20! <https://www.gnu.org/licenses/>.
21
22!> \brief Uncontracted Hamiltonian module
23!> \authors A Al-Refaie
24!> \date 2017
25!>
26!> \note 16/01/2019 - Jakub Benda: Unifom coding style and expanded documentation.
27!>
28module uncontracted_hamiltonian_module
29
30 use const_gbl, only: stdout
31 use consts_mpi_ci, only: mat_dense
32 use precisn, only: wp
33 use basematrix_module, only: basematrix
34 use hamiltonian_module, only: basehamiltonian
35 use parallelization_module, only: grid => process_grid
36 use symbolic_module, only: symbolicelementvector
37 use utility_module, only: compute_total_box, box_index_to_ij
38
39 implicit none
40
41 public uncontracted_hamiltonian
42
43 private
44
45 type, extends(basehamiltonian) :: uncontracted_hamiltonian
46 contains
47 procedure, public :: build_hamiltonian => build_uncontracted_hamiltonian
48 end type uncontracted_hamiltonian
49
50contains
51
52 !> \brief Initialize the type
53 !> \authors A Al-Refaie
54 !> \date 2017
55 !>
56 subroutine build_uncontracted_hamiltonian (this, matrix_elements)
57 class(uncontracted_hamiltonian) :: this
58 class(basematrix), intent(inout) :: matrix_elements
59 type(symbolicelementvector) :: symbolic_elements, ref_symbolic_elements
60 integer :: starting_index, num_csfs, ido, jdo, num_continuum, csf_a_idx, csf_b_idx
61 integer :: num_elements, loop_skip, my_idx, loop_ido, total_vals, m_flag, diag
62 real(wp) :: mat_coeff, element_one
63
64 diag = 0
65 num_continuum = this % options % last_continuum
66 num_csfs = this % options % num_csfs
67 num_elements = 1
68 this % orbitals % MFLG = 0
69
70 call matrix_elements % initialize_matrix_structure(num_csfs, mat_dense, num_csfs)
71 !get the first element and slater it regardless of idiag flag
72 call ref_symbolic_elements % construct
73 call this % slater_rules(this % csfs(1), this % csfs(1), ref_symbolic_elements, diag, .false.)
74
75 element_one = this % evaluate_integrals(ref_symbolic_elements, 0)
76
77 call symbolic_elements % construct
78 call symbolic_elements % clear
79
80 loop_skip = max(1, grid % gprocs)
81 my_idx = max(grid % grank, 0)
82
83 total_vals = compute_total_box(num_csfs, num_csfs)
84 do loop_ido = 1, total_vals, loop_skip
85
86 call matrix_elements % update_pure_L2(.false.)
87 ido = loop_ido + my_idx
88 if (ido > total_vals) cycle
89
90 call box_index_to_ij(ido, num_csfs, csf_a_idx, csf_b_idx)
91
92 if (csf_a_idx == csf_b_idx) then
93 m_flag = 0
94 else
95 m_flag = 1
96 end if
97 if (csf_a_idx < csf_b_idx) cycle
98
99 call this % slater_rules(this % csfs(csf_a_idx), this % csfs(csf_b_idx), symbolic_elements, m_flag, .false.)
100
101 if (m_flag == 0 .and. this % diagonal_flag == 0) then
102 call symbolic_elements % add_symbols(ref_symbolic_elements, -1.0_wp)
103 end if
104
105 mat_coeff = this % evaluate_integrals(symbolic_elements, this % options % phase_correction_flag)
106 if (this % diagonal_flag == 0 .and. this % orbitals % MFLG == 0) mat_coeff = mat_coeff + element_one
107 call matrix_elements % insert_matrix_element(csf_a_idx, csf_b_idx, mat_coeff)
108 num_elements = num_elements + 1
109
110 call symbolic_elements % clear
111
112 end do
113
114 !Flush
115 call matrix_elements % update_pure_L2(.true.)
116
117 !call matrix_elements%print
118 call symbolic_elements % destroy
119 call matrix_elements % finalize_matrix
120
121 write (stdout, "('Num of elements = ',i0)") num_elements
122
123 end subroutine build_uncontracted_hamiltonian
124
125end module uncontracted_hamiltonian_module
MPI SCATCI Constants module.
integer, parameter mat_dense
Matrix is dense.