-
Notifications
You must be signed in to change notification settings - Fork 1
/
ct_methods.F
1751 lines (1530 loc) · 71.9 KB
/
ct_methods.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 Cayley transformation methods
!> \par History
!> 2011.06 created [Rustam Z Khaliullin]
!> \author Rustam Z Khaliullin
! **************************************************************************************************
MODULE ct_methods
USE cp_dbcsr_api, ONLY: &
dbcsr_add, dbcsr_add_on_diag, dbcsr_copy, dbcsr_create, dbcsr_desymmetrize, &
dbcsr_distribution_get, dbcsr_distribution_type, dbcsr_dot, dbcsr_filter, dbcsr_finalize, &
dbcsr_frobenius_norm, dbcsr_func_inverse, dbcsr_function_of_elements, dbcsr_get_diag, &
dbcsr_get_info, dbcsr_get_stored_coordinates, dbcsr_hadamard_product, &
dbcsr_iterator_blocks_left, dbcsr_iterator_next_block, dbcsr_iterator_start, &
dbcsr_iterator_stop, dbcsr_iterator_type, dbcsr_multiply, dbcsr_nblkcols_total, &
dbcsr_nblkrows_total, dbcsr_norm, dbcsr_norm_maxabsnorm, dbcsr_release, &
dbcsr_reserve_block2d, dbcsr_scale, dbcsr_set, dbcsr_set_diag, dbcsr_transposed, &
dbcsr_type, dbcsr_type_no_symmetry, dbcsr_work_create
USE cp_dbcsr_cholesky, ONLY: cp_dbcsr_cholesky_decompose,&
cp_dbcsr_cholesky_invert
USE cp_dbcsr_diag, ONLY: cp_dbcsr_syevd
USE cp_log_handling, ONLY: cp_get_default_logger,&
cp_logger_get_default_unit_nr,&
cp_logger_type
USE ct_types, ONLY: ct_step_env_type
USE input_constants, ONLY: &
cg_dai_yuan, cg_fletcher, cg_fletcher_reeves, cg_hager_zhang, cg_hestenes_stiefel, &
cg_liu_storey, cg_polak_ribiere, cg_zero, tensor_orthogonal, tensor_up_down
USE iterate_matrix, ONLY: matrix_sqrt_Newton_Schulz
USE kinds, ONLY: dp
USE machine, ONLY: m_walltime
#include "./base/base_uses.f90"
IMPLICIT NONE
PRIVATE
CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'ct_methods'
! Public subroutines
PUBLIC :: ct_step_execute, analytic_line_search, diagonalize_diagonal_blocks
CONTAINS
! **************************************************************************************************
!> \brief Performs Cayley transformation
!> \param cts_env ...
!> \par History
!> 2011.06 created [Rustam Z Khaliullin]
!> \author Rustam Z Khaliullin
! **************************************************************************************************
SUBROUTINE ct_step_execute(cts_env)
TYPE(ct_step_env_type) :: cts_env
CHARACTER(len=*), PARAMETER :: routineN = 'ct_step_execute'
INTEGER :: handle, n, preconditioner_type, unit_nr
REAL(KIND=dp) :: gap_estimate, safety_margin
REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: evals
TYPE(cp_logger_type), POINTER :: logger
TYPE(dbcsr_type) :: matrix_pp, matrix_pq, matrix_qp, &
matrix_qp_save, matrix_qq, oo1, &
oo1_sqrt, oo1_sqrt_inv, t_corr, tmp1, &
u_pp, u_qq
!TYPE(dbcsr_type) :: rst_x1, rst_x2
!REAL(KIND=dp) :: ener_tmp
!TYPE(dbcsr_iterator_type) :: iter
!INTEGER :: iblock_row,iblock_col,&
! iblock_row_size,iblock_col_size
!REAL(KIND=dp), DIMENSION(:,:), POINTER :: data_p
CALL timeset(routineN, handle)
logger => cp_get_default_logger()
IF (logger%para_env%is_source()) THEN
unit_nr = cp_logger_get_default_unit_nr(logger, local=.TRUE.)
ELSE
unit_nr = -1
END IF
! check if all input is in place and flags are consistent
IF (cts_env%update_q .AND. (.NOT. cts_env%update_p)) THEN
CPABORT("q-update is possible only with p-update")
END IF
IF (cts_env%tensor_type .EQ. tensor_up_down) THEN
CPABORT("riccati is not implemented for biorthogonal basis")
END IF
IF (.NOT. ASSOCIATED(cts_env%matrix_ks)) THEN
CPABORT("KS matrix is not associated")
END IF
IF (cts_env%use_virt_orbs .AND. (.NOT. cts_env%use_occ_orbs)) THEN
CPABORT("virtual orbs can be used only with occupied orbs")
END IF
IF (cts_env%use_occ_orbs) THEN
IF (.NOT. ASSOCIATED(cts_env%matrix_t)) THEN
CPABORT("T matrix is not associated")
END IF
IF (.NOT. ASSOCIATED(cts_env%matrix_qp_template)) THEN
CPABORT("QP template is not associated")
END IF
IF (.NOT. ASSOCIATED(cts_env%matrix_pq_template)) THEN
CPABORT("PQ template is not associated")
END IF
END IF
IF (cts_env%use_virt_orbs) THEN
IF (.NOT. ASSOCIATED(cts_env%matrix_v)) THEN
CPABORT("V matrix is not associated")
END IF
ELSE
IF (.NOT. ASSOCIATED(cts_env%matrix_p)) THEN
CPABORT("P matrix is not associated")
END IF
END IF
IF (cts_env%tensor_type .NE. tensor_up_down .AND. &
cts_env%tensor_type .NE. tensor_orthogonal) THEN
CPABORT("illegal tensor flag")
END IF
! start real calculations
IF (cts_env%use_occ_orbs) THEN
! create matrices for various ks blocks
CALL dbcsr_create(matrix_pp, &
template=cts_env%p_index_up, &
matrix_type=dbcsr_type_no_symmetry)
CALL dbcsr_create(matrix_qp, &
template=cts_env%matrix_qp_template, &
matrix_type=dbcsr_type_no_symmetry)
CALL dbcsr_create(matrix_qq, &
template=cts_env%q_index_up, &
matrix_type=dbcsr_type_no_symmetry)
CALL dbcsr_create(matrix_pq, &
template=cts_env%matrix_pq_template, &
matrix_type=dbcsr_type_no_symmetry)
! create the residue matrix
CALL dbcsr_create(cts_env%matrix_res, &
template=cts_env%matrix_qp_template)
CALL assemble_ks_qp_blocks(cts_env%matrix_ks, &
cts_env%matrix_p, &
cts_env%matrix_t, &
cts_env%matrix_v, &
cts_env%q_index_down, &
cts_env%p_index_up, &
cts_env%q_index_up, &
matrix_pp, &
matrix_qq, &
matrix_qp, &
matrix_pq, &
cts_env%tensor_type, &
cts_env%use_virt_orbs, &
cts_env%eps_filter)
! create a matrix of single-excitation amplitudes
CALL dbcsr_create(cts_env%matrix_x, &
template=cts_env%matrix_qp_template)
IF (ASSOCIATED(cts_env%matrix_x_guess)) THEN
CALL dbcsr_copy(cts_env%matrix_x, &
cts_env%matrix_x_guess)
IF (cts_env%tensor_type .EQ. tensor_orthogonal) THEN
! bring x from contravariant-covariant representation
! to the orthogonal/cholesky representation
! use res as temporary storage
CALL dbcsr_multiply("N", "N", 1.0_dp, cts_env%q_index_down, &
cts_env%matrix_x, 0.0_dp, cts_env%matrix_res, &
filter_eps=cts_env%eps_filter)
CALL dbcsr_multiply("N", "N", 1.0_dp, cts_env%matrix_res, &
cts_env%p_index_up, 0.0_dp, &
cts_env%matrix_x, &
filter_eps=cts_env%eps_filter)
END IF
ELSE
! set amplitudes to zero
CALL dbcsr_set(cts_env%matrix_x, 0.0_dp)
END IF
!SELECT CASE (cts_env%preconditioner_type)
!CASE (prec_eigenvector_blocks,prec_eigenvector_full)
preconditioner_type = 1
safety_margin = 2.0_dp
gap_estimate = 0.0001_dp
SELECT CASE (preconditioner_type)
CASE (1, 2)
!RZK-warning diagonalization works only with orthogonal tensor!!!
! find a better basis by diagonalizing diagonal blocks
! first pp
CALL dbcsr_create(u_pp, template=matrix_pp, &
matrix_type=dbcsr_type_no_symmetry)
!IF (cts_env%preconditioner_type.eq.prec_eigenvector_full) THEN
IF (.TRUE.) THEN
CALL dbcsr_get_info(matrix_pp, nfullrows_total=n)
ALLOCATE (evals(n))
CALL cp_dbcsr_syevd(matrix_pp, u_pp, evals, &
cts_env%para_env, cts_env%blacs_env)
DEALLOCATE (evals)
ELSE
CALL diagonalize_diagonal_blocks(matrix_pp, u_pp)
END IF
! and now qq
CALL dbcsr_create(u_qq, template=matrix_qq, &
matrix_type=dbcsr_type_no_symmetry)
!IF (cts_env%preconditioner_type.eq.prec_eigenvector_full) THEN
IF (.TRUE.) THEN
CALL dbcsr_get_info(matrix_qq, nfullrows_total=n)
ALLOCATE (evals(n))
CALL cp_dbcsr_syevd(matrix_qq, u_qq, evals, &
cts_env%para_env, cts_env%blacs_env)
DEALLOCATE (evals)
ELSE
CALL diagonalize_diagonal_blocks(matrix_qq, u_qq)
END IF
! apply the transformation to all matrices
CALL matrix_forward_transform(matrix_pp, u_pp, u_pp, &
cts_env%eps_filter)
CALL matrix_forward_transform(matrix_qq, u_qq, u_qq, &
cts_env%eps_filter)
CALL matrix_forward_transform(matrix_qp, u_qq, u_pp, &
cts_env%eps_filter)
CALL matrix_forward_transform(matrix_pq, u_pp, u_qq, &
cts_env%eps_filter)
CALL matrix_forward_transform(cts_env%matrix_x, u_qq, u_pp, &
cts_env%eps_filter)
IF (cts_env%max_iter .GE. 0) THEN
CALL solve_riccati_equation( &
pp=matrix_pp, &
qq=matrix_qq, &
qp=matrix_qp, &
pq=matrix_pq, &
x=cts_env%matrix_x, &
res=cts_env%matrix_res, &
neglect_quadratic_term=cts_env%neglect_quadratic_term, &
conjugator=cts_env%conjugator, &
max_iter=cts_env%max_iter, &
eps_convergence=cts_env%eps_convergence, &
eps_filter=cts_env%eps_filter, &
converged=cts_env%converged)
IF (cts_env%converged) THEN
!IF (unit_nr>0) THEN
! WRITE(unit_nr,*)
! WRITE(unit_nr,'(T6,A)') &
! "RICCATI equations solved"
! CALL m_flush(unit_nr)
!ENDIF
ELSE
CPABORT("RICCATI: CG algorithm has NOT converged")
END IF
END IF
IF (cts_env%calculate_energy_corr) THEN
CALL dbcsr_dot(matrix_qp, cts_env%matrix_x, cts_env%energy_correction)
END IF
CALL dbcsr_release(matrix_pp)
CALL dbcsr_release(matrix_qp)
CALL dbcsr_release(matrix_qq)
CALL dbcsr_release(matrix_pq)
! back-transform to the original basis
CALL matrix_backward_transform(cts_env%matrix_x, u_qq, &
u_pp, cts_env%eps_filter)
CALL dbcsr_release(u_qq)
CALL dbcsr_release(u_pp)
!CASE (prec_cholesky_inverse)
CASE (3)
! RZK-warning implemented only for orthogonal tensors!!!
! generalization to up_down should be easy
CALL dbcsr_create(u_pp, template=matrix_pp, &
matrix_type=dbcsr_type_no_symmetry)
CALL dbcsr_copy(u_pp, matrix_pp)
CALL dbcsr_scale(u_pp, -1.0_dp)
CALL dbcsr_add_on_diag(u_pp, &
ABS(safety_margin*gap_estimate))
CALL cp_dbcsr_cholesky_decompose(u_pp, &
para_env=cts_env%para_env, &
blacs_env=cts_env%blacs_env)
CALL cp_dbcsr_cholesky_invert(u_pp, &
para_env=cts_env%para_env, &
blacs_env=cts_env%blacs_env, &
upper_to_full=.TRUE.)
!CALL dbcsr_scale(u_pp,-1.0_dp)
CALL dbcsr_create(u_qq, template=matrix_qq, &
matrix_type=dbcsr_type_no_symmetry)
CALL dbcsr_copy(u_qq, matrix_qq)
CALL dbcsr_add_on_diag(u_qq, &
ABS(safety_margin*gap_estimate))
CALL cp_dbcsr_cholesky_decompose(u_qq, &
para_env=cts_env%para_env, &
blacs_env=cts_env%blacs_env)
CALL cp_dbcsr_cholesky_invert(u_qq, &
para_env=cts_env%para_env, &
blacs_env=cts_env%blacs_env, &
upper_to_full=.TRUE.)
! transform all riccati matrices (left-right preconditioner)
CALL dbcsr_create(tmp1, template=matrix_qq, &
matrix_type=dbcsr_type_no_symmetry)
CALL dbcsr_multiply("N", "N", 1.0_dp, u_qq, &
matrix_qq, 0.0_dp, tmp1, &
filter_eps=cts_env%eps_filter)
CALL dbcsr_copy(matrix_qq, tmp1)
CALL dbcsr_release(tmp1)
CALL dbcsr_create(tmp1, template=matrix_pp, &
matrix_type=dbcsr_type_no_symmetry)
CALL dbcsr_multiply("N", "N", 1.0_dp, matrix_pp, &
u_pp, 0.0_dp, tmp1, &
filter_eps=cts_env%eps_filter)
CALL dbcsr_copy(matrix_pp, tmp1)
CALL dbcsr_release(tmp1)
CALL dbcsr_create(matrix_qp_save, template=matrix_qp, &
matrix_type=dbcsr_type_no_symmetry)
CALL dbcsr_copy(matrix_qp_save, matrix_qp)
CALL dbcsr_create(tmp1, template=matrix_qp, &
matrix_type=dbcsr_type_no_symmetry)
CALL dbcsr_multiply("N", "N", 1.0_dp, matrix_qp, &
u_pp, 0.0_dp, tmp1, &
filter_eps=cts_env%eps_filter)
CALL dbcsr_multiply("N", "N", 1.0_dp, u_qq, tmp1, &
0.0_dp, matrix_qp, &
filter_eps=cts_env%eps_filter)
CALL dbcsr_release(tmp1)
!CALL dbcsr_print(matrix_qq)
!CALL dbcsr_print(matrix_qp)
!CALL dbcsr_print(matrix_pp)
IF (cts_env%max_iter .GE. 0) THEN
CALL solve_riccati_equation( &
pp=matrix_pp, &
qq=matrix_qq, &
qp=matrix_qp, &
pq=matrix_pq, &
oo=u_pp, &
vv=u_qq, &
x=cts_env%matrix_x, &
res=cts_env%matrix_res, &
neglect_quadratic_term=cts_env%neglect_quadratic_term, &
conjugator=cts_env%conjugator, &
max_iter=cts_env%max_iter, &
eps_convergence=cts_env%eps_convergence, &
eps_filter=cts_env%eps_filter, &
converged=cts_env%converged)
IF (cts_env%converged) THEN
!IF (unit_nr>0) THEN
! WRITE(unit_nr,*)
! WRITE(unit_nr,'(T6,A)') &
! "RICCATI equations solved"
! CALL m_flush(unit_nr)
!ENDIF
ELSE
CPABORT("RICCATI: CG algorithm has NOT converged")
END IF
END IF
IF (cts_env%calculate_energy_corr) THEN
CALL dbcsr_dot(matrix_qp_save, cts_env%matrix_x, cts_env%energy_correction)
END IF
CALL dbcsr_release(matrix_qp_save)
CALL dbcsr_release(matrix_pp)
CALL dbcsr_release(matrix_qp)
CALL dbcsr_release(matrix_qq)
CALL dbcsr_release(matrix_pq)
CALL dbcsr_release(u_qq)
CALL dbcsr_release(u_pp)
CASE DEFAULT
CPABORT("illegal preconditioner type")
END SELECT ! preconditioner type
IF (cts_env%update_p) THEN
IF (cts_env%tensor_type .EQ. tensor_up_down) THEN
CPABORT("orbital update is NYI for this tensor type")
END IF
! transform occupied orbitals
! in a way that preserves the overlap metric
CALL dbcsr_create(oo1, &
template=cts_env%p_index_up, &
matrix_type=dbcsr_type_no_symmetry)
CALL dbcsr_create(oo1_sqrt_inv, &
template=oo1)
CALL dbcsr_create(oo1_sqrt, &
template=oo1)
! Compute (1+tr(X).X)^(-1/2)_up_down
CALL dbcsr_multiply("T", "N", 1.0_dp, cts_env%matrix_x, &
cts_env%matrix_x, 0.0_dp, oo1, &
filter_eps=cts_env%eps_filter)
CALL dbcsr_add_on_diag(oo1, 1.0_dp)
CALL matrix_sqrt_Newton_Schulz(oo1_sqrt, &
oo1_sqrt_inv, &
oo1, &
!if cholesky is used then sqrt
!guess cannot be provided
!matrix_sqrt_inv_guess=cts_env%p_index_up,&
!matrix_sqrt_guess=cts_env%p_index_down,&
threshold=cts_env%eps_filter, &
order=cts_env%order_lanczos, &
eps_lanczos=cts_env%eps_lancsoz, &
max_iter_lanczos=cts_env%max_iter_lanczos)
CALL dbcsr_multiply("N", "N", 1.0_dp, cts_env%p_index_up, &
oo1_sqrt_inv, 0.0_dp, oo1, &
filter_eps=cts_env%eps_filter)
CALL dbcsr_multiply("N", "N", 1.0_dp, oo1, &
cts_env%p_index_down, 0.0_dp, oo1_sqrt, &
filter_eps=cts_env%eps_filter)
CALL dbcsr_release(oo1)
CALL dbcsr_release(oo1_sqrt_inv)
! bring x to contravariant-covariant representation now
CALL dbcsr_create(matrix_qp, &
template=cts_env%matrix_qp_template, &
matrix_type=dbcsr_type_no_symmetry)
CALL dbcsr_multiply("N", "N", 1.0_dp, cts_env%q_index_up, &
cts_env%matrix_x, 0.0_dp, matrix_qp, &
filter_eps=cts_env%eps_filter)
CALL dbcsr_multiply("N", "N", 1.0_dp, matrix_qp, &
cts_env%p_index_down, 0.0_dp, &
cts_env%matrix_x, &
filter_eps=cts_env%eps_filter)
CALL dbcsr_release(matrix_qp)
! update T=T+X or T=T+V.X (whichever is appropriate)
CALL dbcsr_create(t_corr, template=cts_env%matrix_t)
IF (cts_env%use_virt_orbs) THEN
CALL dbcsr_multiply("N", "N", 1.0_dp, cts_env%matrix_v, &
cts_env%matrix_x, 0.0_dp, t_corr, &
filter_eps=cts_env%eps_filter)
CALL dbcsr_add(cts_env%matrix_t, t_corr, &
1.0_dp, 1.0_dp)
ELSE
CALL dbcsr_add(cts_env%matrix_t, cts_env%matrix_x, &
1.0_dp, 1.0_dp)
END IF
! adjust T so the metric is preserved: T=(T+X).(1+tr(X).X)^(-1/2)
CALL dbcsr_multiply("N", "N", 1.0_dp, cts_env%matrix_t, oo1_sqrt, &
0.0_dp, t_corr, filter_eps=cts_env%eps_filter)
CALL dbcsr_copy(cts_env%matrix_t, t_corr)
CALL dbcsr_release(t_corr)
CALL dbcsr_release(oo1_sqrt)
ELSE ! do not update p
IF (cts_env%tensor_type .EQ. tensor_orthogonal) THEN
! bring x to contravariant-covariant representation
CALL dbcsr_create(matrix_qp, &
template=cts_env%matrix_qp_template, &
matrix_type=dbcsr_type_no_symmetry)
CALL dbcsr_multiply("N", "N", 1.0_dp, cts_env%q_index_up, &
cts_env%matrix_x, 0.0_dp, matrix_qp, &
filter_eps=cts_env%eps_filter)
CALL dbcsr_multiply("N", "N", 1.0_dp, matrix_qp, &
cts_env%p_index_down, 0.0_dp, &
cts_env%matrix_x, &
filter_eps=cts_env%eps_filter)
CALL dbcsr_release(matrix_qp)
END IF
END IF
ELSE
CPABORT("illegal occ option")
END IF
CALL timestop(handle)
END SUBROUTINE ct_step_execute
! **************************************************************************************************
!> \brief computes oo, ov, vo, and vv blocks of the ks matrix
!> \param ks ...
!> \param p ...
!> \param t ...
!> \param v ...
!> \param q_index_down ...
!> \param p_index_up ...
!> \param q_index_up ...
!> \param pp ...
!> \param qq ...
!> \param qp ...
!> \param pq ...
!> \param tensor_type ...
!> \param use_virt_orbs ...
!> \param eps_filter ...
!> \par History
!> 2011.06 created [Rustam Z Khaliullin]
!> \author Rustam Z Khaliullin
! **************************************************************************************************
SUBROUTINE assemble_ks_qp_blocks(ks, p, t, v, q_index_down, &
p_index_up, q_index_up, pp, qq, qp, pq, tensor_type, use_virt_orbs, eps_filter)
TYPE(dbcsr_type), INTENT(IN) :: ks, p, t, v, q_index_down, p_index_up, &
q_index_up
TYPE(dbcsr_type), INTENT(OUT) :: pp, qq, qp, pq
INTEGER, INTENT(IN) :: tensor_type
LOGICAL, INTENT(IN) :: use_virt_orbs
REAL(KIND=dp), INTENT(IN) :: eps_filter
CHARACTER(len=*), PARAMETER :: routineN = 'assemble_ks_qp_blocks'
INTEGER :: handle
LOGICAL :: library_fixed
TYPE(dbcsr_type) :: kst, ksv, no, on, oo, q_index_up_nosym, &
sp, spf, t_or, v_or
CALL timeset(routineN, handle)
IF (use_virt_orbs) THEN
! orthogonalize the orbitals
CALL dbcsr_create(t_or, template=t)
CALL dbcsr_create(v_or, template=v)
CALL dbcsr_multiply("N", "N", 1.0_dp, t, p_index_up, &
0.0_dp, t_or, filter_eps=eps_filter)
CALL dbcsr_multiply("N", "N", 1.0_dp, v, q_index_up, &
0.0_dp, v_or, filter_eps=eps_filter)
! KS.T
CALL dbcsr_create(kst, template=t)
CALL dbcsr_multiply("N", "N", 1.0_dp, ks, t_or, &
0.0_dp, kst, filter_eps=eps_filter)
! pp=tr(T)*KS.T
CALL dbcsr_multiply("T", "N", 1.0_dp, t_or, kst, &
0.0_dp, pp, filter_eps=eps_filter)
! qp=tr(V)*KS.T
CALL dbcsr_multiply("T", "N", 1.0_dp, v_or, kst, &
0.0_dp, qp, filter_eps=eps_filter)
CALL dbcsr_release(kst)
! KS.V
CALL dbcsr_create(ksv, template=v)
CALL dbcsr_multiply("N", "N", 1.0_dp, ks, v_or, &
0.0_dp, ksv, filter_eps=eps_filter)
! tr(T)*KS.V
CALL dbcsr_multiply("T", "N", 1.0_dp, t_or, ksv, &
0.0_dp, pq, filter_eps=eps_filter)
! tr(V)*KS.V
CALL dbcsr_multiply("T", "N", 1.0_dp, v_or, ksv, &
0.0_dp, qq, filter_eps=eps_filter)
CALL dbcsr_release(ksv)
CALL dbcsr_release(t_or)
CALL dbcsr_release(v_or)
ELSE ! no virtuals, use projected AOs
! THIS PROCEDURE HAS NOT BEEN UPDATED FOR CHOLESKY p/q_index_up/down
CALL dbcsr_create(sp, template=q_index_down, &
matrix_type=dbcsr_type_no_symmetry)
CALL dbcsr_create(spf, template=q_index_down, &
matrix_type=dbcsr_type_no_symmetry)
! qp=KS*T
CALL dbcsr_multiply("N", "N", 1.0_dp, ks, t, 0.0_dp, qp, &
filter_eps=eps_filter)
! pp=tr(T)*KS.T
CALL dbcsr_multiply("T", "N", 1.0_dp, t, qp, 0.0_dp, pp, &
filter_eps=eps_filter)
! sp=-S_*P
CALL dbcsr_multiply("N", "N", -1.0_dp, q_index_down, p, 0.0_dp, sp, &
filter_eps=eps_filter)
! sp=1/S^-S_.P
SELECT CASE (tensor_type)
CASE (tensor_up_down)
CALL dbcsr_add_on_diag(sp, 1.0_dp)
CASE (tensor_orthogonal)
CALL dbcsr_create(q_index_up_nosym, template=q_index_up, &
matrix_type=dbcsr_type_no_symmetry)
CALL dbcsr_desymmetrize(q_index_up, q_index_up_nosym)
CALL dbcsr_add(sp, q_index_up_nosym, 1.0_dp, 1.0_dp)
CALL dbcsr_release(q_index_up_nosym)
END SELECT
! spf=(1/S^-S_.P)*KS
CALL dbcsr_multiply("N", "N", 1.0_dp, sp, ks, 0.0_dp, spf, &
filter_eps=eps_filter)
! qp=spf*T
CALL dbcsr_multiply("N", "N", 1.0_dp, spf, t, 0.0_dp, qp, &
filter_eps=eps_filter)
SELECT CASE (tensor_type)
CASE (tensor_up_down)
! pq=tr(qp)
CALL dbcsr_transposed(pq, qp, transpose_distribution=.FALSE.)
CASE (tensor_orthogonal)
! pq=sig^.tr(qp)
CALL dbcsr_multiply("N", "T", 1.0_dp, p_index_up, qp, 0.0_dp, pq, &
filter_eps=eps_filter)
library_fixed = .FALSE.
IF (library_fixed) THEN
CALL dbcsr_transposed(qp, pq, transpose_distribution=.FALSE.)
ELSE
CALL dbcsr_create(no, template=qp, &
matrix_type=dbcsr_type_no_symmetry)
CALL dbcsr_multiply("N", "N", 1.0_dp, qp, p_index_up, 0.0_dp, no, &
filter_eps=eps_filter)
CALL dbcsr_copy(qp, no)
CALL dbcsr_release(no)
END IF
END SELECT
! qq=spf*tr(sp)
CALL dbcsr_multiply("N", "T", 1.0_dp, spf, sp, 0.0_dp, qq, &
filter_eps=eps_filter)
SELECT CASE (tensor_type)
CASE (tensor_up_down)
CALL dbcsr_create(oo, template=pp, &
matrix_type=dbcsr_type_no_symmetry)
CALL dbcsr_create(no, template=qp, &
matrix_type=dbcsr_type_no_symmetry)
! first index up
CALL dbcsr_multiply("N", "N", 1.0_dp, q_index_up, qq, 0.0_dp, spf, &
filter_eps=eps_filter)
CALL dbcsr_copy(qq, spf)
CALL dbcsr_multiply("N", "N", 1.0_dp, q_index_up, qp, 0.0_dp, no, &
filter_eps=eps_filter)
CALL dbcsr_copy(qp, no)
CALL dbcsr_multiply("N", "N", 1.0_dp, p_index_up, pp, 0.0_dp, oo, &
filter_eps=eps_filter)
CALL dbcsr_copy(pp, oo)
CALL dbcsr_multiply("N", "N", 1.0_dp, p_index_up, pq, 0.0_dp, on, &
filter_eps=eps_filter)
CALL dbcsr_copy(pq, on)
CALL dbcsr_release(no)
CALL dbcsr_release(oo)
CASE (tensor_orthogonal)
CALL dbcsr_create(oo, template=pp, &
matrix_type=dbcsr_type_no_symmetry)
! both indeces up in the pp block
CALL dbcsr_multiply("N", "N", 1.0_dp, p_index_up, pp, 0.0_dp, oo, &
filter_eps=eps_filter)
CALL dbcsr_multiply("N", "N", 1.0_dp, oo, p_index_up, 0.0_dp, pp, &
filter_eps=eps_filter)
CALL dbcsr_release(oo)
END SELECT
CALL dbcsr_release(sp)
CALL dbcsr_release(spf)
END IF
CALL timestop(handle)
END SUBROUTINE assemble_ks_qp_blocks
! **************************************************************************************************
!> \brief Solves the generalized Riccati or Sylvester eqation
!> using the preconditioned conjugate gradient algorithm
!> qp + qq.x.oo - vv.x.pp - vv.x.pq.x.oo = 0 [oo and vv are optional]
!> qp + qq.x - x.pp - x.pq.x = 0
!> \param pp ...
!> \param qq ...
!> \param qp ...
!> \param pq ...
!> \param oo ...
!> \param vv ...
!> \param x ...
!> \param res ...
!> \param neglect_quadratic_term ...
!> \param conjugator ...
!> \param max_iter ...
!> \param eps_convergence ...
!> \param eps_filter ...
!> \param converged ...
!> \par History
!> 2011.06 created [Rustam Z Khaliullin]
!> 2011.11 generalized [Rustam Z Khaliullin]
!> \author Rustam Z Khaliullin
! **************************************************************************************************
RECURSIVE SUBROUTINE solve_riccati_equation(pp, qq, qp, pq, oo, vv, x, res, &
neglect_quadratic_term, &
conjugator, max_iter, eps_convergence, eps_filter, &
converged)
TYPE(dbcsr_type), INTENT(IN) :: pp, qq
TYPE(dbcsr_type), INTENT(INOUT) :: qp
TYPE(dbcsr_type), INTENT(IN) :: pq
TYPE(dbcsr_type), INTENT(IN), OPTIONAL :: oo, vv
TYPE(dbcsr_type), INTENT(INOUT) :: x
TYPE(dbcsr_type), INTENT(OUT) :: res
LOGICAL, INTENT(IN) :: neglect_quadratic_term
INTEGER, INTENT(IN) :: conjugator, max_iter
REAL(KIND=dp), INTENT(IN) :: eps_convergence, eps_filter
LOGICAL, INTENT(OUT) :: converged
CHARACTER(len=*), PARAMETER :: routineN = 'solve_riccati_equation'
INTEGER :: handle, istep, iteration, nsteps, &
unit_nr, update_prec_freq
LOGICAL :: prepare_to_exit, present_oo, present_vv, &
quadratic_term, restart_conjugator
REAL(KIND=dp) :: best_norm, best_step_size, beta, c0, c1, &
c2, c3, denom, kappa, numer, &
obj_function, t1, t2, tau
REAL(KIND=dp), DIMENSION(3) :: step_size
TYPE(cp_logger_type), POINTER :: logger
TYPE(dbcsr_type) :: aux1, aux2, grad, m, n, oo1, oo2, prec, &
res_trial, step, step_oo, vv_step
!TYPE(dbcsr_type) :: qqqq, pppp, zero_pq, zero_qp
CALL timeset(routineN, handle)
logger => cp_get_default_logger()
IF (logger%para_env%is_source()) THEN
unit_nr = cp_logger_get_default_unit_nr(logger, local=.TRUE.)
ELSE
unit_nr = -1
END IF
t1 = m_walltime()
!IF (level.gt.5) THEN
! CPErrorMessage(cp_failure_level,routineP,"recursion level is too high")
! CPPrecondition(.FALSE.,cp_failure_level,routineP,failure)
!ENDIF
!IF (unit_nr>0) THEN
! WRITE(unit_nr,*) &
! "========== LEVEL ",level,"=========="
!ENDIF
!CALL dbcsr_print(qq)
!CALL dbcsr_print(pp)
!CALL dbcsr_print(qp)
!!CALL dbcsr_print(pq)
!IF (unit_nr>0) THEN
! WRITE(unit_nr,*) &
! "====== END LEVEL ",level,"=========="
!ENDIF
quadratic_term = .NOT. neglect_quadratic_term
present_oo = PRESENT(oo)
present_vv = PRESENT(vv)
! create aux1 matrix and init
CALL dbcsr_create(aux1, template=pp)
CALL dbcsr_copy(aux1, pp)
CALL dbcsr_scale(aux1, -1.0_dp)
! create aux2 matrix and init
CALL dbcsr_create(aux2, template=qq)
CALL dbcsr_copy(aux2, qq)
! create the gradient matrix and init
CALL dbcsr_create(grad, template=x)
CALL dbcsr_set(grad, 0.0_dp)
! create a preconditioner
! RZK-warning how to apply it to up_down tensor?
CALL dbcsr_create(prec, template=x)
!CALL create_preconditioner(prec,aux1,aux2,qp,res,tensor_type,eps_filter)
!CALL dbcsr_set(prec,1.0_dp)
! create the step matrix and init
CALL dbcsr_create(step, template=x)
!CALL dbcsr_hadamard_product(prec,grad,step)
!CALL dbcsr_scale(step,-1.0_dp)
CALL dbcsr_create(n, template=x)
CALL dbcsr_create(m, template=x)
CALL dbcsr_create(oo1, template=pp)
CALL dbcsr_create(oo2, template=pp)
CALL dbcsr_create(res_trial, template=res)
CALL dbcsr_create(vv_step, template=res)
CALL dbcsr_create(step_oo, template=res)
! start conjugate gradient iterations
iteration = 0
converged = .FALSE.
prepare_to_exit = .FALSE.
beta = 0.0_dp
best_step_size = 0.0_dp
best_norm = 1.0E+100_dp
!ecorr=0.0_dp
!change_ecorr=0.0_dp
restart_conjugator = .FALSE.
update_prec_freq = 20
DO
! (re)-compute the residuals
IF (iteration .EQ. 0) THEN
CALL dbcsr_copy(res, qp)
IF (present_oo) THEN
CALL dbcsr_multiply("N", "N", +1.0_dp, qq, x, 0.0_dp, res_trial, &
filter_eps=eps_filter)
CALL dbcsr_multiply("N", "N", +1.0_dp, res_trial, oo, 1.0_dp, res, &
filter_eps=eps_filter)
ELSE
CALL dbcsr_multiply("N", "N", +1.0_dp, qq, x, 1.0_dp, res, &
filter_eps=eps_filter)
END IF
IF (present_vv) THEN
CALL dbcsr_multiply("N", "N", -1.0_dp, x, pp, 0.0_dp, res_trial, &
filter_eps=eps_filter)
CALL dbcsr_multiply("N", "N", +1.0_dp, vv, res_trial, 1.0_dp, res, &
filter_eps=eps_filter)
ELSE
CALL dbcsr_multiply("N", "N", -1.0_dp, x, pp, 1.0_dp, res, &
filter_eps=eps_filter)
END IF
IF (quadratic_term) THEN
IF (present_oo) THEN
CALL dbcsr_multiply("N", "N", +1.0_dp, pq, x, 0.0_dp, oo1, &
filter_eps=eps_filter)
CALL dbcsr_multiply("N", "N", +1.0_dp, oo1, oo, 0.0_dp, oo2, &
filter_eps=eps_filter)
ELSE
CALL dbcsr_multiply("N", "N", +1.0_dp, pq, x, 0.0_dp, oo2, &
filter_eps=eps_filter)
END IF
IF (present_vv) THEN
CALL dbcsr_multiply("N", "N", -1.0_dp, x, oo2, 0.0_dp, res_trial, &
filter_eps=eps_filter)
CALL dbcsr_multiply("N", "N", +1.0_dp, vv, res_trial, 1.0_dp, res, &
filter_eps=eps_filter)
ELSE
CALL dbcsr_multiply("N", "N", -1.0_dp, x, oo2, 1.0_dp, res, &
filter_eps=eps_filter)
END IF
END IF
CALL dbcsr_norm(res, dbcsr_norm_maxabsnorm, norm_scalar=best_norm)
ELSE
CALL dbcsr_add(res, m, 1.0_dp, best_step_size)
CALL dbcsr_add(res, n, 1.0_dp, -best_step_size*best_step_size)
CALL dbcsr_filter(res, eps_filter)
END IF
! check convergence and other exit criteria
converged = (best_norm .LT. eps_convergence)
IF (converged .OR. (iteration .GE. max_iter)) THEN
prepare_to_exit = .TRUE.
END IF
IF (.NOT. prepare_to_exit) THEN
! update aux1=-pp-pq.x.oo and aux2=qq-vv.x.pq
IF (quadratic_term) THEN
IF (iteration .EQ. 0) THEN
IF (present_oo) THEN
CALL dbcsr_multiply("N", "N", -1.0_dp, pq, x, 0.0_dp, oo1, &
filter_eps=eps_filter)
CALL dbcsr_multiply("N", "N", +1.0_dp, oo1, oo, 1.0_dp, aux1, &
filter_eps=eps_filter)
ELSE
CALL dbcsr_multiply("N", "N", -1.0_dp, pq, x, 1.0_dp, aux1, &
filter_eps=eps_filter)
END IF
IF (present_vv) THEN
CALL dbcsr_multiply("N", "N", -1.0_dp, vv, x, 0.0_dp, res_trial, &
filter_eps=eps_filter)
CALL dbcsr_multiply("N", "N", +1.0_dp, res_trial, pq, 1.0_dp, aux2, &
filter_eps=eps_filter)
ELSE
CALL dbcsr_multiply("N", "N", -1.0_dp, x, pq, 1.0_dp, aux2, &
filter_eps=eps_filter)
END IF
ELSE
IF (present_oo) THEN
CALL dbcsr_multiply("N", "N", -best_step_size, pq, step_oo, 1.0_dp, aux1, &
filter_eps=eps_filter)
ELSE
CALL dbcsr_multiply("N", "N", -best_step_size, pq, step, 1.0_dp, aux1, &
filter_eps=eps_filter)
END IF
IF (present_vv) THEN
CALL dbcsr_multiply("N", "N", -best_step_size, vv_step, pq, 1.0_dp, aux2, &
filter_eps=eps_filter)
ELSE
CALL dbcsr_multiply("N", "N", -best_step_size, step, pq, 1.0_dp, aux2, &
filter_eps=eps_filter)
END IF
END IF
END IF
! recompute the gradient, do not update it yet
! use m matrix as a temporary storage
! grad=t(vv).res.t(aux1)+t(aux2).res.t(oo)
IF (present_vv) THEN
CALL dbcsr_multiply("N", "T", 1.0_dp, res, aux1, 0.0_dp, res_trial, &
filter_eps=eps_filter)
CALL dbcsr_multiply("T", "N", 1.0_dp, vv, res_trial, 0.0_dp, m, &
filter_eps=eps_filter)
ELSE
CALL dbcsr_multiply("N", "T", 1.0_dp, res, aux1, 0.0_dp, m, &
filter_eps=eps_filter)
END IF
IF (present_oo) THEN
CALL dbcsr_multiply("T", "N", 1.0_dp, aux1, res, 0.0_dp, res_trial, &
filter_eps=eps_filter)
CALL dbcsr_multiply("N", "T", 1.0_dp, res_trial, oo, 1.0_dp, m, &
filter_eps=eps_filter)
ELSE
CALL dbcsr_multiply("T", "N", 1.0_dp, aux2, res, 1.0_dp, m, &
filter_eps=eps_filter)
END IF
! compute preconditioner
!IF (iteration.eq.0.OR.(mod(iteration,update_prec_freq).eq.0)) THEN
IF (iteration .EQ. 0) THEN
CALL create_preconditioner(prec, aux1, aux2, eps_filter)
!restart_conjugator=.TRUE.
!CALL dbcsr_set(prec,1.0_dp)
!CALL dbcsr_print(prec)
END IF
! compute the conjugation coefficient - beta
IF ((iteration .EQ. 0) .OR. restart_conjugator) THEN
beta = 0.0_dp
ELSE
restart_conjugator = .FALSE.
SELECT CASE (conjugator)
CASE (cg_hestenes_stiefel)
CALL dbcsr_add(grad, m, -1.0_dp, 1.0_dp)
CALL dbcsr_hadamard_product(prec, grad, n)
CALL dbcsr_dot(n, m, numer)
CALL dbcsr_dot(grad, step, denom)
beta = numer/denom
CASE (cg_fletcher_reeves)
CALL dbcsr_hadamard_product(prec, grad, n)
CALL dbcsr_dot(grad, n, denom)
CALL dbcsr_hadamard_product(prec, m, n)
CALL dbcsr_dot(m, n, numer)
beta = numer/denom
CASE (cg_polak_ribiere)
CALL dbcsr_hadamard_product(prec, grad, n)
CALL dbcsr_dot(grad, n, denom)
CALL dbcsr_add(grad, m, -1.0_dp, 1.0_dp)
CALL dbcsr_hadamard_product(prec, grad, n)
CALL dbcsr_dot(n, m, numer)
beta = numer/denom
CASE (cg_fletcher)
CALL dbcsr_hadamard_product(prec, m, n)
CALL dbcsr_dot(m, n, numer)
CALL dbcsr_dot(grad, step, denom)
beta = -1.0_dp*numer/denom
CASE (cg_liu_storey)
CALL dbcsr_dot(grad, step, denom)
CALL dbcsr_add(grad, m, -1.0_dp, 1.0_dp)
CALL dbcsr_hadamard_product(prec, grad, n)
CALL dbcsr_dot(n, m, numer)
beta = -1.0_dp*numer/denom
CASE (cg_dai_yuan)
CALL dbcsr_hadamard_product(prec, m, n)
CALL dbcsr_dot(m, n, numer)
CALL dbcsr_add(grad, m, -1.0_dp, 1.0_dp)
CALL dbcsr_dot(grad, step, denom)
beta = numer/denom
CASE (cg_hager_zhang)
CALL dbcsr_add(grad, m, -1.0_dp, 1.0_dp)
CALL dbcsr_dot(grad, step, denom)
CALL dbcsr_hadamard_product(prec, grad, n)
CALL dbcsr_dot(n, grad, numer)
kappa = 2.0_dp*numer/denom
CALL dbcsr_dot(n, m, numer)