MPI-SCATCI 2.0
An MPI version of SCATCI
Loading...
Searching...
No Matches
BaseIntegral_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 Base integral 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 baseintegral_module
29
31 use mpi_gbl, only: master, myrank
32 use precisn, only: longint, wp
33 use options_module, only: options
34 use parallelization_module, only: grid => process_grid
35 use utility_module, only: unpack_ints
36
37 implicit none
38
39 public baseintegral
40
41 type, abstract :: baseintegral
42 integer, allocatable :: orbital_mapping(:)
43 integer :: num_orbitals
44 real(wp) :: core_energy
45 integer :: positron_flag
46 logical :: quantamoln_flag
47
48 !>Matrix header information
49 integer :: matrix_size
50 integer :: nhe(20), nhd(10)
51 integer :: num_symmetries
52 integer :: number_of_matrix_records
53 integer :: nnuc
54 real(wp) :: dtnuc(41)
55
56 character(NAME_LEN_MAX) :: name
57 contains
58 procedure(generic_finalize), deferred :: finalize_self
59 procedure(generic_initialize), deferred :: initialize_self
60 procedure(generic_load), deferred :: load_integrals
61 procedure(generic_get), deferred :: get_integral_ijklm
62 procedure(generic_geometries), deferred :: write_geometries
63 procedure(generic_destroy), deferred :: destroy_integrals
64
65 procedure, public :: initialize => initialize_base
66 procedure, public :: finalize => finalize_base
67 procedure, public :: write_matrix_header => base_write_header
68 procedure, public :: get_core_energy
69 procedure, public :: get_num_nuclei
70 procedure, public :: get_integralf
71 end type baseintegral
72
73 abstract interface
74 subroutine generic_initialize (this, option)
75 use options_module, only: options
76 import :: baseintegral
77 class(baseintegral) :: this
78 class(options), intent(in) :: option
79 end subroutine generic_initialize
80 end interface
81
82 abstract interface
83 subroutine generic_finalize (this)
84 use options_module, only: options
85 import :: baseintegral
86 class(baseintegral) :: this
87 end subroutine generic_finalize
88 end interface
89
90 abstract interface
91 subroutine generic_load (this, iounit)
92 import :: baseintegral
93 class(baseintegral) :: this
94 integer, intent(in) :: iounit
95 end subroutine generic_load
96 end interface
97
98 abstract interface
99 function generic_get (this, i, j, k, l, m) result(coeff)
100 use precisn, only : wp
101 import :: baseintegral
102 class(baseintegral) :: this
103 integer, intent(in) :: i, j, k, l, m
104 real(wp) :: coeff
105 end function generic_get
106 end interface
107
108 abstract interface
109 subroutine generic_geometries (this, iounit)
110 use precisn, only : wp
111 import :: baseintegral
112 class(baseintegral), intent(in) :: this
113 integer, intent(in) :: iounit
114 end subroutine generic_geometries
115 end interface
116
117 abstract interface
118 subroutine generic_destroy (this)
119 import :: baseintegral
120 class(baseintegral) :: this
121 end subroutine generic_destroy
122 end interface
123
124contains
125
126 !> \brief Initialize the class
127 !> \authors A Al-Refaie
128 !> \date 2017
129 !>
130 subroutine initialize_base (this, option, num_orbitals, mapping)
131
132 class(baseintegral) :: this
133 class(options), intent(in) :: option
134 integer, intent(in) :: num_orbitals
135 integer, intent(in) :: mapping(num_orbitals)
136
137 this % num_orbitals = num_orbitals
138
139 allocate(this % orbital_mapping(num_orbitals))
140
141 this % orbital_mapping(:) = mapping(:)
142 this % positron_flag = option % positron_flag
143 this % quantamoln_flag = option % QuantaMolN
144
145 !>Store needed matrix header information
146 this % num_symmetries = option % num_syms
147
148 this % number_of_matrix_records = option % num_matrix_elements_per_rec
149 this % name = option % name
150
151 call this % initialize_self(option)
152
153 end subroutine initialize_base
154
155
156 !> \brief Finalize the class
157 !> \authors J Benda
158 !> \date 2019
159 !>
160 !> Deallocate memory used by the integrals storage.
161 !>
162 subroutine finalize_base (this)
163
164 class(baseintegral) :: this
165
166 call this % finalize_self
167
168 if (allocated(this % orbital_mapping)) then
169 deallocate (this % orbital_mapping)
170 end if
171
172 end subroutine finalize_base
173
174
175 !> \brief Write data file header
176 !> \authors A Al-Refaie
177 !> \date 2017
178 !>
179 subroutine base_write_header (this, matrix_io, matrix_size)
180
181 class(baseintegral), intent(in) :: this
182 integer, intent(in) :: matrix_io, matrix_size
183
184 if (grid % grank == master) then
185 open (unit = matrix_io, form = 'unformatted')
186 write (matrix_io) matrix_size, this % number_of_matrix_records , 0, matrix_size, &
187 0, this % num_symmetries, 0, 0, 0, 0, this % nnuc, 0, &
188 this % NAME, this % NHE, this % DTNUC
189 end if
190
191 end subroutine base_write_header
192
193
194 !> \brief Retrieve integral value
195 !> \authors A Al-Refaie
196 !> \date 2017
197 !>
198 function get_integralf (this, integral) result(coeff)
199
200 class(baseintegral) :: this
201 integer(longint), intent(in) :: integral(nidx)
202 real(wp) :: coeff
203 integer :: label(8)
204
205 call unpack_ints(integral, label)
206 coeff = this % get_integral_ijklm(label(1), label(2), label(3), label(4), label(5))
207
208 end function get_integralf
209
210
211 !> \brief Retrieve integral value
212 !> \authors A Al-Refaie
213 !> \date 2017
214 !>
215 real(wp) function get_core_energy (this)
216
217 class(baseintegral) :: this
218
219 get_core_energy = this % core_energy
220
221 end function get_core_energy
222
223
224 !> \brief Retrieve integral value
225 !> \authors A Al-Refaie
226 !> \date 2017
227 !>
228 integer function get_num_nuclei (this)
229
230 class(baseintegral) :: this
231
232 get_num_nuclei = this % nnuc
233
234 end function get_num_nuclei
235
236end module baseintegral_module
MPI SCATCI Constants module.
integer, parameter nidx
Number of long integers used to store up to 8 (packed) integral indices.
integer, parameter name_len_max
The maximum length of a name.