30 use const_gbl,
only: stdout
31 use global_utils,
only: indfunc, mprod
32 use integer_packing,
only: pack8ints, unpack8ints
33 use mpi_memory_gbl,
only: mpi_memory_allocate_real, mpi_memory_deallocate_real, mpi_memory_synchronize, local_master
34 use precisn,
only: longint, wp
51 real(wp),
pointer :: one_electron_integral(:)
52 real(wp),
pointer :: two_electron_integral(:)
54 real(wp),
allocatable :: one_electron_integral(:)
55 real(wp),
allocatable :: two_electron_integral(:)
58 integer :: one_electron_window
59 integer :: two_electron_window
61 integer :: num_one_electron_integrals
62 integer :: num_two_electron_integrals
64 real(wp),
allocatable :: xnuc(:),ynuc(:),znuc(:),charge(:)
65 character(len=8),
allocatable :: cname(:)
67 integer :: num_unique_pairs
68 integer(longint),
allocatable :: pair_labels(:,:)
69 integer,
allocatable :: num_orbitals_sym(:)
71 integer :: max_number_pair_sets
72 integer :: num_two_electron_blocks
73 integer,
allocatable :: one_electron_pointer(:)
74 integer,
allocatable :: two_electron_pointer(:)
76 integer :: num_pq, num_rs
77 integer :: num_pair_idx
80 integer,
allocatable :: pair_idx(:)
81 integer,
allocatable :: orbital_idx(:)
82 integer,
allocatable :: symmetry_idx(:)
107 use global_utils,
only: mprod
110 integer :: i, j, k, l
111 integer :: jend, lend
112 integer :: ij, kl, ijkltest
116 do i = 1, this % num_symmetries
120 if (this % num_symmetries == 1)
return
123 do i = 2, this % num_symmetries
132 do i = 2, this % num_symmetries
140 if (this % num_symmetries < 4)
return
142 do i = 2, this % num_symmetries
149 if (i == k .and. j == l) cycle
150 ij = mprod(i, j, 0, stdout)
151 kl = mprod(k, l, 0, stdout)
152 ijkltest = mprod(ij, kl, 0, stdout)
153 if (ijkltest /= 1) cycle
165 integer :: i, j, k, l
166 integer :: jend, lend
167 integer :: ij, kl, ijkltest
168 integer :: pair_number
172 write (stdout, *)
'All IIII blocks of integrals'
174 do i = 1, this % num_symmetries
176 pair_number = pair_number + 1
177 call pack8ints(j, j, j, j, 0, 0, 0, 0, this % pair_labels(1, pair_number))
178 write (stdout, *) pair_number, j, j, j, j
181 if (this % num_symmetries == 1)
return
184 write (stdout, *)
'All IJIJ blocks'
185 do i = 2, this % num_symmetries
189 pair_number = pair_number + 1
190 call pack8ints(jend, lend, jend, lend, 0, 0, 0, 0, this % pair_labels(1, pair_number))
191 write (stdout, *) pair_number, jend, lend, jend, lend
196 write (stdout, *)
'All IIJJ blocks'
197 do i = 2, this % num_symmetries
201 pair_number = pair_number + 1
202 call pack8ints(jend, jend, lend, lend, 0, 0, 0, 0, this % pair_labels(1, pair_number))
203 write (stdout, *) pair_number, jend, jend, lend, lend
207 if (this % num_symmetries < 4)
return
208 write (stdout, *)
'All IJKL blocks'
209 do i = 2, this % num_symmetries
216 if (i == k .and. j == l) cycle
217 ij = mprod(i, j, 0, stdout)
218 kl = mprod(k, l, 0, stdout)
219 ijkltest = mprod(ij, kl, 0, stdout)
220 if (ijkltest /=1 ) cycle
221 pair_number = pair_number + 1
222 call pack8ints(i - 1, j - 1, k - 1, l - 1, 0, 0, 0, 0, this % pair_labels(1, pair_number))
223 write (stdout, *) pair_number, i - 1, j - 1, k - 1, l - 1
235 integer :: label(8), num_orbs(4)
236 integer :: ido, jdo, nam1, nam2, num_pq, num_rs, num_pr
237 integer :: block_number, err
239 this % num_two_electron_blocks = indfunc(this % num_symmetries + 1, 0)
240 this % num_two_electron_blocks = indfunc(this % num_two_electron_blocks + 1, 0)
241 this % num_two_electron_integrals = 1
243 allocate(this % one_electron_pointer(this % num_symmetries + 1), stat = err)
245 call master_memory % track_memory(kind(this % one_electron_pointer), &
246 size(this % one_electron_pointer), err,
'SWEDEN::One electron pointer')
248 this % one_electron_pointer(1) = 1
250 do ido = 1, this % num_symmetries
251 this % one_electron_pointer(ido + 1) = this % one_electron_pointer(ido) + indfunc(this % num_orbitals_sym(ido) + 1, 0)
254 this % num_one_electron_integrals = this % one_electron_pointer(this % num_symmetries + 1) - 1
256 allocate(this % two_electron_pointer(this % num_two_electron_blocks), stat = err)
258 call master_memory % track_memory(kind(this % two_electron_pointer), &
259 size(this % two_electron_pointer), err,
'SWEDEN::Two electron pointer')
261 this % two_electron_pointer(:) = 0
263 write (stdout,
"(//,5X,'Integral Storage Table follows : ',/)")
265 do ido = 1, this % num_unique_pairs
267 call unpack8ints(this % pair_labels(1, ido), label)
270 label(jdo) = label(jdo) + 1
271 num_orbs(jdo) = this % num_orbitals_sym(label(jdo))
275 if (label(1) == label(2))
then
276 num_pq = indfunc(num_orbs(1) + 1, 0)
278 num_pq = num_orbs(1) * num_orbs(2)
280 nam1 = indfunc(max(label(1), label(2)) + 1, 0) - abs(label(1) - label(2))
283 if (label(3) == label(4))
then
284 num_rs = indfunc(num_orbs(3) + 1, 0)
286 num_rs = num_orbs(3) * num_orbs(4)
289 nam2 = indfunc(max(label(3), label(4)) + 1, 0) - abs(label(3)-label(4))
290 block_number = indfunc(max(nam1, nam2) + 1, 0) - abs(nam1 - nam2)
291 this % two_electron_pointer(block_number) = this % num_two_electron_integrals
293 write (stdout,
"(3X,4I3,1X,I5,1X,I10)") label(1), label(2), label(3), label(4), &
294 block_number, this % num_two_electron_integrals
296 if (nam1 == nam2)
then
297 num_pr = indfunc(num_pq + 1, 0)
299 num_pr = num_pq * num_rs
302 this % max_number_pair_sets = max(this % max_number_pair_sets, num_pq, num_rs)
303 if (num_pr > 0) this % num_two_electron_integrals = this % num_two_electron_integrals + num_pr + 1
307 this % num_two_electron_integrals = this % num_two_electron_integrals - 1
309 write (stdout,
"('No of 1 electron integrals = ',i8)") this % num_one_electron_integrals
310 write (stdout,
"('No of 2 electron integrals = ',i8)") this % num_two_electron_integrals
317 integer :: total_num_orbitals, orbital, iorb, isym, max_orbitals, ifail
319 total_num_orbitals = sum(this % num_orbitals_sym(:))
321 allocate(this % orbital_idx(total_num_orbitals), this % symmetry_idx(total_num_orbitals), stat = ifail)
323 call master_memory % track_memory(kind(this % orbital_idx),
size(this % orbital_idx), ifail,
'SWEDEN::pair_idx')
324 call master_memory % track_memory(kind(this % symmetry_idx),
size(this % symmetry_idx), ifail,
'SWEDEN::symmetry_idx')
326 this % orbital_idx(:) = 0
327 this % symmetry_idx(:) = 0
332 do isym = 1, this % num_symmetries
333 max_orbitals = max(max_orbitals, this % num_orbitals_sym(isym)**2)
334 do iorb = 1, this % num_orbitals_sym(isym)
335 orbital = orbital + 1
336 this % orbital_idx(orbital) = iorb
337 this % symmetry_idx(orbital) = isym - 1
341 this % num_pair_idx = max((this % num_symmetries**2 + this % num_symmetries) / 2 + 1, &
342 this % max_number_pair_sets, &
351 integer :: ido, idx, ifail
353 allocate(this % pair_idx(0:this%num_pair_idx), stat = ifail)
354 call master_memory % track_memory(kind(this % pair_idx),
size(this % pair_idx), ifail,
'SWEDEN::pair_idx')
357 this % pair_idx(0) = 0
358 do ido = 1, this % num_pair_idx
359 this % pair_idx(ido) = idx
368 class(
options),
intent(in) :: option
369 integer :: ido, jdo, ifail
371 write (stdout,
"('Using SWEDEN integral')")
373 this % num_unique_pairs = 0
374 this % num_symmetries = option % num_syms
376 allocate(this % num_orbitals_sym(this % num_symmetries), stat = ifail)
377 call master_memory % track_memory(kind(this % num_orbitals_sym), &
378 size(this % num_orbitals_sym), ifail,
'SWEDEN::num_orbitals_sym')
380 this % num_orbitals_sym(:) = option % num_electronic_orbitals_sym(:)
381 this % num_unique_pairs = this % count_num_pairs()
383 allocate(this % pair_labels(2, this % num_unique_pairs), stat = ifail)
384 call master_memory % track_memory(kind(this % pair_labels), &
385 size(this % pair_labels), ifail,
'SWEDEN::Pair labels')
387 write (stdout, *)
'SWEDEN Integral -- Generating pairs'
390 call this % generate_pairs
391 call this % generate_orbital_index
392 call this % generate_pair_index
394 write (stdout,
"('1',/,10x,'D2h Two Electron Integral Box ','Information')")
395 write (stdout,
"( /,10x,'No. of Boxes = ',i4,/)") this % num_unique_pairs
396 write (stdout,
"( 10x,'Box Descripti ons: ',/)")
397 write (stdout, *)
'SWEDEN Integral -- Generating block pointers'
399 call this % generate_pointer_table
411 use params,
only: ctrans1, cpoly
412 use global_utils,
only: intape
413 use scatci_routines,
only: search, chn2e
416 integer,
intent(in) :: iounit
419 character(len=4),
dimension(4) :: LABEL
420 character(len=4),
dimension(33) :: NAM1
421 integer(longint),
dimension(8) :: NAO, NCORE
422 integer(longint),
dimension(20) :: NHE, NOB
423 real(wp),
allocatable :: xtempe(:), xtempp(:)
425 integer(longint) :: I, NSYM, IONEIN
426 integer :: ifail, k, l, ind, ii
431 write (stdout, *)
' Loading SWEDEN Integrals into core'
436 this % one_electron_window = mpi_memory_allocate_real(this % one_electron_integral, 2 * this % num_one_electron_integrals)
442 this % two_electron_window = mpi_memory_allocate_real(this % two_electron_integral, this % num_two_electron_integrals)
444 call mpi_memory_synchronize(this % one_electron_window)
445 call mpi_memory_synchronize(this % two_electron_window)
451 read (iounit) (label(i), i = 1, 4)
453 write (stdout, fmt=
'(/,5X,''Label on Sweden tape = '',4A)') (label(i), i = 1, 4)
455 read (iounit) nsym, en, (nao(i), nob(i), ncore(i), i = 1, nsym)
457 write (stdout,
"(' NSYM =',I5,3X,'NOB =',20I5)") nsym, (nob(i), i = 1, nsym)
458 write (stdout,
"(/,10x,'Nuclear repulsion energy = ',f15.7)") en
459 write (stdout,
"(/10x,' NAO NMO NCORE',/,(10x,2I5,i7))") (nao(i), nob(i), ncore(i), i = 1, nsym)
461 this % core_energy = en
464 if (ncore(i) /= 0)
then
465 write (stdout,
"(10X,'Non-zero CORE in Sweden. Sym no. = ',i5,' Ncore = ',I5,/)") i, ncore(i)
466 write (stdout,
"(10X,'Core Energy = ',F15.7,' Hartrees ',/)") en
474 ionein = ionein + (nob(i) * (nob(i) + 1)) / 2
477 if ((this % two_electron_window == -1 .and. this % one_electron_window == -1) .or. grid % lrank == local_master)
then
478 call search(iounit, ctrans1, ifail)
480 call intape(iounit, this % one_electron_integral, this % num_one_electron_integrals)
483 if (this % positron_flag /= 0)
then
484 allocate(xtempe(this % num_one_electron_integrals), stat = ifail)
485 call master_memory % track_memory(kind(xtempe),
size(xtempe), ifail,
'SWEDEN::xtempe')
487 allocate(xtempp(this % num_one_electron_integrals), stat = ifail)
488 call master_memory % track_memory(kind(xtempp),
size(xtempp), ifail,
'SWEDEN::xtempp')
490 xtempe(1:this % num_one_electron_integrals) = this % one_electron_integral(1:this % num_one_electron_integrals)
491 call intape(iounit, xtempp, this % num_one_electron_integrals)
494 if (this % quantamoln_flag)
then
499 write (413,
"(i6,3x,i6,3x,i6,5x,f12.7)") k - 1, l, l, xtempp(ind)
504 do i = 1, this % num_one_electron_integrals
505 this % one_electron_integral(i + this % num_one_electron_integrals) = 2.0_wp * xtempp(i) - xtempe(i)
514 call chn2e(iounit, stdout, this % two_electron_integral, this % num_two_electron_integrals)
517 1600
format(
' Integrals read successfully:', /, &
518 ' 1-electron integrals, NINT1e =',i10, /, &
519 ' 2-electron integrals, NINT2e =',i10)
521 write (stdout, 1600) this % num_one_electron_integrals, this % num_two_electron_integrals
523 this % nhe(:) = nhe(:)
527 this % nhe(:) = nob(:)
528 this % dtnuc(22:41) = 0
529 this % dtnuc(2:21) = 0
533 call search(iounit, cpoly, ifail)
534 read (iounit) this % nnuc
536 allocate(this % xnuc(this % nnuc), this % ynuc(this % nnuc), this % znuc(this % nnuc), &
537 this % charge(this % nnuc), this % CNAME(this % nnuc))
539 do i = 1, this % nnuc
540 read (iounit) this % cname(i), ii, this % xnuc(i), this % ynuc(i), this % znuc(i), this % charge(i)
547 call mpi_memory_synchronize(this % one_electron_window)
548 call mpi_memory_synchronize(this % two_electron_window)
555 integer,
intent(in) :: i, j, k, l, m
557 integer :: ia, ib, integral_idx
563 integral_idx = this % get_one_electron_index(j, l, m)
564 coeff = this % one_electron_integral(integral_idx)
566 integral_idx = this % get_two_electron_index(i, j, k, l, m)
567 coeff = this % two_electron_integral(integral_idx)
575 integer,
intent(in) :: i, j, pos
576 integer :: m, p, q, ii, jj
578 ii = this % orbital_mapping(i)
579 jj = this % orbital_mapping(j)
581 p = this % orbital_idx(ii)
582 q = this % orbital_idx(jj)
584 m = this % symmetry_idx(ii)
600 integer,
intent(in) :: i, j, k, l, m
602 integer,
dimension(4) :: orb_num_sym, orb_sym, num_orb_sym
605 integer :: mpq, mrs, mpr, aaa, block_pointer, nobrs, nobpq, ipqrs
607 mapped(1) = this % orbital_mapping(i)
608 mapped(2) = this % orbital_mapping(j)
609 mapped(3) = this % orbital_mapping(k)
610 mapped(4) = this % orbital_mapping(l)
613 orb_num_sym(ido) = this % orbital_idx(mapped(ido))
614 symmetry = this % symmetry_idx(mapped(ido)) + 1
615 orb_sym(ido) = symmetry
616 num_orb_sym(ido) = this % num_orbitals_sym(symmetry)
621 mpq = this % pair_idx(orb_sym(1)) + orb_sym(2)
622 mrs = this % pair_idx(orb_sym(3)) + orb_sym(4)
623 mpr = this % pair_idx(mpq) + mrs
625 if (mpq == mrs .and. &
626 ( orb_num_sym(1) < orb_num_sym(3) .or. &
627 (orb_num_sym(1) == orb_num_sym(3) .and. orb_num_sym(2) < orb_num_sym(4)) .or. &
628 (orb_num_sym(1) < orb_num_sym(3) .and. orb_num_sym(2) < orb_num_sym(4)) ))
then
629 aaa = orb_num_sym(1) ; orb_num_sym(1) = orb_num_sym(3) ; orb_num_sym(3) = aaa
630 aaa = orb_sym(1) ; orb_sym(1) = orb_sym(3) ; orb_sym(3) = aaa
631 aaa = orb_num_sym(2) ; orb_num_sym(2) = orb_num_sym(4) ; orb_num_sym(4) = aaa
632 aaa = orb_sym(2) ; orb_sym(2) = orb_sym(4) ; orb_sym(4) = aaa
635 block_pointer = this % two_electron_pointer(mpr) - 1
637 if (orb_sym(1) == orb_sym(2) .and. orb_sym(3) == orb_sym(4) .and. orb_sym(1) == orb_sym(3))
then
638 ipqrs = this % pair_idx(orb_num_sym(1)) + orb_num_sym(2)
640 ipqrs = this % pair_idx(ipqrs - 1) + ipqrs - 1 + this % pair_idx(orb_num_sym(3)) + orb_num_sym(4)
641 else if (orb_sym(1) == orb_sym(3) .and. orb_sym(2) == orb_sym(4))
then
642 ipqrs = (orb_num_sym(1) - 1) * num_orb_sym(2) + orb_num_sym(2) - 1
644 ipqrs = this % pair_idx(ipqrs) + ipqrs + (orb_num_sym(3) - 1) * num_orb_sym(4) + orb_num_sym(4)
645 else if (orb_sym(1) == orb_sym(2) .and. orb_sym(3) == orb_sym(4))
then
646 nobrs = this % pair_idx(num_orb_sym(3)) + num_orb_sym(3)
647 ipqrs = (this % pair_idx(orb_num_sym(1)) + orb_num_sym(2) - 1)*nobrs + this % pair_idx(orb_num_sym(3)) + orb_num_sym(4)
649 nobrs = num_orb_sym(3) * num_orb_sym(4)
650 ipqrs = ((orb_num_sym(1) - 1) * num_orb_sym(2) + orb_num_sym(2) - 1) * nobrs + (orb_num_sym(3) - 1) * num_orb_sym(4) &
662 use params,
only: cpoly
663 use scatci_routines,
only: search
666 integer,
intent(in) :: iounit
667 integer :: ido, ifail, i
669 do ido = 1, this % nnuc
670 write (iounit) this % cname(ido), this % xnuc(ido), this % ynuc(ido), this % znuc(ido), this % charge(ido)
673 write (stdout,
"(/,' Nuclear data X Y Z Charge',/,(3x,a8,2x,4F10.6))") &
674 (this % cname(i), this % xnuc(i), this % ynuc(i), this % znuc(i), this % charge(i), i = 1, this % nnuc)
687 call mpi_memory_deallocate_real(this % one_electron_integral, &
688 size(this % one_electron_integral), this % one_electron_window)
689 call mpi_memory_deallocate_real(this % two_electron_integral, &
690 size(this % two_electron_integral), this % two_electron_window)
697 if (
allocated(this % pair_labels))
then
698 call master_memory % free_memory(kind(this % pair_labels),
size(this % pair_labels))
699 deallocate(this % pair_labels)
703 if (
allocated(this % num_orbitals_sym))
then
704 call master_memory % free_memory(kind(this % num_orbitals_sym),
size(this % num_orbitals_sym))
705 deallocate(this % num_orbitals_sym)
708 if (
allocated(this % one_electron_pointer))
then
709 call master_memory % free_memory(kind(this % one_electron_pointer),
size(this % one_electron_pointer))
710 deallocate(this % one_electron_pointer)
713 if (
allocated(this % two_electron_pointer))
then
714 call master_memory % free_memory(kind(this % two_electron_pointer),
size(this % two_electron_pointer))
715 deallocate(this % two_electron_pointer)
718 if (
allocated(this % pair_idx))
then
719 call master_memory % free_memory(kind(this % pair_idx),
size(this % pair_idx))
720 deallocate(this % pair_idx)
723 if (
allocated(this%orbital_idx))
then
724 call master_memory % free_memory(kind(this % orbital_idx),
size(this % orbital_idx))
725 deallocate(this %orbital_idx)
728 if (
allocated(this % symmetry_idx))
then
729 call master_memory % free_memory(kind(this % symmetry_idx),
size(this % symmetry_idx))
730 deallocate(this % symmetry_idx)