MPI-SCATCI 2.0
An MPI version of SCATCI
Loading...
Searching...
No Matches
UKRMOL_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 UKRmol+ integral module
23!> \authors A Al-Refaie
24!> \date 2017
25!>
26!> The module makes use of the UKRmol+ integral library (GBTOlib).
27!>
28!> \note 16/01/2019 - Jakub Benda: Unifom coding style and expanded documentation.
29!>
30module ukrmol_module
31
32 use precisn, only: longint, wp
33 use const_gbl, only: stdout
34 use baseintegral_module, only: baseintegral
35 use options_module, only: options
36
37 implicit none
38
39 public ukrmolintegral
40
41 private
42
43 type, extends(baseintegral) :: ukrmolintegral
44 integer :: num_csf
45 integer :: num_one_electron_integrals
46 integer :: num_two_electron_integrals
47 integer :: symmetry_type
48 integer :: matrix_io_unit
49 integer :: num_symmetry
50 real(wp) :: exotic_mass
51 integer :: hamiltonian_unit
52 integer, allocatable :: num_orbitals_sym(:)
53 real(wp), allocatable :: xnuc(:), ynuc(:), znuc(:), charge(:)
54 character(len=8), allocatable, dimension(:) :: cname
55 contains
56 procedure, public :: initialize_self => initialize_ukrmol
57 procedure, public :: finalize_self => finalize_ukrmol
58 procedure, public :: load_integrals => load_ukrmol
59 procedure, public :: get_integral_ijklm => get_integral_ukrmol
60 procedure, public :: write_geometries => write_geometries_ukrmol
61 procedure, public :: destroy_integrals => destroy_integral_ukrmol
62 end type
63
64contains
65
66 !> \brief Initialize the type
67 !> \authors A Al-Refaie
68 !> \date 2017
69 !>
70 subroutine initialize_ukrmol (this, option)
71 class(ukrmolintegral) :: this
72 class(options), intent(in) :: option
73
74 write (stdout, "('Using UKRMOL+ integrals')")
75
76 this % num_csf = option % num_csfs
77 this % num_symmetry = option % num_syms
78
79 allocate(this % num_orbitals_sym(this % num_symmetry))
80
81 this % num_orbitals_sym(:) = option % num_electronic_orbitals_sym(:)
82 this % hamiltonian_unit = option % hamiltonian_unit
83 this % exotic_mass = option % exotic_mass
84 this % matrix_io_unit = option % hamiltonian_unit
85 this % name = option % name
86
87 end subroutine initialize_ukrmol
88
89
90 !> \brief Finalize the type
91 !> \authors J Benda
92 !> \date 2019
93 !>
94 !> Deallocate memory used by the integral storage.
95 !>
96 subroutine finalize_ukrmol (this)
97
98 use ukrmol_interface_gbl, only: destroy_ukrmolp
99
100 class(ukrmolintegral) :: this
101
102 if (allocated(this % num_orbitals_sym)) deallocate (this % num_orbitals_sym)
103 if (allocated(this % xnuc)) deallocate (this % xnuc)
104 if (allocated(this % ynuc)) deallocate (this % ynuc)
105 if (allocated(this % znuc)) deallocate (this % znuc)
106 if (allocated(this % charge)) deallocate (this % charge)
107 if (allocated(this % cname)) deallocate (this % cname)
108
109 call destroy_ukrmolp
110
111 end subroutine finalize_ukrmol
112
113
114 !> \brief Read UKRmol+ integrals
115 !> \authors A Al-Refaie
116 !> \date 2017
117 !>
118 !> Uses calls to GBTOlib to perform reading itself.
119 !>
120 subroutine load_ukrmol (this, iounit)
121 use ukrmol_interface_gbl, only: read_ukrmolp_ints, get_geom
122
123 class(ukrmolintegral) :: this
124 integer, intent(in) :: iounit
125 integer :: nalm = 0, lembf, matrix_size = 0,dummy_int
126
127 this % num_one_electron_integrals = -1
128 this % num_two_electron_integrals = -1
129
130 call read_ukrmolp_ints(this % matrix_io_unit, iounit, this % number_of_matrix_records, &
131 this % num_one_electron_integrals, this % num_two_electron_integrals, &
132 this % num_csf, stdout, this % symmetry_type, this % num_symmetry, this % num_orbitals_sym, &
133 this % positron_flag, 1.0_wp, this % name, &
134 nalm, .false.)
135
136 !Get geometries
137 call get_geom(this % nnuc, this % cname, this % xnuc, this % ynuc, this % znuc, this % charge)
138
139 !Read header files locally to store to memory
140 rewind(this % matrix_io_unit)
141 read (this % matrix_io_unit) matrix_size, this % number_of_matrix_records , dummy_int, matrix_size, &
142 dummy_int, this % num_symmetries, dummy_int, dummy_int, dummy_int, dummy_int, this % nnuc, &
143 dummy_int, this % NAME, this % NHE, this % DTNUC
144 this % core_energy = this % DTNUC(1)
145 close (this % matrix_io_unit)
146
147 end subroutine load_ukrmol
148
149
150 !> \brief Get integral value
151 !> \authors A Al-Refaie
152 !> \date 2017
153 !>
154 !> Obtain the integral value from GBTOlib.
155 !>
156 function get_integral_ukrmol (this, i, j, k, l, m) result(coeff)
157 use ukrmol_interface_gbl, only: get_integral
158
159 class(ukrmolintegral) :: this
160 integer, intent(in) :: i, j, k, l, m
161 integer :: a, b, c, d
162 real(wp) :: coeff
163
164 if (i == 0) then
165 a = this % orbital_mapping(j)
166 b = this % orbital_mapping(l)
167 coeff = get_integral(a, b, 0, 0, m)
168 else
169 a = this % orbital_mapping(i)
170 b = this % orbital_mapping(j)
171 c = this % orbital_mapping(k)
172 d = this % orbital_mapping(l)
173 coeff = get_integral(a, b, c, d, m)
174 end if
175
176 end function get_integral_ukrmol
177
178
179 !> \brief Write geometries to file
180 !> \authors A Al-Refaie
181 !> \date 2017
182 !>
183 subroutine write_geometries_ukrmol (this, iounit)
184
185 class(ukrmolintegral), intent(in) :: this
186 integer, intent(in) :: iounit
187 integer :: ido, ifail, i
188
189 do ido = 1, this % nnuc
190 write (iounit) this % cname(ido), this % xnuc(ido), this % ynuc(ido), this % znuc(ido), this % charge(ido)
191 end do
192
193 write (stdout, 120) (this % cname(i), this % xnuc(i), this % ynuc(i), this % znuc(i), this % charge(i), i = 1, this % nnuc)
194 120 FORMAT(/,' Nuclear data X Y Z Charge',/,(3x,a8,2x,4f10.6))
195
196 end subroutine write_geometries_ukrmol
197
198
199 !> \brief Destroy this type
200 !> \authors A Al-Refaie
201 !> \date 2017
202 !>
203 subroutine destroy_integral_ukrmol (this)
204 use ukrmol_interface_gbl, only: destroy_ukrmolp
205
206 class(ukrmolintegral) :: this
207
208 !Where is the destory integral or cleanup function?
209 !MADE ONE!
210 call destroy_ukrmolp()
211
212 end subroutine destroy_integral_ukrmol
213
214end module ukrmol_module