-
Notifications
You must be signed in to change notification settings - Fork 1
/
qs_ot.F
1571 lines (1383 loc) · 66.9 KB
/
qs_ot.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 orbital transformations
!> \par History
!> Added Taylor expansion based computation of the matrix functions (01.2004)
!> added additional rotation variables for non-equivalent occupied orbs (08.2004)
!> \author Joost VandeVondele (06.2002)
! **************************************************************************************************
MODULE qs_ot
USE arnoldi_api, ONLY: arnoldi_extremal
USE cp_dbcsr_api, ONLY: &
dbcsr_add, dbcsr_add_on_diag, dbcsr_copy, dbcsr_distribution_type, dbcsr_filter, &
dbcsr_frobenius_norm, dbcsr_gershgorin_norm, dbcsr_get_block_p, dbcsr_get_info, &
dbcsr_get_occupation, dbcsr_hadamard_product, dbcsr_init_p, dbcsr_iterator_blocks_left, &
dbcsr_iterator_next_block, dbcsr_iterator_start, dbcsr_iterator_stop, dbcsr_iterator_type, &
dbcsr_multiply, dbcsr_release, dbcsr_release_p, dbcsr_scale, dbcsr_scale_by_vector, &
dbcsr_set, dbcsr_transposed, dbcsr_type
USE cp_dbcsr_cholesky, ONLY: cp_dbcsr_cholesky_decompose,&
cp_dbcsr_cholesky_invert,&
cp_dbcsr_cholesky_restore
USE cp_dbcsr_diag, ONLY: cp_dbcsr_heevd,&
cp_dbcsr_syevd
USE kinds, ONLY: dp
USE message_passing, ONLY: mp_comm_type
USE preconditioner, ONLY: apply_preconditioner
USE preconditioner_types, ONLY: preconditioner_type
USE qs_ot_types, ONLY: qs_ot_type
#include "./base/base_uses.f90"
IMPLICIT NONE
PRIVATE
PUBLIC :: qs_ot_get_p
PUBLIC :: qs_ot_get_orbitals
PUBLIC :: qs_ot_get_derivative
PUBLIC :: qs_ot_get_orbitals_ref
PUBLIC :: qs_ot_get_derivative_ref
PUBLIC :: qs_ot_new_preconditioner
PRIVATE :: qs_ot_p2m_diag
PRIVATE :: qs_ot_sinc
PRIVATE :: qs_ot_ref_poly
PRIVATE :: qs_ot_ref_chol
PRIVATE :: qs_ot_ref_lwdn
PRIVATE :: qs_ot_ref_decide
PRIVATE :: qs_ot_ref_update
PRIVATE :: qs_ot_refine
PRIVATE :: qs_ot_on_the_fly_localize
CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_ot'
CONTAINS
! gets ready to use the preconditioner/ or renew the preconditioner
! only keeps a pointer to the preconditioner.
! If you change the preconditioner, you have to call this routine
! you remain responsible of proper deallocate of your preconditioner
! (or you can reuse it on the next step of the computation)
! **************************************************************************************************
!> \brief ...
!> \param qs_ot_env ...
!> \param preconditioner ...
! **************************************************************************************************
SUBROUTINE qs_ot_new_preconditioner(qs_ot_env, preconditioner)
TYPE(qs_ot_type) :: qs_ot_env
TYPE(preconditioner_type), POINTER :: preconditioner
INTEGER :: ncoef
qs_ot_env%preconditioner => preconditioner
qs_ot_env%os_valid = .FALSE.
IF (.NOT. ASSOCIATED(qs_ot_env%matrix_psc0)) THEN
CALL dbcsr_init_p(qs_ot_env%matrix_psc0)
CALL dbcsr_copy(qs_ot_env%matrix_psc0, qs_ot_env%matrix_sc0, 'matrix_psc0')
END IF
IF (.NOT. qs_ot_env%use_dx) THEN
qs_ot_env%use_dx = .TRUE.
CALL dbcsr_init_p(qs_ot_env%matrix_dx)
CALL dbcsr_copy(qs_ot_env%matrix_dx, qs_ot_env%matrix_gx, 'matrix_dx')
IF (qs_ot_env%settings%do_rotation) THEN
CALL dbcsr_init_p(qs_ot_env%rot_mat_dx)
CALL dbcsr_copy(qs_ot_env%rot_mat_dx, qs_ot_env%rot_mat_gx, 'rot_mat_dx')
END IF
IF (qs_ot_env%settings%do_ener) THEN
ncoef = SIZE(qs_ot_env%ener_gx)
ALLOCATE (qs_ot_env%ener_dx(ncoef))
qs_ot_env%ener_dx = 0.0_dp
END IF
END IF
END SUBROUTINE qs_ot_new_preconditioner
! **************************************************************************************************
!> \brief ...
!> \param qs_ot_env ...
!> \param C_NEW ...
!> \param SC ...
!> \param G_OLD ...
!> \param D ...
! **************************************************************************************************
SUBROUTINE qs_ot_on_the_fly_localize(qs_ot_env, C_NEW, SC, G_OLD, D)
!
TYPE(qs_ot_type) :: qs_ot_env
TYPE(dbcsr_type), POINTER :: C_NEW, SC, G_OLD, D
CHARACTER(len=*), PARAMETER :: routineN = 'qs_ot_on_the_fly_localize'
INTEGER, PARAMETER :: taylor_order = 50
REAL(KIND=dp), PARAMETER :: alpha = 0.1_dp, f2_eps = 0.01_dp
INTEGER :: blk, col, col_size, group_handle, &
handle, i, k, n, p, row, row_size
REAL(dp), DIMENSION(:, :), POINTER :: DATA
REAL(KIND=dp) :: expfactor, f2, norm_fro, norm_gct, tmp
TYPE(dbcsr_distribution_type) :: dist
TYPE(dbcsr_iterator_type) :: iter
TYPE(dbcsr_type), POINTER :: C, Gp1, Gp2, GU, U
TYPE(mp_comm_type) :: group
CALL timeset(routineN, handle)
!
!
CALL dbcsr_get_info(C_NEW, nfullrows_total=n, nfullcols_total=k)
!
! C = C*expm(-G)
GU => qs_ot_env%buf1_k_k_nosym ! a buffer
U => qs_ot_env%buf2_k_k_nosym ! a buffer
Gp1 => qs_ot_env%buf3_k_k_nosym ! a buffer
Gp2 => qs_ot_env%buf4_k_k_nosym ! a buffer
C => qs_ot_env%buf1_n_k ! a buffer
!
! compute the derivative of the norm
!-------------------------------------------------------------------
! (x^2+eps)^1/2
f2 = 0.0_dp
CALL dbcsr_copy(C, C_NEW)
CALL dbcsr_iterator_start(iter, C)
DO WHILE (dbcsr_iterator_blocks_left(iter))
CALL dbcsr_iterator_next_block(iter, row, col, DATA, blk, &
row_size=row_size, col_size=col_size)
DO p = 1, col_size ! p
DO i = 1, row_size ! i
tmp = SQRT(DATA(i, p)**2 + f2_eps)
f2 = f2 + tmp
DATA(i, p) = DATA(i, p)/tmp
END DO
END DO
END DO
CALL dbcsr_iterator_stop(iter)
CALL dbcsr_get_info(C, group=group_handle)
CALL group%set_handle(group_handle)
CALL group%sum(f2)
!
!
CALL dbcsr_multiply('T', 'N', 1.0_dp, C, C_NEW, 0.0_dp, GU)
!
! antisymetrize
CALL dbcsr_get_info(GU, distribution=dist)
CALL dbcsr_transposed(U, GU, shallow_data_copy=.FALSE., &
use_distribution=dist, &
transpose_distribution=.FALSE.)
CALL dbcsr_add(GU, U, alpha_scalar=-0.5_dp, beta_scalar=0.5_dp)
!-------------------------------------------------------------------
!
norm_fro = dbcsr_frobenius_norm(GU)
norm_gct = dbcsr_gershgorin_norm(GU)
!write(*,*) 'qs_ot_localize: ||P-I||_f=',norm_fro,' ||P-I||_GCT=',norm_gct
!
!kscale = CEILING(LOG(MIN(norm_fro,norm_gct))/LOG(2.0_dp))
!scale = LOG(MIN(norm_fro,norm_gct))/LOG(2.0_dp)
!write(*,*) 'qs_ot_localize: scale=',scale,' kscale=',kscale
!
! rescale for steepest descent
CALL dbcsr_scale(GU, -alpha)
!
! compute unitary transform
! zeroth and first order
expfactor = 1.0_dp
CALL dbcsr_copy(U, GU)
CALL dbcsr_scale(U, expfactor)
CALL dbcsr_add_on_diag(U, 1.0_dp)
! other orders
CALL dbcsr_copy(Gp1, GU)
DO i = 2, taylor_order
! new power of G
CALL dbcsr_multiply('N', 'N', 1.0_dp, GU, Gp1, 0.0_dp, Gp2)
CALL dbcsr_copy(Gp1, Gp2)
! add to the taylor expansion so far
expfactor = expfactor/REAL(i, KIND=dp)
CALL dbcsr_add(U, Gp1, alpha_scalar=1.0_dp, beta_scalar=expfactor)
norm_fro = dbcsr_frobenius_norm(Gp1)
!write(*,*) 'Taylor expansion i=',i,' norm(X^i)/i!=',norm_fro*expfactor
IF (norm_fro*expfactor .LT. 1.0E-10_dp) EXIT
END DO
!
! rotate MOs
CALL dbcsr_multiply('N', 'N', 1.0_dp, C_NEW, U, 0.0_dp, C)
CALL dbcsr_copy(C_NEW, C)
!
! rotate SC
CALL dbcsr_multiply('N', 'N', 1.0_dp, SC, U, 0.0_dp, C)
CALL dbcsr_copy(SC, C)
!
! rotate D_i
CALL dbcsr_multiply('N', 'N', 1.0_dp, D, U, 0.0_dp, C)
CALL dbcsr_copy(D, C)
!
! rotate G_i-1
IF (ASSOCIATED(G_OLD)) THEN
CALL dbcsr_multiply('N', 'N', 1.0_dp, G_OLD, U, 0.0_dp, C)
CALL dbcsr_copy(G_OLD, C)
END IF
!
CALL timestop(handle)
END SUBROUTINE qs_ot_on_the_fly_localize
! **************************************************************************************************
!> \brief ...
!> \param qs_ot_env ...
!> \param C_OLD ...
!> \param C_TMP ...
!> \param C_NEW ...
!> \param P ...
!> \param SC ...
!> \param update ...
! **************************************************************************************************
SUBROUTINE qs_ot_ref_chol(qs_ot_env, C_OLD, C_TMP, C_NEW, P, SC, update)
!
TYPE(qs_ot_type) :: qs_ot_env
TYPE(dbcsr_type) :: C_OLD, C_TMP, C_NEW, P, SC
LOGICAL, INTENT(IN) :: update
CHARACTER(len=*), PARAMETER :: routineN = 'qs_ot_ref_chol'
INTEGER :: handle, k, n
CALL timeset(routineN, handle)
!
CALL dbcsr_get_info(C_NEW, nfullrows_total=n, nfullcols_total=k)
!
! P = U'*U
CALL cp_dbcsr_cholesky_decompose(P, k, qs_ot_env%para_env, qs_ot_env%blacs_env)
!
! C_NEW = C_OLD*inv(U)
CALL cp_dbcsr_cholesky_restore(C_OLD, k, P, C_NEW, op="SOLVE", pos="RIGHT", &
transa="N", para_env=qs_ot_env%para_env, blacs_env=qs_ot_env%blacs_env)
!
! Update SC if needed
IF (update) THEN
CALL cp_dbcsr_cholesky_restore(SC, k, P, C_TMP, op="SOLVE", pos="RIGHT", &
transa="N", para_env=qs_ot_env%para_env, blacs_env=qs_ot_env%blacs_env)
CALL dbcsr_copy(SC, C_TMP)
END IF
!
CALL timestop(handle)
END SUBROUTINE qs_ot_ref_chol
! **************************************************************************************************
!> \brief ...
!> \param qs_ot_env ...
!> \param C_OLD ...
!> \param C_TMP ...
!> \param C_NEW ...
!> \param P ...
!> \param SC ...
!> \param update ...
! **************************************************************************************************
SUBROUTINE qs_ot_ref_lwdn(qs_ot_env, C_OLD, C_TMP, C_NEW, P, SC, update)
!
TYPE(qs_ot_type) :: qs_ot_env
TYPE(dbcsr_type) :: C_OLD, C_TMP, C_NEW, P, SC
LOGICAL, INTENT(IN) :: update
CHARACTER(len=*), PARAMETER :: routineN = 'qs_ot_ref_lwdn'
INTEGER :: handle, i, k, n
REAL(dp), ALLOCATABLE, DIMENSION(:) :: eig, fun
TYPE(dbcsr_type), POINTER :: V, W
CALL timeset(routineN, handle)
!
CALL dbcsr_get_info(C_NEW, nfullrows_total=n, nfullcols_total=k)
!
V => qs_ot_env%buf1_k_k_nosym ! a buffer
W => qs_ot_env%buf2_k_k_nosym ! a buffer
ALLOCATE (eig(k), fun(k))
!
CALL cp_dbcsr_syevd(P, V, eig, qs_ot_env%para_env, qs_ot_env%blacs_env)
!
! compute the P^(-1/2)
DO i = 1, k
IF (eig(i) .LE. 0.0_dp) &
CPABORT("P not positive definite")
IF (eig(i) .LT. 1.0E-8_dp) THEN
fun(i) = 0.0_dp
ELSE
fun(i) = 1.0_dp/SQRT(eig(i))
END IF
END DO
CALL dbcsr_copy(W, V)
CALL dbcsr_scale_by_vector(V, alpha=fun, side='right')
CALL dbcsr_multiply('N', 'T', 1.0_dp, W, V, 0.0_dp, P)
!
! Update C
CALL dbcsr_multiply('N', 'N', 1.0_dp, C_OLD, P, 0.0_dp, C_NEW)
!
! Update SC if needed
IF (update) THEN
CALL dbcsr_multiply('N', 'N', 1.0_dp, SC, P, 0.0_dp, C_TMP)
CALL dbcsr_copy(SC, C_TMP)
END IF
!
DEALLOCATE (eig, fun)
!
CALL timestop(handle)
END SUBROUTINE qs_ot_ref_lwdn
! **************************************************************************************************
!> \brief ...
!> \param qs_ot_env ...
!> \param C_OLD ...
!> \param C_TMP ...
!> \param C_NEW ...
!> \param P ...
!> \param SC ...
!> \param norm_in ...
!> \param update ...
! **************************************************************************************************
SUBROUTINE qs_ot_ref_poly(qs_ot_env, C_OLD, C_TMP, C_NEW, P, SC, norm_in, update)
!
TYPE(qs_ot_type) :: qs_ot_env
TYPE(dbcsr_type), POINTER :: C_OLD, C_TMP, C_NEW, P
TYPE(dbcsr_type) :: SC
REAL(dp), INTENT(IN) :: norm_in
LOGICAL, INTENT(IN) :: update
CHARACTER(len=*), PARAMETER :: routineN = 'qs_ot_ref_poly'
INTEGER :: handle, irefine, k, n
LOGICAL :: quick_exit
REAL(dp) :: norm, norm_fro, norm_gct, occ_in, &
occ_out, rescale
TYPE(dbcsr_type), POINTER :: BUF1, BUF2, BUF_NOSYM, FT, FY
CALL timeset(routineN, handle)
!
CALL dbcsr_get_info(C_NEW, nfullrows_total=n, nfullcols_total=k)
!
BUF_NOSYM => qs_ot_env%buf1_k_k_nosym ! a buffer
BUF1 => qs_ot_env%buf1_k_k_sym ! a buffer
BUF2 => qs_ot_env%buf2_k_k_sym ! a buffer
FY => qs_ot_env%buf3_k_k_sym ! a buffer
FT => qs_ot_env%buf4_k_k_sym ! a buffer
!
! initialize the norm (already computed in qs_ot_get_orbitals_ref)
norm = norm_in
!
! can we do a quick exit?
quick_exit = .FALSE.
IF (norm .LT. qs_ot_env%settings%eps_irac_quick_exit) quick_exit = .TRUE.
!
! lets refine
rescale = 1.0_dp
DO irefine = 1, qs_ot_env%settings%max_irac
!
! rescaling
IF (norm .GT. 1.0_dp) THEN
CALL dbcsr_scale(P, 1.0_dp/norm)
rescale = rescale/SQRT(norm)
END IF
!
! get the refinement polynomial
CALL qs_ot_refine(P, FY, BUF1, BUF2, qs_ot_env%settings%irac_degree, &
qs_ot_env%settings%eps_irac_filter_matrix)
!
! collect the transformation
IF (irefine .EQ. 1) THEN
CALL dbcsr_copy(FT, FY, name='FT')
ELSE
CALL dbcsr_multiply('N', 'N', 1.0_dp, FT, FY, 0.0_dp, BUF1)
IF (qs_ot_env%settings%eps_irac_filter_matrix .GT. 0.0_dp) THEN
occ_in = dbcsr_get_occupation(buf1)
CALL dbcsr_filter(buf1, qs_ot_env%settings%eps_irac_filter_matrix)
occ_out = dbcsr_get_occupation(buf1)
END IF
CALL dbcsr_copy(FT, BUF1, name='FT')
END IF
!
! quick exit if possible
IF (quick_exit) THEN
EXIT
END IF
!
! P = FY^T * P * FY
CALL dbcsr_multiply('N', 'N', 1.0_dp, P, FY, 0.0_dp, BUF_NOSYM)
IF (qs_ot_env%settings%eps_irac_filter_matrix .GT. 0.0_dp) THEN
occ_in = dbcsr_get_occupation(buf_nosym)
CALL dbcsr_filter(buf_nosym, qs_ot_env%settings%eps_irac_filter_matrix)
occ_out = dbcsr_get_occupation(buf_nosym)
END IF
CALL dbcsr_multiply('N', 'N', 1.0_dp, FY, BUF_NOSYM, 0.0_dp, P)
IF (qs_ot_env%settings%eps_irac_filter_matrix .GT. 0.0_dp) THEN
occ_in = dbcsr_get_occupation(p)
CALL dbcsr_filter(p, qs_ot_env%settings%eps_irac_filter_matrix)
occ_out = dbcsr_get_occupation(p)
END IF
!
! check ||P-1||_gct
CALL dbcsr_add_on_diag(P, -1.0_dp)
norm_fro = dbcsr_frobenius_norm(P)
norm_gct = dbcsr_gershgorin_norm(P)
CALL dbcsr_add_on_diag(P, 1.0_dp)
norm = MIN(norm_gct, norm_fro)
!
! printing
!
! blows up
IF (norm > 1.0E10_dp) THEN
CALL cp_abort(__LOCATION__, &
"Refinement blows up! "// &
"We need you to improve the code, please post your input on "// &
"the forum https://www.cp2k.org/")
END IF
!
! can we do a quick exit next step?
IF (norm .LT. qs_ot_env%settings%eps_irac_quick_exit) quick_exit = .TRUE.
!
! are we done?
IF (norm .LT. qs_ot_env%settings%eps_irac) EXIT
!
END DO
!
! C_NEW = C_NEW * FT * rescale
CALL dbcsr_multiply('N', 'N', rescale, C_OLD, FT, 0.0_dp, C_NEW)
IF (qs_ot_env%settings%eps_irac_filter_matrix .GT. 0.0_dp) THEN
occ_in = dbcsr_get_occupation(c_new)
CALL dbcsr_filter(c_new, qs_ot_env%settings%eps_irac_filter_matrix)
occ_out = dbcsr_get_occupation(c_new)
END IF
!
! update SC = SC * FY * rescale
IF (update) THEN
CALL dbcsr_multiply('N', 'N', rescale, SC, FT, 0.0_dp, C_TMP)
IF (qs_ot_env%settings%eps_irac_filter_matrix .GT. 0.0_dp) THEN
occ_in = dbcsr_get_occupation(c_tmp)
CALL dbcsr_filter(c_tmp, qs_ot_env%settings%eps_irac_filter_matrix)
occ_out = dbcsr_get_occupation(c_tmp)
END IF
CALL dbcsr_copy(SC, C_TMP)
END IF
!
CALL timestop(handle)
END SUBROUTINE qs_ot_ref_poly
! **************************************************************************************************
!> \brief ...
!> \param qs_ot_env1 ...
!> \return ...
! **************************************************************************************************
FUNCTION qs_ot_ref_update(qs_ot_env1) RESULT(update)
!
TYPE(qs_ot_type) :: qs_ot_env1
LOGICAL :: update
update = .FALSE.
SELECT CASE (qs_ot_env1%settings%ot_method)
CASE ("CG")
SELECT CASE (qs_ot_env1%settings%line_search_method)
CASE ("2PNT")
IF (qs_ot_env1%line_search_count .EQ. 2) update = .TRUE.
CASE DEFAULT
CPABORT("NYI")
END SELECT
CASE ("DIIS")
update = .TRUE.
CASE DEFAULT
CPABORT("NYI")
END SELECT
END FUNCTION qs_ot_ref_update
! **************************************************************************************************
!> \brief ...
!> \param qs_ot_env1 ...
!> \param norm_in ...
!> \param ortho_irac ...
! **************************************************************************************************
SUBROUTINE qs_ot_ref_decide(qs_ot_env1, norm_in, ortho_irac)
!
TYPE(qs_ot_type) :: qs_ot_env1
REAL(dp), INTENT(IN) :: norm_in
CHARACTER(LEN=*), INTENT(INOUT) :: ortho_irac
ortho_irac = qs_ot_env1%settings%ortho_irac
IF (norm_in .LT. qs_ot_env1%settings%eps_irac_switch) ortho_irac = "POLY"
END SUBROUTINE qs_ot_ref_decide
! **************************************************************************************************
!> \brief ...
!> \param matrix_c ...
!> \param matrix_s ...
!> \param matrix_x ...
!> \param matrix_sx ...
!> \param matrix_gx_old ...
!> \param matrix_dx ...
!> \param qs_ot_env ...
!> \param qs_ot_env1 ...
! **************************************************************************************************
SUBROUTINE qs_ot_get_orbitals_ref(matrix_c, matrix_s, matrix_x, matrix_sx, &
matrix_gx_old, matrix_dx, qs_ot_env, qs_ot_env1)
!
TYPE(dbcsr_type), POINTER :: matrix_c, matrix_s, matrix_x, matrix_sx, &
matrix_gx_old, matrix_dx
TYPE(qs_ot_type) :: qs_ot_env, qs_ot_env1
CHARACTER(len=*), PARAMETER :: routineN = 'qs_ot_get_orbitals_ref'
CHARACTER(LEN=4) :: ortho_irac
INTEGER :: handle, k, n
LOGICAL :: on_the_fly_loc, update
REAL(dp) :: norm, norm_fro, norm_gct, occ_in, occ_out
TYPE(dbcsr_type), POINTER :: C_NEW, C_OLD, C_TMP, D, G_OLD, P, S, SC
CALL timeset(routineN, handle)
CALL dbcsr_get_info(matrix_c, nfullrows_total=n, nfullcols_total=k)
!
C_NEW => matrix_c
C_OLD => matrix_x ! need to be carefully updated for the gradient !
SC => matrix_sx ! need to be carefully updated for the gradient !
G_OLD => matrix_gx_old ! need to be carefully updated for localization !
D => matrix_dx ! need to be carefully updated for localization !
S => matrix_s
P => qs_ot_env%p_k_k_sym ! a buffer
C_TMP => qs_ot_env%buf1_n_k ! a buffer
!
! do we need to update C_OLD and SC?
update = qs_ot_ref_update(qs_ot_env1)
!
! do we want to on the fly localize?
! for the moment this is set from the input,
! later we might want to localize every n-step or
! when the sparsity increases...
on_the_fly_loc = qs_ot_env1%settings%on_the_fly_loc
!
! compute SC = S*C
IF (ASSOCIATED(S)) THEN
CALL dbcsr_multiply('N', 'N', 1.0_dp, S, C_OLD, 0.0_dp, SC)
IF (qs_ot_env1%settings%eps_irac_filter_matrix .GT. 0.0_dp) THEN
occ_in = dbcsr_get_occupation(sc)
CALL dbcsr_filter(sc, qs_ot_env1%settings%eps_irac_filter_matrix)
occ_out = dbcsr_get_occupation(sc)
END IF
ELSE
CALL dbcsr_copy(SC, C_OLD)
END IF
!
! compute P = C'*SC
CALL dbcsr_multiply('T', 'N', 1.0_dp, C_OLD, SC, 0.0_dp, P)
IF (qs_ot_env1%settings%eps_irac_filter_matrix .GT. 0.0_dp) THEN
occ_in = dbcsr_get_occupation(p)
CALL dbcsr_filter(p, qs_ot_env1%settings%eps_irac_filter_matrix)
occ_out = dbcsr_get_occupation(p)
END IF
!
! check ||P-1||_f and ||P-1||_gct
CALL dbcsr_add_on_diag(P, -1.0_dp)
norm_fro = dbcsr_frobenius_norm(P)
norm_gct = dbcsr_gershgorin_norm(P)
CALL dbcsr_add_on_diag(P, 1.0_dp)
norm = MIN(norm_gct, norm_fro)
CALL qs_ot_ref_decide(qs_ot_env1, norm, ortho_irac)
!
! select the orthogonality method
SELECT CASE (ortho_irac)
CASE ("CHOL")
CALL qs_ot_ref_chol(qs_ot_env, C_OLD, C_TMP, C_NEW, P, SC, update)
CASE ("LWDN")
CALL qs_ot_ref_lwdn(qs_ot_env, C_OLD, C_TMP, C_NEW, P, SC, update)
CASE ("POLY")
CALL qs_ot_ref_poly(qs_ot_env, C_OLD, C_TMP, C_NEW, P, SC, norm, update)
CASE DEFAULT
CPABORT("Wrong argument")
END SELECT
!
! We update the C_i+1 and localization
IF (update) THEN
IF (on_the_fly_loc) THEN
CALL qs_ot_on_the_fly_localize(qs_ot_env, C_NEW, SC, G_OLD, D)
END IF
CALL dbcsr_copy(C_OLD, C_NEW)
END IF
!
CALL timestop(handle)
END SUBROUTINE qs_ot_get_orbitals_ref
! **************************************************************************************************
!> \brief ...
!> \param P ...
!> \param FY ...
!> \param P2 ...
!> \param T ...
!> \param irac_degree ...
!> \param eps_irac_filter_matrix ...
! **************************************************************************************************
SUBROUTINE qs_ot_refine(P, FY, P2, T, irac_degree, eps_irac_filter_matrix)
!----------------------------------------------------------------------
! refinement polynomial of degree 2,3 and 4 (PRB 70, 193102 (2004))
!----------------------------------------------------------------------
TYPE(dbcsr_type), INTENT(inout) :: P, FY, P2, T
INTEGER, INTENT(in) :: irac_degree
REAL(dp), INTENT(in) :: eps_irac_filter_matrix
CHARACTER(len=*), PARAMETER :: routineN = 'qs_ot_refine'
INTEGER :: handle, k
REAL(dp) :: occ_in, occ_out, r
CALL timeset(routineN, handle)
CALL dbcsr_get_info(P, nfullcols_total=k)
SELECT CASE (irac_degree)
CASE (2)
! C_out = C_in * ( 15/8 * I - 10/8 * P + 3/8 * P^2)
r = 3.0_dp/8.0_dp
CALL dbcsr_multiply('N', 'N', r, P, P, 0.0_dp, FY)
IF (eps_irac_filter_matrix .GT. 0.0_dp) THEN
occ_in = dbcsr_get_occupation(fy)
CALL dbcsr_filter(fy, eps_irac_filter_matrix)
occ_out = dbcsr_get_occupation(fy)
END IF
r = -10.0_dp/8.0_dp
CALL dbcsr_add(FY, P, alpha_scalar=1.0_dp, beta_scalar=r)
r = 15.0_dp/8.0_dp
CALL dbcsr_add_on_diag(FY, alpha_scalar=r)
CASE (3)
! C_out = C_in * ( 35/16 * I - 35/16 * P + 21/16 * P^2 - 5/16 P^3)
CALL dbcsr_multiply('N', 'N', 1.0_dp, P, P, 0.0_dp, P2)
IF (eps_irac_filter_matrix .GT. 0.0_dp) THEN
occ_in = dbcsr_get_occupation(p2)
CALL dbcsr_filter(p2, eps_irac_filter_matrix)
occ_out = dbcsr_get_occupation(p2)
END IF
r = -5.0_dp/16.0_dp
CALL dbcsr_multiply('N', 'N', r, P2, P, 0.0_dp, FY)
IF (eps_irac_filter_matrix .GT. 0.0_dp) THEN
occ_in = dbcsr_get_occupation(fy)
CALL dbcsr_filter(fy, eps_irac_filter_matrix)
occ_out = dbcsr_get_occupation(fy)
END IF
r = 21.0_dp/16.0_dp
CALL dbcsr_add(FY, P2, alpha_scalar=1.0_dp, beta_scalar=r)
r = -35.0_dp/16.0_dp
CALL dbcsr_add(FY, P, alpha_scalar=1.0_dp, beta_scalar=r)
r = 35.0_dp/16.0_dp
CALL dbcsr_add_on_diag(FY, alpha_scalar=r)
CASE (4)
! C_out = C_in * ( 315/128 * I - 420/128 * P + 378/128 * P^2 - 180/128 P^3 + 35/128 P^4 )
! = C_in * ( 315/128 * I - 420/128 * P + 378/128 * P^2 + ( - 180/128 * P + 35/128 * P^2 ) * P^2 )
CALL dbcsr_multiply('N', 'N', 1.0_dp, P, P, 0.0_dp, P2) ! P^2
IF (eps_irac_filter_matrix .GT. 0.0_dp) THEN
occ_in = dbcsr_get_occupation(p2)
CALL dbcsr_filter(p2, eps_irac_filter_matrix)
occ_out = dbcsr_get_occupation(p2)
END IF
r = -180.0_dp/128.0_dp
CALL dbcsr_add(T, P, alpha_scalar=0.0_dp, beta_scalar=r) ! T=-180/128*P
r = 35.0_dp/128.0_dp
CALL dbcsr_add(T, P2, alpha_scalar=1.0_dp, beta_scalar=r) ! T=T+35/128*P^2
CALL dbcsr_multiply('N', 'N', 1.0_dp, T, P2, 0.0_dp, FY) ! Y=T*P^2
IF (eps_irac_filter_matrix .GT. 0.0_dp) THEN
occ_in = dbcsr_get_occupation(fy)
CALL dbcsr_filter(fy, eps_irac_filter_matrix)
occ_out = dbcsr_get_occupation(fy)
END IF
r = 378.0_dp/128.0_dp
CALL dbcsr_add(FY, P2, alpha_scalar=1.0_dp, beta_scalar=r) ! Y=Y+378/128*P^2
r = -420.0_dp/128.0_dp
CALL dbcsr_add(FY, P, alpha_scalar=1.0_dp, beta_scalar=r) ! Y=Y-420/128*P
r = 315.0_dp/128.0_dp
CALL dbcsr_add_on_diag(FY, alpha_scalar=r) ! Y=Y+315/128*I
CASE DEFAULT
CPABORT("This irac_order NYI")
END SELECT
CALL timestop(handle)
END SUBROUTINE qs_ot_refine
! **************************************************************************************************
!> \brief ...
!> \param matrix_hc ...
!> \param matrix_x ...
!> \param matrix_sx ...
!> \param matrix_gx ...
!> \param qs_ot_env ...
! **************************************************************************************************
SUBROUTINE qs_ot_get_derivative_ref(matrix_hc, matrix_x, matrix_sx, matrix_gx, &
qs_ot_env)
TYPE(dbcsr_type), POINTER :: matrix_hc, matrix_x, matrix_sx, matrix_gx
TYPE(qs_ot_type) :: qs_ot_env
CHARACTER(len=*), PARAMETER :: routineN = 'qs_ot_get_derivative_ref'
INTEGER :: handle, k, n
REAL(dp) :: occ_in, occ_out
TYPE(dbcsr_type), POINTER :: C, CHC, G, G_dp, HC, SC
CALL timeset(routineN, handle)
CALL dbcsr_get_info(matrix_x, nfullrows_total=n, nfullcols_total=k)
!
C => matrix_x ! NBsf*NOcc
SC => matrix_sx ! NBsf*NOcc need to be up2date
HC => matrix_hc ! NBsf*NOcc
G => matrix_gx ! NBsf*NOcc
CHC => qs_ot_env%buf1_k_k_sym ! buffer
G_dp => qs_ot_env%buf1_n_k_dp ! buffer dp
! C'*(H*C)
CALL dbcsr_multiply('T', 'N', 1.0_dp, C, HC, 0.0_dp, CHC)
IF (qs_ot_env%settings%eps_irac_filter_matrix .GT. 0.0_dp) THEN
occ_in = dbcsr_get_occupation(chc)
CALL dbcsr_filter(chc, qs_ot_env%settings%eps_irac_filter_matrix)
occ_out = dbcsr_get_occupation(chc)
END IF
! (S*C)*(C'*H*C)
CALL dbcsr_multiply('N', 'N', 1.0_dp, SC, CHC, 0.0_dp, G)
IF (qs_ot_env%settings%eps_irac_filter_matrix .GT. 0.0_dp) THEN
occ_in = dbcsr_get_occupation(g)
CALL dbcsr_filter(g, qs_ot_env%settings%eps_irac_filter_matrix)
occ_out = dbcsr_get_occupation(g)
END IF
! G = 2*(1-S*C*C')*H*C
CALL dbcsr_add(G, HC, alpha_scalar=-1.0_dp, beta_scalar=1.0_dp)
!
CALL timestop(handle)
END SUBROUTINE qs_ot_get_derivative_ref
! computes p=x*S*x and the matrix functionals related matrices
! **************************************************************************************************
!> \brief ...
!> \param matrix_x ...
!> \param matrix_sx ...
!> \param qs_ot_env ...
! **************************************************************************************************
SUBROUTINE qs_ot_get_p(matrix_x, matrix_sx, qs_ot_env)
TYPE(dbcsr_type), POINTER :: matrix_x, matrix_sx
TYPE(qs_ot_type) :: qs_ot_env
CHARACTER(len=*), PARAMETER :: routineN = 'qs_ot_get_p'
REAL(KIND=dp), PARAMETER :: rone = 1.0_dp, rzero = 0.0_dp
INTEGER :: handle, k, max_iter, n
LOGICAL :: converged
REAL(KIND=dp) :: max_ev, min_ev, threshold
CALL timeset(routineN, handle)
CALL dbcsr_get_info(matrix_x, nfullrows_total=n, nfullcols_total=k)
! get the overlap
CALL dbcsr_multiply('T', 'N', rone, matrix_x, matrix_sx, rzero, &
qs_ot_env%matrix_p)
! get an upper bound for the largest eigenvalue
! try using lancos first and fall back to gershgorin norm if it fails
max_iter = 30; threshold = 1.0E-03_dp
CALL arnoldi_extremal(qs_ot_env%matrix_p, max_ev, min_ev, converged, threshold, max_iter)
qs_ot_env%largest_eval_upper_bound = MAX(max_ev, ABS(min_ev))
IF (.NOT. converged) qs_ot_env%largest_eval_upper_bound = dbcsr_gershgorin_norm(qs_ot_env%matrix_p)
CALL decide_strategy(qs_ot_env)
IF (qs_ot_env%do_taylor) THEN
CALL qs_ot_p2m_taylor(qs_ot_env)
ELSE
CALL qs_ot_p2m_diag(qs_ot_env)
END IF
IF (qs_ot_env%settings%do_rotation) THEN
CALL qs_ot_generate_rotation(qs_ot_env)
END IF
CALL timestop(handle)
END SUBROUTINE qs_ot_get_p
! **************************************************************************************************
!> \brief computes the rotation matrix rot_mat_u that is associated to a given
!> rot_mat_x using rot_mat_u=exp(rot_mat_x)
!> \param qs_ot_env a valid qs_ot_env
!> \par History
!> 08.2004 created [Joost VandeVondele]
! **************************************************************************************************
SUBROUTINE qs_ot_generate_rotation(qs_ot_env)
TYPE(qs_ot_type) :: qs_ot_env
CHARACTER(len=*), PARAMETER :: routineN = 'qs_ot_generate_rotation'
COMPLEX(KIND=dp), PARAMETER :: cone = (1.0_dp, 0.0_dp), &
czero = (0.0_dp, 0.0_dp)
COMPLEX(KIND=dp), ALLOCATABLE, DIMENSION(:) :: evals_exp
COMPLEX(KIND=dp), DIMENSION(:), POINTER :: data_z
INTEGER :: blk, col, handle, k, row
LOGICAL :: found
REAL(KIND=dp), DIMENSION(:), POINTER :: data_d
TYPE(dbcsr_iterator_type) :: iter
TYPE(dbcsr_type) :: cmat_u, cmat_x
CALL timeset(routineN, handle)
CALL dbcsr_get_info(qs_ot_env%rot_mat_x, nfullrows_total=k)
IF (k /= 0) THEN
CALL dbcsr_copy(cmat_x, qs_ot_env%rot_mat_evec, name='cmat_x')
CALL dbcsr_copy(cmat_u, qs_ot_env%rot_mat_evec, name='cmat_u')
ALLOCATE (evals_exp(k))
! rot_mat_u = exp(rot_mat_x)
! i rot_mat_x is hermitian, so go over the complex variables for diag
!vwCALL cp_cfm_get_info(cmat_x,local_data=local_data_c)
!vwCALL cp_fm_get_info(qs_ot_env%rot_mat_x,local_data=local_data_r)
!vwlocal_data_c=CMPLX(0.0_dp,local_data_r,KIND=dp)
CALL dbcsr_iterator_start(iter, cmat_x)
DO WHILE (dbcsr_iterator_blocks_left(iter))
CALL dbcsr_iterator_next_block(iter, row, col, data_z, blk)
CALL dbcsr_get_block_p(qs_ot_env%rot_mat_x, row, col, data_d, found)
IF (.NOT. found) THEN
WRITE (*, *) routineN//' .NOT.found'
!stop
ELSE
data_z = CMPLX(0.0_dp, data_d, KIND=dp)
END IF
END DO
CALL dbcsr_iterator_stop(iter)
CALL cp_dbcsr_heevd(cmat_x, qs_ot_env%rot_mat_evec, qs_ot_env%rot_mat_evals, &
qs_ot_env%para_env, qs_ot_env%blacs_env)
evals_exp(:) = EXP((0.0_dp, -1.0_dp)*qs_ot_env%rot_mat_evals(:))
CALL dbcsr_copy(cmat_x, qs_ot_env%rot_mat_evec)
CALL dbcsr_scale_by_vector(cmat_x, alpha=evals_exp, side='right')
CALL dbcsr_multiply('N', 'C', cone, cmat_x, qs_ot_env%rot_mat_evec, czero, cmat_u)
CALL dbcsr_copy(qs_ot_env%rot_mat_u, cmat_u, keep_imaginary=.FALSE.)
CALL dbcsr_release(cmat_x)
CALL dbcsr_release(cmat_u)
DEALLOCATE (evals_exp)
END IF
CALL timestop(handle)
END SUBROUTINE qs_ot_generate_rotation
! **************************************************************************************************
!> \brief computes the derivative fields with respect to rot_mat_x
!> \param qs_ot_env valid qs_ot_env. In particular qs_ot_generate_rotation has to be called before
!> and the rot_mat_dedu matrix has to be up to date
!> \par History
!> 08.2004 created [ Joost VandeVondele ]
! **************************************************************************************************
SUBROUTINE qs_ot_rot_mat_derivative(qs_ot_env)
TYPE(qs_ot_type) :: qs_ot_env
CHARACTER(len=*), PARAMETER :: routineN = 'qs_ot_rot_mat_derivative'
COMPLEX(KIND=dp), PARAMETER :: cone = (1.0_dp, 0.0_dp), &
czero = (0.0_dp, 0.0_dp)
INTEGER :: handle, i, j, k
REAL(KIND=dp) :: e1, e2
REAL(KIND=dp), DIMENSION(:, :), POINTER :: data_d
TYPE(dbcsr_type) :: cmat_buf1, cmat_buf2
TYPE(dbcsr_iterator_type) :: iter
COMPLEX(dp), DIMENSION(:, :), POINTER :: data_z
INTEGER::row, col, blk, row_offset, col_offset, row_size, col_size
LOGICAL::found
TYPE(dbcsr_distribution_type) :: dist
CALL timeset(routineN, handle)
CALL dbcsr_get_info(qs_ot_env%rot_mat_u, nfullrows_total=k)
IF (k /= 0) THEN
CALL dbcsr_copy(qs_ot_env%matrix_buf1, qs_ot_env%rot_mat_dedu)
! now we get to the derivative wrt the antisymmetric matrix rot_mat_x
CALL dbcsr_copy(cmat_buf1, qs_ot_env%rot_mat_evec, "cmat_buf1")
CALL dbcsr_copy(cmat_buf2, qs_ot_env%rot_mat_evec, "cmat_buf2")
! init cmat_buf1
!CALL cp_fm_get_info(qs_ot_env%matrix_buf1,matrix_struct=fm_struct, local_data=local_data_r)
!CALL cp_cfm_get_info(cmat_buf1, nrow_local=nrow_local, ncol_local=ncol_local, &
! row_indices=row_indices, col_indices=col_indices, &
! local_data=local_data_c)
!local_data_c=local_data_r
CALL dbcsr_iterator_start(iter, cmat_buf1)
DO WHILE (dbcsr_iterator_blocks_left(iter))
CALL dbcsr_iterator_next_block(iter, row, col, data_z, blk)
CALL dbcsr_get_block_p(qs_ot_env%matrix_buf1, row, col, data_d, found)
data_z = data_d
END DO
CALL dbcsr_iterator_stop(iter)
CALL dbcsr_multiply('T', 'N', cone, cmat_buf1, qs_ot_env%rot_mat_evec, &
czero, cmat_buf2)
CALL dbcsr_multiply('C', 'N', cone, qs_ot_env%rot_mat_evec, cmat_buf2, &
czero, cmat_buf1)
CALL dbcsr_iterator_start(iter, cmat_buf1)
DO WHILE (dbcsr_iterator_blocks_left(iter))
CALL dbcsr_iterator_next_block(iter, row, col, data_z, blk, &
row_size=row_size, col_size=col_size, &
row_offset=row_offset, col_offset=col_offset)
DO j = 1, col_size
DO i = 1, row_size
e1 = qs_ot_env%rot_mat_evals(row_offset + i - 1)
e2 = qs_ot_env%rot_mat_evals(col_offset + j - 1)
data_z(i, j) = data_z(i, j)*cint(e1, e2)
END DO
END DO
END DO
CALL dbcsr_iterator_stop(iter)
CALL dbcsr_multiply('N', 'N', cone, qs_ot_env%rot_mat_evec, cmat_buf1, &
czero, cmat_buf2)
CALL dbcsr_multiply('N', 'C', cone, cmat_buf2, qs_ot_env%rot_mat_evec, &
czero, cmat_buf1)
CALL dbcsr_copy(qs_ot_env%matrix_buf1, cmat_buf1)
CALL dbcsr_get_info(qs_ot_env%matrix_buf3, distribution=dist)
CALL dbcsr_transposed(qs_ot_env%matrix_buf2, qs_ot_env%matrix_buf1, &
shallow_data_copy=.FALSE., use_distribution=dist, &
transpose_distribution=.FALSE.)
CALL dbcsr_add(qs_ot_env%matrix_buf1, qs_ot_env%matrix_buf2, &
alpha_scalar=-1.0_dp, beta_scalar=+1.0_dp)
CALL dbcsr_copy(qs_ot_env%rot_mat_gx, qs_ot_env%matrix_buf1)
CALL dbcsr_release(cmat_buf1)
CALL dbcsr_release(cmat_buf2)
END IF
CALL timestop(handle)
CONTAINS
! **************************************************************************************************
!> \brief ...
!> \param e1 ...
!> \param e2 ...
!> \return ...
! **************************************************************************************************
FUNCTION cint(e1, e2)
REAL(KIND=dp) :: e1, e2
COMPLEX(KIND=dp) :: cint
COMPLEX(KIND=dp) :: l1, l2, x
INTEGER :: I
l1 = (0.0_dp, -1.0_dp)*e1
l2 = (0.0_dp, -1.0_dp)*e2
IF (ABS(l1 - l2) .GT. 0.5_dp) THEN
cint = (EXP(l1) - EXP(l2))/(l1 - l2)
ELSE
x = 1.0_dp
cint = 0.0_dp
DO I = 1, 16
cint = cint + x
x = x*(l1 - l2)/REAL(I + 1, KIND=dp)
END DO
cint = cint*EXP(l2)
END IF
END FUNCTION cint
END SUBROUTINE qs_ot_rot_mat_derivative
!
! decide strategy
! tries to decide if the taylor expansion of cos(sqrt(xsx)) converges rapidly enough
! to make a taylor expansion of the functions cos(sqrt(xsx)) and sin(sqrt(xsx))/sqrt(xsx)
! and their derivatives faster than their computation based on diagonalization
! since xsx can be very small, especially during dynamics, only a few terms might indeed be needed
! we find the necessary order N to have largest_eval_upper_bound**(N+1)/(2(N+1))! < eps_taylor
!
! **************************************************************************************************
!> \brief ...
!> \param qs_ot_env ...
! **************************************************************************************************
SUBROUTINE decide_strategy(qs_ot_env)
TYPE(qs_ot_type) :: qs_ot_env
INTEGER :: N
REAL(KIND=dp) :: num_error
qs_ot_env%do_taylor = .FALSE.
N = 0
num_error = qs_ot_env%largest_eval_upper_bound/(2.0_dp)
DO WHILE (num_error > qs_ot_env%settings%eps_taylor .AND. N <= qs_ot_env%settings%max_taylor)
N = N + 1
num_error = num_error*qs_ot_env%largest_eval_upper_bound/REAL((2*N + 1)*(2*N + 2), KIND=dp)
END DO