MPI-SCATCI  2.0
An MPI version of SCATCI
ALCHEMY_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 
29 
30  use precisn, only: longint, wp
31  use const_gbl, only: stdout
34  use options_module, only: options
35  use global_utils, only: indfunc
36  use integer_packing, only: pack8ints, unpack8ints
37  use scatci_routines, only: rdint, rdints
38 
39  implicit none
40 
41  public alchemyintegral
42 
43  private
44 
45  type, extends(baseintegral) :: alchemyintegral
46 
47  !Our integrals
48  real(wp), pointer :: one_electron_integral(:)
49  real(wp), pointer :: two_electron_integral(:)
50  integer :: num_one_electron_integrals
51  integer :: num_two_electron_integrals
52  integer :: integral_ordering
53  integer :: use_scf = 0
54 
56  integer :: num_unique_pairs
58  integer(longint), allocatable :: pair_labels(:,:)
60  integer, allocatable :: num_orbitals_sym(:)
61 
62  integer :: max_number_pair_sets
63  integer, allocatable :: num_two_electron_blocks(:)
64  integer, allocatable :: num_one_electron_blocks(:)
65  integer, allocatable :: one_electron_pointer(:)
66  integer, allocatable :: two_electron_pointer(:)
67 
68  integer :: num_pq, num_rs
69  integer :: num_pair_idx
70 
72  integer, allocatable :: pair_idx(:)
73  integer, allocatable :: orbital_idx(:)
74  integer, allocatable :: symmetry_idx(:)
75 
76  contains
77 
78  procedure, public :: initialize_self => initialize_alchemy
79  procedure, public :: finalize_self => finalize_alchemy
80  procedure, public :: load_integrals => load_integrals_alchemy
81  procedure, public :: get_integral_ijklm => get_integral_alchemy
82  procedure, public :: destroy_integrals => destroy_integral_alchemy
83  procedure, public :: write_geometries => write_geometries_alchemy
84  procedure, private :: count_num_pairs
85  procedure, private :: generate_pairs
86  procedure, private :: generate_pointer_table
87  procedure, private :: generate_pair_index
88  procedure, private :: generate_orbital_index
89  !procedure, private :: read_one_electron_integrals
90  !procedure, private :: read_two_electron_integrals
91  procedure, private :: get_one_electron_index
92  procedure, private :: get_two_electron_index
93 
94  end type alchemyintegral
95 
96 contains
97 
102  integer function count_num_pairs (this)
103  class(alchemyintegral) :: this
104  integer :: pair_number, mdo, ndo, ido, jdo, nn(5), msym, begin, pass, istart, num_sym
105 
106  num_sym = this % num_symmetries
107  msym = num_sym + num_sym - 1
108  count_num_pairs = 0
109 
110  if (this % use_SCF == 0) then
111  !Standard Alchemy list
112  do mdo = 1, msym
113  begin = mdo / 2 + 1
114  nn(1) = mdo - 1
115  do ndo = begin, num_sym
116  nn(2) = ndo - 1
117  nn(3) = abs(ndo - mdo)
118  do ido = begin, ndo
119  nn(4) = ido - 1
120  nn(5) = abs(ido - mdo)
122  end do
123  end do
124  end do
125  else
126  do pass = 1, 2
127  do mdo = pass, msym
128  begin = mdo / 2 + 1
129  nn(1) = mdo - 1
130  do ndo = mdo / 2 + pass, num_sym
131  nn(2) = ndo - 1
132  nn(3) = abs(ndo - mdo)
133  istart = ndo + 1 - pass
134  if (mdo /= 1 .and. pass == 1) begin = istart
135  do ido = begin, istart
136  nn(4) = ido - 1
137  nn(5) = abs(ido - mdo)
139  end do
140  end do
141  end do
142  msym = msym - 2
143  end do
144  end if
145 
146  end function count_num_pairs
147 
148 
153  subroutine generate_pairs (this)
154  class(alchemyintegral) :: this
155  integer :: pair_number, npqrs, nn(5), mdo, ndo, ido, jdo, I, IBGN, IF, IPASS, J, K, M, MSYM, N, NBGN, begin, pass, num_sym
156  integer(longint) :: packed_label(2)
157  integer,allocatable :: lpqrs(:,:)
158 
159  num_sym = this % num_symmetries
160  msym = num_sym + num_sym - 1
161 
162  allocate(lpqrs(5, this % num_unique_pairs))
163 
164  pair_number = 0
165 
166  if (this % use_SCF == 0) then
167  !Standard Alchemy list
168  do mdo = 1, msym
169  begin = mdo / 2 + 1
170  nn(1) = mdo - 1
171  do ndo = begin, num_sym
172  nn(2) = ndo - 1
173  nn(3) = abs(ndo - mdo)
174  do ido = begin, ndo
175  nn(4) = ido - 1
176  nn(5) = abs(ido - mdo)
177  pair_number = pair_number + 1
178  do jdo = 1, 5
179  lpqrs(jdo,pair_number) = nn(jdo)
180  end do
181  end do
182  end do
183  end do
184  else
185  do ipass = 1, 2
186  do m = ipass, msym
187  ibgn = m / 2 + 1
188  nn(1) = m - 1
189  do n = m / 2 + ipass, num_sym
190  nn(2) = n - 1
191  nn(3) = abs(n - m)
192  IF = n + 1 - ipass
193  if (m /= 1 .and. ipass == 1) ibgn = IF
194  do i = ibgn, IF
195  nn(4) = i - 1
196  nn(5) = abs(i - m)
197  pair_number = pair_number + 1
198  do j = 1, 5
199  lpqrs(j,pair_number) = nn(j)
200  end do
201  end do
202  end do
203  end do
204  msym = msym - 2
205  end do
206  end if
207 
208  if (pair_number /= this % num_unique_pairs) stop "Error in pair count"
209 
210  do ido = 1, this % num_unique_pairs
211  call pack8ints(lpqrs(2,ido), lpqrs(3,ido), lpqrs(4,ido), lpqrs(5,ido), 0, 0, 0, 0, packed_label)
212  if (lpqrs(2,ido) + lpqrs(3,ido) /= lpqrs(1,ido)) packed_label(1) = -packed_label(1)
213  this % pair_labels(1,ido) = packed_label(1)
214  this % pair_labels(2,ido) = packed_label(2)
215  end do
216 
217  write (stdout, "(' NMPRS=',I5//' LPQRS'/(2X,5(2X,I3)))") this % num_unique_pairs, &
218  ((lpqrs(ido,jdo), ido = 1, 5), jdo = 1, this % num_unique_pairs)
219 
220  deallocate(lpqrs)
221 
222  end subroutine generate_pairs
223 
224 
225  subroutine generate_pointer_table (this)
226  class(alchemyintegral) :: this
227  integer :: label(8), num_orbs(4), size_two_pointer
228  integer(longint) :: packed_label(2)
229  integer :: ido, jdo, nam1, nam2, num_pq, num_rs, num_pr
230  integer :: block_number, mvl, md, mm, m, mpa, mpb, mp
231 
232  !Lets allocate our one elctron pointers
233 
234  size_two_pointer = this % num_symmetries * (this % num_symmetries + 1) / 2
235  do ido = 2, this % num_symmetries
236  size_two_pointer = size_two_pointer + (this % num_symmetries - ido + 2) * (this % num_symmetries - ido + 1)
237  end do
238 
239  allocate(this % num_two_electron_blocks(2 * this % num_symmetries + 1))
240 
241  this % num_two_electron_blocks(1) = 0
242  this % num_two_electron_blocks(2) = indfunc(this % num_symmetries + 1, 0)
243 
244  mvl = 2 * this% num_symmetries - 1
245 
246  do ido = 2, mvl
247  md = (ido - 2) / 2
248  mm = this % num_symmetries - md - 1
249  this% num_two_electron_blocks(ido + 1) = this % num_two_electron_blocks(ido) + indfunc(mm + 1, 0)
250  end do
251 
252  allocate(this % one_electron_pointer(this % num_symmetries + 1))
253 
254  this % one_electron_pointer(1) = 1
255 
256  do ido = 1, this % num_symmetries
257  this % one_electron_pointer(ido + 1) = this % one_electron_pointer(ido) + indfunc(this % num_orbitals_sym(ido) + 1, 0)
258  end do
259 
260  allocate(this % two_electron_pointer(this % num_two_electron_blocks(mvl + 1)))
261 
262  do ido =1, this % num_two_electron_blocks(mvl + 1)
263  this % two_electron_pointer(ido) = 0
264  end do
265 
266  this % num_one_electron_integrals = this % one_electron_pointer(this % num_symmetries + 1) - 1
267  this % num_two_electron_integrals = 1
268 
269  write (stdout, "(//,5X,'Integral Storage Table follows : ',/)")
270 
271  do ido = 1, this % num_unique_pairs
273  packed_label(1) = this % pair_labels(1,ido)
274  packed_label(2) = this % pair_labels(2,ido)
275 
276  if (packed_label(1) > 0) then
277  md = 0
278  else
279  md = 1
280  packed_label(1) = -packed_label(1)
281  end if
282 
283  call unpack8ints(packed_label, label)
284 
285  do jdo = 1, 4
286  num_orbs(jdo) = this % num_orbitals_sym(label(jdo) + 1)
287  end do
288 
289  if (md == 0) then
290  m = label(1) + label(2)
291  else
292  m = abs(label(1) - label(2))
293  end if
294 
295  if (m == 0) then
296  md = -1
297  else
298  md = (m - 1) / 2
299  end if
300  mpa = label(1) - md
301  mpb = label(3) - md
302  if (mpa < mpb) then
303  !MP=(MPB*(MPB-1))/2+MPA+this%num_two_electron_blocks(M+1)
304  mp = (mpb * (mpb - 1)) / 2 + mpa + this % num_two_electron_blocks(m + 1)
305  else
306  !MP=(MPA*(MPA-1))/2+MPB+this%num_two_electron_blocks(M+1)
307  mp = (mpa * (mpa - 1)) / 2 + mpb + this % num_two_electron_blocks(m + 1)
308  end if
309 
310  this % two_electron_pointer(mp) = this % num_two_electron_integrals
311  !write(stdout,"('ELEC ',3i8)") mp,this%num_two_electron_integrals,this%two_electron_pointer(mp)
312  if (label(1) == label(2)) then
313  num_pq = indfunc(num_orbs(1) + 1, 0)
314  else
315  num_pq = num_orbs(1) * num_orbs(2)
316  end if
317 
318  if (label(3) == label(4)) then
319  num_rs = indfunc(num_orbs(3) + 1, 0)
320  else
321  num_rs = num_orbs(3) * num_orbs(4)
322  end if
323  this % max_number_pair_sets = max(this % max_number_pair_sets, num_pq, num_rs)
324  if (label(1) == label(3) .and. label(2) == label(4)) then
325  this % num_two_electron_integrals = this % num_two_electron_integrals + indfunc(num_pq + 1, 0)
326  else
327  this % num_two_electron_integrals = this % num_two_electron_integrals + num_pq * num_rs
328  end if
329  end do
330  !Move to the end of the block
331  this % num_two_electron_integrals = this % num_two_electron_integrals - 1
332 
333  write (stdout, "('No of 1 electron integrals = ',i8)") this % num_one_electron_integrals
334  write (stdout, "('No of 2 electron integrals = ',i8)") this % num_two_electron_integrals
335 
336  end subroutine generate_pointer_table
337 
338 
339  subroutine generate_orbital_index (this)
340  class(alchemyintegral) :: this
341  integer :: total_num_orbitals, orbital, iorb, isym, max_orbitals, ifail
342 
343  total_num_orbitals = sum(this % num_orbitals_sym(:))
344 
345  allocate(this % orbital_idx(total_num_orbitals), &
346  this % symmetry_idx(total_num_orbitals), stat = ifail)
347 
348  call master_memory % track_memory(kind(this % orbital_idx), size(this % orbital_idx), ifail, 'SWEDEN::pair_idx')
349  call master_memory % track_memory(kind(this % symmetry_idx), size(this % symmetry_idx), ifail, 'SWEDEN::symmetry_idx')
350 
351  this % orbital_idx(:) = 0
352  this % symmetry_idx(:) = 0
353 
354  max_orbitals = 0
355  orbital = 0
356 
357  do isym = 1,this % num_symmetries
358  max_orbitals = max(max_orbitals, this % num_orbitals_sym(isym)**2)
359  do iorb = 1, this % num_orbitals_sym(isym)
360  orbital = orbital + 1
361  this % orbital_idx(orbital) = iorb
362  this % symmetry_idx(orbital) = isym - 1
363  end do
364  end do
365 
366  this % num_pair_idx = max((this% num_symmetries**2 + this % num_symmetries) / 2 + 1, &
367  this % max_number_pair_sets, &
368  300, &
369  max_orbitals)
370 
371  end subroutine generate_orbital_index
372 
373 
374  subroutine generate_pair_index (this)
375  class(alchemyintegral) :: this
376  integer :: ido, idx, ifail
377 
378  allocate(this % pair_idx(0:this % num_pair_idx), stat = ifail)
379 
380  call master_memory % track_memory(kind(this % pair_idx), size(this % pair_idx), ifail, 'SWEDEN::pair_idx')
381 
382  idx = 0
383  this % pair_idx(0) = 0
384  do ido = 1, this % num_pair_idx
385  this % pair_idx(ido) = idx
386  idx = idx + ido
387  end do
388 
389  end subroutine generate_pair_index
390 
391 
392  subroutine initialize_alchemy (this, option)
393  class(alchemyintegral) :: this
394  class(options), intent(in) :: option
395  integer :: ido, jdo
396 
397  write (stdout, "('Using ALCHEMY integral')")
398 
399  this % num_unique_pairs = 0
400  this % num_symmetries = option % num_syms
401  this % integral_ordering = option % integral_ordering
402  this % use_SCF = option % use_SCF
403 
404  !Lets begin counting
405 
406  allocate(this % num_orbitals_sym(this % num_symmetries))
407 
408  this % num_orbitals_sym(:) = option % num_electronic_orbitals_sym(:)
409  this % num_unique_pairs = this % count_num_pairs()
410 
411  write (stdout, *) 'Number of pairs =',this % num_unique_pairs
412 
413  allocate(this % pair_labels(2, this % num_unique_pairs))
414 
415  write (stdout,*) 'ALCHEMY Integral -- Generating pairs'
416 
418  call this % generate_pairs
419  call this % generate_orbital_index
420  call this % generate_pair_index
421 
422  write (stdout, "('1',/,10x,'D2h Two Electron Integral Box ','Information')")
423  write (stdout, "( /,10x,'No. of Boxes = ',i4,/)") this % num_unique_pairs
424  write (stdout, "( 10x,'Box Descriptions: ',/)")
425  write (stdout, *) 'ALCHEMY Integral -- Generating block pointers'
426 
427  call this % generate_pointer_table
428 
429  end subroutine initialize_alchemy
430 
431 
432  subroutine finalize_alchemy (this)
433  class(alchemyintegral) :: this
434  end subroutine finalize_alchemy
435 
436 
439  subroutine load_integrals_alchemy (this, iounit)
440  use params, only: ctrans1, cpoly
441  use global_utils, only: intape
442  use scatci_routines, only: search, chn2e
443 
444  class(alchemyintegral) :: this
445  integer, intent(in) :: iounit
446  ! ALCHEMY header arrays
447  character(len=4), dimension(4) :: LABEL
448  character(len=4), dimension(33) :: NAM1
449  integer, dimension(8) :: NAO, NCORE
450  integer, dimension(20) :: NHE, NOB
451  real(wp), allocatable :: xtempe(:), xtempp(:)
452 
453  integer :: I, NSYM, IONEIN, NT, ND, LTRI, NNUC, N, IA, IAT, NALM
454  integer :: ifail, k, l, ind, IMAX, IMIN, NEND, IMAXP
455  real(wp) :: EN, SIGN, x1e
456  real(wp) :: vsmall = 1.e-10_wp
457  real(kind=wp), dimension(20) :: charg, xnuc
458 
459  write (stdout, *) ' Loading ALCHEMY Integrals into core'
460 
461  !allocate the integral space, this will be replaced by a caching system
462  !One electron integral
463  allocate(this % one_electron_integral(2 * this % num_one_electron_integrals), stat = ifail)
464  if (ifail /= 0) stop "could not allocate 1 electron integral"
465 
466  !Two electron integral
467  allocate(this % two_electron_integral(this % num_two_electron_integrals), stat = ifail)
468  if (ifail /= 0) stop "could not allocate 2 electron integral"
469 
470  read (iounit) (nam1(i), i = 1, 33), nsym, nt, nnuc, nd, ltri, (nob(i), i = 1, nsym), (nd, i = 1, nt), (nd, i = 1, nt), &
471  (charg(i), i = 1, nnuc), (xnuc(i), i = 1, nnuc)
472 
473  write (stdout, "(/' Transformed integrals read:',/5X,30A4)") (nam1(i), i = 1, 30)
474 
475  en = 0.0_wp
476  do n = 2, nnuc
477  if (abs(charg(n)) < vsmall) cycle
478  do i = 1, n - 1
479  if (abs(charg(i)) < vsmall) cycle
480  en = en + charg(n) * charg(i) / abs(xnuc(n) - xnuc(i))
481  end do
482  end do
483 
484  this % core_energy = en
485 
486  write (stdout, "(' NSYM =',I5,3X,'NOB =',20I5)") nsym, (nob(i), i = 1, nsym)
487  write (stdout, "(/,10x,'Nuclear repulsion energy = ',f15.7)") en
488  write (stdout, "(' NNUC =',I5,3X,'CHARG=',10F10.0/21X,10F10.0)") nnuc, (charg(i), i = 1, nnuc)
489  write (stdout, "(15X,'XNUC =',10F10.5/21X,10F10.5)") (xnuc(i), i = 1, nnuc)
490 
491  ia = 1
492  call rdints(iounit, nsym, nob, ltri, this % num_one_electron_integrals + 1, ia, this % one_electron_integral, nalm)
493  ! print*,this%one_electron_integral(1)
494  if (nalm /= 0) stop "UNABLE TO READ ALCHEMY INTEGRALS"
495 
496  ia = 1
497  call rdints(iounit, nsym, nob, ltri, this % num_one_electron_integrals + 1, ia, this % one_electron_integral, nalm)
498  ! print*,this%one_electron_integral(1)
499  if (nalm /= 0) stop "UNABLE TO READ ALCHEMY INTEGRALS"
500 
501  iat = ia - 1
502  ia = 1
503  allocate(xtempp(this % num_one_electron_integrals))
504  call rdints(iounit, nsym, nob, ltri, this % num_one_electron_integrals + 1, ia, xtempp, nalm)
505  ! print*,xtempp(1)
506 
507  if (nalm == 0) then
508  if (iat /= this % num_one_electron_integrals) stop "INCORRECT NUMBER OF INTEGRALS"
509  sign = 1.0
510  if(abs(this % positron_flag) == 1) sign = -1.0
511  do i = 1, iat
512  x1e = this % one_electron_integral(i) + xtempp(i)
513  this % one_electron_integral(i) = x1e
514  end do
515  !this%one_electron_integral(1:) = this%one_electron_integral(:) + xtempp(:)
516  end if
517 
518  !print *,this%one_electron_integral(1)
519 
520  imax = 0
521  imin = 0
522  nend = -1
523  imaxp = 0
524  ia = 1
525  call rdint(imin, imax, nend, imaxp, iounit, ia, this % num_two_electron_integrals + 1, this % two_electron_integral)
526  if (nend /= 1) stop "UNABLE TO READ ALCHEMY INTEGRALS"
527 
528  deallocate(xtempp)
529 
530  !endif
531 
532  1600 FORMAT(' Integrals read successfully:', /, &
533  ' 1-electron integrals, NINT1e =',i10, /, &
534  ' 2-electron integrals, NINT2e =',i10)
535  write (stdout, 1600) this % num_one_electron_integrals, this % num_two_electron_integrals
536 
537  this % nhe(:) = nhe(:)
538  this % nnuc = nnuc
539  this % dtnuc(:) = 0
540  this % dtnuc(1) = en
541  this % nhe(:) = nob(:)
542  this % dtnuc(22:21+nnuc) = xnuc(1:nnuc)
543  this % dtnuc(2:1+nnuc) = charg(1:nnuc)
544 
545  end subroutine load_integrals_alchemy
546 
547 
548  function get_integral_alchemy (this, i, j, k, l, m) result(coeff)
549  class(alchemyintegral) :: this
550  integer, intent(in) :: i,j,k,l,m
551  real(wp) :: coeff
552  integer :: ia, ib, integral_idx
553 
554  coeff = 0.0
555 
556  !One electron case
557  if (i == 0) then
558  integral_idx = this % get_one_electron_index(j, l, m)
559  coeff = this % one_electron_integral(integral_idx)
560  !print *,coeff
561  else
562  integral_idx = this % get_two_electron_index(i, j, k, l, m)
563  coeff = this % two_electron_integral(integral_idx)
564  end if
565 
566  !WRITE(stdout,*)i, j, k, l, m
567  !write(stdout,*)coeff,integral_idx
568  !write(stdout,*) i,j,k,l,m,coeff,integral_idx
569 
570  end function get_integral_alchemy
571 
572 
573  integer function get_one_electron_index (this, i, j, pos)
574  class(alchemyintegral), intent(in) :: this
575  integer, intent(in) :: i, j, pos
576  integer :: m, p, q, ii, jj
577 
578  ii = this % orbital_mapping(i)
579  jj = this % orbital_mapping(j)
580 
581  p = this % orbital_idx(ii)
582  q = this % orbital_idx(jj)
583 
584  m = this % symmetry_idx(ii)
585 
586  if (p < q) then
587  get_one_electron_index = this % pair_idx(q) + p + this % one_electron_pointer(m + 1) - 1
588  else
589  get_one_electron_index = this % pair_idx(p) + q + this % one_electron_pointer(m + 1) - 1
590  end if
591 
592  if (pos /= 0) get_one_electron_index = get_one_electron_index + this % num_one_electron_integrals
593  !write(stdout,"(10i8)") 0,p,0,q,0,get_one_electron_index,m,this%one_electron_pointer(m+1),this%pair_idx(p)
594 
595  end function get_one_electron_index
596 
597 
598  integer function get_two_electron_index (this, i, j, k, l, m)
599  class(alchemyintegral), intent(in) :: this
600  integer, intent(in) :: i, j, k, l, m
601  integer :: symmetry, mpp(4), nbp(4), npp(4), mapped(4), npq(4), ido, block_pointer, ipqrs, mv, md, mpq, mrs, mpr
602 
603  mapped(1) = this % orbital_mapping(i)
604  mapped(2) = this % orbital_mapping(j)
605  mapped(3) = this % orbital_mapping(k)
606  mapped(4) = this % orbital_mapping(l)
607 
608  do ido = 1, 4
609  npp(ido) = this % orbital_idx(mapped(ido))
610  symmetry = this % symmetry_idx(mapped(ido))
611  mpp(ido) = symmetry
612  nbp(ido) = this % num_orbitals_sym(symmetry + 1)
613  end do
614 
615  if (m /= 0) then
616  mv = mpp(1) + mpp(2)
617  else
618  mv = abs(mpp(1) - mpp(2))
619  end if
620 
621  if (mv /= 0) then
622  md = (mv - 1) / 2
623  else
624  md = -1
625  end if
626 
627  mpq = mpp(1) - md
628  mrs = mpp(3) - md
629  mpr = this % pair_idx(mpq) + mrs + this % num_two_electron_blocks(mv + 1)
630 
631  block_pointer = this % two_electron_pointer(mpr) - 1
632 
633  do ido = 1, 4, 2
634  if (mpp(ido) /= mpp(ido + 1)) then
635  npq(ido) = nbp(ido) * (npp(ido + 1) - 1) + npp(ido)
636  npq(ido + 1) = nbp(ido) * nbp(ido + 1)
637  else
638  if (npp(ido) < npp(ido + 1)) then
639  npq(ido) = this % pair_idx(npp(ido + 1)) + npp(ido)
640  else
641  npq(ido) = this % pair_idx(npp(ido)) + npp(ido + 1)
642  end if
643  npq(ido + 1) = this % pair_idx(nbp(ido) + 1)
644  end if
645  end do
646 
647  if (mpp(1) == mpp(3) .and. mpp(2) == mpp(4)) then
648  if (npq(1) < npq(3)) then
649  ipqrs = this % pair_idx(npq(3)) + npq(1)
650  else
651  ipqrs = this % pair_idx(npq(1)) + npq(3)
652  end if
653  else
654  if (this % integral_ordering == 0) then
655  ipqrs = npq(2) * (npq(3) - 1) + npq(1)
656  else
657  ipqrs = npq(4) * (npq(1) - 1) + npq(3)
658  end if
659  end if
660 
661  get_two_electron_index = ipqrs + block_pointer
662  !print *,'ALC ',block_pointer, IPQRS
663  !write(stdout,"(6i8)") i,j,k,l,m,get_two_electron_index
664 
665  end function get_two_electron_index
666 
667 
668  subroutine write_geometries_alchemy (this, iounit)
669  use params, only : cpoly
670  use scatci_routines, only : search
671 
672  class(alchemyintegral), intent(in) :: this
673  integer, intent(in) :: iounit
674  integer :: ido, ifail
675 
676  !REWIND iounit
677  !CALL SEARCH(iounit,cpoly,ifail)
678  !READ(iounit)nnuc
679  ! ALLOCATE(xnuc(nnuc),ynuc(nnuc),znuc(nnuc),charge(nnuc),CNAME(nnuc))
680  !DO i=1, nnuc
681  ! READ(iounit)cname(i), ii, xnuc(i), ynuc(i), znuc(i), charge(i)
682  !END DO
683 
684  end subroutine write_geometries_alchemy
685 
686 
687  subroutine destroy_integral_alchemy (this)
688  class(alchemyintegral) :: this
689 
690  if (associated(this % one_electron_integral)) deallocate(this % one_electron_integral)
691  if (associated(this % two_electron_integral)) deallocate(this % two_electron_integral)
692 
693  if (allocated(this % pair_labels)) deallocate(this % pair_labels)
694 
696  if (allocated(this % num_orbitals_sym)) deallocate(this % num_orbitals_sym)
697 
698  if (allocated(this % one_electron_pointer)) deallocate(this % one_electron_pointer)
699  if (allocated(this % two_electron_pointer)) deallocate(this % two_electron_pointer)
700 
701  if (allocated(this % pair_idx)) deallocate(this % pair_idx)
702  if (allocated(this % orbital_idx)) deallocate(this % orbital_idx)
703  if (allocated(this % symmetry_idx)) deallocate(this % symmetry_idx)
704 
705  if (allocated(this % num_two_electron_blocks)) deallocate(this % num_two_electron_blocks)
706  if (allocated(this % num_one_electron_blocks)) deallocate(this % num_one_electron_blocks)
707 
708  end subroutine destroy_integral_alchemy
709 
710 end module alchemy_module
alchemy_module
ALCHEMY integral module.
Definition: ALCHEMY_Module.f90:28
memorymanager_module
Memory manager module.
Definition: MemoryManager_module.f90:31
options_module::options
Calculation setup for a single Hamiltonian diagonalization.
Definition: Options_module.f90:74
alchemy_module::finalize_alchemy
subroutine finalize_alchemy(this)
Definition: ALCHEMY_Module.f90:433
alchemy_module::generate_orbital_index
subroutine generate_orbital_index(this)
Definition: ALCHEMY_Module.f90:340
alchemy_module::generate_pointer_table
subroutine generate_pointer_table(this)
Definition: ALCHEMY_Module.f90:226
alchemy_module::get_one_electron_index
integer function get_one_electron_index(this, i, j, pos)
Definition: ALCHEMY_Module.f90:574
alchemy_module::load_integrals_alchemy
subroutine load_integrals_alchemy(this, iounit)
This is just a copy from scatci_routines with only the relevant ALCHEMY parts. This will be replaced ...
Definition: ALCHEMY_Module.f90:440
alchemy_module::generate_pairs
subroutine generate_pairs(this)
?
Definition: ALCHEMY_Module.f90:154
alchemy_module::get_integral_alchemy
real(wp) function get_integral_alchemy(this, i, j, k, l, m)
Definition: ALCHEMY_Module.f90:549
baseintegral_module
Base integral module.
Definition: BaseIntegral_module.f90:28
alchemy_module::generate_pair_index
subroutine generate_pair_index(this)
Definition: ALCHEMY_Module.f90:375
alchemy_module::write_geometries_alchemy
subroutine write_geometries_alchemy(this, iounit)
Definition: ALCHEMY_Module.f90:669
baseintegral_module::baseintegral
Definition: BaseIntegral_module.f90:41
alchemy_module::initialize_alchemy
subroutine initialize_alchemy(this, option)
Definition: ALCHEMY_Module.f90:393
alchemy_module::alchemyintegral
Definition: ALCHEMY_Module.f90:45
options_module
Options module.
Definition: Options_module.f90:41
alchemy_module::get_two_electron_index
integer function get_two_electron_index(this, i, j, k, l, m)
Definition: ALCHEMY_Module.f90:599
alchemy_module::destroy_integral_alchemy
subroutine destroy_integral_alchemy(this)
Definition: ALCHEMY_Module.f90:688
alchemy_module::count_num_pairs
integer function count_num_pairs(this)
?
Definition: ALCHEMY_Module.f90:103
memorymanager_module::master_memory
type(memorymanager), public master_memory
This is the global memory manager.
Definition: MemoryManager_module.f90:69