-
Notifications
You must be signed in to change notification settings - Fork 1
/
rpa_util.F
1296 lines (1079 loc) · 56.8 KB
/
rpa_util.F
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!--------------------------------------------------------------------------------------------------!
! CP2K: A general program to perform molecular dynamics simulations !
! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
! !
! SPDX-License-Identifier: GPL-2.0-or-later !
!--------------------------------------------------------------------------------------------------!
! **************************************************************************************************
!> \brief Utility functions for RPA calculations
!> \par History
!> 06.2019 Moved from rpa_ri_gpw.F [Frederick Stein]
! **************************************************************************************************
MODULE rpa_util
USE cell_types, ONLY: cell_type,&
get_cell
USE cp_blacs_env, ONLY: cp_blacs_env_create,&
cp_blacs_env_release,&
cp_blacs_env_type
USE cp_cfm_types, ONLY: cp_cfm_create,&
cp_cfm_release,&
cp_cfm_set_all,&
cp_cfm_type
USE cp_dbcsr_api, ONLY: &
dbcsr_create, dbcsr_deallocate_matrix, dbcsr_filter, dbcsr_get_info, dbcsr_multiply, &
dbcsr_p_type, dbcsr_release, dbcsr_set, dbcsr_type
USE cp_dbcsr_operations, ONLY: copy_dbcsr_to_fm,&
copy_fm_to_dbcsr,&
dbcsr_allocate_matrix_set,&
dbcsr_deallocate_matrix_set
USE cp_fm_basic_linalg, ONLY: cp_fm_syrk,&
cp_fm_transpose
USE cp_fm_cholesky, ONLY: cp_fm_cholesky_decompose
USE cp_fm_struct, ONLY: cp_fm_struct_create,&
cp_fm_struct_release,&
cp_fm_struct_type
USE cp_fm_types, ONLY: &
cp_fm_copy_general, cp_fm_create, cp_fm_get_info, cp_fm_release, cp_fm_set_all, &
cp_fm_set_element, cp_fm_to_fm, cp_fm_to_fm_submat_general, cp_fm_type
USE dbt_api, ONLY: dbt_destroy,&
dbt_type
USE dgemm_counter_types, ONLY: dgemm_counter_start,&
dgemm_counter_stop,&
dgemm_counter_type
USE hfx_types, ONLY: block_ind_type,&
dealloc_containers,&
hfx_compression_type
USE input_constants, ONLY: wfc_mm_style_gemm,&
wfc_mm_style_syrk
USE kinds, ONLY: dp
USE kpoint_types, ONLY: get_kpoint_info,&
kpoint_release,&
kpoint_type
USE mathconstants, ONLY: z_zero
USE message_passing, ONLY: mp_para_env_release,&
mp_para_env_type
USE mp2_laplace, ONLY: calc_fm_mat_S_laplace
USE parallel_gemm_api, ONLY: parallel_gemm
USE qs_environment_types, ONLY: get_qs_env,&
qs_environment_type
USE rpa_gw_kpoints_util, ONLY: compute_wkp_W
#include "./base/base_uses.f90"
IMPLICIT NONE
PRIVATE
CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'rpa_util'
PUBLIC :: compute_Erpa_by_freq_int, alloc_im_time, calc_mat_Q, Q_trace_and_add_unit_matrix, &
dealloc_im_time, contract_P_omega_with_mat_L, calc_fm_mat_S_rpa, remove_scaling_factor_rpa
CONTAINS
! **************************************************************************************************
!> \brief ...
!> \param qs_env ...
!> \param para_env ...
!> \param dimen_RI ...
!> \param dimen_RI_red ...
!> \param num_integ_points ...
!> \param nspins ...
!> \param fm_mat_Q ...
!> \param fm_mo_coeff_occ ...
!> \param fm_mo_coeff_virt ...
!> \param fm_matrix_Minv_L_kpoints ...
!> \param fm_matrix_L_kpoints ...
!> \param mat_P_global ...
!> \param t_3c_O ...
!> \param matrix_s ...
!> \param kpoints ...
!> \param eps_filter_im_time ...
!> \param cut_memory ...
!> \param nkp ...
!> \param num_cells_dm ...
!> \param num_3c_repl ...
!> \param size_P ...
!> \param ikp_local ...
!> \param index_to_cell_3c ...
!> \param cell_to_index_3c ...
!> \param col_blk_size ...
!> \param do_ic_model ...
!> \param do_kpoints_cubic_RPA ...
!> \param do_kpoints_from_Gamma ...
!> \param do_ri_Sigma_x ...
!> \param my_open_shell ...
!> \param has_mat_P_blocks ...
!> \param wkp_W ...
!> \param cfm_mat_Q ...
!> \param fm_mat_Minv_L_kpoints ...
!> \param fm_mat_L_kpoints ...
!> \param fm_mat_RI_global_work ...
!> \param fm_mat_work ...
!> \param fm_mo_coeff_occ_scaled ...
!> \param fm_mo_coeff_virt_scaled ...
!> \param mat_dm ...
!> \param mat_L ...
!> \param mat_M_P_munu_occ ...
!> \param mat_M_P_munu_virt ...
!> \param mat_MinvVMinv ...
!> \param mat_P_omega ...
!> \param mat_P_omega_kp ...
!> \param mat_work ...
!> \param mo_coeff ...
!> \param fm_scaled_dm_occ_tau ...
!> \param fm_scaled_dm_virt_tau ...
!> \param homo ...
!> \param nmo ...
! **************************************************************************************************
SUBROUTINE alloc_im_time(qs_env, para_env, dimen_RI, dimen_RI_red, num_integ_points, nspins, &
fm_mat_Q, fm_mo_coeff_occ, fm_mo_coeff_virt, &
fm_matrix_Minv_L_kpoints, fm_matrix_L_kpoints, mat_P_global, &
t_3c_O, matrix_s, kpoints, eps_filter_im_time, &
cut_memory, nkp, num_cells_dm, num_3c_repl, &
size_P, ikp_local, &
index_to_cell_3c, &
cell_to_index_3c, &
col_blk_size, &
do_ic_model, do_kpoints_cubic_RPA, &
do_kpoints_from_Gamma, do_ri_Sigma_x, my_open_shell, &
has_mat_P_blocks, wkp_W, &
cfm_mat_Q, fm_mat_Minv_L_kpoints, fm_mat_L_kpoints, &
fm_mat_RI_global_work, fm_mat_work, fm_mo_coeff_occ_scaled, &
fm_mo_coeff_virt_scaled, mat_dm, mat_L, mat_M_P_munu_occ, mat_M_P_munu_virt, &
mat_MinvVMinv, mat_P_omega, mat_P_omega_kp, &
mat_work, mo_coeff, fm_scaled_dm_occ_tau, fm_scaled_dm_virt_tau, homo, nmo)
TYPE(qs_environment_type), POINTER :: qs_env
TYPE(mp_para_env_type), POINTER :: para_env
INTEGER, INTENT(IN) :: dimen_RI, dimen_RI_red, &
num_integ_points, nspins
TYPE(cp_fm_type), INTENT(IN) :: fm_mat_Q
TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:) :: fm_mo_coeff_occ, fm_mo_coeff_virt
TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:, :) :: fm_matrix_Minv_L_kpoints, &
fm_matrix_L_kpoints
TYPE(dbcsr_p_type), INTENT(IN) :: mat_P_global
TYPE(dbt_type), ALLOCATABLE, DIMENSION(:, :), &
INTENT(INOUT) :: t_3c_O
TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_s
TYPE(kpoint_type), POINTER :: kpoints
REAL(KIND=dp), INTENT(IN) :: eps_filter_im_time
INTEGER, INTENT(IN) :: cut_memory
INTEGER, INTENT(OUT) :: nkp, num_cells_dm, num_3c_repl, size_P, &
ikp_local
INTEGER, ALLOCATABLE, DIMENSION(:, :), INTENT(OUT) :: index_to_cell_3c
INTEGER, ALLOCATABLE, DIMENSION(:, :, :), &
INTENT(OUT) :: cell_to_index_3c
INTEGER, DIMENSION(:), POINTER :: col_blk_size
LOGICAL, INTENT(IN) :: do_ic_model, do_kpoints_cubic_RPA, &
do_kpoints_from_Gamma, do_ri_Sigma_x, &
my_open_shell
LOGICAL, ALLOCATABLE, DIMENSION(:, :, :, :, :), &
INTENT(OUT) :: has_mat_P_blocks
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:), &
INTENT(OUT) :: wkp_W
TYPE(cp_cfm_type), INTENT(OUT) :: cfm_mat_Q
TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:, :) :: fm_mat_Minv_L_kpoints, fm_mat_L_kpoints
TYPE(cp_fm_type), INTENT(OUT) :: fm_mat_RI_global_work, fm_mat_work, &
fm_mo_coeff_occ_scaled, &
fm_mo_coeff_virt_scaled
TYPE(dbcsr_p_type), INTENT(OUT) :: mat_dm, mat_L, mat_M_P_munu_occ, &
mat_M_P_munu_virt, mat_MinvVMinv
TYPE(dbcsr_p_type), ALLOCATABLE, &
DIMENSION(:, :, :) :: mat_P_omega
TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: mat_P_omega_kp
TYPE(dbcsr_type), POINTER :: mat_work
TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: mo_coeff
TYPE(cp_fm_type), INTENT(OUT) :: fm_scaled_dm_occ_tau, &
fm_scaled_dm_virt_tau
INTEGER, DIMENSION(:), INTENT(IN) :: homo
INTEGER, INTENT(IN) :: nmo
CHARACTER(LEN=*), PARAMETER :: routineN = 'alloc_im_time'
INTEGER :: cell_grid_dm(3), first_ikp_local, &
handle, i_dim, i_kp, ispin, jquad, &
nspins_P_omega, periodic(3)
INTEGER, DIMENSION(:), POINTER :: row_blk_size
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: wkp_V
TYPE(cell_type), POINTER :: cell
TYPE(cp_fm_struct_type), POINTER :: fm_struct, fm_struct_sub_kp
CALL timeset(routineN, handle)
ALLOCATE (fm_mo_coeff_occ(nspins), fm_mo_coeff_virt(nspins))
CALL cp_fm_create(fm_scaled_dm_occ_tau, mo_coeff(1)%matrix_struct)
CALL cp_fm_set_all(fm_scaled_dm_occ_tau, 0.0_dp)
CALL cp_fm_create(fm_scaled_dm_virt_tau, mo_coeff(1)%matrix_struct)
CALL cp_fm_set_all(fm_scaled_dm_virt_tau, 0.0_dp)
DO ispin = 1, SIZE(mo_coeff)
CALL create_occ_virt_mo_coeffs(fm_mo_coeff_occ(ispin), fm_mo_coeff_virt(ispin), mo_coeff(ispin), &
nmo, homo(ispin))
END DO
num_3c_repl = SIZE(t_3c_O, 2)
IF (do_kpoints_cubic_RPA) THEN
! we always use an odd number of image cells
! CAUTION: also at another point, cell_grid_dm is defined, these definitions have to be identical
DO i_dim = 1, 3
cell_grid_dm(i_dim) = (kpoints%nkp_grid(i_dim)/2)*2 - 1
END DO
num_cells_dm = cell_grid_dm(1)*cell_grid_dm(2)*cell_grid_dm(3)
ALLOCATE (index_to_cell_3c(3, SIZE(kpoints%index_to_cell, 2)))
CPASSERT(SIZE(kpoints%index_to_cell, 1) == 3)
index_to_cell_3c(:, :) = kpoints%index_to_cell(:, :)
ALLOCATE (cell_to_index_3c(LBOUND(kpoints%cell_to_index, 1):UBOUND(kpoints%cell_to_index, 1), &
LBOUND(kpoints%cell_to_index, 2):UBOUND(kpoints%cell_to_index, 2), &
LBOUND(kpoints%cell_to_index, 3):UBOUND(kpoints%cell_to_index, 3)))
cell_to_index_3c(:, :, :) = kpoints%cell_to_index(:, :, :)
ELSE
ALLOCATE (index_to_cell_3c(3, 1))
index_to_cell_3c(:, 1) = 0
ALLOCATE (cell_to_index_3c(0:0, 0:0, 0:0))
cell_to_index_3c(0, 0, 0) = 1
num_cells_dm = 1
END IF
IF (do_kpoints_cubic_RPA .OR. do_kpoints_from_Gamma) THEN
CALL get_sub_para_kp(fm_struct_sub_kp, para_env, dimen_RI, ikp_local, first_ikp_local)
CALL cp_cfm_create(cfm_mat_Q, fm_struct_sub_kp)
CALL cp_cfm_set_all(cfm_mat_Q, z_zero)
ELSE
first_ikp_local = 1
END IF
! if we do kpoints, mat_P has a kpoint and mat_P_omega has the inted
! mat_P(tau, kpoint)
IF (do_kpoints_cubic_RPA .OR. do_kpoints_from_Gamma) THEN
NULLIFY (cell)
CALL get_qs_env(qs_env, cell=cell)
CALL get_cell(cell=cell, periodic=periodic)
CALL get_kpoint_info(kpoints, nkp=nkp)
! compute k-point weights such that functions 1/k^2, 1/k and const function are
! integrated correctly
CALL compute_wkp_W(qs_env, wkp_W, wkp_V, kpoints, cell%h_inv, periodic)
DEALLOCATE (wkp_V)
ELSE
nkp = 1
END IF
IF (do_kpoints_cubic_RPA) THEN
size_P = MAX(num_cells_dm/2 + 1, nkp)
ELSE IF (do_kpoints_from_Gamma) THEN
size_P = MAX(3**(periodic(1) + periodic(2) + periodic(3)), nkp)
ELSE
size_P = 1
END IF
nspins_P_omega = 1
IF (my_open_shell) nspins_P_omega = 2
ALLOCATE (mat_P_omega(num_integ_points, size_P, nspins_P_omega))
DO ispin = 1, nspins_P_omega
DO i_kp = 1, size_P
DO jquad = 1, num_integ_points
NULLIFY (mat_P_omega(jquad, i_kp, ispin)%matrix)
ALLOCATE (mat_P_omega(jquad, i_kp, ispin)%matrix)
CALL dbcsr_create(matrix=mat_P_omega(jquad, i_kp, ispin)%matrix, &
template=mat_P_global%matrix)
CALL dbcsr_set(mat_P_omega(jquad, i_kp, ispin)%matrix, 0.0_dp)
END DO
END DO
END DO
IF (do_kpoints_cubic_RPA .OR. do_kpoints_from_Gamma) THEN
CALL alloc_mat_P_omega(mat_P_omega_kp, 2, size_P, mat_P_global%matrix)
END IF
CALL cp_fm_create(fm_mo_coeff_occ_scaled, fm_mo_coeff_occ(1)%matrix_struct)
CALL cp_fm_to_fm(fm_mo_coeff_occ(1), fm_mo_coeff_occ_scaled)
CALL cp_fm_set_all(matrix=fm_mo_coeff_occ_scaled, alpha=0.0_dp)
CALL cp_fm_create(fm_mo_coeff_virt_scaled, fm_mo_coeff_virt(1)%matrix_struct)
CALL cp_fm_to_fm(fm_mo_coeff_virt(1), fm_mo_coeff_virt_scaled)
CALL cp_fm_set_all(matrix=fm_mo_coeff_virt_scaled, alpha=0.0_dp)
IF (do_kpoints_cubic_RPA .OR. do_kpoints_from_Gamma) THEN
CALL cp_fm_create(fm_mat_RI_global_work, fm_matrix_Minv_L_kpoints(1, 1)%matrix_struct)
CALL cp_fm_to_fm(fm_matrix_Minv_L_kpoints(1, 1), fm_mat_RI_global_work)
CALL cp_fm_set_all(fm_mat_RI_global_work, 0.0_dp)
END IF
ALLOCATE (has_mat_P_blocks(num_cells_dm/2 + 1, cut_memory, cut_memory, num_3c_repl, num_3c_repl))
has_mat_P_blocks = .TRUE.
IF (do_kpoints_cubic_RPA .OR. do_kpoints_from_Gamma) THEN
CALL reorder_mat_L(fm_mat_Minv_L_kpoints, fm_matrix_Minv_L_kpoints, fm_mat_Q%matrix_struct, para_env, mat_L, &
mat_P_global%matrix, dimen_RI, dimen_RI_red, first_ikp_local, ikp_local, fm_struct_sub_kp, &
allocate_mat_L=.FALSE.)
CALL reorder_mat_L(fm_mat_L_kpoints, fm_matrix_L_kpoints, fm_mat_Q%matrix_struct, para_env, mat_L, &
mat_P_global%matrix, dimen_RI, dimen_RI_red, first_ikp_local, ikp_local, fm_struct_sub_kp)
CALL cp_fm_struct_release(fm_struct_sub_kp)
ELSE
CALL reorder_mat_L(fm_mat_Minv_L_kpoints, fm_matrix_Minv_L_kpoints, fm_mat_Q%matrix_struct, para_env, mat_L, &
mat_P_global%matrix, dimen_RI, dimen_RI_red, first_ikp_local)
END IF
! Create Scalapack working matrix for the contraction with the metric
IF (dimen_RI == dimen_RI_red) THEN
CALL cp_fm_create(fm_mat_work, fm_mat_Q%matrix_struct)
CALL cp_fm_to_fm(fm_mat_Q, fm_mat_work)
CALL cp_fm_set_all(matrix=fm_mat_work, alpha=0.0_dp)
ELSE
NULLIFY (fm_struct)
CALL cp_fm_struct_create(fm_struct, template_fmstruct=fm_mat_Q%matrix_struct, &
ncol_global=dimen_RI_red, nrow_global=dimen_RI)
CALL cp_fm_create(fm_mat_work, fm_struct)
CALL cp_fm_set_all(matrix=fm_mat_work, alpha=0.0_dp)
CALL cp_fm_struct_release(fm_struct)
END IF
! Then its DBCSR counter part
IF (.NOT. (do_kpoints_cubic_RPA .OR. do_kpoints_from_Gamma)) THEN
CALL dbcsr_get_info(mat_L%matrix, col_blk_size=col_blk_size, row_blk_size=row_blk_size)
! Create mat_work having the shape of the transposed of mat_L (compare with contract_P_omega_with_mat_L)
NULLIFY (mat_work)
ALLOCATE (mat_work)
CALL dbcsr_create(mat_work, template=mat_L%matrix, row_blk_size=col_blk_size, col_blk_size=row_blk_size)
END IF
IF (do_ri_Sigma_x .OR. do_ic_model) THEN
NULLIFY (mat_MinvVMinv%matrix)
ALLOCATE (mat_MinvVMinv%matrix)
CALL dbcsr_create(mat_MinvVMinv%matrix, template=mat_P_global%matrix)
CALL dbcsr_set(mat_MinvVMinv%matrix, 0.0_dp)
! for kpoints we compute SinvVSinv later with kpoints
IF (.NOT. do_kpoints_from_Gamma) THEN
! get the Coulomb matrix for Sigma_x = G*V
CALL dbcsr_multiply("T", "N", 1.0_dp, mat_L%matrix, mat_L%matrix, &
0.0_dp, mat_MinvVMinv%matrix, filter_eps=eps_filter_im_time)
END IF
END IF
IF (do_ri_Sigma_x) THEN
NULLIFY (mat_dm%matrix)
ALLOCATE (mat_dm%matrix)
CALL dbcsr_create(mat_dm%matrix, template=matrix_s(1)%matrix)
END IF
CALL timestop(handle)
END SUBROUTINE alloc_im_time
! **************************************************************************************************
!> \brief ...
!> \param fm_mo_coeff_occ ...
!> \param fm_mo_coeff_virt ...
!> \param mo_coeff ...
!> \param nmo ...
!> \param homo ...
! **************************************************************************************************
SUBROUTINE create_occ_virt_mo_coeffs(fm_mo_coeff_occ, fm_mo_coeff_virt, mo_coeff, &
nmo, homo)
TYPE(cp_fm_type), INTENT(OUT) :: fm_mo_coeff_occ, fm_mo_coeff_virt
TYPE(cp_fm_type), INTENT(IN) :: mo_coeff
INTEGER, INTENT(IN) :: nmo, homo
CHARACTER(LEN=*), PARAMETER :: routineN = 'create_occ_virt_mo_coeffs'
INTEGER :: handle, icol_global, irow_global
CALL timeset(routineN, handle)
CALL cp_fm_create(fm_mo_coeff_occ, mo_coeff%matrix_struct)
CALL cp_fm_set_all(fm_mo_coeff_occ, 0.0_dp)
CALL cp_fm_to_fm(mo_coeff, fm_mo_coeff_occ)
! set all virtual MO coeffs to zero
DO irow_global = 1, nmo
DO icol_global = homo + 1, nmo
CALL cp_fm_set_element(fm_mo_coeff_occ, irow_global, icol_global, 0.0_dp)
END DO
END DO
CALL cp_fm_create(fm_mo_coeff_virt, mo_coeff%matrix_struct)
CALL cp_fm_set_all(fm_mo_coeff_virt, 0.0_dp)
CALL cp_fm_to_fm(mo_coeff, fm_mo_coeff_virt)
! set all occupied MO coeffs to zero
DO irow_global = 1, nmo
DO icol_global = 1, homo
CALL cp_fm_set_element(fm_mo_coeff_virt, irow_global, icol_global, 0.0_dp)
END DO
END DO
CALL timestop(handle)
END SUBROUTINE create_occ_virt_mo_coeffs
! **************************************************************************************************
!> \brief ...
!> \param mat_P_omega ...
!> \param num_integ_points ...
!> \param size_P ...
!> \param template ...
! **************************************************************************************************
SUBROUTINE alloc_mat_P_omega(mat_P_omega, num_integ_points, size_P, template)
TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: mat_P_omega
INTEGER, INTENT(IN) :: num_integ_points, size_P
TYPE(dbcsr_type), POINTER :: template
CHARACTER(LEN=*), PARAMETER :: routineN = 'alloc_mat_P_omega'
INTEGER :: handle, i_kp, jquad
CALL timeset(routineN, handle)
NULLIFY (mat_P_omega)
CALL dbcsr_allocate_matrix_set(mat_P_omega, num_integ_points, size_P)
DO i_kp = 1, size_P
DO jquad = 1, num_integ_points
ALLOCATE (mat_P_omega(jquad, i_kp)%matrix)
CALL dbcsr_create(matrix=mat_P_omega(jquad, i_kp)%matrix, &
template=template)
CALL dbcsr_set(mat_P_omega(jquad, i_kp)%matrix, 0.0_dp)
END DO
END DO
CALL timestop(handle)
END SUBROUTINE alloc_mat_P_omega
! **************************************************************************************************
!> \brief ...
!> \param fm_mat_L ...
!> \param fm_matrix_Minv_L_kpoints ...
!> \param fm_struct_template ...
!> \param para_env ...
!> \param mat_L ...
!> \param mat_template ...
!> \param dimen_RI ...
!> \param dimen_RI_red ...
!> \param first_ikp_local ...
!> \param ikp_local ...
!> \param fm_struct_sub_kp ...
!> \param allocate_mat_L ...
! **************************************************************************************************
SUBROUTINE reorder_mat_L(fm_mat_L, fm_matrix_Minv_L_kpoints, fm_struct_template, para_env, mat_L, mat_template, &
dimen_RI, dimen_RI_red, first_ikp_local, ikp_local, fm_struct_sub_kp, allocate_mat_L)
TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:, :) :: fm_mat_L, fm_matrix_Minv_L_kpoints
TYPE(cp_fm_struct_type), POINTER :: fm_struct_template
TYPE(mp_para_env_type), POINTER :: para_env
TYPE(dbcsr_p_type), INTENT(OUT) :: mat_L
TYPE(dbcsr_type), INTENT(IN) :: mat_template
INTEGER, INTENT(IN) :: dimen_RI, dimen_RI_red, first_ikp_local
INTEGER, OPTIONAL :: ikp_local
TYPE(cp_fm_struct_type), OPTIONAL, POINTER :: fm_struct_sub_kp
LOGICAL, INTENT(IN), OPTIONAL :: allocate_mat_L
CHARACTER(LEN=*), PARAMETER :: routineN = 'reorder_mat_L'
INTEGER :: handle, ikp, j_size, nblk
INTEGER, DIMENSION(:), POINTER :: col_blk_size, row_blk_size
LOGICAL :: do_kpoints, my_allocate_mat_L
TYPE(cp_blacs_env_type), POINTER :: blacs_env
TYPE(cp_fm_struct_type), POINTER :: fm_struct
TYPE(cp_fm_type) :: fm_mat_L_transposed, fmdummy
CALL timeset(routineN, handle)
do_kpoints = .FALSE.
IF (PRESENT(ikp_local) .AND. PRESENT(fm_struct_sub_kp)) THEN
do_kpoints = .TRUE.
END IF
! Get the fm_struct for fm_mat_L
NULLIFY (fm_struct)
IF (dimen_RI == dimen_RI_red) THEN
fm_struct => fm_struct_template
ELSE
! The template is assumed to be square such that we need a new fm_struct if dimensions are not equal
CALL cp_fm_struct_create(fm_struct, nrow_global=dimen_RI_red, ncol_global=dimen_RI, template_fmstruct=fm_struct_template)
END IF
! Start to allocate the new full matrix
ALLOCATE (fm_mat_L(SIZE(fm_matrix_Minv_L_kpoints, 1), SIZE(fm_matrix_Minv_L_kpoints, 2)))
DO ikp = 1, SIZE(fm_matrix_Minv_L_kpoints, 1)
DO j_size = 1, SIZE(fm_matrix_Minv_L_kpoints, 2)
IF (do_kpoints) THEN
IF (ikp == first_ikp_local .OR. ikp_local == -1) THEN
CALL cp_fm_create(fm_mat_L(ikp, j_size), fm_struct_sub_kp)
CALL cp_fm_set_all(fm_mat_L(ikp, j_size), 0.0_dp)
END IF
ELSE
CALL cp_fm_create(fm_mat_L(ikp, j_size), fm_struct)
CALL cp_fm_set_all(fm_mat_L(ikp, j_size), 0.0_dp)
END IF
END DO
END DO
! For the transposed matric we need a different fm_struct
IF (dimen_RI == dimen_RI_red) THEN
fm_struct => fm_mat_L(first_ikp_local, 1)%matrix_struct
ELSE
CALL cp_fm_struct_release(fm_struct)
! Create a fm_struct with transposed sizes
NULLIFY (fm_struct)
CALL cp_fm_struct_create(fm_struct, nrow_global=dimen_RI, ncol_global=dimen_RI_red, &
template_fmstruct=fm_mat_L(first_ikp_local, 1)%matrix_struct) !, force_block=.TRUE.)
END IF
! Allocate buffer matrix
CALL cp_fm_create(fm_mat_L_transposed, fm_struct)
CALL cp_fm_set_all(matrix=fm_mat_L_transposed, alpha=0.0_dp)
IF (dimen_RI /= dimen_RI_red) CALL cp_fm_struct_release(fm_struct)
CALL cp_fm_get_info(fm_mat_L_transposed, context=blacs_env)
! For k-points copy matrices of your group
! Without kpoints, transpose matrix
! without kpoints, the size of fm_mat_L is 1x1. with kpoints, the size is N_kpoints x 2 (2 for real/complex)
DO ikp = 1, SIZE(fm_matrix_Minv_L_kpoints, 1)
DO j_size = 1, SIZE(fm_matrix_Minv_L_kpoints, 2)
IF (do_kpoints) THEN
IF (ikp_local == ikp .OR. ikp_local == -1) THEN
CALL cp_fm_copy_general(fm_matrix_Minv_L_kpoints(ikp, j_size), fm_mat_L_transposed, para_env)
CALL cp_fm_to_fm(fm_mat_L_transposed, fm_mat_L(ikp, j_size))
ELSE
CALL cp_fm_copy_general(fm_matrix_Minv_L_kpoints(ikp, j_size), fmdummy, para_env)
END IF
ELSE
CALL cp_fm_copy_general(fm_matrix_Minv_L_kpoints(ikp, j_size), fm_mat_L_transposed, blacs_env%para_env)
CALL cp_fm_transpose(fm_mat_L_transposed, fm_mat_L(ikp, j_size))
END IF
END DO
END DO
! Release old matrix
CALL cp_fm_release(fm_matrix_Minv_L_kpoints)
! Release buffer
CALL cp_fm_release(fm_mat_L_transposed)
my_allocate_mat_L = .TRUE.
IF (PRESENT(allocate_mat_L)) my_allocate_mat_L = allocate_mat_L
IF (my_allocate_mat_L) THEN
! Create sparse variant of L
NULLIFY (mat_L%matrix)
ALLOCATE (mat_L%matrix)
IF (dimen_RI == dimen_RI_red) THEN
CALL dbcsr_create(mat_L%matrix, template=mat_template)
ELSE
CALL dbcsr_get_info(mat_template, nblkrows_total=nblk, col_blk_size=col_blk_size)
CALL calculate_equal_blk_size(row_blk_size, dimen_RI_red, nblk)
CALL dbcsr_create(mat_L%matrix, template=mat_template, row_blk_size=row_blk_size, col_blk_size=col_blk_size)
DEALLOCATE (row_blk_size)
END IF
IF (.NOT. (do_kpoints)) THEN
CALL copy_fm_to_dbcsr(fm_mat_L(1, 1), mat_L%matrix)
END IF
END IF
CALL timestop(handle)
END SUBROUTINE reorder_mat_L
! **************************************************************************************************
!> \brief ...
!> \param blk_size_new ...
!> \param dimen_RI_red ...
!> \param nblk ...
! **************************************************************************************************
SUBROUTINE calculate_equal_blk_size(blk_size_new, dimen_RI_red, nblk)
INTEGER, DIMENSION(:), POINTER :: blk_size_new
INTEGER, INTENT(IN) :: dimen_RI_red, nblk
INTEGER :: col_per_blk, remainder
NULLIFY (blk_size_new)
ALLOCATE (blk_size_new(nblk))
remainder = MOD(dimen_RI_red, nblk)
col_per_blk = dimen_RI_red/nblk
! Determine a new distribution for the columns (corresponding to the number of columns)
IF (remainder > 0) blk_size_new(1:remainder) = col_per_blk + 1
blk_size_new(remainder + 1:nblk) = col_per_blk
END SUBROUTINE calculate_equal_blk_size
! **************************************************************************************************
!> \brief ...
!> \param fm_mat_S ...
!> \param do_ri_sos_laplace_mp2 ...
!> \param first_cycle ...
!> \param virtual ...
!> \param Eigenval ...
!> \param homo ...
!> \param omega ...
!> \param omega_old ...
!> \param jquad ...
!> \param mm_style ...
!> \param dimen_RI ...
!> \param dimen_ia ...
!> \param alpha ...
!> \param fm_mat_Q ...
!> \param fm_mat_Q_gemm ...
!> \param do_bse ...
!> \param fm_mat_Q_static_bse_gemm ...
!> \param dgemm_counter ...
!> \param num_integ_points ...
!> \param count_ev_sc_GW ...
! **************************************************************************************************
SUBROUTINE calc_mat_Q(fm_mat_S, do_ri_sos_laplace_mp2, first_cycle, virtual, &
Eigenval, homo, omega, omega_old, jquad, mm_style, dimen_RI, dimen_ia, alpha, fm_mat_Q, fm_mat_Q_gemm, &
do_bse, fm_mat_Q_static_bse_gemm, dgemm_counter, &
num_integ_points, count_ev_sc_GW)
TYPE(cp_fm_type), INTENT(IN) :: fm_mat_S
LOGICAL, INTENT(IN) :: do_ri_sos_laplace_mp2, first_cycle
INTEGER, INTENT(IN) :: virtual
REAL(KIND=dp), DIMENSION(:), INTENT(IN) :: Eigenval
INTEGER, INTENT(IN) :: homo
REAL(KIND=dp), INTENT(IN) :: omega, omega_old
INTEGER, INTENT(IN) :: jquad, mm_style, dimen_RI, dimen_ia
REAL(KIND=dp), INTENT(IN) :: alpha
TYPE(cp_fm_type), INTENT(IN) :: fm_mat_Q, fm_mat_Q_gemm
LOGICAL, INTENT(IN) :: do_bse
TYPE(cp_fm_type), INTENT(IN) :: fm_mat_Q_static_bse_gemm
TYPE(dgemm_counter_type), INTENT(INOUT) :: dgemm_counter
INTEGER, INTENT(IN) :: num_integ_points, count_ev_sc_GW
CHARACTER(LEN=*), PARAMETER :: routineN = 'calc_mat_Q'
INTEGER :: handle
CALL timeset(routineN, handle)
IF (do_ri_sos_laplace_mp2) THEN
! the first index of tau_tj starts with 0 (see mp2_weights)
CALL calc_fm_mat_S_laplace(fm_mat_S, homo, virtual, Eigenval, omega - omega_old)
ELSE
CALL calc_fm_mat_S_rpa(fm_mat_S, first_cycle, virtual, Eigenval, &
homo, omega, omega_old)
END IF
CALL contract_S_to_Q(mm_style, dimen_RI, dimen_ia, alpha, fm_mat_S, fm_mat_Q_gemm, &
fm_mat_Q, dgemm_counter)
! fm_mat_Q_static_bse_gemm does not enter W_ijab (A matrix in TDA), but only full ABBA
! (since only B_ij_bar enters W_ijab)
! Changing jquad, since omega=0 is at last idx
! We enforce W0 for BSE in case of evGW
IF (do_bse .AND. jquad == num_integ_points .AND. count_ev_sc_GW == 1) THEN
CALL cp_fm_to_fm(fm_mat_Q_gemm, fm_mat_Q_static_bse_gemm)
END IF
CALL timestop(handle)
END SUBROUTINE calc_mat_Q
! **************************************************************************************************
!> \brief ...
!> \param fm_mat_S ...
!> \param virtual ...
!> \param Eigenval_last ...
!> \param homo ...
!> \param omega_old ...
! **************************************************************************************************
SUBROUTINE remove_scaling_factor_rpa(fm_mat_S, virtual, Eigenval_last, homo, omega_old)
TYPE(cp_fm_type), INTENT(IN) :: fm_mat_S
INTEGER, INTENT(IN) :: virtual
REAL(KIND=dp), DIMENSION(:), INTENT(IN) :: Eigenval_last
INTEGER, INTENT(IN) :: homo
REAL(KIND=dp), INTENT(IN) :: omega_old
CHARACTER(LEN=*), PARAMETER :: routineN = 'remove_scaling_factor_rpa'
INTEGER :: avirt, handle, i_global, iiB, iocc, &
ncol_local
INTEGER, DIMENSION(:), POINTER :: col_indices
REAL(KIND=dp) :: eigen_diff
CALL timeset(routineN, handle)
! get info of fm_mat_S
CALL cp_fm_get_info(matrix=fm_mat_S, &
ncol_local=ncol_local, &
col_indices=col_indices)
!$OMP PARALLEL DO DEFAULT(NONE) PRIVATE(iiB,iocc,avirt,eigen_diff,i_global) &
!$OMP SHARED(ncol_local,col_indices,Eigenval_last,fm_mat_S,virtual,homo,omega_old)
DO iiB = 1, ncol_local
i_global = col_indices(iiB)
iocc = MAX(1, i_global - 1)/virtual + 1
avirt = i_global - (iocc - 1)*virtual
eigen_diff = Eigenval_last(avirt + homo) - Eigenval_last(iocc)
fm_mat_S%local_data(:, iiB) = fm_mat_S%local_data(:, iiB)/ &
SQRT(eigen_diff/(eigen_diff**2 + omega_old**2))
END DO
CALL timestop(handle)
END SUBROUTINE remove_scaling_factor_rpa
! **************************************************************************************************
!> \brief ...
!> \param fm_mat_S ...
!> \param first_cycle ...
!> \param virtual ...
!> \param Eigenval ...
!> \param homo ...
!> \param omega ...
!> \param omega_old ...
! **************************************************************************************************
SUBROUTINE calc_fm_mat_S_rpa(fm_mat_S, first_cycle, virtual, Eigenval, homo, &
omega, omega_old)
TYPE(cp_fm_type), INTENT(IN) :: fm_mat_S
LOGICAL, INTENT(IN) :: first_cycle
INTEGER, INTENT(IN) :: virtual
REAL(KIND=dp), DIMENSION(:), INTENT(IN) :: Eigenval
INTEGER, INTENT(IN) :: homo
REAL(KIND=dp), INTENT(IN) :: omega, omega_old
CHARACTER(LEN=*), PARAMETER :: routineN = 'calc_fm_mat_S_rpa'
INTEGER :: avirt, handle, i_global, iiB, iocc, &
ncol_local
INTEGER, DIMENSION(:), POINTER :: col_indices
REAL(KIND=dp) :: eigen_diff
CALL timeset(routineN, handle)
! get info of fm_mat_S
CALL cp_fm_get_info(matrix=fm_mat_S, &
ncol_local=ncol_local, &
col_indices=col_indices)
! update G matrix with the new value of omega
IF (first_cycle) THEN
! In this case just update the matrix (symmetric form) with
! SQRT((epsi_a-epsi_i)/((epsi_a-epsi_i)**2+omega**2))
!$OMP PARALLEL DO DEFAULT(NONE) PRIVATE(iiB,iocc,avirt,eigen_diff,i_global) &
!$OMP SHARED(ncol_local,col_indices,Eigenval,fm_mat_S,virtual,homo,omega)
DO iiB = 1, ncol_local
i_global = col_indices(iiB)
iocc = MAX(1, i_global - 1)/virtual + 1
avirt = i_global - (iocc - 1)*virtual
eigen_diff = Eigenval(avirt + homo) - Eigenval(iocc)
fm_mat_S%local_data(:, iiB) = fm_mat_S%local_data(:, iiB)* &
SQRT(eigen_diff/(eigen_diff**2 + omega**2))
END DO
ELSE
! In this case the update has to remove the old omega component thus
! SQRT(((epsi_a-epsi_i)**2+omega_old**2)/((epsi_a-epsi_i)**2+omega**2))
!$OMP PARALLEL DO DEFAULT(NONE) PRIVATE(iiB,iocc,avirt,eigen_diff,i_global) &
!$OMP SHARED(ncol_local,col_indices,Eigenval,fm_mat_S,virtual,homo,omega,omega_old)
DO iiB = 1, ncol_local
i_global = col_indices(iiB)
iocc = MAX(1, i_global - 1)/virtual + 1
avirt = i_global - (iocc - 1)*virtual
eigen_diff = Eigenval(avirt + homo) - Eigenval(iocc)
fm_mat_S%local_data(:, iiB) = fm_mat_S%local_data(:, iiB)* &
SQRT((eigen_diff**2 + omega_old**2)/(eigen_diff**2 + omega**2))
END DO
END IF
CALL timestop(handle)
END SUBROUTINE calc_fm_mat_S_rpa
! **************************************************************************************************
!> \brief ...
!> \param mm_style ...
!> \param dimen_RI ...
!> \param dimen_ia ...
!> \param alpha ...
!> \param fm_mat_S ...
!> \param fm_mat_Q_gemm ...
!> \param fm_mat_Q ...
!> \param dgemm_counter ...
! **************************************************************************************************
SUBROUTINE contract_S_to_Q(mm_style, dimen_RI, dimen_ia, alpha, fm_mat_S, fm_mat_Q_gemm, &
fm_mat_Q, dgemm_counter)
INTEGER, INTENT(IN) :: mm_style, dimen_RI, dimen_ia
REAL(KIND=dp), INTENT(IN) :: alpha
TYPE(cp_fm_type), INTENT(IN) :: fm_mat_S, fm_mat_Q_gemm, fm_mat_Q
TYPE(dgemm_counter_type), INTENT(INOUT) :: dgemm_counter
CHARACTER(LEN=*), PARAMETER :: routineN = 'contract_S_to_Q'
INTEGER :: handle
CALL timeset(routineN, handle)
CALL dgemm_counter_start(dgemm_counter)
SELECT CASE (mm_style)
CASE (wfc_mm_style_gemm)
! waste-fully computes the full symmetrix matrix, but maybe faster than cp_fm_syrk for optimized cp_fm_gemm !!!
CALL parallel_gemm(transa="N", transb="T", m=dimen_RI, n=dimen_RI, k=dimen_ia, alpha=alpha, &
matrix_a=fm_mat_S, matrix_b=fm_mat_S, beta=0.0_dp, &
matrix_c=fm_mat_Q_gemm)
CASE (wfc_mm_style_syrk)
! will only compute the upper half of the matrix, which is fine, since we only use it for cholesky later
CALL cp_fm_syrk(uplo='U', trans='N', k=dimen_ia, alpha=alpha, matrix_a=fm_mat_S, &
ia=1, ja=1, beta=0.0_dp, matrix_c=fm_mat_Q_gemm)
CASE DEFAULT
CPABORT("")
END SELECT
CALL dgemm_counter_stop(dgemm_counter, dimen_RI, dimen_RI, dimen_ia)
! copy/redistribute fm_mat_Q_gemm to fm_mat_Q
CALL cp_fm_set_all(matrix=fm_mat_Q, alpha=0.0_dp)
CALL cp_fm_to_fm_submat_general(fm_mat_Q_gemm, fm_mat_Q, dimen_RI, dimen_RI, 1, 1, 1, 1, &
fm_mat_Q_gemm%matrix_struct%context)
CALL timestop(handle)
END SUBROUTINE contract_S_to_Q
! **************************************************************************************************
!> \brief ...
!> \param dimen_RI ...
!> \param trace_Qomega ...
!> \param fm_mat_Q ...
! **************************************************************************************************
SUBROUTINE Q_trace_and_add_unit_matrix(dimen_RI, trace_Qomega, fm_mat_Q)
INTEGER, INTENT(IN) :: dimen_RI
REAL(KIND=dp), DIMENSION(dimen_RI), INTENT(OUT) :: trace_Qomega
TYPE(cp_fm_type), INTENT(IN) :: fm_mat_Q
CHARACTER(LEN=*), PARAMETER :: routineN = 'Q_trace_and_add_unit_matrix'
INTEGER :: handle, i_global, iiB, j_global, jjB, &
ncol_local, nrow_local
INTEGER, DIMENSION(:), POINTER :: col_indices, row_indices
TYPE(mp_para_env_type), POINTER :: para_env
CALL timeset(routineN, handle)
CALL cp_fm_get_info(matrix=fm_mat_Q, &
nrow_local=nrow_local, &
ncol_local=ncol_local, &
row_indices=row_indices, &
col_indices=col_indices, &
para_env=para_env)
! calculate the trace of Q and add 1 on the diagonal
trace_Qomega = 0.0_dp
!$OMP PARALLEL DO DEFAULT(NONE) PRIVATE(jjB,iiB,i_global,j_global) &
!$OMP SHARED(ncol_local,nrow_local,col_indices,row_indices,trace_Qomega,fm_mat_Q,dimen_RI)
DO jjB = 1, ncol_local
j_global = col_indices(jjB)
DO iiB = 1, nrow_local
i_global = row_indices(iiB)
IF (j_global == i_global .AND. i_global <= dimen_RI) THEN
trace_Qomega(i_global) = fm_mat_Q%local_data(iiB, jjB)
fm_mat_Q%local_data(iiB, jjB) = fm_mat_Q%local_data(iiB, jjB) + 1.0_dp
END IF
END DO
END DO
CALL para_env%sum(trace_Qomega)
CALL timestop(handle)
END SUBROUTINE Q_trace_and_add_unit_matrix
! **************************************************************************************************
!> \brief ...
!> \param dimen_RI ...
!> \param trace_Qomega ...
!> \param fm_mat_Q ...
!> \param para_env_RPA ...
!> \param Erpa ...
!> \param wjquad ...
! **************************************************************************************************
SUBROUTINE compute_Erpa_by_freq_int(dimen_RI, trace_Qomega, fm_mat_Q, para_env_RPA, Erpa, wjquad)
INTEGER, INTENT(IN) :: dimen_RI
REAL(KIND=dp), DIMENSION(dimen_RI), INTENT(IN) :: trace_Qomega
TYPE(cp_fm_type), INTENT(IN) :: fm_mat_Q
TYPE(mp_para_env_type), INTENT(IN) :: para_env_RPA
REAL(KIND=dp), INTENT(INOUT) :: Erpa
REAL(KIND=dp), INTENT(IN) :: wjquad
CHARACTER(LEN=*), PARAMETER :: routineN = 'compute_Erpa_by_freq_int'
INTEGER :: handle, i_global, iiB, info_chol, &
j_global, jjB, ncol_local, nrow_local
INTEGER, DIMENSION(:), POINTER :: col_indices, row_indices
REAL(KIND=dp) :: FComega
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: Q_log
CALL timeset(routineN, handle)
CALL cp_fm_get_info(matrix=fm_mat_Q, &
nrow_local=nrow_local, &
ncol_local=ncol_local, &
row_indices=row_indices, &
col_indices=col_indices)
! calculate Trace(Log(Matrix)) as Log(DET(Matrix)) via cholesky decomposition
CALL cp_fm_cholesky_decompose(matrix=fm_mat_Q, n=dimen_RI, info_out=info_chol)
IF (info_chol .NE. 0) THEN
CALL cp_warn(__LOCATION__, &
"The Cholesky decomposition before inverting the RPA matrix / dielectric "// &
"function failed. "// &
"In case of low-scaling RPA/GW, decreasing EPS_FILTER in the &LOW_SCALING "// &
"section might "// &
"increase the overall accuracy making the matrix positive definite. "// &
"Code will abort.")
END IF
CPASSERT(info_chol == 0)
ALLOCATE (Q_log(dimen_RI))
Q_log = 0.0_dp
!$OMP PARALLEL DO DEFAULT(NONE) PRIVATE(jjB,iiB,i_global,j_global) &
!$OMP SHARED(ncol_local,nrow_local,col_indices,row_indices,Q_log,fm_mat_Q,dimen_RI)
DO jjB = 1, ncol_local
j_global = col_indices(jjB)
DO iiB = 1, nrow_local
i_global = row_indices(iiB)
IF (j_global == i_global .AND. i_global <= dimen_RI) THEN
Q_log(i_global) = 2.0_dp*LOG(fm_mat_Q%local_data(iiB, jjB))
END IF
END DO
END DO
CALL para_env_RPA%sum(Q_log)
! the following frequency integration is Eq. (27) in M. Del Ben et al., JCTC 9, 2654 (2013)
! (https://doi.org/10.1021/ct4002202)
FComega = 0.0_dp
DO iiB = 1, dimen_RI
IF (MODULO(iiB, para_env_RPA%num_pe) /= para_env_RPA%mepos) CYCLE
FComega = FComega + (Q_log(iiB) - trace_Qomega(iiB))/2.0_dp
END DO
Erpa = Erpa + FComega*wjquad
DEALLOCATE (Q_log)
CALL timestop(handle)
END SUBROUTINE compute_Erpa_by_freq_int
! **************************************************************************************************
!> \brief ...
!> \param fm_struct_sub_kp ...
!> \param para_env ...