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
166 integer :: ido, ifail, size_mem
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
174 allocate(this % num_dtr_csf_out(this % num_csfs), stat = ifail)
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')
179 this % num_dtr_csf_out = option % num_dtr_csf_out
180 this % num_electrons = option % num_electrons
182 this % num_MO = orbital % total_num_orbitals
183 this % num_SO = orbital % total_num_spin_orbitals
185 allocate(this % orbitals(this % num_SO), &
186 this % reference_dtrs(this % num_electrons), &
187 this % spin(this % num_SO), stat = ifail)
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')
196 do ido = 1, this % num_SO
197 this % orbitals(ido) = orbital % get_orbital_number(ido)
198 this % spin(ido) = orbital % get_spin(ido)
201 this % reference_dtrs(:) = option % reference_dtrs(:)
202 this % Nints = this % num_SO / 64 + 1
204 allocate(this % ref_int64dtrs(this % Nints), stat = ifail)
206 call master_memory % track_memory(storage_size(this % ref_int64dtrs)/8, &
207 size(this % ref_int64dtrs), ifail,
'CSFMANAGER::ref_int64dtrs')
209 this % ref_int64dtrs(:) = 0
210 write (stdout,
"('Nints = ',i4)") this % Nints
211 this % pt_orb => orbital
213 do ido = 1, this % num_electrons
214 call this % toggle_determinant(this % reference_dtrs(ido), this % ref_int64dtrs, this % Nints)
217 if ((shared_enabled .and. grid % lrank == local_master) .or. .not. shared_enabled)
then
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)
225 call setposwf(option % positron_flag, &
226 option % num_continuum, &
227 option % num_target_sym, &
231 call rdwf(this % num_csfs + 1, &
232 this % index_coeffs, &
234 this % packed_dtrs, &
235 this % length_dtrs, &
237 this % length_coeff, &
250 subroutine finalize_manager (this)
252 class(csfmanager),
intent(inout),
target :: this
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)
259 if (
allocated(this % orbitals))
then
260 call master_memory % free_memory(storage_size(this % orbitals)/8,
size(this % orbitals))
261 deallocate (this % orbitals)
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)
269 if (
allocated(this % spin))
then
270 call master_memory % free_memory(storage_size(this % spin)/8,
size(this % spin))
271 deallocate (this % spin)
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)
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)
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)
306 subroutine create_csfs (this, CSF, orbital_sequence, num_ci_target_sym, continuum_projection)
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(:)
312 integer :: orbital_sequence(this % num_csfs)
313 integer :: csf_det(this % num_electrons)
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
321 this % largest_num_configs = 0
323 write (stdout,
"('Creating our CSFs.........')", advance =
"no")
325 total_number_orbs = 0
327 do ido = 1, this % num_csfs
328 total_number_orbs = total_number_orbs + this % num_dtr_csf_out(ido)
331 write (stdout,
"('Total number of determinants = ',i12)") total_number_orbs
333 allocate(csf(this % num_csfs))
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)
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)
342 do ido = 1, this % num_csfs
343 num_orbs = this % num_dtr_csf_out(ido)
346 csf(ido) % manager => this
347 csf(ido) % num_orbitals = num_orbs
350 csf(ido) % orbital_sequence_number = orbital_sequence(ido)
352 this % largest_num_configs = max(this % largest_num_configs, num_orbs)
354 if (csf(ido) % orbital_sequence_number > 0)
then
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
362 count_number = count_number + 1
364 if (count_number > num_ci_target_sym(symmetry_idx) * 2)
then
365 symmetry_idx = symmetry_idx + 1
371 csf(ido) % is_continuum = .false.
375 orb_counter = orb_counter + 1
376 allocate(csf(ido) % orbital(num_orbs), stat = ifail)
378 csf(ido) % orbital(1) % id_dtr = orb_counter
379 csf(ido) % orbital(1) % csf => csf(ido)
381 if ((shared_enabled .and. grid % lrank == local_master) .or. .not. shared_enabled)
then
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(:)
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)
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
408 orb_counter = orb_counter + 1
409 csf(ido) % orbital(jdo) % id_dtr = orb_counter
410 csf(ido) % orbital(jdo) % csf => csf(ido)
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(:)
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)
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
434 write (stdout,
"('success!')")
435 write (stdout,
"('Largest number of configurations in a single CSF is ',i10)") this % largest_num_configs
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)
443 call master_memory % free_memory(storage_size(this % num_dtr_csf_out)/8,
size(this % num_dtr_csf_out))
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)
449 deallocate(this % num_dtr_csf_out)
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
656 integer :: tz, l, ispin, ishift, nperm, i, j, k, m, n, high, low
658 integer :: particles(1)
659 integer(longint) :: hole, particle, tmp, masklow, maskhigh
660 double precision,
parameter :: phase_dble(0:1) = (/ 1.d0, -1.d0 /)
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)
675 if (hole /= 0_8)
then
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)
690 nperm = nperm + popcnt(iand(det1(j), iand(not(ishft(1_8, n + 1)) + 1, ishft(1_8, m) - 1)))
692 nperm = nperm + popcnt(iand(det1(k), ishft(1_8, m) - 1)) &
693 + popcnt(iand(det1(j), not(ishft(1_8, n + 1)) + 1))
695 nperm = nperm + popcnt(det1(i))
698 phase = phase_dble(iand(nperm, 1))
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 /)
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)
733 idx_particle = idx_particle + 1
735 particles(idx_particle) = tz + ishift
736 particle = iand(particle, particle - 1_8)
738 do while (hole /= 0_8)
741 idx_hole = idx_hole + 1
742 holes(idx_hole) = tz + ishift
743 num_holes = num_holes + 1
745 hole = iand(hole, hole - 1_8)
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)
784 nperm = nperm + popcnt(iand(det1(j), iand(not(ishft(1_8, n + 1)) + 1, ishft(1_8, m) - 1)))
786 nperm = nperm + popcnt(iand(det1(k), ishft(1_8, m) - 1)) &
787 + popcnt(iand(det1(j), not(ishft(1_8, n + 1)) + 1))
789 nperm = nperm + popcnt(det1(l))
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
822 exc(2) = particles(1)
824 exc(4) = particles(2)
826 phase = phase_dble(iand(nperm, 1))