17#ifndef linearAlgebraOperationsInternal_h
18#define linearAlgebraOperationsInternal_h
27#include <unordered_map>
45 const MPI_Comm &mpi_communicator,
46 MPI_Comm & processGridCommunicatorActive,
47 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
48 const unsigned int na,
49 const unsigned int nev,
50 const unsigned int blockSize,
60 const MPI_Comm & mpi_communicator,
62 std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
64 const bool useOnlyThumbRule =
false);
72 const MPI_Comm & mpi_communicator,
73 const unsigned sizeRows,
74 const unsigned sizeColumns,
75 std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
85 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
87 std::unordered_map<unsigned int, unsigned int> & globalToLocalRowIdMap,
88 std::unordered_map<unsigned int, unsigned int>
89 &globalToLocalColumnIdMap);
99 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
101 const MPI_Comm & interComm);
109 template <
typename T>
112 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
113 const std::shared_ptr<
124 template <
typename T>
127 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
129 const MPI_Comm & interComm,
130 const unsigned int broadcastRoot);
141 template <
typename T>
145 const std::shared_ptr<
148 const unsigned int XLocalSize,
149 const unsigned int numberVectors,
150 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
151 const MPI_Comm & interBandGroupComm,
152 const MPI_Comm & mpiComm,
166 template <
typename T,
typename TLowPrec>
170 const std::shared_ptr<
173 const unsigned int XLocalSize,
174 const unsigned int numberVectors,
175 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
176 const MPI_Comm & interBandGroupComm,
177 const MPI_Comm & mpiComm,
191 template <
typename T>
194 T *subspaceVectorsArray,
195 const std::shared_ptr<
198 const unsigned int subspaceVectorsArrayLocalSize,
199 const unsigned int N,
200 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
201 const MPI_Comm & interBandGroupComm,
202 const MPI_Comm & mpiComm,
205 const bool rotationMatTranspose =
false,
206 const bool isRotationMatLowerTria =
false,
207 const bool doCommAfterBandParal =
true);
218 template <
typename T,
typename TLowPrec>
221 T *subspaceVectorsArray,
222 const std::shared_ptr<
225 const unsigned int subspaceVectorsArrayLocalSize,
226 const unsigned int N,
227 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
228 const MPI_Comm & interBandGroupComm,
229 const MPI_Comm & mpiComm,
232 const bool rotationMatTranspose =
false,
233 const bool doCommAfterBandParal =
true);
250 template <
typename T>
254 const std::shared_ptr<
258 const unsigned int subspaceVectorsArrayLocalSize,
259 const unsigned int N,
260 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
261 const unsigned int numberTopVectors,
262 const MPI_Comm & interBandGroupComm,
263 const MPI_Comm & mpiComm,
266 const bool QMatTranspose =
false);
282 template <
typename T,
typename TLowPrec>
286 const std::shared_ptr<
290 const unsigned int subspaceVectorsArrayLocalSize,
291 const unsigned int N,
292 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
293 const unsigned int numberTopVectors,
294 const MPI_Comm & interBandGroupComm,
295 const MPI_Comm & mpiComm,
298 const bool QMatTranspose =
false);
310 template <
typename T,
typename TLowPrec>
313 T *subspaceVectorsArray,
314 const std::shared_ptr<
317 const unsigned int subspaceVectorsArrayLocalSize,
318 const unsigned int N,
319 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
320 const MPI_Comm & interBandGroupComm,
321 const MPI_Comm & mpiComm,
324 const bool rotationMatTranspose =
false,
325 const bool doCommAfterBandParal =
true);
Scalapack wrapper adapted from dealii library and extended implementation to complex datatype.
Definition scalapackWrapper.h:37
Namespace which declares the input parameters and the functions to parse them from the input paramete...
Definition dftParameters.h:35
Definition BLASWrapper.h:35
Contains internal functions used in linearAlgebraOperations.
Definition linearAlgebraOperationsInternal.h:39
void createProcessGridSquareMatrix(const MPI_Comm &mpi_communicator, const unsigned size, std::shared_ptr< const dftfe::ProcessGrid > &processGrid, const dftParameters &dftParams, const bool useOnlyThumbRule=false)
Wrapper function to create a two dimensional processor grid for a square matrix in dftfe::ScaLAPACKMa...
void subspaceRotationSpectrumSplitMixedPrec(const T *X, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, T *Y, const unsigned int subspaceVectorsArrayLocalSize, const unsigned int N, const std::shared_ptr< const dftfe::ProcessGrid > &processGrid, const unsigned int numberTopVectors, const MPI_Comm &interBandGroupComm, const MPI_Comm &mpiComm, const dftfe::ScaLAPACKMatrix< T > &QMat, const dftParameters &dftParams, const bool QMatTranspose=false)
Computes Y^{T}=Q*X^{T}.
void scaleScaLAPACKMat(const std::shared_ptr< const dftfe::ProcessGrid > &processGrid, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, dftfe::ScaLAPACKMatrix< T > &mat, const T scalar)
scale a ScaLAPACKMat with a scalar
void fillParallelOverlapMatrix(const T *X, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, const unsigned int XLocalSize, const unsigned int numberVectors, const std::shared_ptr< const dftfe::ProcessGrid > &processGrid, const MPI_Comm &interBandGroupComm, const MPI_Comm &mpiComm, dftfe::ScaLAPACKMatrix< T > &overlapMatPar, const dftParameters &dftParams)
Computes Sc=X^{T}*Xc and stores in a parallel ScaLAPACK matrix. X^{T} is the subspaceVectorsArray sto...
void broadcastAcrossInterCommScaLAPACKMat(const std::shared_ptr< const dftfe::ProcessGrid > &processGrid, dftfe::ScaLAPACKMatrix< T > &mat, const MPI_Comm &interComm, const unsigned int broadcastRoot)
MPI_Bcast of ScaLAPACKMat across a given inter communicator from a given broadcast root....
void subspaceRotation(T *subspaceVectorsArray, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, const unsigned int subspaceVectorsArrayLocalSize, const unsigned int N, const std::shared_ptr< const dftfe::ProcessGrid > &processGrid, const MPI_Comm &interBandGroupComm, const MPI_Comm &mpiComm, const dftfe::ScaLAPACKMatrix< T > &rotationMatPar, const dftParameters &dftParams, const bool rotationMatTranspose=false, const bool isRotationMatLowerTria=false, const bool doCommAfterBandParal=true)
Computes X^{T}=Q*X^{T} inplace. X^{T} is the subspaceVectorsArray stored in the column major format (...
void subspaceRotationSpectrumSplit(const T *X, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, T *Y, const unsigned int subspaceVectorsArrayLocalSize, const unsigned int N, const std::shared_ptr< const dftfe::ProcessGrid > &processGrid, const unsigned int numberTopVectors, const MPI_Comm &interBandGroupComm, const MPI_Comm &mpiComm, const dftfe::ScaLAPACKMatrix< T > &QMat, const dftParameters &dftParams, const bool QMatTranspose=false)
Computes Y^{T}=Q*X^{T}.
void subspaceRotationCGSMixedPrec(T *subspaceVectorsArray, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, const unsigned int subspaceVectorsArrayLocalSize, const unsigned int N, const std::shared_ptr< const dftfe::ProcessGrid > &processGrid, const MPI_Comm &interBandGroupComm, const MPI_Comm &mpiComm, const dftfe::ScaLAPACKMatrix< T > &rotationMatPar, const dftParameters &dftParams, const bool rotationMatTranspose=false, const bool doCommAfterBandParal=true)
Computes X^{T}=Q*X^{T} inplace. X^{T} is the subspaceVectorsArray stored in the column major format (...
void createProcessGridRectangularMatrix(const MPI_Comm &mpi_communicator, const unsigned sizeRows, const unsigned sizeColumns, std::shared_ptr< const dftfe::ProcessGrid > &processGrid, const dftParameters &dftParams)
Wrapper function to create a two dimensional processor grid for a rectangular matrix in dftfe::ScaLAP...
void fillParallelOverlapMatrixMixedPrec(const T *X, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, const unsigned int XLocalSize, const unsigned int numberVectors, const std::shared_ptr< const dftfe::ProcessGrid > &processGrid, const MPI_Comm &interBandGroupComm, const MPI_Comm &mpiComm, dftfe::ScaLAPACKMatrix< T > &overlapMatPar, const dftParameters &dftParams)
Computes Sc=X^{T}*Xc and stores in a parallel ScaLAPACK matrix. X^{T} is the subspaceVectorsArray sto...
void sumAcrossInterCommScaLAPACKMat(const std::shared_ptr< const dftfe::ProcessGrid > &processGrid, dftfe::ScaLAPACKMatrix< T > &mat, const MPI_Comm &interComm)
Mpi all reduce of ScaLAPACKMat across a given inter communicator. Used for band parallelization.
void subspaceRotationMixedPrec(T *subspaceVectorsArray, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, const unsigned int subspaceVectorsArrayLocalSize, const unsigned int N, const std::shared_ptr< const dftfe::ProcessGrid > &processGrid, const MPI_Comm &interBandGroupComm, const MPI_Comm &mpiComm, const dftfe::ScaLAPACKMatrix< T > &rotationMatPar, const dftParameters &dftParams, const bool rotationMatTranspose=false, const bool doCommAfterBandParal=true)
Computes X^{T}=Q*X^{T} inplace. X^{T} is the subspaceVectorsArray stored in the column major format (...
void setupELPAHandleParameters(const MPI_Comm &mpi_communicator, MPI_Comm &processGridCommunicatorActive, const std::shared_ptr< const dftfe::ProcessGrid > &processGrid, const unsigned int na, const unsigned int nev, const unsigned int blockSize, elpa_t &elpaHandle, const dftParameters &dftParams)
setup ELPA parameters.
void createGlobalToLocalIdMapsScaLAPACKMat(const std::shared_ptr< const dftfe::ProcessGrid > &processGrid, const dftfe::ScaLAPACKMatrix< T > &mat, std::unordered_map< unsigned int, unsigned int > &globalToLocalRowIdMap, std::unordered_map< unsigned int, unsigned int > &globalToLocalColumnIdMap)
Creates global row/column id to local row/column ids for dftfe::ScaLAPACKMatrix.
Contains linear algebra functions used in the implementation of an eigen solver.
Definition linearAlgebraOperations.h:502
Definition pseudoPotentialToDftfeConverter.cc:34