-
Notifications
You must be signed in to change notification settings - Fork 1
/
mp2_grids.F
1526 lines (1223 loc) · 62.8 KB
/
mp2_grids.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 Routines to calculate frequency and time grids (integration points and weights)
!> for correlation methods
!> \par History
!> 05.2019 Refactored from rpa_ri_gpw [Frederick Stein]
! **************************************************************************************************
MODULE mp2_grids
USE cp_fm_types, ONLY: cp_fm_get_info,&
cp_fm_type
USE input_section_types, ONLY: section_vals_type,&
section_vals_val_set
USE kinds, ONLY: dp
USE kpoint_types, ONLY: get_kpoint_info,&
kpoint_env_type,&
kpoint_type
USE machine, ONLY: m_flush
USE mathconstants, ONLY: pi
USE message_passing, ONLY: mp_para_env_release,&
mp_para_env_type
USE minimax_exp, ONLY: get_exp_minimax_coeff
USE minimax_exp_gw, ONLY: get_exp_minimax_coeff_gw
USE minimax_rpa, ONLY: get_rpa_minimax_coeff,&
get_rpa_minimax_coeff_larger_grid
USE qs_environment_types, ONLY: get_qs_env,&
qs_environment_type
USE qs_mo_types, ONLY: get_mo_set,&
mo_set_type
#include "./base/base_uses.f90"
IMPLICIT NONE
PRIVATE
CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'mp2_grids'
PUBLIC :: get_minimax_grid, get_clenshaw_grid, test_least_square_ft, get_l_sq_wghts_cos_tf_t_to_w, &
get_l_sq_wghts_cos_tf_w_to_t, get_l_sq_wghts_sin_tf_t_to_w
CONTAINS
! **************************************************************************************************
!> \brief ...
!> \param para_env ...
!> \param unit_nr ...
!> \param homo ...
!> \param Eigenval ...
!> \param num_integ_points ...
!> \param do_im_time ...
!> \param do_ri_sos_laplace_mp2 ...
!> \param do_print ...
!> \param tau_tj ...
!> \param tau_wj ...
!> \param qs_env ...
!> \param do_gw_im_time ...
!> \param do_kpoints_cubic_RPA ...
!> \param e_fermi ...
!> \param tj ...
!> \param wj ...
!> \param weights_cos_tf_t_to_w ...
!> \param weights_cos_tf_w_to_t ...
!> \param weights_sin_tf_t_to_w ...
!> \param regularization ...
! **************************************************************************************************
SUBROUTINE get_minimax_grid(para_env, unit_nr, homo, Eigenval, num_integ_points, &
do_im_time, do_ri_sos_laplace_mp2, do_print, tau_tj, tau_wj, qs_env, do_gw_im_time, &
do_kpoints_cubic_RPA, e_fermi, tj, wj, weights_cos_tf_t_to_w, &
weights_cos_tf_w_to_t, weights_sin_tf_t_to_w, regularization)
TYPE(mp_para_env_type), INTENT(IN) :: para_env
INTEGER, INTENT(IN) :: unit_nr
INTEGER, DIMENSION(:), INTENT(IN) :: homo
REAL(KIND=dp), DIMENSION(:, :, :), INTENT(IN) :: Eigenval
INTEGER, INTENT(IN) :: num_integ_points
LOGICAL, INTENT(IN) :: do_im_time, do_ri_sos_laplace_mp2, &
do_print
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:), &
INTENT(OUT) :: tau_tj, tau_wj
TYPE(qs_environment_type), POINTER :: qs_env
LOGICAL, INTENT(IN) :: do_gw_im_time, do_kpoints_cubic_RPA
REAL(KIND=dp), INTENT(OUT) :: e_fermi
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:), &
INTENT(INOUT) :: tj, wj
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :), &
INTENT(OUT) :: weights_cos_tf_t_to_w, &
weights_cos_tf_w_to_t, &
weights_sin_tf_t_to_w
REAL(KIND=dp), INTENT(IN), OPTIONAL :: regularization
CHARACTER(LEN=*), PARAMETER :: routineN = 'get_minimax_grid'
INTEGER, PARAMETER :: num_points_per_magnitude = 200
INTEGER :: handle, ierr, ispin, jquad, nspins
LOGICAL :: my_do_kpoints, my_open_shell
REAL(KIND=dp) :: Emax, Emin, max_error_min, my_E_Range, &
my_regularization, scaling
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: x_tw
TYPE(section_vals_type), POINTER :: input
CALL timeset(routineN, handle)
! Test for spin unrestricted
nspins = SIZE(homo)
my_open_shell = (nspins == 2)
! Test whether all necessary variables are available
my_do_kpoints = .FALSE.
IF (.NOT. do_ri_sos_laplace_mp2) THEN
my_do_kpoints = do_kpoints_cubic_RPA
END IF
my_regularization = 0.0_dp
IF (PRESENT(regularization)) THEN
my_regularization = regularization
END IF
IF (my_do_kpoints) THEN
CALL gap_and_max_eig_diff_kpoints(qs_env, para_env, Emin, Emax, e_fermi)
my_E_Range = Emax/Emin
ELSE
IF (qs_env%mp2_env%E_range <= 1.0_dp .OR. qs_env%mp2_env%E_gap <= 0.0_dp) THEN
Emin = HUGE(dp)
Emax = 0.0_dp
DO ispin = 1, nspins
IF (homo(ispin) > 0) THEN
Emin = MIN(Emin, Eigenval(homo(ispin) + 1, 1, ispin) - Eigenval(homo(ispin), 1, ispin))
Emax = MAX(Emax, MAXVAL(Eigenval(:, :, ispin)) - MINVAL(Eigenval(:, :, ispin)))
END IF
END DO
my_E_Range = Emax/Emin
qs_env%mp2_env%e_range = my_e_range
qs_env%mp2_env%e_gap = Emin
CALL get_qs_env(qs_env, input=input)
CALL section_vals_val_set(input, "DFT%XC%WF_CORRELATION%E_RANGE", r_val=my_e_range)
CALL section_vals_val_set(input, "DFT%XC%WF_CORRELATION%E_GAP", r_val=emin)
ELSE
my_E_range = qs_env%mp2_env%E_range
Emin = qs_env%mp2_env%E_gap
Emax = Emin*my_E_range
END IF
END IF
IF (num_integ_points > 20 .AND. my_E_Range < 100.0_dp) THEN
IF (unit_nr > 0) &
CALL cp_warn(__LOCATION__, &
"You requested a large minimax grid (> 20 points) for a small minimax range R (R < 100). "// &
"That may lead to numerical "// &
"instabilities when computing minimax grid weights. You can prevent small ranges by choosing "// &
"a larger basis set with higher angular momenta or alternatively using all-electron calculations.")
END IF
IF (.NOT. do_ri_sos_laplace_mp2) THEN
ALLOCATE (x_tw(2*num_integ_points))
x_tw = 0.0_dp
ierr = 0
IF (num_integ_points .LE. 20) THEN
CALL get_rpa_minimax_coeff(num_integ_points, my_E_Range, x_tw, ierr)
ELSE
CALL get_rpa_minimax_coeff_larger_grid(num_integ_points, my_E_Range, x_tw)
END IF
ALLOCATE (tj(num_integ_points))
tj = 0.0_dp
ALLOCATE (wj(num_integ_points))
wj = 0.0_dp
DO jquad = 1, num_integ_points
tj(jquad) = x_tw(jquad)
wj(jquad) = x_tw(jquad + num_integ_points)
END DO
DEALLOCATE (x_tw)
IF (unit_nr > 0 .AND. do_print) THEN
WRITE (UNIT=unit_nr, FMT="(T3,A,T75,i6)") &
"MINIMAX_INFO| Number of integration points:", num_integ_points
WRITE (UNIT=unit_nr, FMT="(T3,A,T66,F15.4)") &
"MINIMAX_INFO| Gap for the minimax approximation:", Emin
WRITE (UNIT=unit_nr, FMT="(T3,A,T66,F15.4)") &
"MINIMAX_INFO| Range for the minimax approximation:", my_E_Range
WRITE (UNIT=unit_nr, FMT="(T3,A,T54,A,T72,A)") "MINIMAX_INFO| Minimax parameters:", "Weights", "Abscissas"
DO jquad = 1, num_integ_points
WRITE (UNIT=unit_nr, FMT="(T41,F20.10,F20.10)") wj(jquad), tj(jquad)
END DO
CALL m_flush(unit_nr)
END IF
! scale the minimax parameters
tj(:) = tj(:)*Emin
wj(:) = wj(:)*Emin
ELSE
! When we perform SOS-MP2, we need an additional factor of 2 for the energies (compare with mp2_laplace.F)
! We do not need weights etc. for the cosine transform
! We do not scale Emax because it is not needed for SOS-MP2
Emin = Emin*2.0_dp
END IF
! set up the minimax time grid
IF (do_im_time .OR. do_ri_sos_laplace_mp2) THEN
ALLOCATE (x_tw(2*num_integ_points))
x_tw = 0.0_dp
IF (num_integ_points .LE. 20) THEN
CALL get_exp_minimax_coeff(num_integ_points, my_E_Range, x_tw)
ELSE
CALL get_exp_minimax_coeff_gw(num_integ_points, my_E_Range, x_tw)
END IF
! For RPA we include already a factor of two (see later steps)
scaling = 2.0_dp
IF (do_ri_sos_laplace_mp2) scaling = 1.0_dp
ALLOCATE (tau_tj(0:num_integ_points))
tau_tj = 0.0_dp
ALLOCATE (tau_wj(num_integ_points))
tau_wj = 0.0_dp
DO jquad = 1, num_integ_points
tau_tj(jquad) = x_tw(jquad)/scaling
tau_wj(jquad) = x_tw(jquad + num_integ_points)/scaling
END DO
DEALLOCATE (x_tw)
IF (unit_nr > 0 .AND. do_print) THEN
WRITE (UNIT=unit_nr, FMT="(T3,A,T66,F15.4)") &
"MINIMAX_INFO| Range for the minimax approximation:", my_E_Range
! For testing the gap
WRITE (UNIT=unit_nr, FMT="(T3,A,T66,F15.4)") &
"MINIMAX_INFO| Gap:", Emin
WRITE (UNIT=unit_nr, FMT="(T3,A,T54,A,T72,A)") &
"MINIMAX_INFO| Minimax parameters of the time grid:", "Weights", "Abscissas"
DO jquad = 1, num_integ_points
WRITE (UNIT=unit_nr, FMT="(T41,F20.10,F20.10)") tau_wj(jquad), tau_tj(jquad)
END DO
CALL m_flush(unit_nr)
END IF
! scale grid from [1,R] to [Emin,Emax]
tau_tj(:) = tau_tj(:)/Emin
tau_wj(:) = tau_wj(:)/Emin
IF (.NOT. do_ri_sos_laplace_mp2) THEN
ALLOCATE (weights_cos_tf_t_to_w(num_integ_points, num_integ_points))
weights_cos_tf_t_to_w = 0.0_dp
CALL get_l_sq_wghts_cos_tf_t_to_w(num_integ_points, tau_tj, weights_cos_tf_t_to_w, tj, &
Emin, Emax, max_error_min, num_points_per_magnitude, &
my_regularization)
! get the weights for the cosine transform W^c(iw) -> W^c(it)
ALLOCATE (weights_cos_tf_w_to_t(num_integ_points, num_integ_points))
weights_cos_tf_w_to_t = 0.0_dp
CALL get_l_sq_wghts_cos_tf_w_to_t(num_integ_points, tau_tj, weights_cos_tf_w_to_t, tj, &
Emin, Emax, max_error_min, num_points_per_magnitude, &
my_regularization)
IF (do_gw_im_time) THEN
! get the weights for the sine transform Sigma^sin(it) -> Sigma^sin(iw) (PRB 94, 165109 (2016), Eq. 71)
ALLOCATE (weights_sin_tf_t_to_w(num_integ_points, num_integ_points))
weights_sin_tf_t_to_w = 0.0_dp
CALL get_l_sq_wghts_sin_tf_t_to_w(num_integ_points, tau_tj, weights_sin_tf_t_to_w, tj, &
Emin, Emax, max_error_min, num_points_per_magnitude, &
my_regularization)
IF (unit_nr > 0) THEN
WRITE (UNIT=unit_nr, FMT="(T3,A,T66,ES15.2)") &
"MINIMAX_INFO| Maximum deviation of the imag. time fit:", max_error_min
END IF
END IF
END IF
END IF
CALL timestop(handle)
END SUBROUTINE get_minimax_grid
! **************************************************************************************************
!> \brief ...
!> \param para_env ...
!> \param para_env_RPA ...
!> \param unit_nr ...
!> \param homo ...
!> \param virtual ...
!> \param Eigenval ...
!> \param num_integ_points ...
!> \param num_integ_group ...
!> \param color_rpa_group ...
!> \param fm_mat_S ...
!> \param my_do_gw ...
!> \param ext_scaling ...
!> \param a_scaling ...
!> \param tj ...
!> \param wj ...
! **************************************************************************************************
SUBROUTINE get_clenshaw_grid(para_env, para_env_RPA, unit_nr, homo, virtual, Eigenval, num_integ_points, &
num_integ_group, color_rpa_group, fm_mat_S, my_do_gw, &
ext_scaling, a_scaling, tj, wj)
TYPE(mp_para_env_type), INTENT(IN) :: para_env, para_env_RPA
INTEGER, INTENT(IN) :: unit_nr
INTEGER, DIMENSION(:), INTENT(IN) :: homo, virtual
REAL(KIND=dp), DIMENSION(:, :, :), INTENT(IN) :: Eigenval
INTEGER, INTENT(IN) :: num_integ_points, num_integ_group, &
color_rpa_group
TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: fm_mat_S
LOGICAL, INTENT(IN) :: my_do_gw
REAL(KIND=dp), INTENT(IN) :: ext_scaling
REAL(KIND=dp), INTENT(OUT) :: a_scaling
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:), &
INTENT(OUT) :: tj, wj
CHARACTER(LEN=*), PARAMETER :: routineN = 'get_clenshaw_grid'
INTEGER :: handle, jquad, nspins
LOGICAL :: my_open_shell
CALL timeset(routineN, handle)
nspins = SIZE(homo)
my_open_shell = (nspins == 2)
! Now, start to prepare the different grid
ALLOCATE (tj(num_integ_points))
tj = 0.0_dp
ALLOCATE (wj(num_integ_points))
wj = 0.0_dp
DO jquad = 1, num_integ_points - 1
tj(jquad) = jquad*pi/(2.0_dp*num_integ_points)
wj(jquad) = pi/(num_integ_points*SIN(tj(jquad))**2)
END DO
tj(num_integ_points) = pi/2.0_dp
wj(num_integ_points) = pi/(2.0_dp*num_integ_points*SIN(tj(num_integ_points))**2)
IF (my_do_gw .AND. ext_scaling > 0.0_dp) THEN
a_scaling = ext_scaling
ELSE
CALL calc_scaling_factor(a_scaling, para_env, para_env_RPA, homo, virtual, Eigenval, &
num_integ_points, num_integ_group, color_rpa_group, &
tj, wj, fm_mat_S)
END IF
IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T56,F25.5)') 'INTEG_INFO| Scaling parameter:', a_scaling
wj(:) = wj(:)*a_scaling
CALL timestop(handle)
END SUBROUTINE get_clenshaw_grid
! **************************************************************************************************
!> \brief ...
!> \param a_scaling_ext ...
!> \param para_env ...
!> \param para_env_RPA ...
!> \param homo ...
!> \param virtual ...
!> \param Eigenval ...
!> \param num_integ_points ...
!> \param num_integ_group ...
!> \param color_rpa_group ...
!> \param tj_ext ...
!> \param wj_ext ...
!> \param fm_mat_S ...
! **************************************************************************************************
SUBROUTINE calc_scaling_factor(a_scaling_ext, para_env, para_env_RPA, homo, virtual, Eigenval, &
num_integ_points, num_integ_group, color_rpa_group, &
tj_ext, wj_ext, fm_mat_S)
REAL(KIND=dp), INTENT(OUT) :: a_scaling_ext
TYPE(mp_para_env_type), INTENT(IN) :: para_env, para_env_RPA
INTEGER, DIMENSION(:), INTENT(IN) :: homo, virtual
REAL(KIND=dp), DIMENSION(:, :, :), INTENT(IN) :: Eigenval
INTEGER, INTENT(IN) :: num_integ_points, num_integ_group, &
color_rpa_group
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:), &
INTENT(IN) :: tj_ext, wj_ext
TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: fm_mat_S
CHARACTER(LEN=*), PARAMETER :: routineN = 'calc_scaling_factor'
INTEGER :: handle, icycle, jquad, ncol_local, &
ncol_local_beta, nspins
LOGICAL :: my_open_shell
REAL(KIND=dp) :: a_high, a_low, a_scaling, conv_param, eps, first_deriv, left_term, &
right_term, right_term_ref, right_term_ref_beta, step
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: cottj, D_ia, D_ia_beta, iaia_RI, &
iaia_RI_beta, M_ia, M_ia_beta
TYPE(mp_para_env_type), POINTER :: para_env_col, para_env_col_beta
CALL timeset(routineN, handle)
nspins = SIZE(homo)
my_open_shell = (nspins == 2)
eps = 1.0E-10_dp
ALLOCATE (cottj(num_integ_points))
! calculate the cotangent of the abscissa tj
DO jquad = 1, num_integ_points
cottj(jquad) = 1.0_dp/TAN(tj_ext(jquad))
END DO
CALL calc_ia_ia_integrals(para_env_RPA, homo(1), virtual(1), ncol_local, right_term_ref, Eigenval(:, 1, 1), &
D_ia, iaia_RI, M_ia, fm_mat_S(1), para_env_col)
! In the open shell case do point 1-2-3 for the beta spin
IF (my_open_shell) THEN
CALL calc_ia_ia_integrals(para_env_RPA, homo(2), virtual(2), ncol_local_beta, right_term_ref_beta, Eigenval(:, 1, 2), &
D_ia_beta, iaia_RI_beta, M_ia_beta, fm_mat_S(2), para_env_col_beta)
right_term_ref = right_term_ref + right_term_ref_beta
END IF
! bcast the result
IF (para_env%mepos == 0) THEN
CALL para_env%bcast(right_term_ref, 0)
ELSE
right_term_ref = 0.0_dp
CALL para_env%bcast(right_term_ref, 0)
END IF
! 5) start iteration for solving the non-linear equation by bisection
! find limit, here step=0.5 seems a good compromise
conv_param = 100.0_dp*EPSILON(right_term_ref)
step = 0.5_dp
a_low = 0.0_dp
a_high = step
right_term = -right_term_ref
DO icycle = 1, num_integ_points*2
a_scaling = a_high
CALL calculate_objfunc(a_scaling, left_term, first_deriv, num_integ_points, my_open_shell, &
M_ia, cottj, wj_ext, D_ia, D_ia_beta, M_ia_beta, &
ncol_local, ncol_local_beta, num_integ_group, color_rpa_group, &
para_env, para_env_col, para_env_col_beta)
left_term = left_term/4.0_dp/pi*a_scaling
IF (ABS(left_term) > ABS(right_term) .OR. ABS(left_term + right_term) <= conv_param) EXIT
a_low = a_high
a_high = a_high + step
END DO
IF (ABS(left_term + right_term) >= conv_param) THEN
IF (a_scaling >= 2*num_integ_points*step) THEN
a_scaling = 1.0_dp
ELSE
DO icycle = 1, num_integ_points*2
a_scaling = (a_low + a_high)/2.0_dp
CALL calculate_objfunc(a_scaling, left_term, first_deriv, num_integ_points, my_open_shell, &
M_ia, cottj, wj_ext, D_ia, D_ia_beta, M_ia_beta, &
ncol_local, ncol_local_beta, num_integ_group, color_rpa_group, &
para_env, para_env_col, para_env_col_beta)
left_term = left_term/4.0_dp/pi*a_scaling
IF (ABS(left_term) > ABS(right_term)) THEN
a_high = a_scaling
ELSE
a_low = a_scaling
END IF
IF (ABS(a_high - a_low) < 1.0e-5_dp) EXIT
END DO
END IF
END IF
a_scaling_ext = a_scaling
CALL para_env%bcast(a_scaling_ext, 0)
DEALLOCATE (cottj)
DEALLOCATE (iaia_RI)
DEALLOCATE (D_ia)
DEALLOCATE (M_ia)
CALL mp_para_env_release(para_env_col)
IF (my_open_shell) THEN
DEALLOCATE (iaia_RI_beta)
DEALLOCATE (D_ia_beta)
DEALLOCATE (M_ia_beta)
CALL mp_para_env_release(para_env_col_beta)
END IF
CALL timestop(handle)
END SUBROUTINE calc_scaling_factor
! **************************************************************************************************
!> \brief ...
!> \param para_env_RPA ...
!> \param homo ...
!> \param virtual ...
!> \param ncol_local ...
!> \param right_term_ref ...
!> \param Eigenval ...
!> \param D_ia ...
!> \param iaia_RI ...
!> \param M_ia ...
!> \param fm_mat_S ...
!> \param para_env_col ...
! **************************************************************************************************
SUBROUTINE calc_ia_ia_integrals(para_env_RPA, homo, virtual, ncol_local, right_term_ref, Eigenval, &
D_ia, iaia_RI, M_ia, fm_mat_S, para_env_col)
TYPE(mp_para_env_type), INTENT(IN) :: para_env_RPA
INTEGER, INTENT(IN) :: homo, virtual
INTEGER, INTENT(OUT) :: ncol_local
REAL(KIND=dp), INTENT(OUT) :: right_term_ref
REAL(KIND=dp), DIMENSION(:), INTENT(IN) :: Eigenval
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:), &
INTENT(OUT) :: D_ia, iaia_RI, M_ia
TYPE(cp_fm_type), INTENT(IN) :: fm_mat_S
TYPE(mp_para_env_type), POINTER :: para_env_col
CHARACTER(LEN=*), PARAMETER :: routineN = 'calc_ia_ia_integrals'
INTEGER :: avirt, color_col, color_row, handle, &
i_global, iiB, iocc, nrow_local
INTEGER, DIMENSION(:), POINTER :: col_indices, row_indices
REAL(KIND=dp) :: eigen_diff
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: iaia_RI_dp
TYPE(mp_para_env_type), POINTER :: para_env_row
CALL timeset(routineN, handle)
! calculate the (ia|ia) RI integrals
! ----------------------------------
! 1) get info fm_mat_S
CALL cp_fm_get_info(matrix=fm_mat_S, &
nrow_local=nrow_local, &
ncol_local=ncol_local, &
row_indices=row_indices, &
col_indices=col_indices)
! allocate the local buffer of iaia_RI integrals (dp kind)
ALLOCATE (iaia_RI_dp(ncol_local))
iaia_RI_dp = 0.0_dp
! 2) perform the local multiplication SUM_K (ia|K)*(ia|K)
DO iiB = 1, ncol_local
iaia_RI_dp(iiB) = iaia_RI_dp(iiB) + DOT_PRODUCT(fm_mat_S%local_data(:, iiB), fm_mat_S%local_data(:, iiB))
END DO
! 3) sum the result with the processes of the RPA_group having the same columns
! _______ia______ _
! | | | | | | |
! --> | 1 | 5 | 9 | 13| SUM --> | |
! |___|__ |___|___| |_|
! | | | | | | |
! --> | 2 | 6 | 10| 14| SUM --> | |
! K |___|___|___|___| |_| (ia|ia)_RI
! | | | | | | |
! --> | 3 | 7 | 11| 15| SUM --> | |
! |___|___|___|___| |_|
! | | | | | | |
! --> | 4 | 8 | 12| 16| SUM --> | |
! |___|___|___|___| |_|
!
color_col = fm_mat_S%matrix_struct%context%mepos(2)
ALLOCATE (para_env_col)
CALL para_env_col%from_split(para_env_RPA, color_col)
CALL para_env_col%sum(iaia_RI_dp)
! convert the iaia_RI_dp into double-double precision
ALLOCATE (iaia_RI(ncol_local))
DO iiB = 1, ncol_local
iaia_RI(iiB) = iaia_RI_dp(iiB)
END DO
DEALLOCATE (iaia_RI_dp)
! 4) calculate the right hand term, D_ia is the matrix containing the
! orbital energy differences, M_ia is the diagonal of the full RPA 'excitation'
! matrix
ALLOCATE (D_ia(ncol_local))
ALLOCATE (M_ia(ncol_local))
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)
D_ia(iiB) = eigen_diff
END DO
DO iiB = 1, ncol_local
M_ia(iiB) = D_ia(iiB)*D_ia(iiB) + 2.0_dp*D_ia(iiB)*iaia_RI(iiB)
END DO
right_term_ref = 0.0_dp
DO iiB = 1, ncol_local
right_term_ref = right_term_ref + (SQRT(M_ia(iiB)) - D_ia(iiB) - iaia_RI(iiB))
END DO
right_term_ref = right_term_ref/2.0_dp
! sum the result with the processes of the RPA_group having the same row
color_row = fm_mat_S%matrix_struct%context%mepos(1)
ALLOCATE (para_env_row)
CALL para_env_row%from_split(para_env_RPA, color_row)
! allocate communication array for rows
CALL para_env_row%sum(right_term_ref)
CALL mp_para_env_release(para_env_row)
CALL timestop(handle)
END SUBROUTINE calc_ia_ia_integrals
! **************************************************************************************************
!> \brief ...
!> \param a_scaling ...
!> \param left_term ...
!> \param first_deriv ...
!> \param num_integ_points ...
!> \param my_open_shell ...
!> \param M_ia ...
!> \param cottj ...
!> \param wj ...
!> \param D_ia ...
!> \param D_ia_beta ...
!> \param M_ia_beta ...
!> \param ncol_local ...
!> \param ncol_local_beta ...
!> \param num_integ_group ...
!> \param color_rpa_group ...
!> \param para_env ...
!> \param para_env_col ...
!> \param para_env_col_beta ...
! **************************************************************************************************
SUBROUTINE calculate_objfunc(a_scaling, left_term, first_deriv, num_integ_points, my_open_shell, &
M_ia, cottj, wj, D_ia, D_ia_beta, M_ia_beta, &
ncol_local, ncol_local_beta, num_integ_group, color_rpa_group, &
para_env, para_env_col, para_env_col_beta)
REAL(KIND=dp), INTENT(IN) :: a_scaling
REAL(KIND=dp), INTENT(INOUT) :: left_term, first_deriv
INTEGER, INTENT(IN) :: num_integ_points
LOGICAL, INTENT(IN) :: my_open_shell
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:), &
INTENT(IN) :: M_ia, cottj, wj, D_ia, D_ia_beta, &
M_ia_beta
INTEGER, INTENT(IN) :: ncol_local, ncol_local_beta, &
num_integ_group, color_rpa_group
TYPE(mp_para_env_type), INTENT(IN) :: para_env, para_env_col
TYPE(mp_para_env_type), POINTER :: para_env_col_beta
INTEGER :: iiB, jquad
REAL(KIND=dp) :: first_deriv_beta, left_term_beta, omega
left_term = 0.0_dp
first_deriv = 0.0_dp
left_term_beta = 0.0_dp
first_deriv_beta = 0.0_dp
DO jquad = 1, num_integ_points
! parallelize over integration points
IF (MODULO(jquad, num_integ_group) /= color_rpa_group) CYCLE
omega = a_scaling*cottj(jquad)
DO iiB = 1, ncol_local
! parallelize over ia elements in the para_env_row group
IF (MODULO(iiB, para_env_col%num_pe) /= para_env_col%mepos) CYCLE
! calculate left_term
left_term = left_term + wj(jquad)* &
(LOG(1.0_dp + (M_ia(iiB) - D_ia(iiB)**2)/(omega**2 + D_ia(iiB)**2)) - &
(M_ia(iiB) - D_ia(iiB)**2)/(omega**2 + D_ia(iiB)**2))
first_deriv = first_deriv + wj(jquad)*cottj(jquad)**2* &
((-M_ia(iiB) + D_ia(iiB)**2)**2/((omega**2 + D_ia(iiB)**2)**2*(omega**2 + M_ia(iiB))))
END DO
IF (my_open_shell) THEN
DO iiB = 1, ncol_local_beta
! parallelize over ia elements in the para_env_row group
IF (MODULO(iiB, para_env_col_beta%num_pe) /= para_env_col_beta%mepos) CYCLE
! calculate left_term
left_term_beta = left_term_beta + wj(jquad)* &
(LOG(1.0_dp + (M_ia_beta(iiB) - D_ia_beta(iiB)**2)/(omega**2 + D_ia_beta(iiB)**2)) - &
(M_ia_beta(iiB) - D_ia_beta(iiB)**2)/(omega**2 + D_ia_beta(iiB)**2))
first_deriv_beta = &
first_deriv_beta + wj(jquad)*cottj(jquad)**2* &
((-M_ia_beta(iiB) + D_ia_beta(iiB)**2)**2/((omega**2 + D_ia_beta(iiB)**2)**2*(omega**2 + M_ia_beta(iiB))))
END DO
END IF
END DO
! sum the contribution from all proc, starting form the row group
CALL para_env%sum(left_term)
CALL para_env%sum(first_deriv)
IF (my_open_shell) THEN
CALL para_env%sum(left_term_beta)
CALL para_env%sum(first_deriv_beta)
left_term = left_term + left_term_beta
first_deriv = first_deriv + first_deriv_beta
END IF
END SUBROUTINE calculate_objfunc
! **************************************************************************************************
!> \brief Calculate integration weights for the tau grid (in dependency of the omega node)
!> \param num_integ_points ...
!> \param tau_tj ...
!> \param weights_cos_tf_t_to_w ...
!> \param omega_tj ...
!> \param E_min ...
!> \param E_max ...
!> \param max_error ...
!> \param num_points_per_magnitude ...
!> \param regularization ...
! **************************************************************************************************
SUBROUTINE get_l_sq_wghts_cos_tf_t_to_w(num_integ_points, tau_tj, weights_cos_tf_t_to_w, omega_tj, &
E_min, E_max, max_error, num_points_per_magnitude, &
regularization)
INTEGER, INTENT(IN) :: num_integ_points
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:), &
INTENT(IN) :: tau_tj
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :), &
INTENT(INOUT) :: weights_cos_tf_t_to_w
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:), &
INTENT(IN) :: omega_tj
REAL(KIND=dp), INTENT(IN) :: E_min, E_max
REAL(KIND=dp), INTENT(INOUT) :: max_error
INTEGER, INTENT(IN) :: num_points_per_magnitude
REAL(KIND=dp), INTENT(IN) :: regularization
CHARACTER(LEN=*), PARAMETER :: routineN = 'get_l_sq_wghts_cos_tf_t_to_w'
INTEGER :: handle, iii, info, jjj, jquad, lwork, &
num_x_nodes
INTEGER, ALLOCATABLE, DIMENSION(:) :: iwork
REAL(KIND=dp) :: multiplicator, omega
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: sing_values, tau_wj_work, vec_UTy, work, &
x_values, y_values
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: mat_A, mat_SinvVSinvSigma, &
mat_SinvVSinvT, mat_U
CALL timeset(routineN, handle)
! take num_points_per_magnitude points per 10-interval
num_x_nodes = (INT(LOG10(E_max/E_min)) + 1)*num_points_per_magnitude
! take at least as many x points as integration points to have clear
! input for the singular value decomposition
num_x_nodes = MAX(num_x_nodes, num_integ_points)
ALLOCATE (x_values(num_x_nodes))
x_values = 0.0_dp
ALLOCATE (y_values(num_x_nodes))
y_values = 0.0_dp
ALLOCATE (mat_A(num_x_nodes, num_integ_points))
mat_A = 0.0_dp
ALLOCATE (tau_wj_work(num_integ_points))
tau_wj_work = 0.0_dp
ALLOCATE (sing_values(num_integ_points))
sing_values = 0.0_dp
ALLOCATE (mat_U(num_x_nodes, num_x_nodes))
mat_U = 0.0_dp
ALLOCATE (mat_SinvVSinvT(num_x_nodes, num_integ_points))
mat_SinvVSinvT = 0.0_dp
! double the value nessary for 'A' to achieve good performance
lwork = 8*num_integ_points*num_integ_points + 12*num_integ_points + 2*num_x_nodes
ALLOCATE (work(lwork))
work = 0.0_dp
ALLOCATE (iwork(8*num_integ_points))
iwork = 0
ALLOCATE (mat_SinvVSinvSigma(num_integ_points, num_x_nodes))
mat_SinvVSinvSigma = 0.0_dp
ALLOCATE (vec_UTy(num_x_nodes))
vec_UTy = 0.0_dp
max_error = 0.0_dp
! loop over all omega frequency points
DO jquad = 1, num_integ_points
! set the x-values logarithmically in the interval [Emin,Emax]
multiplicator = (E_max/E_min)**(1.0_dp/(REAL(num_x_nodes, KIND=dp) - 1.0_dp))
DO iii = 1, num_x_nodes
x_values(iii) = E_min*multiplicator**(iii - 1)
END DO
omega = omega_tj(jquad)
! y=2x/(x^2+omega_k^2)
DO iii = 1, num_x_nodes
y_values(iii) = 2.0_dp*x_values(iii)/((x_values(iii))**2 + omega**2)
END DO
! calculate mat_A
DO jjj = 1, num_integ_points
DO iii = 1, num_x_nodes
mat_A(iii, jjj) = COS(omega*tau_tj(jjj))*EXP(-x_values(iii)*tau_tj(jjj))
END DO
END DO
! Singular value decomposition of mat_A
CALL DGESDD('A', num_x_nodes, num_integ_points, mat_A, num_x_nodes, sing_values, mat_U, num_x_nodes, &
mat_SinvVSinvT, num_x_nodes, work, lwork, iwork, info)
CPASSERT(info == 0)
! integration weights = V Sigma U^T y
! 1) V*Sigma
DO jjj = 1, num_integ_points
DO iii = 1, num_integ_points
! mat_SinvVSinvSigma(iii, jjj) = mat_SinvVSinvT(jjj, iii)/sing_values(jjj)
mat_SinvVSinvSigma(iii, jjj) = mat_SinvVSinvT(jjj, iii)*sing_values(jjj) &
/(regularization**2 + sing_values(jjj)**2)
END DO
END DO
! 2) U^T y
CALL DGEMM('T', 'N', num_x_nodes, 1, num_x_nodes, 1.0_dp, mat_U, num_x_nodes, y_values, num_x_nodes, &
0.0_dp, vec_UTy, num_x_nodes)
! 3) (V*Sigma) * (U^T y)
CALL DGEMM('N', 'N', num_integ_points, 1, num_x_nodes, 1.0_dp, mat_SinvVSinvSigma, num_integ_points, vec_UTy, &
num_x_nodes, 0.0_dp, tau_wj_work, num_integ_points)
weights_cos_tf_t_to_w(jquad, :) = tau_wj_work(:)
CALL calc_max_error_fit_tau_grid_with_cosine(max_error, omega, tau_tj, tau_wj_work, x_values, &
y_values, num_integ_points, num_x_nodes)
END DO ! jquad
DEALLOCATE (x_values, y_values, mat_A, tau_wj_work, sing_values, mat_U, mat_SinvVSinvT, &
work, iwork, mat_SinvVSinvSigma, vec_UTy)
CALL timestop(handle)
END SUBROUTINE get_l_sq_wghts_cos_tf_t_to_w
! **************************************************************************************************
!> \brief Calculate integration weights for the tau grid (in dependency of the omega node)
!> \param num_integ_points ...
!> \param tau_tj ...
!> \param weights_sin_tf_t_to_w ...
!> \param omega_tj ...
!> \param E_min ...
!> \param E_max ...
!> \param max_error ...
!> \param num_points_per_magnitude ...
!> \param regularization ...
! **************************************************************************************************
SUBROUTINE get_l_sq_wghts_sin_tf_t_to_w(num_integ_points, tau_tj, weights_sin_tf_t_to_w, omega_tj, &
E_min, E_max, max_error, num_points_per_magnitude, regularization)
INTEGER, INTENT(IN) :: num_integ_points
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:), &
INTENT(IN) :: tau_tj
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :), &
INTENT(INOUT) :: weights_sin_tf_t_to_w
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:), &
INTENT(IN) :: omega_tj
REAL(KIND=dp), INTENT(IN) :: E_min, E_max
REAL(KIND=dp), INTENT(OUT) :: max_error
INTEGER, INTENT(IN) :: num_points_per_magnitude
REAL(KIND=dp), INTENT(IN) :: regularization
CHARACTER(LEN=*), PARAMETER :: routineN = 'get_l_sq_wghts_sin_tf_t_to_w'
INTEGER :: handle, iii, info, jjj, jquad, lwork, &
num_x_nodes
INTEGER, ALLOCATABLE, DIMENSION(:) :: iwork
REAL(KIND=dp) :: chi2_min_jquad, multiplicator, omega
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: sing_values, tau_wj_work, vec_UTy, work, &
work_array, x_values, y_values
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: mat_A, mat_SinvVSinvSigma, &
mat_SinvVSinvT, mat_U
CALL timeset(routineN, handle)
! take num_points_per_magnitude points per 10-interval
num_x_nodes = (INT(LOG10(E_max/E_min)) + 1)*num_points_per_magnitude
! take at least as many x points as integration points to have clear
! input for the singular value decomposition
num_x_nodes = MAX(num_x_nodes, num_integ_points)
ALLOCATE (x_values(num_x_nodes))
x_values = 0.0_dp
ALLOCATE (y_values(num_x_nodes))
y_values = 0.0_dp
ALLOCATE (mat_A(num_x_nodes, num_integ_points))
mat_A = 0.0_dp
ALLOCATE (tau_wj_work(num_integ_points))
tau_wj_work = 0.0_dp
ALLOCATE (work_array(2*num_integ_points))
work_array = 0.0_dp
ALLOCATE (sing_values(num_integ_points))
sing_values = 0.0_dp
ALLOCATE (mat_U(num_x_nodes, num_x_nodes))
mat_U = 0.0_dp
ALLOCATE (mat_SinvVSinvT(num_x_nodes, num_integ_points))
mat_SinvVSinvT = 0.0_dp
! double the value nessary for 'A' to achieve good performance
lwork = 8*num_integ_points*num_integ_points + 12*num_integ_points + 2*num_x_nodes
ALLOCATE (work(lwork))
work = 0.0_dp
ALLOCATE (iwork(8*num_integ_points))
iwork = 0
ALLOCATE (mat_SinvVSinvSigma(num_integ_points, num_x_nodes))
mat_SinvVSinvSigma = 0.0_dp
ALLOCATE (vec_UTy(num_x_nodes))
vec_UTy = 0.0_dp
max_error = 0.0_dp
! loop over all omega frequency points
DO jquad = 1, num_integ_points
chi2_min_jquad = 100.0_dp
! set the x-values logarithmically in the interval [Emin,Emax]
multiplicator = (E_max/E_min)**(1.0_dp/(REAL(num_x_nodes, KIND=dp) - 1.0_dp))
DO iii = 1, num_x_nodes
x_values(iii) = E_min*multiplicator**(iii - 1)
END DO
omega = omega_tj(jquad)
! y=2x/(x^2+omega_k^2)
DO iii = 1, num_x_nodes
! y_values(iii) = 2.0_dp*x_values(iii)/((x_values(iii))**2+omega**2)
y_values(iii) = 2.0_dp*omega/((x_values(iii))**2 + omega**2)
END DO
! calculate mat_A
DO jjj = 1, num_integ_points
DO iii = 1, num_x_nodes
mat_A(iii, jjj) = SIN(omega*tau_tj(jjj))*EXP(-x_values(iii)*tau_tj(jjj))
END DO
END DO
! Singular value decomposition of mat_A
CALL DGESDD('A', num_x_nodes, num_integ_points, mat_A, num_x_nodes, sing_values, mat_U, num_x_nodes, &
mat_SinvVSinvT, num_x_nodes, work, lwork, iwork, info)
CPASSERT(info == 0)
! integration weights = V Sigma U^T y
! 1) V*Sigma
DO jjj = 1, num_integ_points
DO iii = 1, num_integ_points
! mat_SinvVSinvSigma(iii, jjj) = mat_SinvVSinvT(jjj, iii)/sing_values(jjj)
mat_SinvVSinvSigma(iii, jjj) = mat_SinvVSinvT(jjj, iii)*sing_values(jjj) &
/(regularization**2 + sing_values(jjj)**2)
END DO
END DO
! 2) U^T y
CALL DGEMM('T', 'N', num_x_nodes, 1, num_x_nodes, 1.0_dp, mat_U, num_x_nodes, y_values, num_x_nodes, &
0.0_dp, vec_UTy, num_x_nodes)
! 3) (V*Sigma) * (U^T y)
CALL DGEMM('N', 'N', num_integ_points, 1, num_x_nodes, 1.0_dp, mat_SinvVSinvSigma, num_integ_points, vec_UTy, &
num_x_nodes, 0.0_dp, tau_wj_work, num_integ_points)
weights_sin_tf_t_to_w(jquad, :) = tau_wj_work(:)
CALL calc_max_error_fit_tau_grid_with_sine(max_error, omega, tau_tj, tau_wj_work, x_values, &
y_values, num_integ_points, num_x_nodes)
END DO ! jquad
DEALLOCATE (x_values, y_values, mat_A, tau_wj_work, work_array, sing_values, mat_U, mat_SinvVSinvT, &
work, iwork, mat_SinvVSinvSigma, vec_UTy)
CALL timestop(handle)