MPI-SCATCI 2.0
An MPI version of SCATCI
Loading...
Searching...
No Matches
Contracted_Hamiltonian_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 Contracted Hamiltonian module
23!> \authors A Al-Refaie
24!> \date 2017
25!>
26!> This module handles the calculation of the contracted scattering hamiltonian
27!> through the usage of the Contracted_Hamiltonian class. Eq. references refer to the paper by
28!> J. Tennyson: J. Phys B.: At. Mol. Opt. Phys. 29 (1996) 1817-1828
29!>
30!> \note 30/01/2017 - Ahmed Al-Refaie: Initial Revision
31!> \note 16/01/2019 - Jakub Benda: Unifom coding style and expanded documentation.
32!>
33module contracted_hamiltonian_module
34
35 use const_gbl, only: stdout
37 use mpi_gbl, only: mpi_reduceall_max
38 use precisn, only: longint, wp
39 use basematrix_module, only: basematrix
40 use contracted_symbolic_module, only: contractedsymbolicelementvector
41 use hamiltonian_module, only: basehamiltonian
42 use memorymanager_module, only: master_memory
43 use parallelization_module, only: grid => process_grid
44 use symbolic_module, only: symbolicelementvector
45 use target_rmat_ci_module, only: target_rmat_ci
46 use timing_module, only: master_timer
47 use utility_module, only: triangular_index_to_ij, compute_total_triangular, compute_total_box, box_index_to_ij, &
48 pack_ints, unpack_ints
49
50 implicit none
51
52 public contracted_hamiltonian
53
54 private
55
56 !> \brief Computation of Hamiltonian
57 !> \authors A Al-Refaie
58 !> \date 2017
59 !>
60 !> This class computes the hamitonian using the contracted prototype scheme
61 !> given by J. Tennyson: J. Phys B.: At. Mol. Opt. Phys. 29 (1996) 1817-1828.
62 !> The notation i1,n1,j1,i2,n2,j2 referes to the variables i,n,j,i',n',j' given in table 1.
63 !> The class refer to the ones desribed in the Paper in Section 3.
64 !>
65 type, extends(basehamiltonian) :: contracted_hamiltonian
66
67 class(target_rmat_ci), pointer :: rmat_ci(:) !< Our Target coefficients
68 integer :: start_l2 !< Start of the L2 CSFS
69 integer :: total_num_csf !< Total number of CSFs we're dealing with
70 integer :: num_l2_functions
71 integer :: l2_mat_offset !< Matrix starting position of the L2 portion
72 integer :: csf_skip !< Used to move to the next CSF (most continuum ones come in pairs so the next one is skipped usually)
73 integer :: non_zero_prototype = 0 !< Counter for calculated prototypes that are non_zero
74
75 contains
76
77 procedure, public :: initialize => initialize_contracted_hamiltonian
78 procedure, public :: build_hamiltonian => build_contracted_hamiltonian
79
80 !----------class routines----------------------!
81 procedure, private :: evaluate_class_1
82 procedure, private :: evaluate_class_2
83 procedure, private :: evaluate_class_3
84 procedure, private :: evaluate_class_4
85 procedure, private :: evaluate_class_5
86 procedure, private :: evaluate_class_6
87 procedure, private :: evaluate_class_7
88 procedure, private :: evaluate_class_8
89 procedure, private :: evaluate_class_2_and_8
90 !--------class routines------------------------!
91 !procedure, private :: evaluate_class_7
92 !---Contraction routines
93 procedure, private :: fast_contract_class_1_3_diag
94 procedure, private :: fast_contract_class_1_3_offdiag
95 procedure, private :: fast_contract_class_567
96
97 !---------Expansion routines-------------------!
98 procedure, private :: expand_diagonal
99 procedure, private :: expand_off_diagonal
100 procedure, private :: expand_continuum_l2
101 procedure, private :: expand_off_diagonal_general
102
103 procedure, private :: expand_continuum_l2_eval
104 procedure, private :: expand_diagonal_eval
105 procedure, private :: contr_expand_diagonal_eval
106 procedure, private :: contr_expand_off_diagonal_eval
107 procedure, private :: contr_expand_continuum_l2_eval
108 procedure, private :: expand_off_diagonal_eval
109 procedure, private :: expand_off_diagonal_gen_eval
110
111 !-------auxilary--routines-----------------!
112 procedure, private :: compute_matrix_ij
113 procedure, private :: get_target_coeff
114 procedure, private :: reduce_prototypes
115 procedure, private :: get_starting_index
116
117 end type contracted_hamiltonian
118
119contains
120
121 !> \brief Initializes the contracted hamiltonian and supplies target coefficients
122 !> \authors A Al-Refaie
123 !> \date 2017
124 !>
125 !> \param[inout] this Hamiltonian object to update.
126 !> \param[in] rmat_ci Our Target_RMat_CI array containing the cofficients.
127 !>
128 subroutine initialize_contracted_hamiltonian (this, rmat_ci)
129 class(contracted_hamiltonian) :: this
130 class(target_rmat_ci), target, intent(in) :: rmat_ci(:)
131
132 !Assign the rmat pointer
133 this % rmat_ci => rmat_ci
134 !Our skip value is 2 for the moment
135 this % csf_skip = 2
136 !Where the L2 starts
137 this % start_L2 = this % options % last_continuum + 1
138 !total number of CSFS
139 this % total_num_csf = this % options % num_csfs
140 !Calculate where each L2 CSF belongs in the matrix
141 this % L2_mat_offset = this % options % contracted_mat_size - this % total_num_csf
142 this % num_L2_functions = this % options % num_L2_CSF
143
144 end subroutine initialize_contracted_hamiltonian
145
146
147 !> \brief Builds the Contracted Hamiltonian class by class
148 !> \authors A Al-Refaie
149 !> \date 2017
150 !>
151 !> \param[inout] this Hamiltonian object to update.
152 !> \param[out] matrix_elements A BaseMatrix object, stores the calculated elements.
153 !>
154 subroutine build_contracted_hamiltonian (this, matrix_elements)
155 class(contracted_hamiltonian) :: this
156 class(basematrix), intent(inout) :: matrix_elements
157 integer :: num_target_sym, target_sym, num_targets_per_sym, i1, i2, n1, n2, l1, loop_skip, loop_ido, my_idx, ido
158 integer :: num_csfs, m, n, num_continuum, total_l1_cont_elems, matrix_block_size, estimated_mat_size
159
160 !The block size is the size of the continuum and continuum+L2 section of the matrix
161 matrix_block_size = this % options % contracted_mat_size - this % total_num_csf + this % options % last_continuum
162
163 write (stdout, "('continuum block_size =',i8)") matrix_block_size
164 !We tell the matrix storage class the contracted size, the block size and the fact that it is sparse
165 call matrix_elements % initialize_matrix_structure(this % options % contracted_mat_size, mat_sparse, matrix_block_size)
166
167 !Get the number of target symmetries
168 num_target_sym = this % options % num_target_sym
169
170 !------------------------------------------------------------------------------------------------------------------!
171 !-----------------------------CONTINUUM CALCULATIONS---------------------------------------------------------------!
172 !------------------------------------------------------------------------------------------------------------------!
173
174 this % orbitals % MFLG = 0
175
176 !-----------------------------------------------CLASS 1---------------------------------------------------------
177 call master_timer % start_timer("Class 1")
178
179 !Perform a class one ('diagonal same target symmetry') calculation on all target symmetries
180 do i1 = 1, num_target_sym
181 call this % evaluate_class_1(i1, this % options % num_target_state_sym(i1), matrix_elements)
182 end do
183 call master_timer % stop_timer("Class 1")
184 write (stdout, "('Class 1 complete')")
185 flush(stdout)
186 call master_timer % report_timers
187 !stop
188 !-----------------------------------------------CLASS 3---------------------------------------------------------
189
190
191 !-----------------------------------------------CLASSES 567---------------------------------------------------------
192 !These classes deal with differing target symmetries
193 !call master_timer%start_timer("Class 56")
194 !Loop through every combination of target symmetry
195 !do i1 = 1, num_target_sym-1
196 ! do i2=i1+1,num_target_sym
197 ! if(i1==i2) cycle !We've already done same target symmetries so ignore
198 !Check if we have the same continuum symmetry
199 ! if(this%options%lambda_continuum_orbitals_target(i1) == this%options%lambda_continuum_orbitals_target(i2)) then
200 ! if(this%options%gerade_sym_continuum_orbital_target(i1) == this%options%gerade_sym_continuum_orbital_target(i2)) then
201
202 ! cycle
203 ! endif
204 ! endif
205 ! enddo
206 !enddo
207 !call master_timer%stop_timer("Class 56")
208
209 call master_memory % print_memory_report
210
211 !Before moving to the pure L2 calculations we need to clear out any other continuum calculations from the matrix format
212 !call matrix_elements%update_continuum(.true.)
213
214 !--------------------------------------------------------------------------------------------------------------------------!
215 !----------------------------------------------PURE L2 Calculations--------------------------------------------------------!
216 !--------------------------------------------------------------------------------------------------------------------------!
217
218 this % orbitals % MFLG = 1
219
220 call master_timer % start_timer("Class 3")
221 !Perform a class three ('off-diagonal same target symmetry') calculation on all target symmetries
222 do i1 = 1, num_target_sym
223 call this % evaluate_class_3(i1, this % options % num_target_state_sym(i1), matrix_elements)
224 end do
225 call master_timer % stop_timer("Class 3")
226
227 !call master_timer%report_timers
228
229 write (stdout, "('Class 3 complete')")
230 flush(stdout)
231
232 call master_timer % report_timers
233 call master_memory % print_memory_report
234
235 !-----------------------------------------------CLASS 4---------------------------------------------------------
236 !Class four is a difficult one as we have to loop through the entire L2 functions which are quite large
237 loop_skip = max(1, grid % gprocs)
238 my_idx = max(grid % grank, 0)
239
240 call master_timer % start_timer("Class 4")
241 do i1 = 1, num_target_sym
242 total_l1_cont_elems = this % options % num_target_state_sym(i1) * this % options % num_continuum_orbitals_target(i1)
243 do loop_ido = this % start_L2, this % total_num_csf, loop_skip
244 l1 = my_idx + loop_ido
245
246 if (l1 > this % total_num_csf) then
247 !This is dummy to ensure synchornization with other procs
248 do ido = 1, total_l1_cont_elems
249 call matrix_elements % update_pure_L2(.false.)
250 end do
251 else
252 !Class four (continuum-L2 calculation)
253 call this % evaluate_class_4(i1, this % options % num_target_state_sym(i1), l1, matrix_elements)
254 end if
255 end do
256 end do
257
258 call master_timer % stop_timer("Class 4")
259 write (stdout, "('Class 4 complete')")
260 flush(stdout)
261
262 call master_timer % report_timers
263 !call matrix_elements%update_pure_L2(.true.)
264 !call master_memory%print_memory_report
265
266 !-----------------------------------------------CLASSES 2+8---------------------------------------------------------
267 call master_timer % start_timer("Class 567")
268
269 do i1 = 1, num_target_sym - 1
270 do i2 = i1 + 1, num_target_sym
271
272 if (i1 == i2) cycle !We've already done same target symmetries so ignore
273 !write(stdout,"('class 567: ',2i8)") i1,i2
274 !Check if we have the same continuum symmetry
275 if (this % options % lambda_continuum_orbitals_target(i1) == &
276 this % options % lambda_continuum_orbitals_target(i2)) then
277 if (this % options % gerade_sym_continuum_orbital_target(i1) == &
278 this % options % gerade_sym_continuum_orbital_target(i2)) then
279 !If we do then do class five ('diagonal same continuum symmetry') calculation
280 ! write(stdout,"('class 5')")
281 call this % evaluate_class_5(i1, this % options % num_target_state_sym(i1), i2, &
282 this % options % num_target_state_sym(i2), matrix_elements)
283 !Then a class six ('off-diagonal same continuum symmetry') calculation
284 ! write(stdout,"('class 6')")
285 call this % evaluate_class_6 (i1, this % options % num_target_state_sym(i1), i2, &
286 this % options % num_target_state_sym(i2), matrix_elements)
287 !call master_timer%report_timers
288 !call master_memory%print_memory_report
289 cycle
290 end if
291 end if
292
293 !write(stdout,"('class 7')")
294 !Otherwise do a class seven ('differing continuum symmetry') calculation
295 call this % evaluate_class_7(i1, this % options % num_target_state_sym(i1), i2, &
296 this % options % num_target_state_sym(i2), matrix_elements)
297 !call master_timer%report_timers
298 !call master_memory%print_memory_report
299
300 end do
301 end do
302
303 call master_timer % stop_timer("Class 567")
304 write (stdout, "('Class 567 complete')")
305 flush(stdout)
306
307 call master_timer % report_timers
308 call master_memory % print_memory_report
309 call master_timer % start_timer("Class 2+8")
310 !Here I've combined both class two ('L2 diagonal') and class eight ('L2 offdiagonal') calculations
311 !call this%evaluate_class_2(matrix_elements)
312 call this % evaluate_class_2_and_8(matrix_elements)
313 !call this%evaluate_class_8(matrix_elements)
314 call master_timer % stop_timer("Class 2+8")
315
316 write (stdout, "('Class 2+8 complete')")
317 flush(stdout)
318
319 call master_timer % report_timers
320 call master_memory % print_memory_report
321
322 !Clear up remaining L2 elements in a matrix formats stash
323 call matrix_elements % update_pure_L2(.true.)
324
325 !The matrix may have some last step thats required
326 call matrix_elements % finalize_matrix
327
328 !Print out interesting information
329 n = (this % options % num_csfs * (this % options % num_csfs + 1)) / 2
330 m = (this % options % contracted_mat_size * (this % options % contracted_mat_size + 1)) / 2
331
332 write (stdout, 6000) m, matrix_elements % get_size()
333 6000 format(//,' TOTAL H Matrix ELEMENTS =',i15,/,' NON-ZERO ELEMENTS evaluated =',i15)
334 write (stdout, 6010) n, this % non_zero_prototype
335 6010 format(' Total prototype ELEMENTS =',i15,/,' Non-zero prototype ELEMENTS =',i15)
336 ! WRITE(stdout,6020)nint0, nint00
337 6020 format(/,' Number (prototype) integrals =',i15,/,' Compressed (prototype) integrals =',i15)
338 write (stdout, 6030) this % number_of_integrals
339 6030 format(' Number integrals evaluated =',i15)
340
341 end subroutine build_contracted_hamiltonian
342
343
344 !> \brief Builds the class 1 part of the matrix
345 !> \authors A Al-Refaie
346 !> \date 2017
347 !>
348 !> This relates to the building of diagonal elements of the target states with the same target symmetries.
349 !>
350 !> \param[inout] this Hamiltonian object to query.
351 !> \param[in] i1 The target symmetry.
352 !> \param[in] num_targets The number of target states withing the symmetry.
353 !> \param[out] matrix_elements A BaseMatrix object, stores the calculated elements.
354 !>
355 subroutine evaluate_class_1 (this, i1, num_targets, matrix_elements)
356 class(contracted_hamiltonian) :: this
357 class(basematrix), intent(inout) :: matrix_elements
358 integer, intent(in) :: i1
359 integer, intent(in) :: num_targets
360
361 !>These are the contracted prototypes for each combination of target states
362 type(contractedsymbolicelementvector) :: master_prototype
363 !>This stores te resultant slater rule calculation of prototypes before contracting for each state
364 type(symbolicelementvector) :: temp_prototype
365 !>The starting matrix index for each target state combination
366 integer, allocatable :: matrix_index(:,:)
367
368 !>The number of target state combinations
369 integer :: num_ci
370 !------------------------------target state variables
371 integer :: n1, n2
372 !----------------------------Looping variables------------------------------!
373 integer :: starting_index, ido, jdo, ii, jj, ij, err
374 !----------------------------MPI looping variables--------------------------
375 integer :: loop_idx, loop_ido, my_idx, loop_skip
376 !>The total number of configurations with this target symmetry
377 integer :: num_configurations
378 !>CSF indicies
379 integer :: config_idx_a, config_idx_b
380 !The target coefficient
381 real(wp), allocatable :: alpha_coeffs(:,:)
382 !The resultant matrix coefficient
383 real(wp), allocatable :: mat_coeffs(:,:)
384 !Number of continuum orbitals
385 integer :: continuum_orbital, num_continuum_orbitals
386
387 allocate(matrix_index(2,num_targets*(num_targets+1)/2), stat = err)
388 call master_memory % track_memory(storage_size(matrix_index)/8, size(matrix_index), err, 'CLASS1::MATRIX_INDEX')
389 allocate(alpha_coeffs(num_targets*(num_targets+1)/2,1), stat = err)
390 call master_memory % track_memory(storage_size(alpha_coeffs)/8, size(alpha_coeffs), err, 'CLASS1::alpha_coeffs')
391 allocate(mat_coeffs(num_targets*(num_targets+1)/2,1), stat = err)
392 call master_memory % track_memory(storage_size(mat_coeffs)/8, size(mat_coeffs), err, 'CLASS1::mat_coeffs')
393 !Get number of continuum orbitals
394 num_continuum_orbitals = this % options % num_continuum_orbitals_target(i1)
395
396 !number of target state combinations (triangular form)
397 num_ci = num_targets * (num_targets + 1) / 2
398
399 ii = 1
400 do n1 = 1, this % options % num_target_state_sym(i1)
401 do n2 = n1, this % options % num_target_state_sym(i1)
402 !Construct a master prototype for each target state combo
403 !And get the starting matrix index
404 call this % compute_matrix_ij(i1, n1, 1, i1, n2, 1, matrix_index(1,ii), matrix_index(2,ii))
405 ii = ii + 1
406 end do
407 end do
408
409 call master_prototype % construct(num_ci, 1)
410
411 !Construct the temporary prototype
412 call temp_prototype % construct
413
414 !get the starting index of the configurations
415 starting_index = this % get_starting_index(i1)
416
417 !Get the number of configurations
418 num_configurations = this % options % num_ci_target_sym(i1)
419
420 !Set up our mpi variables
421 loop_skip = max(1, grid % lprocs)
422 my_idx = max(grid % lrank, 0)
423 call master_timer % start_timer("Class 1 prototype")
424 !This loop is commonly found in a lot of the classes. Basically it involves striding the loop across mpi processes
425 !For example if nprocs is 5, process 0 will do indicies 1,6,11..etc while process 1 will do 2,7,12.. and so on
426 do loop_ido = 1, num_configurations, loop_skip
427 this % orbitals % MFLG = 0
428 !Calculate our real index
429 ido = loop_ido + my_idx
430 !If we've gon past (which is possible in this setup then skip
431 if (ido > num_configurations) cycle
432 !Calculate our diagonal configuration index
433 config_idx_a = starting_index + (ido - 1) * this % csf_skip
434 !Slater it to get our temporary prototypes
435 call this % slater_rules(this % csfs(config_idx_a), this % csfs(config_idx_a), temp_prototype, 0, .false.)
436 alpha_coeffs = 0.0_wp
437 !Loop through each target state
438 ii = 1
439 do n1 = 1, this % options % num_target_state_sym(i1)
440 do n2 = n1, this % options % num_target_state_sym(i1)
441 ! Apply our coefficeints for each combination to the master, Left hand of Eq. 6
442 alpha_coeffs(ii,1) = this % get_target_coeff(i1, ido, n1) * this % get_target_coeff(i1, ido, n2)
443 ii = ii + 1
444 end do
445 end do
446 call master_prototype % add_symbols(temp_prototype, alpha_coeffs)
447 !call this%fast_contract_class_1_3_diag(i1,num_targets,ido,temp_prototype,master_prototype)
448
449 !Clear the temporary prototype for usage again
450 call temp_prototype % clear
451
452 !Now the off diagonal
453 do jdo = 1, ido - 1
454
455 !Get the index
456 config_idx_b = starting_index + (jdo - 1) * this % csf_skip
457 !Slater
458 call this % slater_rules(this % csfs(config_idx_b), this % csfs(config_idx_a), temp_prototype, 1, .false.)
459
460 !No prototypes? then we leave
461 if (temp_prototype % is_empty()) cycle
462
463 !Now contract the 'off diagonal' diagonal
464 ii = 1
465 do n1 = 1, this % options % num_target_state_sym(i1)
466 do n2 = n1, this % options % num_target_state_sym(i1)
467 ! (cimn*cim'n' + cim'n*cimn')H_imj,imj' Right hand of Eq. 6
468 alpha_coeffs(ii,1) = this % get_target_coeff(i1, ido, n1) * this % get_target_coeff(i1, jdo, n2) &
469 + this % get_target_coeff(i1, jdo, n1) * this % get_target_coeff(i1, ido, n2)
470 ! call master_prototype(ii)%add_symbols(temp_prototype,alpha)
471 ii = ii + 1
472 end do
473 end do
474 !print *,temp_prototype%get_size()
475 !print *,jdo
476 call master_prototype % add_symbols(temp_prototype, alpha_coeffs)
477 !call this%fast_contract_class_1_3_offdiag(i1,num_targets,ido,jdo,temp_prototype,master_prototype)
478 !Clear the prototype for the next round
479 call temp_prototype % clear
480 end do
481 end do
482 call master_timer % stop_timer("Class 1 prototype")
483 ij = 1
484
485 call master_prototype % synchronize_symbols()
486
487 !Without expansion, evaluate the first elements
488 !do ii=1,num_targets
489 ! do jj=ii,num_targets
490 ! mat_coeff = this%evaluate_integrals(master_prototype(ij),NORMAL_PHASE)
491 !Insert
492 !these are special as they may add the eigenvalues after
493 !We just let the master rank add them since they'll be combined later on
494 ! if(myrank == master .and. this%diagonal_flag == NO_PURE_TARGET_INTEGRAL_DIFF) mat_coeff = mat_coeff + this%rmat_ci(i1)%eigenvalues(ii)
495 ! call matrix_elements%insert_matrix_element(matrix_index(2,ij),matrix_index(1,ij),mat_coeff,1)
496 !Possibly perform update (who knows really?)
497 ! call matrix_elements%update_continuum(.false.)
498 ! ij=ij+1
499 ! enddo
500 !enddo
501
502 !now we have our prototype we can start the expansion
503 config_idx_a = starting_index + (num_configurations - 1) * this % csf_skip
504 !Get which continuum orbital this belongs to
505 continuum_orbital = this % csfs(config_idx_a) % orbital_sequence_number
506
507 !Set up our mpi variables
508 loop_skip = max(1, grid % gprocs)
509 my_idx = max(grid % grank, 0)
510
511 call master_timer % start_timer("Class 1 Expansion")
512 !Now loop through our J
513 do loop_ido = 1, num_continuum_orbitals, loop_skip
514 ido = loop_ido + my_idx
515 call matrix_elements % update_pure_L2(.false., num_targets * (num_targets - 1))
516
517 if (ido > num_continuum_orbitals) cycle
518 !Loop through each target combo
519 mat_coeffs = 0
520 !Immediately expand the correct prototype and then immediately evaluate the matrix element
521 call this % contr_expand_diagonal_eval(continuum_orbital, master_prototype, ido, mat_coeffs)
522 ij = 1
523 do ii = 1, num_targets
524 do jj = ii, num_targets
525 !Insert into our matrix format
526 if (this % diagonal_flag == no_pure_target_integral_diff .and. ii == jj) then
527 mat_coeffs(ij,1) = mat_coeffs(ij,1) + this % rmat_ci(i1) % eigenvalues(ii)
528 end if
529 call matrix_elements % insert_matrix_element(matrix_index(2,ij) + ido - 1, &
530 matrix_index(1,ij) + ido - 1, mat_coeffs(ij,1), 8)
531 !An update may occur
532 ij = ij + 1
533 end do
534 end do
535 end do
536
537 call master_timer % stop_timer("Class 1 Expansion")
538 !Destroy our temporary prototype
539 call temp_prototype %destroy
540 !Destroy all of the masters
541 !do ii=1,num_ci
542 call master_prototype % destroy
543 !enddo
544 call master_memory % free_memory(storage_size(matrix_index)/8, size(matrix_index))
545 call master_memory % free_memory(storage_size(alpha_coeffs)/8, size(alpha_coeffs))
546 call master_memory % free_memory(storage_size(mat_coeffs)/8, size(mat_coeffs))
547 deallocate(matrix_index, alpha_coeffs, mat_coeffs)
548
549 this % non_zero_prototype = this % non_zero_prototype + 1
550
551 !call master_timer%report_timers
552 end subroutine evaluate_class_1
553
554
555 !> \brief Builds the class 2 part of the matrix
556 !> \authors A Al-Refaie
557 !> \date 2017
558 !>
559 !> This has be deprecated, instead class 2 and 8 have been combined into a faster method.
560 !>
561 !> \param[inout] this Hamiltonian object to query.
562 !> \param[out] matrix_elements A BaseMatrix object, stores the calculated elements.
563 !>
564 subroutine evaluate_class_2 (this, matrix_elements)
565 class(contracted_hamiltonian) :: this
566 class(basematrix), intent(inout) :: matrix_elements
567 type(symbolicelementvector) :: symbolic_elements
568 real(wp) :: mat_coeff
569 integer :: l1, mat_offset, loop_ido, loop_skip, my_idx
570
571 this % diagonal_flag = 0
572
573 !Construct our symbolic elements
574 call symbolic_elements % construct
575
576 loop_skip = max(1, grid % gprocs)
577 my_idx = max(1, grid % grank)
578
579 !Loop through the L2 functions
580 do l1 = this % start_L2, this % total_num_csf
581 call matrix_elements % update_pure_L2(.false.)
582 !Slater them for the prototypes
583 call this % slater_rules(this % csfs(l1), this % csfs(l1), symbolic_elements, 0)
584 !Empty? then leave!
585 if (symbolic_elements % is_empty()) cycle
586 !Evaluate the integral
587 mat_coeff = this % evaluate_integrals(symbolic_elements, normal_phase)
588 !Insert using an offset
589 call matrix_elements % insert_matrix_element(l1 + this % L2_mat_offset, l1 + this % L2_mat_offset, mat_coeff, 2)
590 !Clear the symbols for the next step
591 call symbolic_elements % clear
592
593 this % non_zero_prototype = this % non_zero_prototype + 1
594 end do
595
596 !Destroy the symbols
597 call symbolic_elements % destroy
598
599 end subroutine evaluate_class_2
600
601
602 !> \brief Builds the class 3 part of the matrix
603 !> \authors A Al-Refaie
604 !> \date 2017
605 !>
606 !> This relates to building the off-diagonal elements of the target states with the same target symmetries.
607 !>
608 !> \param[inout] this Hamiltonian object to query.
609 !> \param[in] i1 The target symmetry.
610 !> \param[in] num_targets The number of target states withing the symmetry.
611 !> \param[out] matrix_elements A BaseMatrix object, stores the calculated elements.
612 !>
613 subroutine evaluate_class_3 (this, i1, num_targets, matrix_elements)
614 class(contracted_hamiltonian) :: this
615 class(basematrix), intent(inout) :: matrix_elements
616 integer, intent(in) :: i1, num_targets
617
618 !>These are the contracted prototypes for each combination of target states
619 type(contractedsymbolicelementvector) :: master_prototype
620 !>This stores te resultant slater rule calculation of prototypes before contracting for each state
621 type(symbolicelementvector) :: temp_prototype
622 !>The starting matrix index for each target state combination
623 integer, allocatable :: matrix_index(:,:)
624 !>The number of target state combinations
625 integer :: num_ci
626 !------------------------------target state variables
627 integer :: n1, n2
628 !----------------------------Looping variables------------------------------!
629 integer :: starting_index, ido, jdo, ii
630 !----------------------------MPI looping variables--------------------------
631 integer :: loop_idx, loop_ido, my_idx, loop_skip, err
632 !>The total number of configurations with this target symmetry
633 integer :: num_configurations
634 !>CSF indicies
635 integer :: config_idx_a, config_idx_b
636 !The target coefficient
637 real(wp), allocatable :: alpha_coeffs(:,:)
638 !The resultant matrix coefficient
639 real(wp), allocatable :: mat_coeffs(:,:)
640 !Number of continuum orbitals
641 integer :: continuum_orbital_a, continuum_orbital_b, num_continuum_orbitals, total_orbs
642
643 allocate(matrix_index(2,num_targets*(num_targets+1)/2), stat = err)
644 call master_memory % track_memory(storage_size(matrix_index)/8, size(matrix_index), err, 'CLASS3::MATRIX_INDEX')
645
646 allocate(alpha_coeffs(num_targets*(num_targets+1)/2,1), stat = err)
647 call master_memory % track_memory(storage_size(alpha_coeffs)/8, size(alpha_coeffs), err, 'CLASS3::alpha_coeffs')
648
649 allocate(mat_coeffs(num_targets*(num_targets+1)/2,1), stat = err)
650 call master_memory % track_memory(storage_size(mat_coeffs)/8, size(mat_coeffs), err, 'CLASS3::mat_coeffs')
651
652 !Get number of continuum orbitals
653 num_continuum_orbitals = this % options % num_continuum_orbitals_target(i1)
654
655 !number of target state combinations (triangular form
656 num_ci = num_targets * (num_targets + 1) / 2
657
658 ii = 1
659 do n1 = 1, this % options % num_target_state_sym(i1)
660 do n2 = n1, this % options % num_target_state_sym(i1)
661 !Construct a master prototype for each target state combo
662 !And get the starting matrix index
663 call this % compute_matrix_ij(i1, n1, 2, i1, n2, 1, matrix_index(1,ii), matrix_index(2,ii))
664 ii = ii + 1
665 end do
666 end do
667
668 call master_prototype % construct(num_targets * (num_targets + 1) / 2, 1)
669
670 !Construct the temporary prototype
671 call temp_prototype % construct
672
673 !get the starting index of the configurations
674 starting_index = this % get_starting_index(i1)
675
676 !Get the number of configurations
677 num_configurations = this % options % num_ci_target_sym(i1)
678
679 !Set up our mpi variables
680 loop_skip = max(1, grid % lprocs)
681 my_idx = max(grid % lrank, 0)
682 call master_timer % start_timer("Class 3 Prototype")
683
684 !Do the diagonal
685 do loop_ido = 1, num_configurations, loop_skip
686 !Calculate our real index
687 ido = loop_ido + my_idx
688 !If we've gon past (which is possible in this setup then skip
689 if (ido > num_configurations) cycle
690 !Calculate our diagonal configuration index
691 config_idx_a = starting_index + (ido - 1) * this % csf_skip
692 !The second configuration is simply offset by one
693 config_idx_b = config_idx_a + 1
694 !Slater it to get our temporary prototypes
695 call this % slater_rules(this % csfs(config_idx_b), this % csfs(config_idx_a), temp_prototype, 0, .false.)
696
697 alpha_coeffs = 0.0_wp
698 !Similar contraction with class 1 Left hand of Eq. 6
699 ii = 1
700 do n1 = 1, this % options % num_target_state_sym(i1)
701 do n2 = n1, this % options % num_target_state_sym(i1)
702 alpha_coeffs(ii,1) = this % get_target_coeff(i1, ido, n1) * this % get_target_coeff(i1, ido, n2)
703 ii = ii + 1
704 end do
705 end do
706 call master_prototype % add_symbols(temp_prototype, alpha_coeffs)
707 !call this%fast_contract_class_1_3_diag(i1,num_targets,ido,temp_prototype,master_prototype)
708 !Clear temporary prototype
709 call temp_prototype % clear
710
711 !Now loop through offdiagonal
712 do jdo = 1, num_configurations
713 !Skip diagonal (already done above)
714 if (jdo == ido) cycle
715 !Get the second CSF index offset by one
716 config_idx_b = starting_index + 1 + (jdo - 1) * this % csf_skip
717 !Slater
718 call this % slater_rules(this % csfs(config_idx_b), this % csfs(config_idx_a), temp_prototype, 1, .false.)
719
720 !No prototypes? then we leave
721 if (temp_prototype % is_empty()) cycle
722 alpha_coeffs = 0.0_wp
723 !Contract similary to class 1
724 ! (cimn*cim'n' + cim'n*cimn')H_imj,imj' Right hand of Eq. 6
725 ii = 1
726 do n1 = 1, this % options % num_target_state_sym(i1)
727 do n2 = n1, this % options % num_target_state_sym(i1)
728 alpha_coeffs(ii,1) = this % get_target_coeff(i1, ido, n1) * this % get_target_coeff(i1, jdo, n2)
729 ii = ii + 1
730 end do
731 end do
732 call master_prototype % add_symbols(temp_prototype, alpha_coeffs)
733 !call this%fast_contract_class_1_3_offdiag(i1,num_targets,ido,jdo,temp_prototype,master_prototype)
734 call temp_prototype % clear
735 end do
736 end do
737 call master_timer % stop_timer("Class 3 Prototype")
738 !call master_timer % report_timers
739 call master_prototype % synchronize_symbols()
740
741 !Again without expansion, evaluate the first elements
742 !ii=1
743 !do n1=1,num_targets
744 ! do n2=n1,num_targets
745 ! cimn*cimn'
746
747 ! call master_prototype(ii)%synchronize_symbols()
748
749 ! if(myrank == master) then
750 ! mat_coeffs(1) = this%evaluate_integrals(master_prototype(ii),NORMAL_PHASE)
751 !
752 ! call matrix_elements%insert_matrix_element(matrix_index(1,ii),matrix_index(2,ii),mat_coeffs(1),8)
753 ! write(2014,*) matrix_index(1,ii),matrix_index(2,ii),mat_coeffs(1)
754 ! endif
755 ! call matrix_elements%update_pure_L2(.false.)
756 !this time it is possible that if the target states are not the same then we 'flip' the index over and insert it
757 ! if (n1 /= n2) then
758 ! if(myrank == master) then
759 ! call matrix_elements%insert_matrix_element(matrix_index(1,ii)-1,matrix_index(2,ii)+1,mat_coeffs(1),8)
760 ! write(2014,*) matrix_index(1,ii)-1,matrix_index(2,ii)+1,mat_coeffs(1)
761 ! endif
762 ! call matrix_elements%update_pure_L2(.false.)
763 ! endif
764 !Calaculate the next index
765 ! matrix_index(1,ii) = matrix_index(1,ii) - 1
766
767 ! ii=ii+1
768
769 ! enddo
770 !enddo
771
772 !Now get expansion variable variables
773 config_idx_a = starting_index + (num_configurations - 1) * this % csf_skip
774 config_idx_b = config_idx_a + 1
775
776 continuum_orbital_a = this % csfs(config_idx_a) % orbital_sequence_number
777 continuum_orbital_b = this % csfs(config_idx_b) % orbital_sequence_number
778 call master_timer % start_timer("Class 3 Expansion")
779
780 !Set up our mpi variables
781 loop_skip = max(1, grid % gprocs)
782 my_idx = max(grid % grank, 0)
783
784 total_orbs = compute_total_box(num_continuum_orbitals, num_continuum_orbitals)
785
786 do loop_ido = 1, total_orbs, loop_skip
787
788 loop_idx = loop_ido + my_idx
789 call matrix_elements % update_pure_L2(.false., num_targets * (num_targets - 1))
790 !Skip if beyond total
791 if (loop_idx > total_orbs) cycle
792 !Convert the 1-D index into 2-D
793 call box_index_to_ij(loop_idx, num_continuum_orbitals, ido, jdo)
794
795 if (ido > num_continuum_orbitals) cycle
796 if (jdo > num_continuum_orbitals) cycle
797 if (jdo < ido + 1) cycle
798 if (ido == jdo) cycle
799
800 !Loop through each continuum orbital
801 !do ido=1,num_continuum_orbitals-1
802 !do jdo=max(3,ido+1),num_continuum_orbitals
803 ii = 1
804 !Loop through each target state
805 !Expand and Evaluate the matrix element
806
807 call this % contr_expand_off_diagonal_eval(continuum_orbital_a, continuum_orbital_b, &
808 master_prototype, ido, jdo, mat_coeffs)
809
810 do n1 = 1, num_targets
811 do n2 = n1, num_targets
812 !Insert into the matrix
813 call matrix_elements % insert_matrix_element(matrix_index(2,ii) + jdo - 1, &
814 matrix_index(1,ii) + ido - 2, mat_coeffs(ii,1), 8)
815 !print *,matrix_index(2,ii)+jdo-2,matrix_index(1,ii)+ido-1
816 !write(2014,*) matrix_index(1,ii)+jdo-2,matrix_index(2,ii)+ido-1,mat_coeffs(ii)
817 !call matrix_elements%update_continuum(.false.)
818 !this time it is possible that if the target states are not the same then we 'flip' the index over and insert it
819 if (n1 /= n2) then
820 call matrix_elements % insert_matrix_element (matrix_index(2,ii) + ido - 1, &
821 matrix_index(1,ii) + jdo - 2, mat_coeffs(ii,1), 8)
822 !write(2014,*) matrix_index(1,ii)+ido-2,matrix_index(2,ii)+jdo-1,mat_coeffs(ii)
823 !call matrix_elements%update_continuum(.false.)
824 end if
825 ii = ii + 1
826 end do
827 end do
828 !enddo
829 end do
830
831 call master_timer % stop_timer("Class 3 Expansion")
832 !call master_timer%report_timers
833
834 !Clean up master
835 call master_prototype % destroy
836
837 !Clean up temporary
838 call temp_prototype % destroy
839
840 this % non_zero_prototype = this % non_zero_prototype + 1
841
842 call master_memory % free_memory(storage_size(matrix_index)/8, size(matrix_index))
843 call master_memory % free_memory(storage_size(alpha_coeffs)/8, size(alpha_coeffs))
844 call master_memory % free_memory(storage_size(mat_coeffs)/8, size(mat_coeffs))
845 deallocate(matrix_index, alpha_coeffs, mat_coeffs)
846
847 end subroutine evaluate_class_3
848
849
850 !> \brief Builds the class 4 part of the matrix
851 !> \authors A Al-Refaie
852 !> \date 2017
853 !>
854 !> This relates to building the target-L2 off-diagonal elements of the hamiltonian.
855 !>
856 !> \param[inout] this Hamiltonian object to query.
857 !> \param[in] i1 The target symmetry
858 !> \param[in] num_states The number of target states withing the symmetry
859 !> \param[in] l1 The index to the L2 function
860 !> \param[out] matrix_elements A BaseMatrix object, stores the calculated elements
861 !>
862 subroutine evaluate_class_4 (this, i1, num_states, l1, matrix_elements)
863
864 class(contracted_hamiltonian) :: this
865 class(basematrix), intent(inout) :: matrix_elements
866 integer, intent(in) :: i1, num_states, l1
867
868 integer :: config_idx_a ! CSF index
869 integer :: num_continuum_orbitals ! Number of continuum orbitals
870 integer :: num_configurations ! Number of configurations within the target symmetry
871 integer :: n1 ! Target state index
872 integer :: continuum_orbital ! Continuum orbital index
873
874 integer :: ido, jdo, starting_index, i, j ! Loop variables
875 integer :: loop_skip, loop_idx, my_idx, err ! MPI loop variables
876 integer, allocatable :: matrix_index(:) ! The starting matrix index for each state
877
878 type(contractedsymbolicelementvector) :: master_prototype ! The contracted prototypes for each target state
879 type(symbolicelementvector) :: temp_prototype ! Resultant slater rule calculation of symbols before contracting for each state
880 real(wp), allocatable :: alpha_coeffs(:,:), mat_coeffs(:) ! Target coefficient and matrix element coefficient
881
882 allocate(matrix_index(num_states), stat = err)
883 call master_memory % track_memory(storage_size(matrix_index)/8, size(matrix_index), err, 'CLASS4::MATRIX_INDEX')
884
885 allocate(alpha_coeffs(num_states,1), stat = err)
886 call master_memory % track_memory(storage_size(alpha_coeffs)/8, size(alpha_coeffs), err, 'CLASS4::alpha_coeffs')
887
888 allocate(mat_coeffs(num_states), stat = err)
889 call master_memory % track_memory(storage_size(mat_coeffs)/8, size(mat_coeffs), err, 'CLASS4::mat_coeffs')
890
891 !Get number of continuum orbitals
892 num_continuum_orbitals = this % options % num_continuum_orbitals_target(i1)
893 !Get number of configuration state functions for this symmtery
894 num_configurations = this % options % num_ci_target_sym(i1)
895
896 !get the starting index of the configurations
897 starting_index = this % get_starting_index(i1)
898
899 do n1 = 1, num_states
900 !Construct and clear a symbol list for each state
901 !Compute starting index for each state
902 call this % compute_matrix_ij(i1, n1, 1, 1, 1, 1, matrix_index(n1), j)
903 end do
904
905 call master_prototype % construct(num_states, 1)
906
907 !Construct our temporary symbols
908 call temp_prototype % construct
909 !Clear them just in case
910 call temp_prototype % clear
911
912 !Setup our loop variables
913 !loop_skip = max(1,nprocs)
914 !my_idx = max(myrank,0)
915
916 !Again like class 1 and 3 we loop in a strided fashion across mpi nodes (or in a normal fashion without mpi)
917 do ido = 1, num_configurations
918 !Calculate our real index
919 !loop_idx = ido + my_idx
920 !If we've gone past it then leave
921 !if(loop_idx > num_configurations) cycle
922 !Caculate our CSf index
923 config_idx_a = starting_index + (ido - 1) * this % csf_skip
924
925 !Slater it with the L2 functions
926 call this % slater_rules(this % csfs(l1), this % csfs(config_idx_a), temp_prototype, 1, .false.)
927 !Nothing? Skip
928 if (temp_prototype % is_empty()) cycle
929
930 !Contract of the form cimn*Himj,l (Eq. 7)
931 do n1 = 1, num_states
932 alpha_coeffs(n1,1) = this % get_target_coeff(i1, ido, n1)
933 end do
934
935 call master_prototype % add_symbols(temp_prototype, alpha_coeffs)
936 call temp_prototype % clear
937 end do
938
939 call temp_prototype % clear
940
941 !Get last CSF in symmetry
942 config_idx_a = starting_index + (num_configurations - 1) * this % csf_skip
943 !Get the associated continuum orbital
944 continuum_orbital = this % csfs(config_idx_a) % orbital_sequence_number
945
946 !This is the L2 matrix position
947 j = l1 + this % L2_mat_offset
948
949 !do n1=1, num_states
950
951 ! call master_prototype(n1)%synchronize_symbols
952 !enddo
953
954 !For each state, evaluate the integrals without expansion and store
955 !do n1=1,num_states
956 ! mat_coeffs(n1) = this%evaluate_integrals(master_prototype(n1),NORMAL_PHASE)
957 ! call matrix_elements%insert_matrix_element(j,matrix_index(n1),mat_coeffs(n1),4)
958 ! call matrix_elements%update_continuum(.false.)
959 !enddo
960
961 !Start of continuum orbital expansion
962 do ido = 1, num_continuum_orbitals
963 !ido = loop_idx + my_idx
964 !call matrix_elements%update_pure_L2(.false.,num_states)
965 !if(ido > num_continuum_orbitals) cycle
966
967 call this % contr_expand_continuum_L2_eval(continuum_orbital, master_prototype, ido, mat_coeffs)
968
969 do n1 = 1, num_states
970 !Expand and evaluate
971 !Store coefficient into matrix object
972 call matrix_elements % insert_matrix_element(j, matrix_index(n1) + ido - 1, mat_coeffs(n1), 8, 0.0_wp)
973 !Perform update if needed
974 call matrix_elements % update_pure_L2(.false.)
975 end do
976 end do
977
978 !Clean up prototype symbols
979 call master_prototype % destroy
980 call temp_prototype % destroy
981
982 this % non_zero_prototype = this % non_zero_prototype + 1
983
984 call master_memory % free_memory(storage_size(matrix_index)/8, size(matrix_index))
985 call master_memory % free_memory(storage_size(alpha_coeffs)/8, size(alpha_coeffs))
986 call master_memory % free_memory(storage_size(mat_coeffs)/8, size(mat_coeffs))
987
988 deallocate(matrix_index, alpha_coeffs, mat_coeffs)
989
990 end subroutine evaluate_class_4
991
992
993 !> \brief Builds the class 5 part of the matrix
994 !> \authors A Al-Refaie
995 !> \date 2017
996 !>
997 !> This relates to building the 'diagonal' off-diagonal elements of the hamiltonian of differing target symmetries
998 !> but same continua symmetry.
999 !>
1000 !> \param[inout] this Hamiltonian object to query.
1001 !> \param[in] i1 Target symmetry 1
1002 !> \param[in] num_target_1 The number of target states within target symmetry 1
1003 !> \param[in] i2 Target symmetry 2
1004 !> \param[in] num_target_2 The number of target states within target symmetry 2
1005 !> \param[out] matrix_elements A BaseMatrix object, stores the calculated elements
1006 !>
1007 subroutine evaluate_class_5 (this, i1, num_target_1, i2, num_target_2, matrix_elements)
1008 class(contracted_hamiltonian) :: this
1009 class(basematrix), intent(inout) :: matrix_elements
1010 integer, intent(in) :: i1, i2, num_target_1, num_target_2
1011
1012 !>Contracted prototypes for each target state combination
1013 !type(SymbolicElementVector),target :: master_prototype(num_target_1,num_target_2)
1014 !type(SymbolicElementVector),pointer :: mst_pntr(:)
1015 type(contractedsymbolicelementvector) :: master_prototype
1016 !>Starting index for each
1017 integer, allocatable :: matrix_index(:,:,:)
1018 real(wp), allocatable :: alpha_coeffs(:,:)
1019 !>Temporary storage of Slater rule calculations
1020 type(symbolicelementvector) :: temp_prototype
1021 !>Configuration and target state variables
1022 integer :: config_idx_a,config_idx_b,n1,n2
1023 !> Number of continuum orbitals
1024 integer :: num_continuum_orbitals
1025 !>Counts for number of configurations of target symmetry 1,2 and in total TS1*TS2
1026 integer :: num_configurations_a,num_configurations_b,total_configurations
1027 !------------------------------MPI LOOP VARIABLES----------------------------------!
1028 integer :: loop_idx,my_idx,loop_skip,loop_ido
1029 !---------------------------------LOOP VARIABLES--------------------------------!
1030 integer :: ido,jdo,starting_index_a,starting_index_b,i,j,err
1031 !> Current continuum orbital
1032 integer :: continuum_orbital
1033 !
1034
1035 !> Target coefficient and matrix coefficient
1036 real(wp) :: alpha
1037 real(wp), allocatable :: mat_coeffs(:,:)
1038
1039 allocate(matrix_index(2, num_target_1, num_target_2), stat = err)
1040 call master_memory % track_memory(storage_size(matrix_index)/8, size(matrix_index), err, 'CLASS5::MATRIX_INDEX')
1041
1042 allocate(alpha_coeffs(num_target_1, num_target_2), stat = err)
1043 call master_memory % track_memory(storage_size(alpha_coeffs)/8, size(alpha_coeffs), err, 'CLASS5::alpha_coeffs')
1044
1045 allocate(mat_coeffs(num_target_1, num_target_2), stat = err)
1046 call master_memory % track_memory(storage_size(mat_coeffs)/8, size(mat_coeffs), err, 'CLASS5::mat_coeffs')
1047
1048 !Get number of continuum orbitals, which is a minimum of the two
1049 num_continuum_orbitals = min(this % options % num_continuum_orbitals_target(i1), &
1050 this % options % num_continuum_orbitals_target(i2))
1051
1052 !get the starting index of the configuration of target symmetry 1
1053 starting_index_a = this % get_starting_index(i1)
1054
1055 !get the starting index of the configurations of target symmetry 2
1056 starting_index_b = this % get_starting_index(i2)
1057
1058 !total number of configurations of target symmetry 1
1059 num_configurations_a = this % options % num_ci_target_sym(i1)
1060
1061 !total number of configurations of target symmetry 2
1062 num_configurations_b = this % options % num_ci_target_sym(i2)
1063
1064 !Compute starting matrix index and construct contracted symbolic prototypes for each combination of target states
1065 do n1 = 1, num_target_1
1066 do n2 = 1, num_target_2
1067 call this % compute_matrix_ij(i1, n1, 1, i2, n2, 1, matrix_index(1,n1,n2), matrix_index(2,n1,n2))
1068 end do
1069 end do
1070 call master_prototype % construct(num_target_1, num_target_2)
1071
1072 !Construct
1073 call temp_prototype % construct
1074 call temp_prototype % clear
1075
1076 !mst_pntr(1:num_target_1*num_target_2) => master_prototype(:,:)
1077 !mat_ptr(1:num_target_1*num_target_2) => mat_coeffs(:,:)
1078 !Compute the total number of configurations to loop through (num_config_a*num_config_b)
1079 total_configurations = compute_total_box(num_configurations_a, num_configurations_b)
1080
1081 !Set up our mpi variables
1082 loop_skip = max(1, grid % lprocs)
1083 my_idx = max(grid % lrank, 0)
1084
1085 call master_timer % start_timer("Class 5 prototype")
1086 !This loop differs from classes 1,3 and 4. Instead the 2 dimensional loop is collapsed into one dimension
1087 !And the expected CSFs are computed from the single index. This was done to improve load balancing significantly
1088 !and ensure better OpenMP looping when it is added eventually. If this wasn;t done then the last MPI process will have the largest
1089 !chunk of configurations to deal compared to the others.
1090 !Again this loop is strided by MPI processes
1091 do loop_ido = 1, total_configurations, loop_skip
1092 alpha_coeffs = 0.0_wp
1093
1094 !Calculate the real loop index
1095 loop_idx = loop_ido + my_idx
1096 !Skip if beyond total
1097 if (loop_idx > total_configurations) cycle
1098 !Convert the 1-D index into 2-D
1099 call box_index_to_ij(loop_idx, num_configurations_a, ido, jdo)
1100
1101 alpha_coeffs = 0.0_wp
1102
1103 !Use the two dimensional id to calculate the CSF index
1104 config_idx_a = starting_index_a + (ido - 1) * this % csf_skip
1105 config_idx_b = starting_index_b + (jdo - 1) * this % csf_skip
1106
1107 !Slater for prototype symbols
1108 call this % slater_rules(this % csfs(config_idx_b), this % csfs(config_idx_a), temp_prototype, 1, .false.)
1109
1110 if (temp_prototype % is_empty()) cycle
1111 !Perform contraction on each combination of target state of the form cimn*ci'm'n'*Himj,i'm'j' seen in Eq. 8
1112 !do n1=1,num_target_1
1113 ! do n2=1,num_target_2
1114 ! alpha = this%get_target_coeff(i1,ido,n1)*this%get_target_coeff(i2,jdo,n2)
1115 ! call master_prototype(n1,n2)%add_symbols(temp_prototype,alpha)
1116 ! enddo
1117 !enddo
1118 !call this%fast_contract_class_567(i1,i2,num_target_1,num_target_2,ido,jdo,temp_prototype,master_prototype)
1119 do n1 = 1, num_target_1
1120 do n2 = 1, num_target_2
1121 alpha_coeffs(n1,n2) = this % get_target_coeff(i1,ido,n1) * this % get_target_coeff(i2,jdo,n2)
1122 !call master_prototype(n1,n2)%add_symbols(temp_prototype,alpha)
1123 end do
1124 end do
1125 call master_prototype % add_symbols(temp_prototype, alpha_coeffs)
1126
1127 !clear the temporary symbols
1128 call temp_prototype % clear
1129 end do
1130
1131 call master_timer % stop_timer("Class 5 prototype")
1132 call temp_prototype % destroy
1133
1134 !call master_memory%print_memory_report()
1135 !write(stdout,"('Class 5 num symbols ',4i16)") master_prototype%get_size(),num_target_1,num_target_2,master_prototype%get_size()*num_target_1*num_target_2
1136 !Get the last CSf and find the associated orbital
1137 config_idx_a = starting_index_a + (num_configurations_a - 1) * this % csf_skip
1138 continuum_orbital = this % csfs(config_idx_a) % orbital_sequence_number
1139
1140 !Evaluate integrals for the first element without expansion. store and possible update the matrix.
1141 !do n1=1,num_target_1
1142 ! do n2=1,num_target_2
1143 ! mat_coeffs(n1,n2) = this%evaluate_integrals(master_prototype(n1,n2),NORMAL_PHASE)
1144 ! call matrix_elements%insert_matrix_element(matrix_index(2,n1,n2),matrix_index(1,n1,n2),mat_coeff,5)
1145 ! call matrix_elements%update_continuum(.false.)
1146 ! enddo
1147 !enddo
1148
1149 call master_timer % start_timer("Class 5 Synchonize")
1150 !write(stdout,"('Synch-Start ')")
1151 !call master_prototype(n1,n2)%print
1152 call master_prototype % synchronize_symbols
1153 !write(stdout,"('Synch-End ')")
1154 call master_timer % stop_timer("Class 5 Synchonize")
1155
1156 !Set up our mpi variables
1157 loop_skip = max(1, grid % gprocs)
1158 my_idx = max(grid % grank, 0)
1159
1160 !write(stdout,"('Expand ')")
1161 do loop_ido = 1, num_continuum_orbitals, loop_skip
1162 !Calculate the real loop index
1163 ido = loop_ido + my_idx
1164 call matrix_elements % update_pure_L2(.false., num_target_1 * num_target_2)
1165 !Skip if beyond total
1166 if (ido > num_continuum_orbitals) cycle
1167 mat_coeffs = 0
1168 call this % contr_expand_diagonal_eval(continuum_orbital, master_prototype, ido, mat_coeffs)
1169 !write(stdout,"('Orb ',i8)") ido
1170 do n1 = 1, num_target_1
1171 do n2 = 1, num_target_2
1172 !Expand and evaluate
1173 !Insert into the matrix
1174 call matrix_elements % insert_matrix_element(matrix_index(2,n1,n2) + ido - 1, &
1175 matrix_index(1,n1,n2) + ido - 1, mat_coeffs(n1,n2), 8)
1176 !Update matrix if needed.
1177 end do
1178 end do
1179 end do
1180
1181 !---Cleanup----!
1182
1183 call master_prototype % destroy
1184
1185 this % non_zero_prototype = this % non_zero_prototype + 1
1186
1187 call master_memory % free_memory(storage_size(matrix_index)/8, size(matrix_index))
1188 call master_memory % free_memory(storage_size(alpha_coeffs)/8, size(alpha_coeffs))
1189 call master_memory % free_memory(storage_size(mat_coeffs)/8, size(mat_coeffs))
1190
1191 deallocate(matrix_index, alpha_coeffs, mat_coeffs)
1192
1193 end subroutine evaluate_class_5
1194
1195
1196 !> \brief Builds the class 6 part of the matrix
1197 !> \authors A Al-Refaie
1198 !> \date 2017
1199 !>
1200 !> This relates to building the off-diagonal elements of the hamiltonian of differing target symmetries
1201 !> but same continua symmetry.
1202 !>
1203 !> \param[inout] this Hamiltonian object to query.
1204 !> \param[in] i1 Target symmetry 1
1205 !> \param[in] num_target_1 The number of target states within target symmetry 1
1206 !> \param[in] i2 Target symmetry 2
1207 !> \param[in] num_target_2 The number of target states within target symmetry 2
1208 !> \param[out] matrix_elements A BaseMatrix object, stores the calculated elements.
1209 !>
1210 subroutine evaluate_class_6 (this, i1, num_target_1, i2, num_target_2, matrix_elements)
1211 class(contracted_hamiltonian) :: this
1212 integer, intent(in) :: i1, i2, num_target_1, num_target_2
1213 class(basematrix), intent(inout) :: matrix_elements
1214
1215 !>Contracted prototypes for each target state combination
1216 !type(SymbolicElementVector),target :: master_prototype(num_target_1,num_target_2)
1217 !type(SymbolicElementVector),pointer :: mst_pntr(:)
1218 type(contractedsymbolicelementvector) :: master_prototype
1219
1220 !>Starting index for each
1221 integer, allocatable :: matrix_index(:,:,:)
1222 real(wp), allocatable :: alpha_coeffs(:,:)
1223
1224 !>Temporary storage of Slater rule calculations
1225 type(symbolicelementvector) :: temp_prototype
1226
1227 !>Configuration and target state variables
1228 integer :: config_idx_a, config_idx_b, n1, n2
1229
1230 !> Number of continuum orbitals
1231 integer :: num_continuum_orbitals_a, num_continuum_orbitals_b
1232
1233 !>Counts for number of configurations of target symmetry 1,2 and in total TS1*TS2
1234 integer :: num_configurations_a, num_configurations_b, total_configurations
1235
1236 !------------------------------MPI LOOP VARIABLES----------------------------------!
1237 integer :: loop_idx, my_idx, loop_skip, loop_ido
1238
1239 !---------------------------------LOOP VARIABLES--------------------------------!
1240 integer :: ido, jdo, starting_index_a, starting_index_b, i, j, err
1241
1242 !> Current continuum orbital for target symmetry 1 and 2
1243 integer :: continuum_orbital_a, continuum_orbital_b, total_orbs
1244
1245 !> Target coefficient and matrix coefficient
1246 real(wp) :: alpha
1247 real(wp), allocatable :: mat_coeffs(:,:)
1248
1249 allocate(matrix_index(2, num_target_1, num_target_2), stat = err)
1250 call master_memory % track_memory(storage_size(matrix_index)/8, size(matrix_index), err, 'CLASS6::MATRIX_INDEX')
1251
1252 allocate(alpha_coeffs(num_target_1, num_target_2), stat = err)
1253 call master_memory % track_memory(storage_size(alpha_coeffs)/8, size(alpha_coeffs), err, 'CLASS6::alpha_coeffs')
1254
1255 allocate(mat_coeffs(num_target_1, num_target_2), stat = err)
1256 call master_memory % track_memory(storage_size(mat_coeffs)/8, size(mat_coeffs), err, 'CLASS6::mat_coeffs')
1257
1258 !Get number of continuum orbitals for target symmetry 1
1259 num_continuum_orbitals_a = this % options % num_continuum_orbitals_target(i1)
1260
1261 !Get number of continuum orbitals for target symmetry 2
1262 num_continuum_orbitals_b = this % options % num_continuum_orbitals_target(i2)
1263
1264 !get the starting index of the configurations for target symmetry 1
1265 starting_index_a = this % get_starting_index(i1)
1266
1267 !get the starting index of the configurations for target symmetry 2
1268 starting_index_b = this % get_starting_index(i2)
1269
1270 !get the number of CSFS for target symmetry 1
1271 num_configurations_a = this % options % num_ci_target_sym(i1)
1272
1273 !get the number of CSFS for target symmetry 2
1274 num_configurations_b = this % options % num_ci_target_sym(i2)
1275
1276 !Compute starting matrix index and construct contracted symbolic prototypes for each combination of target states
1277 do n1 = 1, num_target_1
1278 do n2 = 1, num_target_2
1279 call this % compute_matrix_ij(i1, n1, 1, i2, n2, 1, matrix_index(1,n1,n2), matrix_index(2,n1,n2))
1280 end do
1281 end do
1282 call master_prototype % construct(num_target_1, num_target_2)
1283
1284 !Construct
1285 call temp_prototype % construct
1286 call temp_prototype % clear
1287 !mst_pntr(1:num_target_1*num_target_2) => master_prototype(:,:)
1288 !mat_ptr(1:num_target_1*num_target_2) => mat_coeffs(:,:)
1289 !Compute the total number of configurations to loop through (num_config_a*num_config_b)
1290 total_configurations = compute_total_box(num_configurations_a, num_configurations_b)
1291
1292 !Set up our mpi variables
1293 loop_skip = max(1, grid % lprocs)
1294 my_idx = max(grid % lrank, 0)
1295
1296 config_idx_a = 0
1297 config_idx_b = 0
1298
1299 call master_timer % start_timer("Class 6 prototype")
1300 !This loop differs from classes 1,3 and 4. Instead the 2 dimensional loop is collapsed into one dimension
1301 !And the expected CSFs are computed from the single index. This was done to improve load balancing significantly
1302 !and ensure better OpenMP looping when it is added eventually. If this wasn't done then the last MPI process will have the largest
1303 !chunk of configurations to deal compared to the others.
1304 !Again this loop is strided by MPI processes
1305 do loop_ido = 1, total_configurations, loop_skip
1306 alpha_coeffs = 0.0_wp
1307
1308 !Calculate the real loop index
1309 loop_idx = loop_ido + my_idx
1310 !Skip if beyond total
1311 if (loop_idx > total_configurations) cycle
1312 !Convert the 1-D index into 2-D
1313 call box_index_to_ij(loop_idx, num_configurations_a, ido, jdo)
1314
1315 !Use the two dimensional id to calculate the CSF index
1316 config_idx_a = starting_index_a + (ido - 1) * this % csf_skip
1317 config_idx_b = starting_index_b + 1 + (jdo - 1) * this % csf_skip
1318
1319 !Slater for prototype symbols
1320 call this % slater_rules(this % csfs(config_idx_b), this % csfs(config_idx_a), temp_prototype, 1, .false.)
1321
1322 if (temp_prototype % is_empty()) cycle
1323 !Perform contraction on each combination of target state of the form cimn*ci'm'n'*Himj,i'm'j' seen in Eq. 8
1324 !do n1=1,num_target_1
1325 ! do n2=1,num_target_2
1326 ! alpha = this%get_target_coeff(i1,ido,n1)*this%get_target_coeff(i2,jdo,n2)
1327 ! call master_prototype(n1,n2)%add_symbols(temp_prototype,alpha)
1328 ! enddo
1329 !enddo
1330 !call this%fast_contract_class_567(i1,i2,num_target_1,num_target_2,ido,jdo,temp_prototype,master_prototype)
1331 do n1 = 1, num_target_1
1332 do n2 = 1, num_target_2
1333 alpha_coeffs(n1,n2) = this % get_target_coeff(i1,ido,n1) * this % get_target_coeff(i2,jdo,n2)
1334 !call master_prototype(n1,n2)%add_symbols(temp_prototype,alpha)
1335 end do
1336 end do
1337 call master_prototype % add_symbols(temp_prototype, alpha_coeffs)
1338
1339 !clear the temporary symbols
1340 call temp_prototype % clear
1341 end do
1342
1343 call master_timer % stop_timer("Class 6 prototype")
1344 call temp_prototype % clear
1345
1346 ! When there are more ranks in the shared memory communicator than pairs of target CSFs,
1347 ! some ranks exit the previous DO cycle without setting config_idx_a and config_idx_b. These are needed
1348 ! just to get number of
1349 n1 = config_idx_a; call mpi_reduceall_max(n1, config_idx_a, grid % lcomm)
1350 n2 = config_idx_b; call mpi_reduceall_max(n2, config_idx_b, grid % lcomm)
1351
1352 !Get the continuum orbital numbers for both of the last CSFS
1353 continuum_orbital_a = this % csfs(config_idx_a) % orbital_sequence_number
1354 continuum_orbital_b = this % csfs(config_idx_b) % orbital_sequence_number
1355
1356 call master_timer % start_timer("Class 6 Synchonize")
1357
1358 !call master_prototype(n1,n2)%print
1359 call master_prototype % synchronize_symbols
1360 call master_timer % stop_timer("Class 6 Synchonize")
1361
1362 total_orbs = compute_total_box(num_continuum_orbitals_a, num_continuum_orbitals_b)
1363
1364 !Set up our mpi variables
1365 loop_skip = max(1, grid % gprocs)
1366 my_idx = max(grid % grank, 0)
1367
1368 !Loop through each orbital in target symmetry 1
1369 do loop_ido = 1, total_orbs, loop_skip
1370 !Loop through each orbital in target symmetry 2
1371 !do jdo=1,num_continuum_orbitals_b
1372 loop_idx = loop_ido + my_idx
1373 call matrix_elements%update_pure_L2(.false.,num_target_1*num_target_2)
1374
1375 !Skip if beyond total
1376 if(loop_idx > total_orbs) cycle
1377
1378 !Convert the 1-D index into 2-D
1379 call box_index_to_ij(loop_idx,num_continuum_orbitals_a,ido,jdo)
1380 if (ido == jdo) cycle !Class 5 already handled the diagonal so ignore
1381
1382 !Loop through each state in target symmetry 1
1383 !Expand and evaluate immediately
1384 call this % contr_expand_off_diagonal_eval (continuum_orbital_a, continuum_orbital_b, &
1385 master_prototype, ido, jdo, mat_coeffs)
1386 do n1 = 1, num_target_1
1387 !Loop through each state in target symmetry 2
1388 do n2 = 1, num_target_2
1389 !Insert into the matrix
1390 call matrix_elements % insert_matrix_element (matrix_index(2,n1,n2) + jdo - 1, &
1391 matrix_index(1,n1,n2) + ido - 1, &
1392 mat_coeffs(n1,n2), 8)
1393 !Update if neccessary
1394 end do
1395 end do
1396 end do
1397
1398 !Cleanup
1399
1400 call master_prototype % destroy
1401 call temp_prototype % destroy
1402
1403 this % non_zero_prototype = this % non_zero_prototype + 1
1404
1405 call master_memory % free_memory(storage_size(matrix_index)/8, size(matrix_index))
1406 call master_memory % free_memory(storage_size(alpha_coeffs)/8, size(alpha_coeffs))
1407 call master_memory % free_memory(storage_size(mat_coeffs)/8, size(mat_coeffs))
1408
1409 deallocate(matrix_index, alpha_coeffs, mat_coeffs)
1410
1411 end subroutine evaluate_class_6
1412
1413
1414 !> \brief Builds the class 7 part of the matrix
1415 !> \authors A Al-Refaie
1416 !> \date 2017
1417 !>
1418 !> This relates to building the off-diagonal elements of the hamiltonian of differing target symmetries and continua symmetry
1419 !>
1420 !> \param[inout] this Hamiltonian object to query.
1421 !> \param[in] i1 Target symmetry 1
1422 !> \param[in] num_target_1 The number of target states within target symmetry 1
1423 !> \param[in] i2 Target symmetry 2
1424 !> \param[in] num_target_2 The number of target states within target symmetry 2
1425 !> \param[out] matrix_elements A BaseMatrix object, stores the calculated elements
1426 !>
1427 subroutine evaluate_class_7 (this, i1, num_target_1, i2, num_target_2, matrix_elements)
1428 class(contracted_hamiltonian) :: this
1429 integer, intent(in) :: i1, i2, num_target_1, num_target_2
1430 class(basematrix), intent(inout) :: matrix_elements
1431
1432 ! Contracted prototypes for each target state combination
1433 !type(SymbolicElementVector), target :: master_prototype(num_target_1, num_target_2)
1434 !type(SymbolicElementVector), pointer :: mst_pntr(:)
1435 type(contractedsymbolicelementvector) :: master_prototype
1436
1437 ! Starting index for each
1438 integer, allocatable :: matrix_index(:,:,:)
1439 real(wp), allocatable :: alpha_coeffs(:,:)
1440
1441 ! Temporary storage of Slater rule calculations
1442 type(symbolicelementvector) :: temp_prototype
1443
1444 ! Configuration and target state variables
1445 integer :: config_idx_a, config_idx_b, n1, n2
1446
1447 ! Number of continuum orbitals
1448 integer :: num_continuum_orbitals_a, num_continuum_orbitals_b
1449
1450 ! Counts for number of configurations of target symmetry 1,2 and in total TS1*TS2
1451 integer :: num_configurations_a, num_configurations_b, total_configurations
1452
1453 !------------------------------MPI LOOP VARIABLES----------------------------------!
1454 integer :: loop_idx, my_idx, loop_skip, loop_ido
1455
1456 !---------------------------------LOOP VARIABLES--------------------------------!
1457 integer :: ido, jdo, starting_index_a, starting_index_b, i, j, err
1458
1459 ! Current continuum orbital for target symmetry 1 and 2
1460 integer :: continuum_orbital_a, continuum_orbital_b, total_orbs
1461
1462 ! Target coefficient and matrix coefficient
1463 real(wp) :: alpha
1464 real(wp), allocatable :: mat_coeffs(:,:)
1465
1466 allocate(matrix_index(2,num_target_1,num_target_2), stat = err)
1467 call master_memory % track_memory(storage_size(matrix_index)/8, size(matrix_index), err, 'CLASS7::MATRIX_INDEX')
1468
1469 allocate(alpha_coeffs(num_target_1,num_target_2), stat = err)
1470 call master_memory % track_memory(storage_size(alpha_coeffs)/8, size(alpha_coeffs), err, 'CLASS7::alpha_coeffs')
1471
1472 allocate(mat_coeffs(num_target_1,num_target_2), stat = err)
1473 call master_memory % track_memory(storage_size(mat_coeffs)/8, size(mat_coeffs), err, 'CLASS7::mat_coeffs')
1474
1475 config_idx_a = 0
1476 config_idx_b = 0
1477
1478 !Get number of continuum orbitals for target symmetry 1
1479 num_continuum_orbitals_a = this % options % num_continuum_orbitals_target(i1)
1480 !Get number of continuum orbitals for target symmetry 2
1481 num_continuum_orbitals_b = this % options % num_continuum_orbitals_target(i2)
1482
1483 !get the starting index of the configurations for target symmetry 1
1484 starting_index_a = this % get_starting_index(i1)
1485 !get the starting index of the configurations for target symmetry 2
1486 starting_index_b = this % get_starting_index(i2)
1487
1488 !get the number of CSFS for target symmetry 1
1489 num_configurations_a = this % options % num_ci_target_sym(i1)
1490 !get the number of CSFS for target symmetry 2
1491 num_configurations_b = this % options % num_ci_target_sym(i2)
1492
1493 !Compute starting matrix index and construct contracted symbolic prototypes for each combination of target states
1494 do n1 = 1, num_target_1
1495 do n2= 1, num_target_2
1496 call this % compute_matrix_ij(i1, n1, 1, i2, n2, 1, matrix_index(1,n1,n2), matrix_index(2,n1,n2))
1497 end do
1498 end do
1499
1500 call master_prototype % construct(num_target_1, num_target_2)
1501
1502 !Construct
1503 call temp_prototype % construct
1504 call temp_prototype % clear
1505 !mst_pntr(1:num_target_1*num_target_2) => master_prototype(:,:)
1506 !mat_ptr(1:num_target_1*num_target_2) => mat_coeffs(:,:)
1507 !Compute the total number of configurations to loop through (num_config_a*num_config_b)
1508 total_configurations = compute_total_box(num_configurations_a, num_configurations_b)
1509
1510 !Set up our mpi variables
1511 loop_skip = max(1, grid % lprocs)
1512 my_idx = max(grid % lrank, 0)
1513
1514 call master_timer % start_timer("Class 7 prototype")
1515 !This loop differs from classes 1,3 and 4. Instead the 2 dimensional loop is collapsed into one dimension
1516 !And the expected CSFs are computed from the single index. This was done to improve load balancing significantly
1517 !and ensure better OpenMP looping when it is added eventually. If this wasn't done then the last MPI process will have the largest
1518 !chunk of configurations to deal compared to the others.
1519 !Again this loop is strided by MPI processes
1520 do loop_ido = 1, total_configurations, loop_skip
1521 alpha_coeffs = 0.0_wp
1522
1523 !Calculate the real loop index
1524 loop_idx = loop_ido + my_idx
1525 !Skip if beyond total
1526 if (loop_idx > total_configurations) cycle
1527 !Convert the 1-D index into 2-D
1528 call box_index_to_ij(loop_idx, num_configurations_a, ido, jdo)
1529
1530 !Use the two dimensional id to calculate the CSF index
1531 config_idx_a = starting_index_a + (ido - 1) * this % csf_skip
1532 config_idx_b = starting_index_b + (jdo - 1) * this % csf_skip
1533
1534 !Slater for prototype symbols
1535 call this % slater_rules(this % csfs(config_idx_b), this % csfs(config_idx_a), temp_prototype, 1, .false.)
1536
1537 !Cycle if empty
1538 if (temp_prototype % is_empty()) cycle
1539 !Perform contraction on each combination of target state of the form cimn*ci'm'n'*Himj,i'm'j' seen in Eq. 8
1540 !do n1=1,num_target_1
1541 ! do n2=1,num_target_2
1542 ! alpha = this%get_target_coeff(i1,ido,n1)*this%get_target_coeff(i2,jdo,n2)
1543
1544 ! call master_prototype(n1,n2)%add_symbols(temp_prototype,alpha)
1545 ! enddo
1546 !enddo
1547 !call this%fast_contract_class_567(i1,i2,num_target_1,num_target_2,ido,jdo,temp_prototype,master_prototype)
1548 do n1 = 1, num_target_1
1549 do n2 = 1, num_target_2
1550 alpha_coeffs(n1,n2) = this % get_target_coeff(i1, ido, n1) * this % get_target_coeff(i2, jdo, n2)
1551 !call master_prototype(n1,n2)%add_symbols(temp_prototype,alpha)
1552 end do
1553 end do
1554 call master_prototype % add_symbols(temp_prototype, alpha_coeffs)
1555
1556 !clear the temporary symbols
1557 call temp_prototype % clear
1558 end do
1559
1560 call master_timer % stop_timer("Class 7 prototype")
1561 call temp_prototype % clear
1562
1563 ! When there are more ranks in the shared memory communicator than pairs of target CSFs,
1564 ! some ranks exit the previous DO cycle without setting config_idx_a and config_idx_b. These are needed
1565 ! just to get number of
1566 n1 = config_idx_a; call mpi_reduceall_max(n1, config_idx_a, grid % lcomm)
1567 n2 = config_idx_b; call mpi_reduceall_max(n2, config_idx_b, grid % lcomm)
1568
1569 !Get the continuum orbital numbers for both of the last CSFS
1570 continuum_orbital_a = this % csfs(config_idx_a) % orbital_sequence_number
1571 continuum_orbital_b = this % csfs(config_idx_b) % orbital_sequence_number
1572
1573 call master_timer % start_timer("Class 7 Synchonize")
1574 call master_prototype % synchronize_symbols
1575 call master_timer % stop_timer("Class 7 Synchonize")
1576
1577 total_orbs = compute_total_box(num_continuum_orbitals_a, num_continuum_orbitals_b)
1578
1579 !Set up our mpi variables
1580 loop_skip = max(1, grid % gprocs)
1581 my_idx = max(grid % grank, 0)
1582
1583 call master_timer % start_timer("Class 7 expansion")
1584
1585 !Loop through each orbital in target symmetry 1
1586 do loop_ido = 1, total_orbs, loop_skip
1587 !Loop through each orbital in target symmetry 2
1588 !do jdo=1,num_continuum_orbitals_b
1589 loop_idx = loop_ido + my_idx
1590 call matrix_elements % update_pure_L2(.false., num_target_1 * num_target_2)
1591 !Skip if beyond total
1592 if (loop_idx > total_orbs) cycle
1593 !Convert the 1-D index into 2-D
1594 call box_index_to_ij(loop_idx, num_continuum_orbitals_a, ido, jdo)
1595 !Loop through each state in target symmetry 1
1596 !Expand and evaluate immediately
1597 call this % expand_off_diagonal_gen_eval(continuum_orbital_a, continuum_orbital_b, master_prototype, &
1598 ido, jdo, mat_coeffs)
1599
1600 do n1 = 1, num_target_1
1601 !Loop through each state in target symmetry 2
1602 do n2 = 1, num_target_2
1603 !Insert into the matrix
1604 call matrix_elements % insert_matrix_element(matrix_index(2,n1,n2) + jdo - 1, &
1605 matrix_index(1,n1,n2) + ido - 1, &
1606 mat_coeffs(n1,n2), 8)
1607 !Update if neccessary
1608 end do
1609 end do
1610 end do
1611
1612 call master_timer % stop_timer("Class 7 expansion")
1613
1614 !-------------Cleanup---------------!
1615 call temp_prototype % destroy
1616 call master_prototype % destroy
1617
1618 this % non_zero_prototype = this % non_zero_prototype + 1
1619
1620 call master_memory % free_memory(storage_size(matrix_index)/8, size(matrix_index))
1621 call master_memory % free_memory(storage_size(alpha_coeffs)/8, size(alpha_coeffs))
1622 call master_memory % free_memory(storage_size(mat_coeffs)/8, size(mat_coeffs))
1623
1624 deallocate(matrix_index, alpha_coeffs, mat_coeffs)
1625
1626 end subroutine evaluate_class_7
1627
1628
1629 !> \brief Builds the class 8 part of the matrix
1630 !> \authors A Al-Refaie
1631 !> \date 2017
1632 !>
1633 !> \deprecated Now replaced with the fast class_2_and_8 combined subroutine.
1634 !>
1635 !> \param[inout] this Hamiltonian object to query.
1636 !> \param[out] matrix_elements A BaseMatrix object, stores the calculated elements
1637 !>
1638 subroutine evaluate_class_8 (this, matrix_elements)
1639 class(contracted_hamiltonian) :: this
1640 class(basematrix), intent(inout) :: matrix_elements
1641
1642 type(symbolicelementvector) :: symbolic_elements ! The smbolic elements
1643 integer :: l1, l2 ! L2 CSf indicies
1644 real(wp) :: mat_coeff ! The calculated matrix coefficient
1645 integer :: loop_ido, my_idx, loop_skip
1646
1647 !Construct the symbolic element object
1648 call symbolic_elements % construct
1649
1650 my_idx = max(1, grid % grank)
1651 loop_skip = max(1, grid % gprocs)
1652
1653 !Loop through each L2 function
1654 do loop_ido = this % start_L2, this % total_num_csf - 1, loop_skip
1655 l1 = loop_ido + my_idx
1656 call matrix_elements % update_pure_L2(.false.)
1657
1658 !Loop through the lower triangular of L2
1659 do l2 = this % start_L2, l1 - 1
1660
1661 !Check for an update here to ensure synchonization with all processes
1662 if (l1 > this % total_num_csf - 1) exit
1663
1664 !Slater the L2 functions
1665 call this % slater_rules(this % csfs(l1), this % csfs(l2), symbolic_elements, 1, .false.)
1666 !If empty then cycle
1667 if (symbolic_elements % is_empty()) cycle
1668 !Evaluate the integral
1669 mat_coeff = this % evaluate_integrals(symbolic_elements, normal_phase)
1670 !Insert the coffieicnt into the correct matrix position
1671 call matrix_elements % insert_matrix_element(l2 + this % L2_mat_offset, l1 + this % L2_mat_offset, mat_coeff, 8)
1672 !Clear for the next step
1673 call symbolic_elements % clear
1674
1675 this % non_zero_prototype = this % non_zero_prototype + 1
1676 end do
1677 end do
1678
1679 !Cleanup the symbols
1680 call symbolic_elements % destroy
1681
1682 end subroutine evaluate_class_8
1683
1684
1685 !> \brief Builds both the class 2 and class 8 (Pure L2) parts of the matrix elements
1686 !> \authors A Al-Refaie
1687 !> \date 2017
1688 !>
1689 !> This suborutine builds both the diagonal and off-diagonal elements of the L2 part of the matrix.
1690 !> It collapses the 2-dimensional loop into a one dimensional loop and then moves through the entire lower
1691 !> triangular portion o the matrix to calculate very quicky and in parallel (if mpi enabled) the matrix elements.
1692 !>
1693 !> \param[inout] this Hamiltonian object to query.
1694 !> \param[out] matrix_elements A BaseMatrix object, stores the calculated elements.
1695 !>
1696 subroutine evaluate_class_2_and_8 (this, matrix_elements)
1697 class(contracted_hamiltonian) :: this
1698 class(basematrix), intent(inout) :: matrix_elements
1699
1700 integer :: l1, l2 ! L2 indicies to be calculated
1701 type(symbolicelementvector) :: symbolic_elements ! The symbolic elements evaluated from Slaters rules
1702 real(wp) :: mat_coeff ! The calculated matrix coefficients
1703 integer :: loop_idx, ido ! Looping indicies
1704 integer :: total_vals ! Total number of CSFS to loop through
1705 integer :: loop_skip, my_idx ! MPI variables
1706 integer :: trig_n ! The total number of L2 CSFS
1707 integer :: diag ! Wheter we are diagonal or not
1708
1709 call symbolic_elements % construct ! Construct our symbol object
1710 trig_n = this % num_L2_functions ! Calculate the number of L2 functions
1711 total_vals = compute_total_box(trig_n, trig_n) ! Compute the total number of CSFS to loop through
1712
1713 this % diagonal_flag = 0
1714
1715 !Our MPI looping variasbles
1716 loop_skip = max(1, grid % gprocs)
1717 my_idx= max(grid % grank, 0)
1718
1719 write (stdout, "('start,total: ',2i12)") this % start_L2, total_vals
1720
1721 ! Whilst similar in concept to classes 5,6,7. Instead of collapsing an N^2 to one dimension it collapses an
1722 ! N(N+1)/2 loop into one dimension. The index calculation is a little mroe complex but offers the same performance
1723 ! benefits as classes 5,6 and 7
1724 do ido = 1, total_vals, loop_skip
1725
1726 !Update if needed
1727 call matrix_elements % update_pure_L2(.false.)
1728 call symbolic_elements % clear
1729
1730 !Calculate the real index
1731 loop_idx = ido + my_idx
1732 !If we've gone over then cycle
1733 if (loop_idx > total_vals) cycle
1734 !Convert the one dimensional index to two dimensions
1735 call box_index_to_ij(loop_idx, trig_n, l1, l2)
1736 !write (stdout, "('idx = ',i12,' l1 = ',i12,' l2 ',i12)") loop_idx, l1, l2
1737 !Get the correct L2 offsets
1738 l1 = l1 + this % start_L2 - 1
1739 l2 = l2 + this % start_L2 - 1
1740 !write (stdout, "('idx = ',i12,' d l1 = ',i12,'d l2 ',i12,'mat ',2i8)") loop_idx, l1, l2, l2 + this % L2_mat_offset, l1 + this % L2_mat_offset
1741
1742 if ((l2 + this % L2_mat_offset) < (l1 + this % L2_mat_offset)) cycle
1743
1744 !As a precaution, if we've gone past either, skip
1745 if (l1 > this % total_num_csf .or. l2 > this % total_num_csf) cycle
1746 if (l1 < this % start_L2 .or. l2 < this % start_L2 ) cycle
1747 !write(stdout,"('l1,2 = ',2i8)") l1,l2
1748 !stop "Error in L2 indexing"
1749 !endif
1750
1751 diag = 1
1752 if (l1 == l2) diag = 0
1753
1754 !Slater the CSFS
1755 call this % slater_rules(this % csfs(l1), this % csfs(l2), symbolic_elements, diag, .false.)
1756 !If empty, skip
1757 if (symbolic_elements % is_empty()) cycle
1758 !Evaluate normally
1759 mat_coeff = this % evaluate_integrals(symbolic_elements, normal_phase)
1760
1761 !if(l2/=l1 .and. abs(mat_coeff) >= DEFAULT_MATELEM_THRESHOLD) write(1235,"(3i8,es16.5)") diag,l2+this%L2_mat_offset,l1+this%L2_mat_offset,mat_coeff
1762 !if(l1==l2) call this%rmat_ci(1)%modify_L2_diagonal(mat_coeff)
1763
1764 !Insert the matrix element into the correct place
1765 call matrix_elements % insert_matrix_element(l2 + this % L2_mat_offset, l1 + this % L2_mat_offset, mat_coeff, 8)
1766
1767 !Clear for the next cycle
1768 this % non_zero_prototype = this % non_zero_prototype + 1
1769
1770 end do
1771
1772 !Cleanup
1773 call symbolic_elements % destroy
1774
1775 end subroutine evaluate_class_2_and_8
1776
1777
1778 !> \brief Finds the starting matrix index for paticular target symmetries and target states
1779 !> \authors A Al-Refaie
1780 !> \date 2017
1781 !>
1782 !> \param[inout] this Hamiltonian object to query.
1783 !> \param[in] i1 Target symmetry 1
1784 !> \param[in] n1 Target state of target symmetry 1
1785 !> \param[in] j1 Continuum orbital of Target state of target symmetry 1
1786 !> \param[in] i2 Target symmetry 2
1787 !> \param[in] n2 Target state of target symmetry 2
1788 !> \param[in] j2 Continuum orbital of Target state of target symmetry 2
1789 !> \param[out] i Resultant matrix co-ordinate i
1790 !> \param[out] j Resultant matrix co-ordinate j
1791 !>
1792 subroutine compute_matrix_ij (this, i1, n1, j1, i2, n2, j2, i, j)
1793 class(contracted_hamiltonian) :: this
1794 integer, intent(in) :: i1, n1, j1, i2, n2, j2
1795 integer, intent(out) :: i, j
1796 integer :: ido
1797
1798 i = 0
1799
1800 do ido = 2, i1
1801 i = i + this % options % num_target_state_sym(ido - 1) &
1802 * this % options % num_continuum_orbitals_target(ido - 1)
1803 end do
1804
1805 do ido = 2, n1
1806 i = i + this % options % num_continuum_orbitals_target(i1)
1807 end do
1808
1809 i = i + j1
1810 j = 0
1811
1812 do ido = 2, i2
1813 j = j + this % options % num_target_state_sym(ido - 1) &
1814 * this % options % num_continuum_orbitals_target(ido - 1)
1815 enddo
1816
1817 do ido = 2, n2
1818 j = j + this % options % num_continuum_orbitals_target(i2)
1819 end do
1820
1821 j = j + j2
1822
1823 end subroutine compute_matrix_ij
1824
1825!----------------------------------EXPANSION-----------------------------------------------------!
1826
1827 !> Deprecated
1828 subroutine expand_diagonal (this, continuum_orbital, master_prototype, j1, expanded_prototype)
1829 class(contracted_hamiltonian) :: this
1830 class(symbolicelementvector), intent(in) :: master_prototype, expanded_prototype
1831 integer, intent(in) :: continuum_orbital, j1
1832
1833 integer(longint) :: integral(NIDX)
1834 real(wp) :: integral_coeff
1835 integer :: lwd(8), lw1, lw2, ido, i, num_integrals
1836
1837 num_integrals = master_prototype % get_size()
1838
1839 do ido = 1, num_integrals
1840
1841 !Get the integral
1842 call master_prototype % get_coeff_and_integral(ido, integral_coeff, integral)
1843 call unpack_ints(integral, lwd)
1844
1845 lw1 = 0
1846 lw2 = 0
1847
1848 do i = 1, 4
1849 if (lwd(i) == continuum_orbital) then
1850 if(lw1 == 0) then
1851 lw1 = i
1852 else
1853 lw2 = i
1854 end if
1855 end if
1856 end do
1857
1858 if (lw1 == 0) then
1859 call expanded_prototype % insert_symbol(integral, integral_coeff, .false.)
1860 cycle
1861 else if (lw2 == 0) then
1862 lwd(lw1) = continuum_orbital + j1 - 1
1863 else
1864 lwd(lw1) = continuum_orbital + j1 - 1
1865 lwd(lw2) = continuum_orbital + j1 - 1
1866 end if
1867
1868 call expanded_prototype % insert_ijklm_symbol(lwd(1), lwd(2), lwd(3), lwd(4), lwd(5), integral_coeff, .false.)
1869
1870 end do
1871
1872 end subroutine expand_diagonal
1873
1874 !> Deprecated
1875 subroutine expand_off_diagonal (this, continuum_orbital_a, continuum_orbital_b, master_prototype, ja, jb, expanded_prototype)
1876 class(contracted_hamiltonian) :: this
1877 class(symbolicelementvector), intent(in) :: master_prototype, expanded_prototype
1878 integer, intent(in) :: continuum_orbital_a, continuum_orbital_b, ja, jb
1879
1880 integer(longint) :: integral(NIDX)
1881 real(wp) :: integral_coeff
1882 integer :: lwd(8), lwa, lwb, ido, i, num_integrals
1883
1884 num_integrals = master_prototype % get_size()
1885
1886 do ido = 1, num_integrals
1887
1888 !Get the integral
1889 call master_prototype % get_coeff_and_integral(ido, integral_coeff, integral)
1890 call unpack_ints(integral, lwd)
1891
1892 lwa = 0
1893 lwb = 0
1894
1895 do i = 1, 4
1896 if (lwd(i) == continuum_orbital_a) lwa = i
1897 if (lwd(i) == continuum_orbital_b) lwb = i
1898 end do
1899
1900 !No occurances
1901 if (max(lwa, lwb) == 0) then
1902 call expanded_prototype % insert_symbol(integral, integral_coeff, .false.)
1903 cycle
1904 else if (lwb == 0) then
1905 lwd(lwa) = continuum_orbital_a + ja - 1
1906 else if (lwa == 0) then
1907 lwd(lwb) = continuum_orbital_b + jb - 2
1908 else
1909 lwd(lwa) = continuum_orbital_a + ja - 1
1910 lwd(lwb) = continuum_orbital_b + jb - 2
1911 endif
1912
1913 !One occurance
1914
1915 call expanded_prototype % insert_ijklm_symbol(lwd(1), lwd(2), lwd(3), lwd(4), lwd(5), integral_coeff, .false.)
1916
1917 end do
1918
1919 end subroutine expand_off_diagonal
1920
1921
1922 !> Deprecated
1923 subroutine expand_continuum_l2 (this, continuum_orbital, master_prototype, j1, expanded_prototype)
1924 class(contracted_hamiltonian) :: this
1925 class(symbolicelementvector), intent(in) :: master_prototype, expanded_prototype
1926 integer, intent(in) :: continuum_orbital, j1
1927
1928 integer(longint) :: integral(NIDX)
1929 real(wp) :: integral_coeff
1930 integer :: lwd(8), lw, ido, i, num_integrals
1931
1932 num_integrals = master_prototype % get_size()
1933
1934 do ido = 1, num_integrals
1935
1936 !Get the integral
1937 call master_prototype % get_coeff_and_integral(ido, integral_coeff, integral)
1938 call unpack_ints(integral, lwd)
1939
1940 lw = 0
1941
1942 do i = 1, 4
1943 if (lwd(i) == continuum_orbital) lw=i
1944 end do
1945
1946 if (lw == 0) then
1947 call expanded_prototype % insert_symbol(integral, integral_coeff, .false.)
1948 cycle
1949 else
1950 lwd(lw) = continuum_orbital + j1 - 1
1951 end if
1952
1953 call expanded_prototype % insert_ijklm_symbol(lwd(1), lwd(2), lwd(3), lwd(4), lwd(5), integral_coeff, .false.)
1954
1955 end do
1956
1957 end subroutine expand_continuum_l2
1958
1959
1960 !> Deprecated
1961 subroutine expand_off_diagonal_general (this, continuum_orbital_a, continuum_orbital_b, master_prototype, &
1962 ja, jb, expanded_prototype)
1963 class(contracted_hamiltonian) :: this
1964 class(symbolicelementvector), intent(in) :: master_prototype, expanded_prototype
1965 integer, intent(in) :: continuum_orbital_a, continuum_orbital_b, ja, jb
1966
1967 integer(longint) :: integral(NIDX)
1968 real(wp) :: integral_coeff
1969 integer :: lwd(8), lwa, lwb, ido, i, num_integrals
1970
1971 num_integrals = master_prototype % get_size()
1972
1973 do ido = 1, num_integrals
1974
1975 !Get the integral
1976 call master_prototype % get_coeff_and_integral(ido, integral_coeff, integral)
1977 call unpack_ints(integral, lwd)
1978
1979 lwa = 0
1980 lwb = 0
1981
1982 do i = 1, 4
1983 if (lwd(i) == continuum_orbital_a) lwa = i
1984 if (lwd(i) == continuum_orbital_b) lwb = i
1985 end do
1986
1987 !No occurances
1988 if (max(lwa, lwb) == 0) then
1989 call expanded_prototype % insert_symbol(integral, integral_coeff, .false.)
1990 cycle
1991 else if(lwb == 0)then
1992 lwd(lwa) = continuum_orbital_a + ja - 1
1993 else if (lwa == 0) then
1994 lwd(lwb) = continuum_orbital_b + jb - 1
1995 else
1996 lwd(lwa) = continuum_orbital_a + ja - 1
1997 lwd(lwb) = continuum_orbital_b + jb - 1
1998 end if
1999
2000 !IF(ipair(lwd(1))+lwd(2) > ipair(lwd(3))+lwd(4)) then
2001 call expanded_prototype % insert_ijklm_symbol(lwd(1), lwd(2), lwd(3), lwd(4), lwd(5), integral_coeff, .false.)
2002 !else
2003 !call expanded_prototype%insert_ijklm_symbol(lwd(3),lwd(4),lwd(1),lwd(2),lwd(5),integral_coeff)
2004 !endif
2005 end do
2006
2007 end subroutine expand_off_diagonal_general
2008
2009
2010 !> \brief Expands the diagonal prototype symbols and evaluates for a single matrix element for co-ordiate (j1)
2011 !> \authors A Al-Refaie
2012 !> \date 2017
2013 !>
2014 !> \param[inout] this Hamiltonian object to query.
2015 !> \param[in] continuum_orbital Prototype continuum orbital
2016 !> \param[in] master_prototype Prototype symbolic elements to be expanded
2017 !> \param[in] j1 Desired continuum orbital
2018 !> \param[out] mat_coeffs Output array for the matrix element.
2019 !>
2020 subroutine expand_diagonal_eval (this, continuum_orbital, master_prototype, j1, mat_coeffs)
2021 class(contracted_hamiltonian) :: this
2022 integer, intent(in) :: continuum_orbital, j1
2023 type(symbolicelementvector), intent(in) :: master_prototype(:)
2024 real(wp), intent(out) :: mat_coeffs(:)
2025
2026 integer(longint) :: integral(NIDX)
2027 real(wp) :: integral_coeff
2028 integer :: lwd(8), lw1, lw2, ido, i, prototypes, num_prototypes, num_integrals
2029
2030 num_prototypes = size(master_prototype)
2031 num_integrals = master_prototype(1) % get_size()
2032 mat_coeffs = 0
2033
2034 if (num_integrals == 0) return
2035
2036 do ido = 1, num_integrals
2037
2038 !Get the integral
2039 call master_prototype(1) % get_coeff_and_integral(ido, integral_coeff, integral)
2040 call unpack_ints(integral, lwd)
2041
2042 integral_coeff = 1.0_wp
2043
2044 lw1 = 0
2045 lw2 = 0
2046
2047 do i = 1, 4
2048 if (lwd(i) == continuum_orbital) then
2049 if (lw1 == 0) then
2050 lw1 = i
2051 else
2052 lw2 = i
2053 end if
2054 end if
2055 end do
2056
2057 if (lw1 == 0) then
2058 integral_coeff = this % evaluate_integrals_singular(integral, 1.0_wp, normal_phase)
2059 do prototypes = 1, num_prototypes
2060 mat_coeffs(prototypes) = mat_coeffs(prototypes) &
2061 + master_prototype(prototypes) % get_coefficient(ido) * integral_coeff
2062 end do
2063 cycle
2064 else if (lw2 == 0) then
2065 lwd(lw1) = continuum_orbital + j1 - 1
2066 else
2067 lwd(lw1) = continuum_orbital + j1 - 1
2068 lwd(lw2) = continuum_orbital + j1 - 1
2069 end if
2070
2071 call pack_ints(lwd(1), lwd(2), lwd(3), lwd(4), lwd(5), 0, 0, 0, integral)
2072 integral_coeff = this % evaluate_integrals_singular(integral, 1.0_wp, normal_phase)
2073
2074 do prototypes = 1, num_prototypes
2075 mat_coeffs(prototypes) = mat_coeffs(prototypes) + master_prototype(prototypes)%get_coefficient(ido) * integral_coeff
2076 end do
2077
2078 end do
2079
2080 end subroutine expand_diagonal_eval
2081
2082
2083 !> \brief Expands the diagonal prototype symbols and evaluates for a single matrix element for co-ordiate (j1)
2084 !> \authors A Al-Refaie
2085 !> \date 2017
2086 !>
2087 !> \param[inout] this Hamiltonian object to query.
2088 !> \param[in] continuum_orbital Prototype continuum orbital
2089 !> \param[in] master_prototype Prototype symbolic elements to be expanded
2090 !> \param[in] j1 Desired continuum orbital
2091 !> \param[out] mat_coeffs Output array for the matrix element.
2092 !>
2093 subroutine contr_expand_diagonal_eval (this, continuum_orbital, master_prototype, j1, mat_coeffs)
2094 class(contracted_hamiltonian) :: this
2095 type(contractedsymbolicelementvector), intent(in) :: master_prototype
2096 integer, intent(in) :: continuum_orbital, j1
2097 real(wp), intent(out) :: mat_coeffs(:,:)
2098
2099 integer(longint) :: integral(NIDX)
2100 real(wp) :: integral_coeff
2101
2102 integer :: lwd(8), lw1, lw2, ido, i, prototypes, num_prototypes, num_integrals, num_targets_1, num_targets_2, n1, n2
2103
2104 num_integrals = master_prototype % get_size()
2105 num_targets_1 = master_prototype % get_num_targets_sym1()
2106 num_targets_2 = master_prototype % get_num_targets_sym2()
2107 mat_coeffs = 0
2108
2109 if (num_integrals == 0) return
2110
2111 do ido = 1, num_integrals
2112
2113 !Get the integral
2114 integral = master_prototype % get_integral_label(ido)
2115 call unpack_ints(integral, lwd)
2116
2117 integral_coeff = 1.0_wp
2118 lw1 = 0
2119 lw2 = 0
2120
2121 do i = 1, 4
2122 if (lwd(i) == continuum_orbital) then
2123 if (lw1 == 0) then
2124 lw1 = i
2125 else
2126 lw2 = i
2127 end if
2128 end if
2129 end do
2130
2131 if (lw1 == 0) then
2132 integral_coeff = this % evaluate_integrals_singular(integral, 1.0_wp, normal_phase)
2133 do n1 = 1,num_targets_1
2134 do n2 = 1,num_targets_2
2135 mat_coeffs(n1,n2) = mat_coeffs(n1,n2) + master_prototype % get_coefficient(ido, n1, n2) * integral_coeff
2136 end do
2137 end do
2138 cycle
2139 else if (lw2 == 0) then
2140 lwd(lw1) = continuum_orbital + j1 - 1
2141 else
2142 lwd(lw1) = continuum_orbital + j1 - 1
2143 lwd(lw2) = continuum_orbital + j1 - 1
2144 end if
2145
2146 call pack_ints(lwd(1), lwd(2), lwd(3), lwd(4), lwd(5), 0, 0, 0, integral)
2147 integral_coeff = this % evaluate_integrals_singular(integral, 1.0_wp, normal_phase)
2148
2149 do n1 = 1,num_targets_1
2150 do n2 = 1,num_targets_2
2151 mat_coeffs(n1,n2) = mat_coeffs(n1,n2) + master_prototype % get_coefficient(ido, n1, n2) * integral_coeff
2152 end do
2153 end do
2154
2155 end do
2156
2157 end subroutine contr_expand_diagonal_eval
2158
2159
2160 !> \brief Expands the continuum-L2 prototype symbols and evaluates for a single matrix element for co-ordinate (j1),L2
2161 !> \authors A Al-Refaie
2162 !> \date 2017
2163 !>
2164 !> \param[inout] this Hamiltonian object to query.
2165 !> \param[in] continuum_orbital Prototype continuum orbital
2166 !> \param[in] master_prototype Prototype symbolic elements to be expanded
2167 !> \param[in] j1 Desired continuum orbital
2168 !> \param[out] mat_coeffs Output array for the matrix element.
2169 !>
2170 subroutine expand_continuum_l2_eval (this, continuum_orbital, master_prototype, j1, mat_coeffs)
2171 class(contracted_hamiltonian) :: this
2172 class(symbolicelementvector), intent(in) :: master_prototype(:)
2173 integer, intent(in) :: continuum_orbital, j1
2174 real(wp), intent(out) :: mat_coeffs(:)
2175
2176 integer(longint) :: integral(NIDX)
2177 real(wp) :: integral_coeff
2178 integer :: lwd(8), lw, ido, i, prototypes, num_prototypes, num_integrals
2179
2180 num_prototypes = size(master_prototype)
2181 num_integrals = master_prototype(1) % get_size()
2182
2183 mat_coeffs = 0.0_wp
2184
2185 if (num_integrals == 0) return
2186
2187 do ido = 1, num_integrals
2188
2189 !Get the integral
2190 call master_prototype(1) % get_coeff_and_integral(ido, integral_coeff, integral)
2191 call unpack_ints(integral, lwd)
2192
2193 integral_coeff = 1.0_wp
2194 lw = 0
2195
2196 do i = 1, 4
2197 if (lwd(i) == continuum_orbital) lw = i
2198 end do
2199
2200 if (lw == 0) then
2201 integral_coeff = this % evaluate_integrals_singular(integral, 1.0_wp, normal_phase)
2202 do prototypes = 1, num_prototypes
2203 mat_coeffs(prototypes) = mat_coeffs(prototypes) &
2204 + master_prototype(prototypes) % get_coefficient(ido) * integral_coeff
2205 end do
2206 else
2207 lwd(lw) = continuum_orbital + j1 - 1
2208 end if
2209
2210 call pack_ints(lwd(1), lwd(2), lwd(3), lwd(4), lwd(5), 0, 0, 0, integral)
2211
2212 integral_coeff = this % evaluate_integrals_singular(integral, 1.0_wp, normal_phase)
2213
2214 do prototypes = 1, num_prototypes
2215 mat_coeffs(prototypes) = mat_coeffs(prototypes) + master_prototype(prototypes)%get_coefficient(ido) * integral_coeff
2216 end do
2217
2218 end do
2219
2220 end subroutine expand_continuum_l2_eval
2221
2222
2223 !> \brief Expands the continuum-L2 prototype symbols and evaluates for a single matrix element for co-ordinate (j1),L2
2224 !> \authors A Al-Refaie
2225 !> \date 2017
2226 !>
2227 !> \param[inout] this Hamiltonian object to query.
2228 !> \param[in] continuum_orbital Prototype continuum orbital
2229 !> \param[in] master_prototype Prototype symbolic elements to be expanded
2230 !> \param[in] j1 Desired continuum orbital
2231 !> \param[out] mat_coeffs Output array for the matrix element.
2232 !>
2233 subroutine contr_expand_continuum_l2_eval (this, continuum_orbital, master_prototype, j1, mat_coeffs)
2234 class(contracted_hamiltonian) :: this
2235 class(contractedsymbolicelementvector), intent(in) :: master_prototype
2236 integer, intent(in) :: continuum_orbital, j1
2237 real(wp), intent(out) :: mat_coeffs(:)
2238
2239 integer(longint) :: integral(NIDX)
2240 real(wp) :: integral_coeff
2241 integer :: lwd(8), lw, ido, i, num_targets_1, num_targets_2, n1, n2, num_integrals
2242
2243 num_targets_1 = master_prototype % get_num_targets_sym1()
2244 num_integrals = master_prototype % get_size()
2245 mat_coeffs = 0.0_wp
2246
2247 if (num_integrals == 0) return
2248
2249 do ido = 1, num_integrals
2250
2251 !Get the integral
2252 integral = master_prototype % get_integral_label(ido)
2253 call unpack_ints(integral, lwd)
2254
2255 integral_coeff = 1.0_wp
2256
2257 lw = 0
2258
2259 do i = 1, 4
2260 if (lwd(i) == continuum_orbital) lw = i
2261 end do
2262
2263 if (lw == 0) then
2264 integral_coeff = this % evaluate_integrals_singular(integral, 1.0_wp, normal_phase)
2265 do n1 = 1, num_targets_1
2266 mat_coeffs(n1) = mat_coeffs(n1) + master_prototype % get_coefficient(ido, n1, 1) * integral_coeff
2267 end do
2268 cycle
2269 else
2270 lwd(lw) = continuum_orbital + j1 - 1
2271 end if
2272
2273 call pack_ints(lwd(1), lwd(2), lwd(3), lwd(4), lwd(5), 0, 0, 0, integral)
2274
2275 integral_coeff = this % evaluate_integrals_singular(integral, 1.0_wp, normal_phase)
2276
2277 do n1 = 1, num_targets_1
2278 mat_coeffs(n1) = mat_coeffs(n1) + master_prototype % get_coefficient(ido, n1, 1) * integral_coeff
2279 end do
2280
2281 end do
2282
2283 end subroutine contr_expand_continuum_l2_eval
2284
2285
2286 !> \brief Expands the off-diagonal same/similar symmetry prototype symbols and evaluates for a single matrix element for co-ordinate offset (ja,jb)
2287 !> \authors A Al-Refaie
2288 !> \date 2017
2289 !>
2290 !> \param[inout] this Hamiltonian object to query.
2291 !> \param[in] continuum_orbital_a Prototype continuum orbital for target symmetry a
2292 !> \param[in] continuum_orbital_b Prototype continuum orbital for target symmetry b
2293 !> \param[in] master_prototype Prototype symbolic elements to be expanded
2294 !> \param[in] ja Desired continuum orbital of target symmetry a
2295 !> \param[in] jb Desired continuum orbital of target symmetry b
2296 !> \param[in] mat_coeffs Output array for the matrix element.
2297 !>
2298 subroutine expand_off_diagonal_eval (this, continuum_orbital_a, continuum_orbital_b, master_prototype, ja, jb, mat_coeffs)
2299 class(contracted_hamiltonian) :: this
2300 class(symbolicelementvector), intent(in) :: master_prototype(:)
2301 integer, intent(in) :: continuum_orbital_a, continuum_orbital_b, ja, jb
2302 real(wp), intent(out) :: mat_coeffs(:)
2303
2304 integer(longint) :: integral(NIDX)
2305 real(wp) :: integral_coeff
2306 integer :: lwd(8), lwa, lwb, ido, i, prototypes, num_prototypes, num_integrals
2307
2308 num_prototypes = size(master_prototype)
2309 num_integrals = master_prototype(1) % get_size()
2310 mat_coeffs = 0
2311
2312 if (num_integrals == 0) return
2313
2314 do ido = 1, num_integrals
2315
2316 !Get the integral
2317 call master_prototype(1) % get_coeff_and_integral(ido, integral_coeff, integral)
2318 call unpack_ints(integral, lwd)
2319
2320 integral_coeff = 1.0_wp
2321
2322 lwa = 0
2323 lwb = 0
2324
2325 do i = 1, 4
2326 if (lwd(i) == continuum_orbital_a) lwa = i
2327 if (lwd(i) == continuum_orbital_b) lwb = i
2328 end do
2329
2330 !No occurances
2331 if (max(lwa, lwb) == 0) then
2332 integral_coeff = this % evaluate_integrals_singular(integral, 1.0_wp, normal_phase)
2333
2334 do prototypes = 1,num_prototypes
2335 mat_coeffs(prototypes) = mat_coeffs(prototypes) &
2336 + master_prototype(prototypes) % get_coefficient(ido) * integral_coeff
2337 end do
2338 !expand_off_diagonal_eval = expand_off_diagonal_eval+this%evaluate_integrals_singular(integral,integral_coeff,NORMAL_PHASE)
2339 cycle
2340 else if (lwb == 0) then
2341 lwd(lwa) = continuum_orbital_a + ja - 1
2342 else if (lwa == 0) then
2343 lwd(lwb) = continuum_orbital_b + jb - 2
2344 else
2345 lwd(lwa) = continuum_orbital_a + ja - 1
2346 lwd(lwb) = continuum_orbital_b + jb - 2
2347 end if
2348
2349 call pack_ints(lwd(1), lwd(2), lwd(3), lwd(4), lwd(5), 0, 0, 0, integral)
2350
2351 integral_coeff = this % evaluate_integrals_singular(integral, 1.0_wp, normal_phase)
2352
2353 do prototypes = 1,num_prototypes
2354 mat_coeffs(prototypes) = mat_coeffs(prototypes) + master_prototype(prototypes)%get_coefficient(ido) * integral_coeff
2355 end do
2356
2357 end do
2358
2359 end subroutine expand_off_diagonal_eval
2360
2361
2362 !> \brief Expands the off-diagonal same/similar symmetry prototype symbols and evaluates for a single matrix element for co-ordinate offset (ja,jb)
2363 !> \authors A Al-Refaie
2364 !> \date 2017
2365 !>
2366 !> \param[inout] this Hamiltonian object to query.
2367 !> \param[in] continuum_orbital_a Prototype continuum orbital for target symmetry a
2368 !> \param[in] continuum_orbital_b Prototype continuum orbital for target symmetry b
2369 !> \param[in] master_prototype Prototype symbolic elements to be expanded
2370 !> \param[in] ja Desired continuum orbital of target symmetry a
2371 !> \param[in] jb Desired continuum orbital of target symmetry b
2372 !> \param[out] mat_coeffs Output array for the matrix element.
2373 !>
2374 subroutine contr_expand_off_diagonal_eval (this, continuum_orbital_a, continuum_orbital_b, master_prototype, ja, jb, mat_coeffs)
2375 class(contracted_hamiltonian) :: this
2376 class(contractedsymbolicelementvector), intent(in) :: master_prototype
2377 integer, intent(in) :: continuum_orbital_a, continuum_orbital_b, ja, jb
2378 real(wp), intent(out) :: mat_coeffs(:,:)
2379
2380 integer(longint) :: integral(NIDX)
2381 real(wp) :: integral_coeff
2382
2383 integer :: lwd(8), lwa, lwb, ido, i, prototypes, num_prototypes, num_integrals, num_targets_1, num_targets_2, n1, n2
2384
2385
2386 num_integrals = master_prototype % get_size()
2387 num_targets_1 = master_prototype % get_num_targets_sym1()
2388 num_targets_2 = master_prototype % get_num_targets_sym2()
2389 mat_coeffs = 0
2390
2391 if (num_integrals == 0) return
2392
2393 do ido = 1, num_integrals
2394
2395 !Get the integral
2396 integral = master_prototype % get_integral_label(ido)
2397 call unpack_ints(integral, lwd)
2398
2399 integral_coeff = 1.0_wp
2400
2401 lwa = 0
2402 lwb = 0
2403
2404 do i = 1, 4
2405 if (lwd(i) == continuum_orbital_a) lwa = i
2406 if (lwd(i) == continuum_orbital_b) lwb = i
2407 end do
2408
2409 !No occurances
2410 if (max(lwa, lwb) == 0) then
2411 integral_coeff = this % evaluate_integrals_singular(integral, 1.0_wp, normal_phase)
2412
2413 do n1 = 1, num_targets_1
2414 do n2 = 1, num_targets_2
2415 mat_coeffs(n1,n2) = mat_coeffs(n1,n2) + master_prototype % get_coefficient(ido, n1, n2) * integral_coeff
2416 end do
2417 end do
2418 !expand_off_diagonal_eval = expand_off_diagonal_eval+this%evaluate_integrals_singular(integral,integral_coeff,NORMAL_PHASE)
2419 cycle
2420 else if (lwb == 0) then
2421 lwd(lwa) = continuum_orbital_a + ja - 1
2422 else if (lwa == 0) then
2423 lwd(lwb) = continuum_orbital_b + jb - 2
2424 else
2425 lwd(lwa) = continuum_orbital_a + ja - 1
2426 lwd(lwb) = continuum_orbital_b + jb - 2
2427 end if
2428
2429 call pack_ints(lwd(1), lwd(2), lwd(3), lwd(4), lwd(5), 0, 0, 0, integral)
2430
2431 integral_coeff = this % evaluate_integrals_singular(integral, 1.0_wp, normal_phase)
2432
2433 do n1 = 1, num_targets_1
2434 do n2 = 1, num_targets_2
2435 mat_coeffs(n1,n2) = mat_coeffs(n1,n2) + master_prototype % get_coefficient(ido, n1, n2) * integral_coeff
2436 end do
2437 end do
2438
2439 end do
2440
2441 end subroutine contr_expand_off_diagonal_eval
2442
2443
2444 !> \brief Expands the off-diagonal generic prototype symbols and evaluates for a single matrix element for co-ordinate (ja,jb)
2445 !> \authors A Al-Refaie
2446 !> \date 2017
2447 !>
2448 !> \param[inout] this Hamiltonian object to query.
2449 !> \param[in] continuum_orbital_a Prototype continuum orbital for target symmetry a
2450 !> \param[in] continuum_orbital_b Prototype continuum orbital for target symmetry b
2451 !> \param[in] master_prototype Prototype symbolic elements to be expanded
2452 !> \param[in] ja Desired continuum orbital of target symmetry a
2453 !> \param[in] jb Desired continuum orbital of target symmetry b
2454 !> \param[out] mat_coeffs Output array for the matrix element.
2455 !>
2456 subroutine expand_off_diagonal_gen_eval (this, continuum_orbital_a, continuum_orbital_b, master_prototype, ja, jb, mat_coeffs)
2457 class(contracted_hamiltonian) :: this
2458 class(contractedsymbolicelementvector), intent(in) :: master_prototype
2459 integer, intent(in) :: continuum_orbital_a, continuum_orbital_b, ja, jb
2460 real(wp), intent(out) :: mat_coeffs(:,:)
2461
2462 integer(longint) :: integral(NIDX)
2463 real(wp) :: integral_coeff
2464 integer :: lwd(8), lwa, lwb, ido, i, n1, n2, num_targets_1, num_targets_2, ii, num_integrals
2465
2466 num_targets_1 = master_prototype % get_num_targets_sym1()
2467 num_targets_2 = master_prototype % get_num_targets_sym2()
2468 num_integrals = master_prototype % get_size()
2469 mat_coeffs = 0
2470
2471 if (num_integrals == 0) return
2472
2473 do ido = 1, num_integrals
2474
2475 !Get the integral
2476 integral = master_prototype % get_integral_label(ido)
2477 call unpack_ints(integral, lwd)
2478
2479 integral_coeff = 1.0_wp
2480
2481 lwa = 0
2482 lwb = 0
2483
2484 do i = 1, 4
2485 if (lwd(i) == continuum_orbital_a) lwa = i
2486 if (lwd(i) == continuum_orbital_b) lwb = i
2487 end do
2488
2489 !No occurances
2490 if (max(lwa, lwb) == 0) then
2491 integral_coeff = this % evaluate_integrals_singular(integral, 1.0_wp, normal_phase)
2492 do n1 = 1,num_targets_1
2493 do n2 = 1,num_targets_2
2494 mat_coeffs(n1,n2) = mat_coeffs(n1,n2) + master_prototype % get_coefficient(ido, n1, n2) * integral_coeff
2495 end do
2496 end do
2497 cycle
2498 else if(lwb == 0)then
2499 lwd(lwa) = continuum_orbital_a + ja - 1
2500 else if (lwa == 0) then
2501 lwd(lwb) = continuum_orbital_b + jb - 1
2502 else
2503 lwd(lwa) = continuum_orbital_a + ja - 1
2504 lwd(lwb) = continuum_orbital_b + jb - 1
2505 end if
2506
2507 call pack_ints(lwd(1), lwd(2), lwd(3), lwd(4), lwd(5), 0, 0, 0, integral)
2508
2509 integral_coeff = this % evaluate_integrals_singular(integral, 1.0_wp, normal_phase)
2510
2511 do n1 = 1, num_targets_1
2512 do n2 = 1, num_targets_2
2513 mat_coeffs(n1,n2) = mat_coeffs(n1,n2) + master_prototype % get_coefficient(ido, n1, n2) * integral_coeff
2514 end do
2515 end do
2516
2517 end do
2518
2519 end subroutine expand_off_diagonal_gen_eval
2520
2521
2522 !> \brief (Not used) Compresses several prototypes into a single prototype at position 1
2523 !> \authors A Al-Refaie
2524 !> \date 2017
2525 !>
2526 !> To be used for the OpenMP implementation by combining multiple symbols.
2527 !>
2528 !> \param[inout] this Hamiltonian object to query.
2529 !> \param[in] prototypes A symbolic prototype array to be compressed.
2530 !>
2531 subroutine reduce_prototypes (this, prototypes)
2532 class(contracted_hamiltonian), intent(in) :: this
2533 class(symbolicelementvector), intent(in) :: prototypes(:)
2534
2535 integer :: num_prototypes, ido
2536
2537 num_prototypes = size(prototypes)
2538
2539 if (num_prototypes <= 1) return
2540
2541 do ido = 2, num_prototypes
2542 call prototypes(1) % add_symbols(prototypes(ido))
2543 end do
2544
2545 end subroutine reduce_prototypes
2546
2547
2548 !> \brief Gets the starting index configuration state functions at a specific symmetry
2549 !> \authors A Al-Refaie
2550 !> \date 2017
2551 !>
2552 !> \param[inout] this Hamiltonian object to query.
2553 !> \param[in] symmetry Target symmetry.
2554 !>
2555 !> \result The starting CSF index for a given target symmetry
2556 !>
2557 integer function get_starting_index (this, symmetry)
2558 class(contracted_hamiltonian), intent(in) :: this
2559 integer, intent(in) :: symmetry
2560
2561 integer :: ido
2562
2563 get_starting_index = 1
2564 do ido = 2, symmetry
2565 get_starting_index = get_starting_index + this % options % num_ci_target_sym(ido - 1) * this % csf_skip
2566 end do
2567
2568 end function get_starting_index
2569
2570
2571 !> \brief Simply a wrapper to aid in clarity
2572 !> \authors A Al-Refaie
2573 !> \date 2017
2574 !>
2575 !> \param[inout] this Hamiltonian object to query.
2576 !> \param[in] i Target symmetry.
2577 !> \param[in] m Target state.
2578 !> \param[in] n Target configuration.
2579 !>
2580 !> \result target coefficient
2581 !>
2582 real(wp) function get_target_coeff (this, i, m, n)
2583 class(contracted_hamiltonian), intent(in) :: this
2584
2585 integer :: i, m, n
2586
2587 get_target_coeff = this % rmat_ci(i) % eigenvectors(m, n)
2588
2589 end function get_target_coeff
2590
2591
2592 subroutine fast_contract_class_1_3_diag (this, i1, num_targets, ido, temp_prototype, master_prototypes)
2593 class(contracted_hamiltonian) :: this
2594 integer, intent(in) :: i1, ido, num_targets
2595 class(symbolicelementvector), intent(in) :: temp_prototype, master_prototypes(num_targets * (num_targets + 1) / 2)
2596
2597 integer(longint) :: label(NIDX)
2598
2599 integer :: found_idx, ii, n1, n2, num_symbols, jdo
2600 real(wp) :: alpha(num_targets * (num_targets + 1) / 2), coeff
2601
2602 ii = 1
2603
2604 !get our coeffs
2605 do n1 = 1, num_targets
2606 do n2 = n1, num_targets
2607 ! Apply our coefficeints for each combination to the master, Left hand of Eq. 6
2608 alpha(ii) = this % get_target_coeff(i1, ido, n1) * this % get_target_coeff(i1, ido, n2)
2609 !call master_prototype(ii)%add_symbols(temp_prototype,alpha)
2610 ii = ii + 1
2611 end do
2612 end do
2613
2614 num_symbols = temp_prototype % get_size()
2615
2616 do jdo = 1, num_symbols
2617 call temp_prototype % get_coeff_and_integral(jdo, coeff, label)
2618 !We only need to perform the expensive check once per insertion
2619 if (master_prototypes(1) % check_same_integral(label, found_idx)) then
2620 ii = 1
2621 !get our coeffs
2622 do n1 = 1, num_targets
2623 do n2 = n1, num_targets
2624 call master_prototypes(ii) % modify_coeff(found_idx, alpha(ii) * coeff)
2625 ii = ii + 1
2626 end do
2627 end do
2628 else
2629 ii = 1
2630 !get our coeffs
2631 do n1 = 1, num_targets
2632 do n2 = n1, num_targets
2633 !Insert without a check
2634 call master_prototypes(ii) % insert_symbol(label, alpha(ii) * coeff, .false.)
2635 ii = ii + 1
2636 end do
2637 end do
2638 end if
2639 end do
2640
2641 end subroutine fast_contract_class_1_3_diag
2642
2643
2644 subroutine fast_contract_class_1_3_offdiag (this, i1, num_targets, s1, s2, temp_prototype, master_prototypes)
2645 class(contracted_hamiltonian) :: this
2646 integer, intent(in) :: i1, s1, s2, num_targets
2647 class(symbolicelementvector), intent(in) :: temp_prototype, master_prototypes(num_targets * (num_targets + 1) / 2)
2648
2649 integer(longint) :: label(NIDX)
2650
2651 integer :: found_idx, ii, n1, n2, num_symbols, jdo
2652 real(wp) :: alpha(num_targets * (num_targets + 1) / 2), coeff
2653
2654 ii = 1
2655
2656 !get our coeffs
2657 do n1 = 1, num_targets
2658 do n2 = n1, num_targets
2659 ! Apply our coefficeints for each combination to the master, Left hand of Eq. 6
2660 alpha(ii) = this % get_target_coeff(i1, s1, n1) * this % get_target_coeff(i1, s2, n2) &
2661 + this % get_target_coeff(i1, s2, n1) * this % get_target_coeff(i1, s1, n2)
2662 !call master_prototype(ii)%add_symbols(temp_prototype,alpha)
2663 ii = ii + 1
2664 end do
2665 end do
2666
2667 num_symbols = temp_prototype % get_size()
2668
2669 do jdo = 1, num_symbols
2670 call temp_prototype % get_coeff_and_integral(jdo, coeff, label)
2671 !We only need to perform the expensive check once per insertion
2672 if (master_prototypes(1) % check_same_integral(label, found_idx)) then
2673 ii = 1
2674 !get our coeffs
2675 do n1 = 1, num_targets
2676 do n2 = n1, num_targets
2677 call master_prototypes(ii) % modify_coeff(found_idx, alpha(ii) * coeff)
2678 ii = ii + 1
2679 end do
2680 end do
2681 else
2682 ii = 1
2683 !get our coeffs
2684 do n1 = 1, num_targets
2685 do n2 = n1, num_targets
2686 !Insert without a check
2687 call master_prototypes(ii) % insert_symbol(label, alpha(ii) * coeff, .false.)
2688 ii = ii + 1
2689 end do
2690 end do
2691 end if
2692 end do
2693
2694 end subroutine fast_contract_class_1_3_offdiag
2695
2696
2697 subroutine fast_contract_class_567 (this, i1, i2, num_targets1, num_targets2, s1, s2, temp_prototype, master_prototypes)
2698 class(contracted_hamiltonian) :: this
2699 integer, intent(in) :: i1, i2, s1, s2, num_targets1, num_targets2
2700 class(symbolicelementvector), intent(in) :: temp_prototype, master_prototypes(num_targets1, num_targets2)
2701
2702 integer(longint) :: label(NIDX)
2703
2704 integer :: found_idx, ii, n1, n2, num_symbols, jdo
2705 real(wp) :: alpha(num_targets1, num_targets2), coeff
2706
2707 alpha = 0
2708 !get our coeffs
2709 do n1 = 1, num_targets1
2710 do n2 = 1, num_targets2
2711 ! Apply our coefficeints for each combination to the master, Left hand of Eq. 6
2712 alpha(n1,n2) = this % get_target_coeff(i1, s1, n1) * this % get_target_coeff(i2, s2, n2)
2713 !call master_prototype(ii)%add_symbols(temp_prototype,alpha)
2714 end do
2715 end do
2716
2717 num_symbols = temp_prototype % get_size()
2718
2719 do jdo = 1, num_symbols
2720 call temp_prototype % get_coeff_and_integral(jdo, coeff, label)
2721 !We only need to perform the expensive check once per insertion
2722 if (master_prototypes(1,1) % check_same_integral(label, found_idx)) then
2723 !get our coeffs
2724 do n1 = 1, num_targets1
2725 do n2 = 1, num_targets2
2726 call master_prototypes(n1,n2) % modify_coeff(found_idx, alpha(n1,n2) * coeff)
2727 end do
2728 end do
2729 else
2730 !get our coeffs
2731 do n1 = 1, num_targets1
2732 do n2 = 1, num_targets2
2733 !Insert without a check
2734 call master_prototypes(n1,n2) % insert_symbol(label, alpha(n1,n2) * coeff, .false.)
2735 end do
2736 end do
2737 end if
2738 end do
2739
2740 end subroutine fast_contract_class_567
2741
2742end module contracted_hamiltonian_module
MPI SCATCI Constants module.
integer, parameter nidx
Number of long integers used to store up to 8 (packed) integral indices.
integer, parameter normal_phase
No phase correction.
integer, parameter mat_sparse
Matrix is sparse.
integer, parameter no_pure_target_integral_diff
Same as diff but dont evalute integrals belonging to the target state.