-
Notifications
You must be signed in to change notification settings - Fork 34
/
linearAlgebraOperations.h
949 lines (900 loc) · 37.1 KB
/
linearAlgebraOperations.h
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
// ---------------------------------------------------------------------
//
// Copyright (c) 2017-2022 The Regents of the University of Michigan and DFT-FE
// authors.
//
// This file is part of the DFT-FE code.
//
// The DFT-FE code is free software; you can use it, redistribute
// it, and/or modify it under the terms of the GNU Lesser General
// Public License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// The full text of the license can be found in the file LICENSE at
// the top level of the DFT-FE distribution.
//
// ---------------------------------------------------------------------
#ifndef linearAlgebraOperations_h
#define linearAlgebraOperations_h
#include <elpaScalaManager.h>
#include <headers.h>
#include <operator.h>
#include "process_grid.h"
#include "scalapackWrapper.h"
#include "dftParameters.h"
#include <BLASWrapper.h>
namespace dftfe
{
//
// extern declarations for blas-lapack routines
//
#ifndef DOXYGEN_SHOULD_SKIP_THIS
extern "C"
{
void
dgemv_(const char * TRANS,
const unsigned int *M,
const unsigned int *N,
const double * alpha,
const double * A,
const unsigned int *LDA,
const double * X,
const unsigned int *INCX,
const double * beta,
double * C,
const unsigned int *INCY);
void
sgemv_(const char * TRANS,
const unsigned int *M,
const unsigned int *N,
const float * alpha,
const float * A,
const unsigned int *LDA,
const float * X,
const unsigned int *INCX,
const float * beta,
float * C,
const unsigned int *INCY);
void
zgemv_(const char * TRANS,
const unsigned int * M,
const unsigned int * N,
const std::complex<double> *alpha,
const std::complex<double> *A,
const unsigned int * LDA,
const std::complex<double> *X,
const unsigned int * INCX,
const std::complex<double> *beta,
std::complex<double> * C,
const unsigned int * INCY);
void
cgemv_(const char * TRANS,
const unsigned int * M,
const unsigned int * N,
const std::complex<float> *alpha,
const std::complex<float> *A,
const unsigned int * LDA,
const std::complex<float> *X,
const unsigned int * INCX,
const std::complex<float> *beta,
std::complex<float> * C,
const unsigned int * INCY);
void
dsymv_(const char * UPLO,
const unsigned int *N,
const double * alpha,
const double * A,
const unsigned int *LDA,
const double * X,
const unsigned int *INCX,
const double * beta,
double * C,
const unsigned int *INCY);
void
dgesv_(int * n,
int * nrhs,
double *a,
int * lda,
int * ipiv,
double *b,
int * ldb,
int * info);
void
dsysv_(const char *UPLO,
const int * n,
const int * nrhs,
double * a,
const int * lda,
int * ipiv,
double * b,
const int * ldb,
double * work,
const int * lwork,
int * info);
void
dscal_(const unsigned int *n,
const double * alpha,
double * x,
const unsigned int *inc);
void
sscal_(const unsigned int *n,
const float * alpha,
float * x,
const unsigned int *inc);
void
zscal_(const unsigned int * n,
const std::complex<double> *alpha,
std::complex<double> * x,
const unsigned int * inc);
void
zdscal_(const unsigned int * n,
const double * alpha,
std::complex<double> *x,
const unsigned int * inc);
void
daxpy_(const unsigned int *n,
const double * alpha,
const double * x,
const unsigned int *incx,
double * y,
const unsigned int *incy);
void
saxpy_(const unsigned int *n,
const float * alpha,
const float * x,
const unsigned int *incx,
float * y,
const unsigned int *incy);
void
dgemm_(const char * transA,
const char * transB,
const unsigned int *m,
const unsigned int *n,
const unsigned int *k,
const double * alpha,
const double * A,
const unsigned int *lda,
const double * B,
const unsigned int *ldb,
const double * beta,
double * C,
const unsigned int *ldc);
void
sgemm_(const char * transA,
const char * transB,
const unsigned int *m,
const unsigned int *n,
const unsigned int *k,
const float * alpha,
const float * A,
const unsigned int *lda,
const float * B,
const unsigned int *ldb,
const float * beta,
float * C,
const unsigned int *ldc);
void
dsyevd_(const char * jobz,
const char * uplo,
const unsigned int *n,
double * A,
const unsigned int *lda,
double * w,
double * work,
const unsigned int *lwork,
int * iwork,
const unsigned int *liwork,
int * info);
void
dsygvx_(const int * itype,
const char * jobz,
const char * range,
const char * uplo,
const int * n,
double * a,
const int * lda,
double * b,
const int * ldb,
const double *vl,
const double *vu,
const int * il,
const int * iu,
const double *abstol,
int * m,
double * w,
double * z,
const int * ldz,
double * work,
const int * lwork,
int * iwork,
int * ifail,
int * info);
void
dsyevx_(const char * jobz,
const char * range,
const char * uplo,
const int * n,
double * a,
const int * lda,
const double *vl,
const double *vu,
const int * il,
const int * iu,
const double *abstol,
int * m,
double * w,
double * z,
const int * ldz,
double * work,
const int * lwork,
int * iwork,
int * ifail,
int * info);
double
dlamch_(const char *cmach);
void
dsyevr_(const char * jobz,
const char * range,
const char * uplo,
const unsigned int *n,
double * A,
const unsigned int *lda,
const double * vl,
const double * vu,
const unsigned int *il,
const unsigned int *iu,
const double * abstol,
const unsigned int *m,
double * w,
double * Z,
const unsigned int *ldz,
unsigned int * isuppz,
double * work,
const int * lwork,
int * iwork,
const int * liwork,
int * info);
void
dsyrk_(const char * uplo,
const char * trans,
const unsigned int *n,
const unsigned int *k,
const double * alpha,
const double * A,
const unsigned int *lda,
const double * beta,
double * C,
const unsigned int *ldc);
void
dsyr_(const char * uplo,
const unsigned int *n,
const double * alpha,
const double * X,
const unsigned int *incx,
double * A,
const unsigned int *lda);
void
dsyr2_(const char * uplo,
const unsigned int *n,
const double * alpha,
const double * x,
const unsigned int *incx,
const double * y,
const unsigned int *incy,
double * a,
const unsigned int *lda);
void
dcopy_(const unsigned int *n,
const double * x,
const unsigned int *incx,
double * y,
const unsigned int *incy);
void
scopy_(const unsigned int *n,
const float * x,
const unsigned int *incx,
float * y,
const unsigned int *incy);
void
zgemm_(const char * transA,
const char * transB,
const unsigned int * m,
const unsigned int * n,
const unsigned int * k,
const std::complex<double> *alpha,
const std::complex<double> *A,
const unsigned int * lda,
const std::complex<double> *B,
const unsigned int * ldb,
const std::complex<double> *beta,
std::complex<double> * C,
const unsigned int * ldc);
void
cgemm_(const char * transA,
const char * transB,
const unsigned int * m,
const unsigned int * n,
const unsigned int * k,
const std::complex<float> *alpha,
const std::complex<float> *A,
const unsigned int * lda,
const std::complex<float> *B,
const unsigned int * ldb,
const std::complex<float> *beta,
std::complex<float> * C,
const unsigned int * ldc);
void
zheevd_(const char * jobz,
const char * uplo,
const unsigned int * n,
std::complex<double> *A,
const unsigned int * lda,
double * w,
std::complex<double> *work,
const unsigned int * lwork,
double * rwork,
const unsigned int * lrwork,
int * iwork,
const unsigned int * liwork,
int * info);
void
zheevr_(const char * jobz,
const char * range,
const char * uplo,
const unsigned int * n,
std::complex<double> *A,
const unsigned int * lda,
const double * vl,
const double * vu,
const unsigned int * il,
const unsigned int * iu,
const double * abstol,
const unsigned int * m,
double * w,
std::complex<double> *Z,
const unsigned int * ldz,
unsigned int * isuppz,
std::complex<double> *work,
const int * lwork,
double * rwork,
const int * lrwork,
int * iwork,
const int * liwork,
int * info);
void
zherk_(const char * uplo,
const char * trans,
const unsigned int * n,
const unsigned int * k,
const double * alpha,
const std::complex<double> *A,
const unsigned int * lda,
const double * beta,
std::complex<double> * C,
const unsigned int * ldc);
void
zcopy_(const unsigned int * n,
const std::complex<double> *x,
const unsigned int * incx,
std::complex<double> * y,
const unsigned int * incy);
void
ccopy_(const unsigned int * n,
const std::complex<float> *x,
const unsigned int * incx,
std::complex<float> * y,
const unsigned int * incy);
std::complex<double>
zdotc_(const unsigned int * N,
const std::complex<double> *X,
const unsigned int * INCX,
const std::complex<double> *Y,
const unsigned int * INCY);
double
ddot_(const unsigned int *N,
const double * X,
const unsigned int *INCX,
const double * Y,
const unsigned int *INCY);
double
dnrm2_(const unsigned int *n, const double *x, const unsigned int *incx);
double
dznrm2_(const unsigned int * n,
const std::complex<double> *x,
const unsigned int * incx);
void
zaxpy_(const unsigned int * n,
const std::complex<double> *alpha,
const std::complex<double> *x,
const unsigned int * incx,
std::complex<double> * y,
const unsigned int * incy);
void
caxpy_(const unsigned int * n,
const std::complex<float> *alpha,
const std::complex<float> *x,
const unsigned int * incx,
std::complex<float> * y,
const unsigned int * incy);
void
dpotrf_(const char * uplo,
const unsigned int *n,
double * a,
const unsigned int *lda,
int * info);
void
dpotri_(const char * uplo,
const unsigned int *n,
double * A,
const unsigned int *lda,
int * info);
void
zpotrf_(const char * uplo,
const unsigned int * n,
std::complex<double> *a,
const unsigned int * lda,
int * info);
void
dtrtri_(const char * uplo,
const char * diag,
const unsigned int *n,
double * a,
const unsigned int *lda,
int * info);
void
ztrtri_(const char * uplo,
const char * diag,
const unsigned int * n,
std::complex<double> *a,
const unsigned int * lda,
int * info);
// LU decomoposition of a general matrix
void
dgetrf_(int *M, int *N, double *A, int *lda, int *IPIV, int *INFO);
// generate inverse of a matrix given its LU decomposition
void
dgetri_(int * N,
double *A,
int * lda,
int * IPIV,
double *WORK,
int * lwork,
int * INFO);
}
#endif
inline void
xgemm(const char * transA,
const char * transB,
const unsigned int *m,
const unsigned int *n,
const unsigned int *k,
const double * alpha,
const double * A,
const unsigned int *lda,
const double * B,
const unsigned int *ldb,
const double * beta,
double * C,
const unsigned int *ldc)
{
dgemm_(transA, transB, m, n, k, alpha, A, lda, B, ldb, beta, C, ldc);
}
inline void
xgemm(const char * transA,
const char * transB,
const unsigned int *m,
const unsigned int *n,
const unsigned int *k,
const float * alpha,
const float * A,
const unsigned int *lda,
const float * B,
const unsigned int *ldb,
const float * beta,
float * C,
const unsigned int *ldc)
{
sgemm_(transA, transB, m, n, k, alpha, A, lda, B, ldb, beta, C, ldc);
}
inline void
xgemm(const char * transA,
const char * transB,
const unsigned int * m,
const unsigned int * n,
const unsigned int * k,
const std::complex<double> *alpha,
const std::complex<double> *A,
const unsigned int * lda,
const std::complex<double> *B,
const unsigned int * ldb,
const std::complex<double> *beta,
std::complex<double> * C,
const unsigned int * ldc)
{
zgemm_(transA, transB, m, n, k, alpha, A, lda, B, ldb, beta, C, ldc);
}
inline void
xgemm(const char * transA,
const char * transB,
const unsigned int * m,
const unsigned int * n,
const unsigned int * k,
const std::complex<float> *alpha,
const std::complex<float> *A,
const unsigned int * lda,
const std::complex<float> *B,
const unsigned int * ldb,
const std::complex<float> *beta,
std::complex<float> * C,
const unsigned int * ldc)
{
cgemm_(transA, transB, m, n, k, alpha, A, lda, B, ldb, beta, C, ldc);
}
inline void
xscal(const unsigned int *n,
const double * alpha,
double * x,
const unsigned int *inc)
{
dscal_(n, alpha, x, inc);
}
inline void
xscal(const unsigned int * n,
const std::complex<double> *alpha,
std::complex<double> * x,
const unsigned int * inc)
{
zscal_(n, alpha, x, inc);
}
inline void
xcopy(const unsigned int *n,
const double * x,
const unsigned int *incx,
double * y,
const unsigned int *incy)
{
dcopy_(n, x, incx, y, incy);
}
inline void
xcopy(const unsigned int * n,
const std::complex<double> *x,
const unsigned int * incx,
std::complex<double> * y,
const unsigned int * incy)
{
zcopy_(n, x, incx, y, incy);
}
/**
* @brief Contains linear algebra functions used in the implementation of an eigen solver
*
* @author Phani Motamarri, Sambit Das
*/
namespace linearAlgebraOperations
{
/** @brief Compute inverse of serial matrix using LAPACK LU factorization
*/
void
inverse(double *A, int N);
/** @brief Calculates an estimate of lower and upper bounds of a matrix using
* k-step Lanczos method.
*
* @param operatorMatrix An object which has access to the given matrix
* @param vect A dummy vector
* @return std::pair<double,double> An estimate of the lower and upper bound of the given matrix
*/
template <typename T, dftfe::utils::MemorySpace memorySpace>
std::pair<double, double>
lanczosLowerUpperBoundEigenSpectrum(
const std::shared_ptr<dftfe::linearAlgebra::BLASWrapper<memorySpace>>
& BLASWrapperPtr,
operatorDFTClass<memorySpace> & operatorMatrix,
dftfe::linearAlgebra::MultiVector<T, memorySpace> &X,
dftfe::linearAlgebra::MultiVector<T, memorySpace> &Y,
dftfe::linearAlgebra::MultiVector<T, memorySpace> &Z,
const dftParameters & dftParams);
/** @brief Apply Chebyshev filter to a given subspace
*
* @param[in] operatorMatrix An object which has access to the given matrix
* @param[in,out] X Given subspace as a dealii array representing multiple
* fields as a flattened array. In-place update of the given subspace.
* @param[in] numberComponents Number of multiple-fields
* @param[in] m Chebyshev polynomial degree
* @param[in] a lower bound of unwanted spectrum
* @param[in] b upper bound of unwanted spectrum
* @param[in] a0 lower bound of wanted spectrum
*/
template <typename T, dftfe::utils::MemorySpace memorySpace>
void
chebyshevFilter(operatorDFTClass<memorySpace> &operatorMatrix,
dftfe::linearAlgebra::MultiVector<T, memorySpace> &X,
dftfe::linearAlgebra::MultiVector<T, memorySpace> &Y,
const unsigned int m,
const double a,
const double b,
const double a0);
template <typename T, typename TFP32, dftfe::utils::MemorySpace memorySpace>
void
chebyshevFilterSinglePrec(
const std::shared_ptr<dftfe::linearAlgebra::BLASWrapper<memorySpace>>
& BLASWrapperPtr,
operatorDFTClass<memorySpace> & operatorMatrix,
dftfe::linearAlgebra::MultiVector<T, memorySpace> & X,
dftfe::linearAlgebra::MultiVector<T, memorySpace> & Y,
dftfe::linearAlgebra::MultiVector<TFP32, memorySpace> &X_SP,
dftfe::linearAlgebra::MultiVector<TFP32, memorySpace> &Y_SP,
std::vector<double> eigenvalues,
const unsigned int m,
const double a,
const double b,
const double a0);
/** @brief Orthogonalize given subspace using GramSchmidt orthogonalization
*
* @param[in,out] X Given subspace as flattened array of multi-vectors.
* In-place update of the given subspace
* @param[in] numberComponents Number of multiple-fields
* @param[in] mpiComm global communicator
*/
template <typename T>
void
gramSchmidtOrthogonalization(T * X,
const unsigned int numberComponents,
const unsigned int numberDofs,
const MPI_Comm & mpiComm);
/** @brief Orthogonalize given subspace using Lowden orthogonalization for double data-type
* (serial version using LAPACK)
*
* @param[in,out] X Given subspace as flattened array of multi-vectors.
* In-place update of the given subspace
* @param[in] numberComponents Number of multiple-fields
* @param[in] mpiComm global communicator
* @return flag indicating success/failure. 1 for failure, 0 for success
*/
unsigned int
lowdenOrthogonalization(std::vector<dataTypes::number> &X,
const unsigned int numberComponents,
const MPI_Comm & mpiComm,
const dftParameters & dftParams);
/** @brief Orthogonalize given subspace using Pseudo-Gram-Schmidt orthogonalization
* (serial version using LAPACK, parallel version using ScaLAPACK)
*
* @param[in,out] X Given subspace as flattened array of multi-vectors.
* In-place update of the given subspace
* @param[in] numberComponents Number of multiple-fields
* @param[in] mpiCommParent parent communicator
* @param[in] interBandGroupComm interpool communicator for parallelization
* over band groups
* @param[in] mpiComm domain decomposition communicator
*
* @return flag indicating success/failure. 1 for failure, 0 for success
*/
template <typename T>
unsigned int
pseudoGramSchmidtOrthogonalization(elpaScalaManager & elpaScala,
T * X,
const unsigned int numberComponents,
const unsigned int numberDofs,
const MPI_Comm & mpiCommParent,
const MPI_Comm & interBandGroupComm,
const MPI_Comm & mpiCommDomain,
const bool useMixedPrec,
const dftParameters &dftParams);
/** @brief Compute Rayleigh-Ritz projection
* (serial version using LAPACK, parallel version using ScaLAPACK)
*
* @param[in] operatorMatrix An object which has access to the given matrix
* @param[in,out] X Given subspace as flattened array of multi-vectors.
* In-place rotated subspace
* @param[in] numberComponents Number of vectors
* @param[in] mpiCommDomain parent communicator
* @param[in] interBandGroupComm interpool communicator for parallelization
* over band groups
* @param[in] mpiCommDomain domain decomposition communicator
* @param[out] eigenValues of the Projected Hamiltonian
*/
template <typename T>
void
rayleighRitzGEP(
operatorDFTClass<dftfe::utils::MemorySpace::HOST> &operatorMatrix,
elpaScalaManager & elpaScala,
T * X,
const unsigned int numberComponents,
const unsigned int numberDofs,
const MPI_Comm & mpiCommParent,
const MPI_Comm & interBandGroupComm,
const MPI_Comm & mpiCommDomain,
std::vector<double> & eigenValues,
const bool useMixedPrec,
const dftParameters & dftParams);
/** @brief Compute Rayleigh-Ritz projection
* (serial version using LAPACK, parallel version using ScaLAPACK)
*
* @param[in] operatorMatrix An object which has access to the given matrix
* @param[in,out] X Given subspace as flattened array of multi-vectors.
* In-place rotated subspace
* @param[in] numberComponents Number of vectors
* @param[in] mpiCommParent parent mpi communicator
* @param[in] interBandGroupComm interpool communicator for parallelization
* over band groups
* @param[in] mpiCommDomain domain decomposition communicator
* @param[out] eigenValues of the Projected Hamiltonian
*/
template <typename T>
void
rayleighRitz(
operatorDFTClass<dftfe::utils::MemorySpace::HOST> &operatorMatrix,
elpaScalaManager & elpaScala,
T * X,
const unsigned int numberComponents,
const unsigned int numberDofs,
const MPI_Comm & mpiCommParent,
const MPI_Comm & interBandGroupComm,
const MPI_Comm & mpiCommDomain,
std::vector<double> & eigenValues,
const dftParameters & dftParams,
const bool doCommAfterBandParal = true);
/** @brief Compute Rayleigh-Ritz projection in case of spectrum split using direct diagonalization
* (serial version using LAPACK, parallel version using ScaLAPACK)
*
* @param[in] operatorMatrix An object which has access to the given matrix
* @param[in] X Given subspace as flattened array of multi-vectors.
* @param[out] Y rotated subspace of top states
* @param[in] numberComponents Number of vectors
* @param[in] numberCoreStates Number of core states to be used for
* spectrum splitting
* @param[in] mpiCommParent parent mpi communicator
* @param[in] interBandGroupComm interpool communicator for parallelization
* over band groups
* @param[in] mpiCommDomain domain decomposition communicator
* @param[out] eigenValues of the Projected Hamiltonian
*/
template <typename T>
void
rayleighRitzGEPSpectrumSplitDirect(
operatorDFTClass<dftfe::utils::MemorySpace::HOST> &operatorMatrix,
elpaScalaManager & elpaScala,
T * X,
T * Y,
const unsigned int numberComponents,
const unsigned int numberDofs,
const unsigned int numberCoreStates,
const MPI_Comm & mpiCommParent,
const MPI_Comm & interBandGroupComm,
const MPI_Comm & mpiCommDomain,
const bool useMixedPrec,
std::vector<double> & eigenValues,
const dftParameters & dftParams);
/** @brief Compute Rayleigh-Ritz projection in case of spectrum split using direct diagonalization
* (serial version using LAPACK, parallel version using ScaLAPACK)
*
* @param[in] operatorMatrix An object which has access to the given matrix
* @param[in] X Given subspace as flattened array of multi-vectors.
* @param[out] Y rotated subspace of top states
* @param[in] numberComponents Number of vectors
* @param[in] numberCoreStates Number of core states to be used for
* spectrum splitting
* @param[in] mpiCommParent parent mpi communicator
* @param[in] interBandGroupComm interpool communicator for parallelization
* over band groups
* @param[in] mpiCommDomain domain decomposition communicator
* @param[out] eigenValues of the Projected Hamiltonian
*/
template <typename T>
void
rayleighRitzSpectrumSplitDirect(
operatorDFTClass<dftfe::utils::MemorySpace::HOST> &operatorMatrix,
elpaScalaManager & elpaScala,
const T * X,
T * Y,
const unsigned int numberComponents,
const unsigned int numberDofs,
const unsigned int numberCoreStates,
const MPI_Comm & mpiCommParent,
const MPI_Comm & interBandGroupComm,
const MPI_Comm & mpiCommDomain,
const bool useMixedPrec,
std::vector<double> & eigenValues,
const dftParameters & dftParams);
/** @brief Compute residual norm associated with eigenValue problem of the given operator
*
* @param[in] operatorMatrix An object which has access to the given matrix
* @param[in] X Given subspace as STL vector of dealii vectors
* @param[in] eigenValues eigenValues of the operator
* @param[in] mpiCommParent parent mpi communicator
* @param[in] mpiCommDomain domain decomposition communicator
* @param[out] residualNorms of the eigen Value problem
*/
template <typename T>
void
computeEigenResidualNorm(
operatorDFTClass<dftfe::utils::MemorySpace::HOST> &operatorMatrix,
T * X,
const std::vector<double> & eigenValues,
const unsigned int numberComponents,
const unsigned int numberDofs,
const MPI_Comm & mpiCommParent,
const MPI_Comm & mpiCommDomain,
const MPI_Comm & interBandGroupComm,
std::vector<double> & residualNorm,
const dftParameters & dftParams);
/** @brief Compute first order response in density matrix with respect to perturbation in the Hamiltonian.
* Perturbation is computed in the eigenbasis.
*/
template <typename T>
void
densityMatrixEigenBasisFirstOrderResponse(
operatorDFTClass<dftfe::utils::MemorySpace::HOST> &operatorMatrix,
T * X,
const unsigned int N,
const unsigned int numberLocalDofs,
const MPI_Comm & mpiCommParent,
const MPI_Comm & mpiCommDomain,
const MPI_Comm & interBandGroupComm,
const std::vector<double> & eigenValues,
const double fermiEnergy,
std::vector<double> &densityMatDerFermiEnergy,
elpaScalaManager & elpaScala,
const dftParameters &dftParams);
/**
* @brief Compute projection of the operator into a subspace spanned by a given orthogonal basis HProjConj=X^{T}*HConj*XConj
*
* @param X Vector of Vectors containing multi-wavefunction fields
* @param numberComponents number of wavefunctions associated with a given node
* @param ProjMatrix projected small matrix
*/
void
XtHX(operatorDFTClass<dftfe::utils::MemorySpace::HOST> &operatorMatrix,
const dataTypes::number * X,
const unsigned int numberComponents,
const unsigned int numberLocalDofs,
const MPI_Comm & mpiCommDomain,
const MPI_Comm & interBandGroupComm,
const dftParameters & dftParams,
std::vector<dataTypes::number> & ProjHam);
/**
* @brief Compute projection of the operator into a subspace spanned by a given orthogonal basis HProjConj=X^{T}*HConj*XConj
*
* @param X Vector of Vectors containing multi-wavefunction fields
* @param numberComponents number of wavefunctions associated with a given node
* @param processGrid two-dimensional processor grid corresponding to the parallel projHamPar
* @param projHamPar parallel ScaLAPACKMatrix which stores the computed projection
* of the operation into the given subspace
*/
void
XtHX(operatorDFTClass<dftfe::utils::MemorySpace::HOST> &operatorMatrix,
const dataTypes::number * X,
const unsigned int numberComponents,
const unsigned int numberLocalDofs,
const std::shared_ptr<const dftfe::ProcessGrid> & processGrid,
const MPI_Comm & mpiCommDomain,
const MPI_Comm & interBandGroupComm,
const dftParameters & dftParams,
dftfe::ScaLAPACKMatrix<dataTypes::number> & projHamPar,
const bool onlyHPrimePartForFirstOrderDensityMatResponse = false);
/**
* @brief Compute projection of the operator into a subspace spanned by a given orthogonal basis HProjConj=X^{T}*HConj*XConj
*
* @param X Vector of Vectors containing multi-wavefunction fields
* @param totalNumberComponents number of wavefunctions associated with a given node
* @param singlePrecComponents number of wavecfuntions starting from the first for
* which the project Hamiltionian block will be computed in single
* procession. However the cross blocks will still be computed in double
* precision.
* @param processGrid two-dimensional processor grid corresponding to the parallel projHamPar
* @param projHamPar parallel ScaLAPACKMatrix which stores the computed projection
* of the operation into the given subspace
*/
void
XtHXMixedPrec(
operatorDFTClass<dftfe::utils::MemorySpace::HOST> &operatorMatrix,
const dataTypes::number * X,
const unsigned int totalNumberComponents,
const unsigned int singlePrecComponents,
const unsigned int numberLocalDofs,
const std::shared_ptr<const dftfe::ProcessGrid> & processGrid,
const MPI_Comm & mpiCommDomain,
const MPI_Comm & interBandGroupComm,
const dftParameters & dftParams,
dftfe::ScaLAPACKMatrix<dataTypes::number> & projHamPar,
const bool onlyHPrimePartForFirstOrderDensityMatResponse = false);
} // namespace linearAlgebraOperations
} // namespace dftfe
#endif