MPI-SCATCI  2.0
An MPI version of SCATCI
SWEDEN_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 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
37  use options_module, only: options
38  use parallelization_module, only: grid => process_grid
39  use timing_module, only: master_timer
40 
41  implicit none
42 
43  public swedenintegral
44 
45  private
46 
47  type, extends(baseintegral) :: swedenintegral
48 
49  !Our integrals
50 #ifdef mpithree
51  real(wp), pointer :: one_electron_integral(:)
52  real(wp), pointer :: two_electron_integral(:)
53 #else
54  real(wp), allocatable :: one_electron_integral(:)
55  real(wp), allocatable :: two_electron_integral(:)
56 #endif
57 
58  integer :: one_electron_window
59  integer :: two_electron_window
60 
61  integer :: num_one_electron_integrals
62  integer :: num_two_electron_integrals
63 
64  real(wp), allocatable :: xnuc(:),ynuc(:),znuc(:),charge(:)
65  character(len=8), allocatable :: cname(:)
66 
67  integer :: num_unique_pairs
68  integer(longint), allocatable :: pair_labels(:,:)
69  integer, allocatable :: num_orbitals_sym(:)
70 
71  integer :: max_number_pair_sets
72  integer :: num_two_electron_blocks
73  integer, allocatable :: one_electron_pointer(:)
74  integer, allocatable :: two_electron_pointer(:)
75 
76  integer :: num_pq, num_rs
77  integer :: num_pair_idx
78 
80  integer, allocatable :: pair_idx(:)
81  integer, allocatable :: orbital_idx(:)
82  integer, allocatable :: symmetry_idx(:)
83 
84  contains
85 
86  procedure, public :: initialize_self => initialize_sweden
87  procedure, public :: finalize_self => finalize_sweden
88  procedure, public :: load_integrals => load_integrals_sweden
89  procedure, public :: get_integral_ijklm => get_integral_sweden
90  procedure, public :: write_geometries => write_geometries_sweden
91  procedure, public :: destroy_integrals => destroy_integral_sweden
92  procedure, private :: count_num_pairs
93  procedure, private :: generate_pairs
94  procedure, private :: generate_pointer_table
95  procedure, private :: generate_pair_index
96  procedure, private :: generate_orbital_index
97  !procedure, private :: read_one_electron_integrals
98  !procedure, private :: read_two_electron_integrals
99  procedure, private :: get_one_electron_index
100  procedure, private :: get_two_electron_index
101 
102  end type swedenintegral
103 
104 contains
105 
106  integer function count_num_pairs (this)
107  use global_utils, only: mprod
108 
109  class(swedenintegral) :: this
110  integer :: i, j, k, l
111  integer :: jend, lend
112  integer :: ij, kl, ijkltest
113 
114  count_num_pairs = 0
115 
116  do i = 1, this % num_symmetries
118  end do
119 
120  if (this % num_symmetries == 1) return
121 
122  !Generate IJIJ pairs
123  do i = 2, this % num_symmetries
124  jend = i - 1
125  do j = 1, jend
126  lend = j - 1
128  end do
129  end do
130 
131  !Generate IIJJ pairs
132  do i = 2, this % num_symmetries
133  jend = i - 1
134  do j = 1, jend
135  lend = j - 1
137  end do
138  end do
139 
140  if (this % num_symmetries < 4) return
141 
142  do i = 2, this % num_symmetries
143  jend = i - 1
144  do j = 1, jend
145  do k = 2, i
146  lend = k - 1
147  if (i == k) lend = j
148  do l = 1, lend
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
155  end do
156  end do
157  end do
158  end do
159 
160  end function count_num_pairs
161 
162 
163  subroutine generate_pairs (this)
164  class(swedenintegral) :: this
165  integer :: i, j, k, l
166  integer :: jend, lend
167  integer :: ij, kl, ijkltest
168  integer :: pair_number
169 
170  pair_number = 0
171 
172  write (stdout, *) 'All IIII blocks of integrals'
173 
174  do i = 1, this % num_symmetries
175  j = i - 1
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
179  end do
180 
181  if (this % num_symmetries == 1) return
182 
183  !Generate IJIJ pairs
184  write (stdout, *) 'All IJIJ blocks'
185  do i = 2, this % num_symmetries
186  jend = i - 1
187  do j = 1, jend
188  lend = j - 1
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
192  end do
193  end do
194 
195  !Generate IIJJ pairs
196  write (stdout, *) 'All IIJJ blocks'
197  do i = 2, this % num_symmetries
198  jend = i - 1
199  do j = 1, jend
200  lend = j - 1
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
204  end do
205  end do
206 
207  if (this % num_symmetries < 4) return
208  write (stdout, *)'All IJKL blocks'
209  do i = 2, this % num_symmetries
210  jend = i - 1
211  do j = 1, jend
212  do k = 2, i
213  lend = k - 1
214  if (i == k) lend = j
215  do l = 1, lend
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
224  end do
225  end do
226  end do
227  end do
228 
229  end subroutine generate_pairs
230 
231 
232  subroutine generate_pointer_table (this)
233  class(swedenintegral) :: this
234 
235  integer :: label(8), num_orbs(4)
236  integer :: ido, jdo, nam1, nam2, num_pq, num_rs, num_pr
237  integer :: block_number, err
238 
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
242 
243  allocate(this % one_electron_pointer(this % num_symmetries + 1), stat = err)
244 
245  call master_memory % track_memory(kind(this % one_electron_pointer), &
246  size(this % one_electron_pointer), err, 'SWEDEN::One electron pointer')
247 
248  this % one_electron_pointer(1) = 1
249 
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)
252  end do
253 
254  this % num_one_electron_integrals = this % one_electron_pointer(this % num_symmetries + 1) - 1
255 
256  allocate(this % two_electron_pointer(this % num_two_electron_blocks), stat = err)
257 
258  call master_memory % track_memory(kind(this % two_electron_pointer), &
259  size(this % two_electron_pointer), err, 'SWEDEN::Two electron pointer')
260 
261  this % two_electron_pointer(:) = 0
262 
263  write (stdout, "(//,5X,'Integral Storage Table follows : ',/)")
264 
265  do ido = 1, this % num_unique_pairs
267  call unpack8ints(this % pair_labels(1, ido), label)
268 
269  do jdo = 1, 4
270  label(jdo) = label(jdo) + 1
271  num_orbs(jdo) = this % num_orbitals_sym(label(jdo))
272  end do
273 
274  ! if IIJJ
275  if (label(1) == label(2)) then
276  num_pq = indfunc(num_orbs(1) + 1, 0)
277  else
278  num_pq = num_orbs(1) * num_orbs(2)
279  end if
280  nam1 = indfunc(max(label(1), label(2)) + 1, 0) - abs(label(1) - label(2))
281 
282  ! if JJ
283  if (label(3) == label(4)) then
284  num_rs = indfunc(num_orbs(3) + 1, 0)
285  else
286  num_rs = num_orbs(3) * num_orbs(4)
287  end if
288 
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
292 
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
295 
296  if (nam1 == nam2) then
297  num_pr = indfunc(num_pq + 1, 0)
298  else
299  num_pr = num_pq * num_rs
300  end if
301 
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
304  end do
305 
306  !Move to the end of the block
307  this % num_two_electron_integrals = this % num_two_electron_integrals - 1
308 
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
311 
312  end subroutine generate_pointer_table
313 
314 
315  subroutine generate_orbital_index (this)
316  class(swedenintegral) :: this
317  integer :: total_num_orbitals, orbital, iorb, isym, max_orbitals, ifail
318 
319  total_num_orbitals = sum(this % num_orbitals_sym(:))
320 
321  allocate(this % orbital_idx(total_num_orbitals), this % symmetry_idx(total_num_orbitals), stat = ifail)
322 
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')
325 
326  this % orbital_idx(:) = 0
327  this % symmetry_idx(:) = 0
328 
329  max_orbitals = 0
330  orbital = 0
331 
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
338  end do
339  end do
340 
341  this % num_pair_idx = max((this % num_symmetries**2 + this % num_symmetries) / 2 + 1, &
342  this % max_number_pair_sets, &
343  300, &
344  max_orbitals)
345 
346  end subroutine generate_orbital_index
347 
348 
349  subroutine generate_pair_index (this)
350  class(swedenintegral) :: this
351  integer :: ido, idx, ifail
352 
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')
355 
356  idx = 0
357  this % pair_idx(0) = 0
358  do ido = 1, this % num_pair_idx
359  this % pair_idx(ido) = idx
360  idx = idx + ido
361  end do
362 
363  end subroutine generate_pair_index
364 
365 
366  subroutine initialize_sweden (this, option)
367  class(swedenintegral) :: this
368  class(options), intent(in) :: option
369  integer :: ido, jdo, ifail
370 
371  write (stdout, "('Using SWEDEN integral')")
372 
373  this % num_unique_pairs = 0
374  this % num_symmetries = option % num_syms
375 
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')
379 
380  this % num_orbitals_sym(:) = option % num_electronic_orbitals_sym(:)
381  this % num_unique_pairs = this % count_num_pairs()
382 
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')
386 
387  write (stdout, *) 'SWEDEN Integral -- Generating pairs'
388 
390  call this % generate_pairs
391  call this % generate_orbital_index
392  call this % generate_pair_index
393 
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'
398 
399  call this % generate_pointer_table
400 
401  end subroutine initialize_sweden
402 
403 
404  subroutine finalize_sweden (this)
405  class(swedenintegral) :: this
406  end subroutine finalize_sweden
407 
408 
410  subroutine load_integrals_sweden (this, iounit)
411  use params, only: ctrans1, cpoly
412  use global_utils, only: intape
413  use scatci_routines, only: search, chn2e
414 
415  class(swedenintegral) :: this
416  integer, intent(in) :: iounit
417 
418  ! Sweden header arrays
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(:)
424 
425  integer(longint) :: I, NSYM, IONEIN
426  integer :: ifail, k, l, ind, ii
427  real(wp) :: EN
428 
429  call master_timer % start_timer('SWEDEN load')
430 
431  write (stdout, *) ' Loading SWEDEN Integrals into core'
432 
433  !allocate the integral space, this will be replaced by a caching system
434  !One electron integral
435  !allocate(this%one_electron_integral(2*this%num_one_electron_integrals),stat=ifail)
436  this % one_electron_window = mpi_memory_allocate_real(this % one_electron_integral, 2 * this % num_one_electron_integrals)
437 
438  !call master_memory%track_memory(kind(this%one_electron_integral),size(this%one_electron_integral),ifail,'One electron integral')
439  !if(ifail /=0) stop "could not allocate 1 electron integral"
440 
441  !Two electron integral
442  this % two_electron_window = mpi_memory_allocate_real(this % two_electron_integral, this % num_two_electron_integrals)
443 
444  call mpi_memory_synchronize(this % one_electron_window)
445  call mpi_memory_synchronize(this % two_electron_window)
446 
447  !if( (this%two_electron_window == -1 .and. this%one_electron_window == -1) .or. local_rank == local_master) then
448  !call master_memory%track_memory(kind(this%two_electron_integral),size(this%two_electron_integral),ifail,'One electron integral')
449  !if(ifail /=0) stop "could not allocate 2 electron integral"
450 
451  read (iounit) (label(i), i = 1, 4)
452 
453  write (stdout, fmt='(/,5X,''Label on Sweden tape = '',4A)') (label(i), i = 1, 4)
454 
455  read (iounit) nsym, en, (nao(i), nob(i), ncore(i), i = 1, nsym)
456 
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)
460 
461  this % core_energy = en
462 
463  do i = 1, nsym
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
467  end if
468  end do
469 
470  !Do Erro checkin here will ignore for now
471 
472  ionein = 0
473  do i = 1, nsym
474  ionein = ionein + (nob(i) * (nob(i) + 1)) / 2
475  end do
476 
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)
479  read (iounit)
480  call intape(iounit, this % one_electron_integral, this % num_one_electron_integrals)
481 
482  !Positron case
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')
486 
487  allocate(xtempp(this % num_one_electron_integrals), stat = ifail)
488  call master_memory % track_memory(kind(xtempp), size(xtempp), ifail, 'SWEDEN::xtempp')
489 
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)
492 
493  !Are we using Quantemol-N?
494  if (this % quantamoln_flag) then
495  ind = 0
496  do k = 1, nsym
497  DO l = 1, nob(k)
498  ind = ind + l
499  write (413, "(i6,3x,i6,3x,i6,5x,f12.7)") k - 1, l, l, xtempp(ind)
500  end do
501  end do
502  end if
503 
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)
506  end do
507 
508  call master_memory % free_memory(kind(xtempp), size(xtempp))
509  deallocate(xtempp)
510  call master_memory % free_memory(kind(xtempe), size(xtempe))
511  deallocate(xtempe)
512  end if
513 
514  call chn2e(iounit, stdout, this % two_electron_integral, this % num_two_electron_integrals)
515  end if
516 
517  1600 format(' Integrals read successfully:', /, &
518  ' 1-electron integrals, NINT1e =',i10, /, &
519  ' 2-electron integrals, NINT2e =',i10)
520 
521  write (stdout, 1600) this % num_one_electron_integrals, this % num_two_electron_integrals
522 
523  this % nhe(:) = nhe(:)
524  this % nnuc = 0
525  this % dtnuc(:) = 0
526  this % dtnuc(1) = en
527  this % nhe(:) = nob(:)
528  this % dtnuc(22:41) = 0
529  this % dtnuc(2:21) = 0
530 
531  !Reading Geometries
532  rewind iounit
533  call search(iounit, cpoly, ifail)
534  read (iounit) this % nnuc
535 
536  allocate(this % xnuc(this % nnuc), this % ynuc(this % nnuc), this % znuc(this % nnuc), &
537  this % charge(this % nnuc), this % CNAME(this % nnuc))
538 
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)
541  end do
542 
543  call master_timer % stop_timer('SWEDEN load')
544 
545  !endif
546 
547  call mpi_memory_synchronize(this % one_electron_window)
548  call mpi_memory_synchronize(this % two_electron_window)
549 
550  end subroutine load_integrals_sweden
551 
552 
553  function get_integral_sweden (this, i, j, k, l, m) result(coeff)
554  class(swedenintegral) :: this
555  integer, intent(in) :: i, j, k, l, m
556  real(wp) :: coeff
557  integer :: ia, ib, integral_idx
558 
559  coeff = 0.0
560 
561  !One electron case
562  if (i == 0) then
563  integral_idx = this % get_one_electron_index(j, l, m)
564  coeff = this % one_electron_integral(integral_idx)
565  else
566  integral_idx = this % get_two_electron_index(i, j, k, l, m)
567  coeff = this % two_electron_integral(integral_idx)
568  end if
569 
570  end function get_integral_sweden
571 
572 
573  integer function get_one_electron_index (this, i, j, pos)
574  class(swedenintegral), 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(swedenintegral), intent(in) :: this
600  integer, intent(in) :: i, j, k, l, m
601  integer :: symmetry
602  integer, dimension(4) :: orb_num_sym, orb_sym, num_orb_sym
603  integer :: mapped(4)
604  integer :: ido
605  integer :: mpq, mrs, mpr, aaa, block_pointer, nobrs, nobpq, ipqrs
606 
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)
611 
612  do ido = 1, 4
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)
617  end do
618 
619  !write(stdout,"(4i4,2x,4i4,2x,4i4)") orb_num_sym(1:4),orb_sym(1:4),num_orb_sym(1:4)
620 
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
624 
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
633  end if
634 
635  block_pointer = this % two_electron_pointer(mpr) - 1
636 
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)
639  !IF(ipqrs.GT.nobtest)WRITE(6,910)ipqrs, nobtest
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
643  ! IF(ipqrs.GT.nobtest)WRITE(6,910)ipqrs, nobtest
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)
648  else
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) &
651  + orb_num_sym(4)
652  end if
653 
654  get_two_electron_index = ipqrs + block_pointer + 1
655 
656  !write(stdout,"(6i8)") i,j,k,l,m,get_two_electron_index
657 
658  end function get_two_electron_index
659 
660 
661  subroutine write_geometries_sweden (this, iounit)
662  use params, only: cpoly
663  use scatci_routines, only: search
664 
665  class(swedenintegral), intent(in) :: this
666  integer, intent(in) :: iounit
667  integer :: ido, ifail, i
668 
669  do ido = 1, this % nnuc
670  write (iounit) this % cname(ido), this % xnuc(ido), this % ynuc(ido), this % znuc(ido), this % charge(ido)
671  end do
672 
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)
675 
676  end subroutine write_geometries_sweden
677 
678 
679  subroutine destroy_integral_sweden (this)
680  class(swedenintegral) :: this
681 
682  !if(allocated(this%one_electron_integral)) then
683  !call master_memory%free_memory(kind(this%one_electron_integral),size(this%one_electron_integral))
684  !deallocate(this%one_electron_integral)
685  !endif
686 
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)
691 
692  !if(allocated(this%two_electron_integral)) then
693  !call master_memory%free_memory(kind(this%two_electron_integral),size(this%two_electron_integral))
694  !deallocate(this%two_electron_integral)
695  !endif
696 
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)
700  end if
701 
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)
706  end if
707 
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)
711  end if
712 
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)
716  end if
717 
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)
721  end if
722 
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)
726  end if
727 
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)
731  end if
732 
733  end subroutine destroy_integral_sweden
734 
735 end module sweden_module
sweden_module::get_two_electron_index
integer function get_two_electron_index(this, i, j, k, l, m)
Definition: SWEDEN_Module.F90:599
memorymanager_module
Memory manager module.
Definition: MemoryManager_module.f90:31
timing_module
Timing module.
Definition: Timing_Module.f90:28
sweden_module::swedenintegral
Definition: SWEDEN_Module.F90:47
options_module::options
Calculation setup for a single Hamiltonian diagonalization.
Definition: Options_module.f90:74
sweden_module::generate_pointer_table
subroutine generate_pointer_table(this)
Definition: SWEDEN_Module.F90:233
sweden_module::get_one_electron_index
integer function get_one_electron_index(this, i, j, pos)
Definition: SWEDEN_Module.F90:574
parallelization_module::process_grid
type(processgrid) process_grid
Definition: Parallelization_module.F90:80
sweden_module::load_integrals_sweden
subroutine load_integrals_sweden(this, iounit)
This is just a copy from scatci_routines with only the relevant SWEDEN parts.
Definition: SWEDEN_Module.F90:411
baseintegral_module
Base integral module.
Definition: BaseIntegral_module.f90:28
parallelization_module
Distribution of processes into a grid.
Definition: Parallelization_module.F90:29
sweden_module::generate_orbital_index
subroutine generate_orbital_index(this)
Definition: SWEDEN_Module.F90:316
sweden_module::destroy_integral_sweden
subroutine destroy_integral_sweden(this)
Definition: SWEDEN_Module.F90:680
sweden_module::get_integral_sweden
real(wp) function get_integral_sweden(this, i, j, k, l, m)
Definition: SWEDEN_Module.F90:554
sweden_module::generate_pair_index
subroutine generate_pair_index(this)
Definition: SWEDEN_Module.F90:350
sweden_module::write_geometries_sweden
subroutine write_geometries_sweden(this, iounit)
Definition: SWEDEN_Module.F90:662
baseintegral_module::baseintegral
Definition: BaseIntegral_module.f90:41
sweden_module
SWEDEN integral module.
Definition: SWEDEN_Module.F90:28
sweden_module::count_num_pairs
integer function count_num_pairs(this)
Definition: SWEDEN_Module.F90:107
options_module
Options module.
Definition: Options_module.f90:41
sweden_module::initialize_sweden
subroutine initialize_sweden(this, option)
Definition: SWEDEN_Module.F90:367
timing_module::master_timer
type(timer), public master_timer
Definition: Timing_Module.f90:74
sweden_module::generate_pairs
subroutine generate_pairs(this)
Definition: SWEDEN_Module.F90:164
memorymanager_module::master_memory
type(memorymanager), public master_memory
This is the global memory manager.
Definition: MemoryManager_module.f90:69
sweden_module::finalize_sweden
subroutine finalize_sweden(this)
Definition: SWEDEN_Module.F90:405