MPI-SCATCI 2.0
An MPI version of SCATCI
Loading...
Searching...
No Matches
CSF_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 CSF module
23!> \authors A Al-Refaie, J Benda
24!> \date 2017 - 2019
25!>
26!> This module handles reading and storing of configuration state functions from CONGEN.
27!> It provides a clean way of handling CSFs than multiple arrays! Hopefully!
28!>
29!> When MPI-3 shared memory is available, the configurations (\ref CSFManager::int64dtrs) and
30!> coefficients (\ref CSFManager::coeffs) are stored in an array shared across all processes.
31!>
32!> \note 16/01/2019 - Jakub Benda: Unifom coding style and expanded documentation.
33!>
34module csf_module
35
36 use const_gbl, only: stdout
37 use mpi_memory_gbl, only: mpi_memory_allocate_integer_2dim, mpi_memory_allocate_real, mpi_memory_synchronize, &
38 mpi_memory_deallocate_integer_2dim, mpi_memory_deallocate_real, local_master
39 use mpi_gbl, only: mpi_mod_barrier, shared_enabled
40 use precisn, only: longint, wp
41 use scatci_routines, only: rdnfto, rdwf, setposwf
42 use memorymanager_module, only: master_memory
43 use options_module, only: options
44 use orbital_module, only: orbitaltable
45 use parallelization_module, only: grid => process_grid
46
47 implicit none
48
49 private
50
51 public csfobject, csfmanager, csforbital
52
53 !> \brief Configuration state function factory
54 !> \authors A Al-Refaie, J Benda
55 !> \date 2017 - 2019
56 !>
57 !> This is a central object that takes care of reading and storing all configurations state functions (CSFs)
58 !> generated by CONGEN. The determinants forming the CSFs are compressed to a long bit array, which can be
59 !> shared by all processes in a MPI group. The manager class creates instances of \ref CSFObject (individual
60 !> CSFs) and \ref CSFOrbital (individual determinants). Both of them just effectively address into shared arrays
61 !> manager by this object.
62 !>
63 type :: csfmanager
64#ifdef mpithree
65 integer(longint), pointer :: int64dtrs(:,:) !< Unpacked determinants as a bit array (spin-orbital occupation numbers).
66 real(wp), pointer :: coeffdtrs(:) !< Coefficients of determinants stored in int64dtrs.
67#else
68 integer(longint), allocatable :: int64dtrs(:,:) !< Unpacked determinants as a bit array (spin-orbital occupation numbers).
69 real(wp), allocatable :: coeffdtrs(:) !< Coefficients of determinants stored in int64dtrs.
70#endif
71
72 integer :: int64_win = -1 !< MPI-3 shared memory window ID for int64dtrs.
73 integer :: coeffs_win = -1 !< MPI-3 shared memory window ID for coeffs.
74
75 integer :: megul !< CONGEN output unit with generated configuration state functions.
76 integer :: num_csfs !< Old name: NOCSF - # of CSFs
77 integer :: length_coeff
78 integer :: length_dtrs
79 integer, pointer :: num_dtr_csf_out(:) !< Old Name: NODO - Number of determinants in CSF in output
80
81 !---------CSF body------------------------------!
82 integer, allocatable :: index_dtrs(:) !< Old Name: INDO - Index of determinants
83 integer, allocatable :: index_coeffs(:) !< Old Name: ICDO - index of coefficients
84 integer, allocatable :: packed_dtrs(:) !< Old Name: NDO - Packed slater determinants
85 real(wp), allocatable :: coeffs(:) !< Old Name: NCO - coeffcients
86
87 class(options), pointer :: option_ptr
88 integer :: largest_num_configs
89
90 !This is for the fast slater rules
91 integer :: num_mo
92 integer :: num_so
93 integer :: num_electrons
94 integer :: nints !< Number of integers forming the shared determinant storage (bit array).
95
96 integer, allocatable :: orbitals(:)
97 integer, allocatable :: spin(:)
98 integer, allocatable :: reference_dtrs(:)
99 integer(longint), allocatable :: ref_int64dtrs(:)
100 class(orbitaltable), pointer :: pt_orb
101 contains
102 procedure, public :: initialize => read_csf_body
103 procedure, public :: print => print_all_csfs
104 procedure, public :: create_csfs
105 procedure, public :: toggle_determinant
106 procedure, public :: finalize => finalize_manager
107 end type csfmanager
108
109
110 !> \brief Single determinant
111 !> \authors A Al-Refaie, J Benda
112 !> \date 2017 - 2019
113 !>
114 !> Class holding data for a single determinant that forms a part of a configuration state function (CSF).
115 !> At the moment, it holds just offset in the shared contiguous determinant bit array, which is managed
116 !> by \ref CSFManager type.
117 !>
118 type csforbital
119 type(csfobject), pointer :: csf !< Pointer to the CSFObject with CSF this determinant belongs to.
120 integer :: id_dtr !< Offset in the global determinant array.
121 contains
122 procedure, public :: compare_excitations_fast
123 procedure, public :: get_determinants
124 final :: destroy_orbitals
125 end type csforbital
126
127
128 !> \brief Single configuration state function
129 !> \authors A Al-Refaie, J Benda
130 !> \date 2017 - 2019
131 !>
132 !> Set of determinants comprising a single configuration state function (CSF). The determinants are
133 !> represented by instances of the \ref CSFOrbital type.
134 !>
135 type csfobject
136 type(csfmanager), pointer :: manager !< Pointer to the creator CSFManager object.
137 type(csforbital), allocatable :: orbital(:) !< List of determinants.
138 integer :: num_orbitals
139 integer :: id
140 integer :: orbital_sequence_number !< Pointer to mapping of energies.
141 logical :: is_continuum !< Whether this is a continuum function or not.
142 integer :: target_symmetry !< The i discussed in the paper.
143 integer :: continuum_symmetry !< The gamma discussed in the paper.
144 integer :: continuum_orbital_idx !< The j = 1 or 2 discussed in the paper.
145 integer :: target_configuration !< The 'm' discussed in the paper.
146 contains
147 procedure, public :: print => print_csf
148 procedure, public :: check_slater => does_it_obey_slater
149 end type csfobject
150
151contains
152
153 !> \brief Reads the CSFs
154 !> \authors A Al-Refaie
155 !> \date 2017
156 !>
157 !> Reads all configurations from CONGEN output using the SCATCI's RDWF routine
158 !> and stores them in the CONGEN format. The translation from CONGEN format to MPI-SCATCI
159 !> internal (bitset) representation is done in \ref create_csfs.
160 !>
161 subroutine read_csf_body (this, option, orbital)
162 class(csfmanager), intent(inout) :: this
163 class(options), intent(in) :: option
164 class(orbitaltable), intent(in), target :: orbital
165
166 integer :: ido, ifail, size_mem
167
168 ! Store important values for the CSF generation
169 this % num_csfs = option % num_csfs
170 this % length_dtrs = option % length_dtrs
171 this % length_coeff = option % length_coeff
172 this % megul = option % megul
173
174 allocate(this % num_dtr_csf_out(this % num_csfs), stat = ifail)
175
176 call master_memory % track_memory(storage_size(this % num_dtr_csf_out)/8, &
177 size(this % num_dtr_csf_out), ifail, 'CSFMANAGER::num_dtr_csr_out')
178
179 this % num_dtr_csf_out = option % num_dtr_csf_out
180 this % num_electrons = option % num_electrons
181
182 this % num_MO = orbital % total_num_orbitals
183 this % num_SO = orbital % total_num_spin_orbitals
184
185 allocate(this % orbitals(this % num_SO), &
186 this % reference_dtrs(this % num_electrons), &
187 this % spin(this % num_SO), stat = ifail)
188
189 call master_memory % track_memory(storage_size(this % orbitals)/8, &
190 size(this % orbitals), ifail, 'CSFMANAGER::orbitals')
191 call master_memory % track_memory(storage_size(this % reference_dtrs)/8, &
192 size(this % reference_dtrs), ifail, 'CSFMANAGER::reference_dtrs')
193 call master_memory % track_memory(storage_size(this % spin)/8, &
194 size(this % spin), ifail, 'CSFMANAGER::spin')
195
196 do ido = 1, this % num_SO
197 this % orbitals(ido) = orbital % get_orbital_number(ido)
198 this % spin(ido) = orbital % get_spin(ido)
199 end do
200
201 this % reference_dtrs(:) = option % reference_dtrs(:)
202 this % Nints = this % num_SO / 64 + 1
203
204 allocate(this % ref_int64dtrs(this % Nints), stat = ifail)
205
206 call master_memory % track_memory(storage_size(this % ref_int64dtrs)/8, &
207 size(this % ref_int64dtrs), ifail, 'CSFMANAGER::ref_int64dtrs')
208
209 this % ref_int64dtrs(:) = 0
210 write (stdout, "('Nints = ',i4)") this % Nints
211 this % pt_orb => orbital
212
213 do ido = 1, this % num_electrons
214 call this % toggle_determinant(this % reference_dtrs(ido), this % ref_int64dtrs, this % Nints)
215 end do
216
217 if ((shared_enabled .and. grid % lrank == local_master) .or. .not. shared_enabled) then
218 ! Allocate arrays for reading from congen output
219 allocate(this % index_dtrs(this % num_csfs + 1), &
220 this % index_coeffs(this % num_csfs + 1), &
221 this % packed_dtrs(this % length_dtrs), &
222 this % coeffs(this % length_coeff), stat = ifail)
223
224 ! Rewind and reposition CONGEN file to the beginning of the wavefunction data
225 call setposwf(option % positron_flag, &
226 option % num_continuum, &
227 option % num_target_sym, &
228 this % megul)
229
230 ! Read from congen
231 call rdwf(this % num_csfs + 1, &
232 this % index_coeffs, &
233 this % index_dtrs, &
234 this % packed_dtrs, &
235 this % length_dtrs, &
236 this % coeffs, &
237 this % length_coeff, &
238 this % megul)
239 end if
240
241 end subroutine read_csf_body
242
243
244 !> \brief Release all memory used by the instance
245 !> \authors J Benda
246 !> \date 2019 - 2022
247 !>
248 !> \param[inout] this Manager object to update.
249 !>
250 subroutine finalize_manager (this)
251
252 class(csfmanager), intent(inout), target :: this
253
254 if (associated(this % num_dtr_csf_out)) then
255 call master_memory % free_memory(storage_size(this % num_dtr_csf_out)/8, size(this % num_dtr_csf_out))
256 deallocate (this % num_dtr_csf_out)
257 end if
258
259 if (allocated(this % orbitals)) then
260 call master_memory % free_memory(storage_size(this % orbitals)/8, size(this % orbitals))
261 deallocate (this % orbitals)
262 end if
263
264 if (allocated(this % reference_dtrs)) then
265 call master_memory % free_memory(storage_size(this % reference_dtrs)/8, size(this % reference_dtrs))
266 deallocate (this % reference_dtrs)
267 end if
268
269 if (allocated(this % spin)) then
270 call master_memory % free_memory(storage_size(this % spin)/8, size(this % spin))
271 deallocate (this % spin)
272 end if
273
274 if (allocated(this % ref_int64dtrs)) then
275 call master_memory % free_memory(storage_size(this % ref_int64dtrs)/8, size(this % ref_int64dtrs))
276 deallocate (this % ref_int64dtrs)
277 end if
278
279 if (allocated(this % index_dtrs)) deallocate (this % index_dtrs)
280 if (allocated(this % index_coeffs)) deallocate (this % index_coeffs)
281 if (allocated(this % packed_dtrs)) deallocate (this % packed_dtrs)
282 if (allocated(this % coeffs)) deallocate (this % coeffs)
283
284 call mpi_memory_deallocate_integer_2dim(this % int64dtrs, size(this % int64dtrs), this % int64_win, grid % lcomm)
285 call mpi_memory_deallocate_real(this % coeffdtrs, size(this % coeffdtrs), this % coeffs_win, grid % lcomm)
286
287 end subroutine finalize_manager
288
289
290 !> \brief This subroutine seperates out the CSFs and creates individual objects out of them
291 !> \authors A Al-Refaie, J Benda
292 !> \date 2017 - 2019
293 !>
294 !> The unpacked determinants forming each CSF are stored as bit fields in \ref CSFManager::int64dtrs (potentially
295 !> shared across the MPI processes participating in the diagonalization) of length equal to the number
296 !> of spin-orbitals and values 0 and 1 representing absence or presence of an electron in the spin-orbital.
297 !> The coefficients of these determinants are then adjusted to be compatible with ascending order
298 !> of spin-orbitals in the determinant (may result in sign flip).
299 !>
300 !> \param[inout] this Manager object to update.
301 !> \param[inout] CSF Configuration set to construct (will be allocated).
302 !> \param[in] orbital_sequence Pointer to orbital sequence (KPT).
303 !> \param[in] num_ci_target_sym Number of CI components of each target symmetry (NCTGT).
304 !> \param[in] continuum_projection Lambda value of the continuum orbitals associated with each target state (MCONT).
305 !>
306 subroutine create_csfs (this, CSF, orbital_sequence, num_ci_target_sym, continuum_projection)
307
308 class(csfmanager), intent(inout), target :: this
309 type(csfobject), intent(inout), allocatable, target :: CSF(:)
310 integer, intent(in) :: num_ci_target_sym(:), continuum_projection(:)
311
312 integer :: orbital_sequence(this % num_csfs)
313 integer :: csf_det(this % num_electrons)
314
315 real(wp) :: phase_fold
316 integer :: ido, jdo, kdo, ifail, num_orbs, start_index_dtr, start_index_coeff, num_spin_orb
317 integer :: dtr_size, elec_idx, total_number_orbs, orb_counter, count_number, symmetry_idx
318
319 symmetry_idx = 1
320 count_number = 1
321 this % largest_num_configs = 0
322
323 write (stdout, "('Creating our CSFs.........')", advance = "no")
324
325 total_number_orbs = 0
326 orb_counter = 0
327 do ido = 1, this % num_csfs
328 total_number_orbs = total_number_orbs + this % num_dtr_csf_out(ido)
329 end do
330
331 write (stdout, "('Total number of determinants = ',i12)") total_number_orbs
332
333 allocate(csf(this % num_csfs))
334
335 this % int64_win = mpi_memory_allocate_integer_2dim(this % int64dtrs, this % Nints, total_number_orbs, grid % lcomm)
336 this % coeffs_win = mpi_memory_allocate_real(this % coeffdtrs, total_number_orbs, grid % lcomm)
337
338 call mpi_mod_barrier(ifail, grid % gcomm)
339 call mpi_memory_synchronize(this % int64_win, grid % lcomm)
340 call mpi_memory_synchronize(this % coeffs_win, grid % lcomm)
341
342 do ido = 1, this % num_csfs
343 num_orbs = this % num_dtr_csf_out(ido)
344
345 csf(ido) % id = ido
346 csf(ido) % manager => this
347 csf(ido) % num_orbitals = num_orbs
348
349 ! Store the orbital number from mkpt
350 csf(ido) % orbital_sequence_number = orbital_sequence(ido)
351
352 this % largest_num_configs = max(this % largest_num_configs, num_orbs)
353
354 if (csf(ido) % orbital_sequence_number > 0) then
355
356 csf(ido) % is_continuum = .true.
357 csf(ido) % continuum_orbital_idx = mod(ido - 1, 2)
358 csf(ido) % target_symmetry = symmetry_idx
359 csf(ido) % target_configuration = ((count_number - 1) / 2) + 1
360 csf(ido) % continuum_symmetry = continuum_projection(symmetry_idx) + 1
361
362 count_number = count_number + 1
363
364 if (count_number > num_ci_target_sym(symmetry_idx) * 2) then
365 symmetry_idx = symmetry_idx + 1
366 count_number = 1
367 end if
368
369 else
370
371 csf(ido) % is_continuum = .false.
372
373 end if
374
375 orb_counter = orb_counter + 1
376 allocate(csf(ido) % orbital(num_orbs), stat = ifail) ! Allocate our orbitals
377
378 csf(ido) % orbital(1) % id_dtr = orb_counter
379 csf(ido) % orbital(1) % csf => csf(ido)
380
381 if ((shared_enabled .and. grid % lrank == local_master) .or. .not. shared_enabled) then
382 !Assign the coefficient to the first
383 start_index_dtr = this % index_dtrs(ido)
384 start_index_coeff = this % index_coeffs(ido) - 1
385 num_spin_orb = this % packed_dtrs(start_index_dtr)
386 this % int64dtrs(:,orb_counter) = this % ref_int64dtrs(:)
387 csf_det(:) = this % reference_dtrs(:)
388
389 do kdo = 1, num_spin_orb
390 call this % toggle_determinant(this % packed_dtrs(start_index_dtr + kdo), &
391 this % int64dtrs(:,orb_counter), this % Nints)
392 call this % toggle_determinant(this % packed_dtrs(start_index_dtr + num_spin_orb + kdo), &
393 this % int64dtrs(:,orb_counter), this % Nints)
394 elec_idx = this % pt_orb % get_electron_number(this % packed_dtrs(start_index_dtr + kdo))
395 csf_det(elec_idx) = this % packed_dtrs(start_index_dtr + num_spin_orb + kdo)
396 end do
397
398 phase_fold = 1_wp
399 call qsortdets(csf_det, phase_fold)
400 this % coeffdtrs(orb_counter) = this % coeffs(start_index_coeff + 1) * phase_fold
401 start_index_dtr = start_index_dtr + 2 * num_spin_orb + 1
402 end if
403
404 !Move to the next orbital
405
406 !Do for the rest
407 do jdo = 2, num_orbs
408 orb_counter = orb_counter + 1
409 csf(ido) % orbital(jdo) % id_dtr = orb_counter
410 csf(ido) % orbital(jdo) % csf => csf(ido)
411
412 if ((shared_enabled .and. grid % lrank == local_master) .or. .not. shared_enabled) then
413 num_spin_orb = this % packed_dtrs(start_index_dtr)
414 this % int64dtrs(:,orb_counter) = this % ref_int64dtrs(:)
415 csf_det(:) = this % reference_dtrs(:)
416
417 do kdo = 1, num_spin_orb
418 call this % toggle_determinant(this % packed_dtrs(start_index_dtr + kdo), &
419 this % int64dtrs(:,orb_counter), this % Nints)
420 call this % toggle_determinant(this % packed_dtrs(start_index_dtr + num_spin_orb + kdo), &
421 this % int64dtrs(:,orb_counter), this % Nints)
422 elec_idx = this % pt_orb % get_electron_number(this % packed_dtrs(start_index_dtr + kdo))
423 csf_det(elec_idx) = this % packed_dtrs(start_index_dtr + num_spin_orb + kdo)
424 end do
425
426 phase_fold = 1_wp
427 call qsortdets(csf_det, phase_fold)
428 this % coeffdtrs(orb_counter) = this % coeffs(start_index_coeff + jdo) * phase_fold
429 start_index_dtr = start_index_dtr + 2 * num_spin_orb + 1
430 end if
431 end do
432 end do
433
434 write (stdout, "('success!')")
435 write (stdout, "('Largest number of configurations in a single CSF is ',i10)") this % largest_num_configs
436
437 call mpi_mod_barrier(ifail, grid % gcomm)
438 call mpi_memory_synchronize(this % int64_win, grid % lcomm)
439 call mpi_memory_synchronize(this % coeffs_win, grid % lcomm)
440 !call master_memory % free_memory(storage_size(this % index_dtrs)/8, size(this % index_dtrs))
441 !call master_memory % free_memory(storage_size(this % index_coeffs)/8, size(this % index_coeffs))
442 !call master_memory % free_memory(storage_size(this % packed_dtrs)/8, size(this % packed_dtrs))
443 call master_memory % free_memory(storage_size(this % num_dtr_csf_out)/8,size(this % num_dtr_csf_out))
444
445 ! We don't need them anymore
446 if ((shared_enabled .and. grid % lrank == local_master) .or. .not. shared_enabled) then
447 deallocate(this % index_dtrs, this % index_coeffs, this % packed_dtrs, this % coeffs)
448 end if
449 deallocate(this % num_dtr_csf_out)
450
451 end subroutine create_csfs
452
453
454 subroutine print_csf (this)
455 class(csfobject), intent(in) :: this
456 integer(longint) :: mask
457 integer :: sporbs_per_int, idet, isporb, iint
458
459 sporbs_per_int = bit_size(this % manager % int64dtrs)
460 do idet = 1, this % num_orbitals
461 write (stdout, '(6x,a,i3,a,e23.15e3,a)', advance = 'no') 'det', idet, ': ', &
462 this % manager % coeffdtrs(this % orbital(idet) % id_dtr), ' :'
463 do isporb = 1, this % manager % nints * sporbs_per_int
464 iint = (isporb - 1) / sporbs_per_int + 1
465 mask = 2**mod(isporb - 1, sporbs_per_int)
466 if (iand(mask, this % manager % int64dtrs(iint, this % orbital(idet) % id_dtr)) /= 0) then
467 write (stdout, '(1x,i0)', advance = 'no') isporb
468 end if
469 end do
470 write (stdout, '()')
471 end do
472
473 end subroutine print_csf
474
475
476 subroutine print_all_csfs (this, csf)
477 class(csfmanager), intent(in) :: this
478 class(csfobject), intent(in) :: csf(:)
479 integer :: ido
480
481 do ido = 1, this % num_csfs
482 call csf(ido) % print
483 end do
484
485 end subroutine print_all_csfs
486
487
488 subroutine destroy_orbitals (this)
489 type(csforbital) :: this
490 end subroutine destroy_orbitals
491
492
493 subroutine compare_excitations_fast (this, that, num_electrons, coeff, num_excitations, output_dtrs)
494 class(csforbital), intent(in) :: this, that
495 class(csfmanager), pointer :: manager
496 integer, intent(in) :: num_electrons
497 integer, intent(out) :: num_excitations
498 integer, intent(out) :: output_dtrs(4)
499 real(wp), intent(out) :: coeff
500 real(wp) :: phase
501
502 phase = 1
503 coeff = 0.0_wp
504 output_dtrs = 0
505
506 manager => this % csf % manager
507
508 !DIR$ FORCEINLINE
509 call get_excitation(manager % int64dtrs(:, this % id_dtr), &
510 manager % int64dtrs(:, that % id_dtr), &
511 output_dtrs, num_excitations, phase, manager % Nints)
512
513 coeff = phase * manager % coeffdtrs(this % id_dtr) &
514 * manager % coeffdtrs(that % id_dtr)
515
516 end subroutine compare_excitations_fast
517
518
519 integer function does_it_obey_slater (this, that)
520 class(csfobject), intent(in) :: this, that
521 class(csfmanager), pointer :: manager
522 integer :: min_num_excitation, ido
523
524 manager => this % manager
525
526 does_it_obey_slater = pzero(manager % int64dtrs(:, this % orbital(1) % id_dtr), &
527 manager % int64dtrs(:, that % orbital(1) % id_dtr), &
528 manager % Nints)
529
530 end function does_it_obey_slater
531
532
533 subroutine toggle_determinant (this, det0, det1, Nints)
534 class(csfmanager), intent(in) :: this
535 integer, intent(in) :: Nints
536 integer :: det0
537 integer(longint) :: det1(Nints)
538 integer(longint) :: bit_set
539 integer :: spin_position, array_position, local_orbital
540
541 bit_set = 0
542 array_position = (det0 - 1) / 64 + 1
543 ! odd is alpha even is beta
544
545 local_orbital = det0 - (array_position - 1) * 64 - 1
546 bit_set = ibset(bit_set, local_orbital)
547
548 det1(array_position) = ieor(det1(array_position), bit_set)
549
550 end subroutine toggle_determinant
551
552
553 integer function pzero (det1, det2, Nint)
554 integer, intent(in) :: nint
555 integer(longint), intent(in) :: det1(nint), det2(nint)
556 integer(longint) :: xor1
557 integer :: l, nexcitations = 0
558
559 xor1 = ieor(det1(1), det2(1))
560 pzero = popcnt(ieor(get_beta(xor1), get_alpha(xor1)))
561
562 do l = 2, nint
563 xor1 = ieor(det1(l), det2(l))
564 pzero = pzero + popcnt(ieor(get_beta(xor1), get_alpha(xor1)))
565 end do
566
567 end function pzero
568
569
570 ! \brief Select only beta orbitals
571 ! \date 2017 - 2020
572 ! \author A Al-Refaie, J Benda
573 !
574 ! Given a section of a bit array representing the population of spin-orbitals, mask out the the alpha
575 ! orbitals (odd bits), keeping only the even bits. On return, the even bits are shifted one bit back,
576 ! to odd positions.
577 !
578 ! Note that rather than masking out odd bits and then shifting the result, we do it the other way round (first
579 ! shift, then mask out *even* bits). This is because in some circumstances the leading non-zero bit is reserved
580 ! for sign and constructing a mask that explicitly sets it may result in "Arithmetic Overflow" errors.
581 !
582 elemental integer function get_beta (det)
583 integer(longint), intent(in) :: det
584 integer(longint), parameter :: mask = int(z'5555555555555555', longint) ! = 0101...0101
585
586 get_beta = iand(ishft(det, -1), mask)
587
588 end function get_beta
589
590
591 ! \brief Select only alpha orbitals
592 ! \date 2017 - 2020
593 ! \author A Al-Refaie, J Benda
594 !
595 ! Given a section of a bit array representing the population of spin-orbitals, mask out the the beta
596 ! orbitals (even bits), keeping only the odd bits.
597 !
598 elemental integer function get_alpha (det)
599 integer(longint), intent(in) :: det
600 integer(longint), parameter :: mask = int(z'5555555555555555', longint) ! = 0101...0101
601
602 get_alpha = iand(det, mask)
603
604 end function get_alpha
605
606
607 integer function n_excitations (det1, det2, Nint)
608 integer, intent(in) :: nint
609 integer(longint), intent(in) :: det1(nint), det2(nint)
610
611 integer :: l
612
613 n_excitations = popcnt(ieor(det1(1), det2(1)))
614
615 do l = 2, nint
616 n_excitations = n_excitations + popcnt(ieor(det1(l), det2(l)))
617 end do
618
619 n_excitations = ishft(n_excitations, -1)
620
621 end function n_excitations
622
623
624 subroutine get_excitation (det1, det2, exc, degree, phase, Nint)
625 integer, intent(in) :: Nint
626 integer(longint), intent(in) :: det1(Nint,2), det2(Nint,2)
627 integer, intent(out) :: exc(4)
628 integer, intent(out) :: degree
629 double precision, intent(inout) :: phase
630
631 degree = n_excitations(det1, det2, nint)
632
633 select case (degree)
634 case (3:)
635 return
636 case (2)
637 call get_double_excitation(det1, det2, exc, phase, nint)
638 return
639 case (1)
640 call get_single_excitation(det1, det2, exc(3:4), phase, nint)
641 return
642 case (0)
643 return
644 end select
645
646 end subroutine get_excitation
647
648
649 subroutine get_single_excitation (det1, det2, exc, phase, Nint)
650 integer, intent(in) :: Nint
651 integer(longint), intent(in) :: det1(Nint)
652 integer(longint), intent(in) :: det2(Nint)
653 integer, intent(inout) :: exc(2)
654 double precision, intent(out) :: phase
655
656 integer :: tz, l, ispin, ishift, nperm, i, j, k, m, n, high, low
657 integer :: holes(1)
658 integer :: particles(1)
659 integer(longint) :: hole, particle, tmp, masklow, maskhigh
660 double precision, parameter :: phase_dble(0:1) = (/ 1.d0, -1.d0 /)
661
662 nperm = 0
663 exc(:) = 0
664 ishift = -63
665 do l = 1, nint
666 ishift = ishift + 64
667 if (det1(l) == det2(l)) cycle
668 tmp = ieor(det1(l), det2(l))
669 particle = iand(tmp, det2(l))
670 hole = iand(tmp, det1(l))
671 if (particle /= 0_8) then
672 tz = trailz(particle)
673 exc(2) = tz + ishift
674 end if
675 if (hole /= 0_8) then
676 tz = trailz(hole)
677 !exc(0,1,ispin) = 1
678 exc(1) = tz + ishift
679 end if
680 if (exc(1) /= 0 .and. exc(2) /=0 ) then
681 low = min(exc(1), exc(2))
682 high = max(exc(1), exc(2))
683 j = ishft(low - 1, -6) + 1
684 n = iand(low - 1, 63)
685 k = ishft(high - 1, -6) + 1
686 m = iand(high - 1, 63)
687 ! masklow = not(ishft(1_8,n+1))+1
688 ! maskhigh = ishft(1_8,m)-1
689 if (j == k) then
690 nperm = nperm + popcnt(iand(det1(j), iand(not(ishft(1_8, n + 1)) + 1, ishft(1_8, m) - 1)))
691 else
692 nperm = nperm + popcnt(iand(det1(k), ishft(1_8, m) - 1)) &
693 + popcnt(iand(det1(j), not(ishft(1_8, n + 1)) + 1))
694 do i = j + 1, k - 1
695 nperm = nperm + popcnt(det1(i))
696 end do
697 end if
698 phase = phase_dble(iand(nperm, 1))
699 return
700 end if
701 end do
702
703 end subroutine get_single_excitation
704
705
706 subroutine get_double_excitation (det1, det2, exc, phase, Nint)
707 integer, intent(in) :: Nint
708 integer(longint), intent(in) :: det1(Nint), det2(Nint)
709 integer, intent(out) :: exc(4)
710 double precision, intent(out) :: phase
711 integer :: l, ispin, idx_hole, idx_particle, ishift
712 integer :: i, j, k, m, n, high, low, a, b, c, d, nperm, tz, nexc, num_holes
713 integer(longint) :: hole, particle, tmp, spin_total, masklow, maskhigh
714 integer :: holes(2), particles(2), comp(2)
715 double precision, parameter :: phase_dble(0:1) = (/ 1.d0, -1.d0 /)
716
717 exc(:) = 0
718 nexc = 0
719 nperm = 0
720 idx_particle = 0
721 idx_hole = 0
722 ishift = -63
723 num_holes = 0
724 do l = 1, nint
725 ishift = ishift + 64
726 if (det1(l) == det2(l)) cycle
727 tmp = ieor(det1(l), det2(l))
728 particle = iand(tmp, det2(l))
729 hole = iand(tmp, det1(l))
730 do while (particle /= 0_8)
731 tz = trailz(particle)
732 nexc = nexc + 1
733 idx_particle = idx_particle + 1
734 !exc(2*idx_particle) = tz+ishift
735 particles(idx_particle) = tz + ishift
736 particle = iand(particle, particle - 1_8)
737 end do
738 do while (hole /= 0_8)
739 tz = trailz(hole)
740 nexc = nexc + 1
741 idx_hole = idx_hole + 1
742 holes(idx_hole) = tz + ishift
743 num_holes = num_holes + 1
744 !exc( idx_hole*2 -1) = tz+ishift
745 hole = iand(hole, hole - 1_8)
746 end do
747 if (nexc == 4) exit
748 end do
749
750 !if(exc(4) < exc(2)) then
751 ! tmp = exc(4)
752 ! exc(4) = exc(2)
753 ! exc(2) = tmp
754 ! endif
755 ! if(exc(3) < exc(1)) then
756 ! tmp = exc(3)
757 ! exc(1) = exc(3)
758 ! exc(3) = tmp
759 ! endif
760 !if(holes(2) < holes(1)) then
761 ! tmp = holes(1)
762 ! holes(1) = holes(2)
763 ! holes(2) = tmp
764 !endif
765
766 !write(stdout,*) particles
767
768 !if(particles(2) < particles(1)) then
769 ! tmp = particles(1)
770 ! particles(1) = particles(2)
771 ! particles(2) = tmp
772 !endif
773 ! write(stdout,*) particles
774
775 do i = 1, num_holes
776 low = min(particles(i), holes(i))
777 high = max(particles(i), holes(i))
778 j = ishft(low - 1, -6) + 1
779 n = iand(low - 1, 63)
780 k = ishft(high - 1, -6) + 1
781 m = iand(high - 1, 63)
782
783 if (j == k) then
784 nperm = nperm + popcnt(iand(det1(j), iand(not(ishft(1_8, n + 1)) + 1, ishft(1_8, m) - 1)))
785 else
786 nperm = nperm + popcnt(iand(det1(k), ishft(1_8, m) - 1)) &
787 + popcnt(iand(det1(j), not(ishft(1_8, n + 1)) + 1))
788 do l = j + 1, k - 1
789 nperm = nperm + popcnt(det1(l))
790 end do
791 end if
792 end do
793
794 ! spin_total = iand(particles(1),1) + iand(particles(2),1) + iand(holes(1),1) + iand(holes(2),1)
795 ! if (spin_total == 0 .or. spin_total == 4) then
796
797 if (num_holes == 2) then
798 a = min(holes(1), particles(1))
799 b = max(holes(1), particles(1))
800 c = min(holes(2), particles(2))
801 d = max(holes(2), particles(2))
802 if (c > a .and. c < b .and. d > b) nperm = nperm + 1
803 end if
804
805 ! end do
806 !exc(1) = holes(1)
807 !exc(2) = particles(1)
808 !exc(3) = holes(2)
809 !exc(4) = particles(2)
810 !if(holes(1) > holes(2)) then
811 ! comp(1) = particles(1)
812 ! comp(2) = particles(2)
813 ! else
814 ! comp(1) = particles(2)
815 ! comp(2) = particles(1)
816 ! endif
817 ! if(comp(1) > comp(2)) nperm = nperm + 1
818 ! exit
819 ! end if
820 ! end do
821 exc(1) = holes(1)
822 exc(2) = particles(1)
823 exc(3) = holes(2)
824 exc(4) = particles(2)
825
826 phase = phase_dble(iand(nperm, 1))
827
828 end subroutine get_double_excitation
829
830
831 subroutine assign_pointer (targ, point, Nints)
832 integer, intent(in) :: Nints
833 integer(longint), target, intent(in) :: targ(Nints)
834 integer(longint), pointer, intent(out) :: point(:)
835
836 point(1:nints) => targ(1:nints)
837
838 end subroutine assign_pointer
839
840
841 subroutine get_determinants (this, dtrs, nelec)
842 integer, intent(in) :: nelec
843 class(csforbital), intent(in) :: this
844 integer, intent(out) :: dtrs(nelec)
845 integer(longint) :: current_so
846 integer :: tz, ishift = -63, ints, idx_dtr, nexec
847
848 idx_dtr = 0
849 dtrs = 0
850 ishift = -63
851
852 do ints = 1, this % csf % manager % Nints
853 ishift = ishift + 64
854 current_so = this % csf % manager % int64dtrs(ints, this % id_dtr)
855 do while (current_so /= 0)
856 tz = trailz(current_so)
857 idx_dtr = idx_dtr + 1
858 dtrs(idx_dtr) = tz + ishift
859 current_so = iand(current_so, current_so - 1)
860 end do
861 end do
862
863 end subroutine get_determinants
864
865
866 recursive subroutine qsortdets(A,phase)
867 integer, intent(inout), dimension(:) :: a
868 real(wp), intent(inout) :: phase
869 integer :: iq
870
871 if (size(a) > 1) then
872 call partition(a, iq, phase)
873 call qsortdets(a(:iq-1), phase)
874 call qsortdets(a(iq:), phase)
875 end if
876
877 end subroutine qsortdets
878
879
880 subroutine partition (A, marker, phase)
881 integer, intent(inout) :: A(:)
882 integer, intent(out) :: marker
883 real(wp), intent(inout) :: phase
884 integer :: i, j
885 integer :: temp
886 integer :: x ! pivot point
887
888 x = a(1)
889 i = 0
890 j = size(a) + 1
891
892 do
893 j = j - 1
894 do
895 if (a(j) <= x) exit
896 j = j-1
897 end do
898 i = i + 1
899 do
900 if (a(i) >= x) exit
901 i = i + 1
902 end do
903 if (i < j) then
904 ! exchange A(i) and A(j)
905 temp = a(i)
906 a(i) = a(j)
907 a(j) = temp
908 phase = -phase
909 elseif (i == j) then
910 marker = i + 1
911 return
912 else
913 marker = i
914 return
915 end if
916 end do
917
918 end subroutine partition
919
920end module csf_module