DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
linearAlgebraOperationsDevice.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (c) 2017-2025 The Regents of the University of Michigan and DFT-FE
4// authors.
5//
6// This file is part of the DFT-FE code.
7//
8// The DFT-FE code is free software; you can use it, redistribute
9// it, and/or modify it under the terms of the GNU Lesser General
10// Public License as published by the Free Software Foundation; either
11// version 2.1 of the License, or (at your option) any later version.
12// The full text of the license can be found in the file LICENSE at
13// the top level of the DFT-FE distribution.
14//
15// ---------------------------------------------------------------------
16
17#if defined(DFTFE_WITH_DEVICE)
18# ifndef linearAlgebraOperationsDevice_h
19# define linearAlgebraOperationsDevice_h
20
21# include <headers.h>
22# include <operator.h>
23# include "process_grid.h"
24# include "scalapackWrapper.h"
25# include "elpaScalaManager.h"
27# include "dftParameters.h"
28# include <BLASWrapper.h>
29
30namespace dftfe
31{
32 extern "C"
33 {
34 void
35 dsyevd_(const char *jobz,
36 const char *uplo,
37 const unsigned int *n,
38 double *A,
39 const unsigned int *lda,
40 double *w,
41 double *work,
42 const unsigned int *lwork,
43 int *iwork,
44 const unsigned int *liwork,
45 int *info);
46
47 void
48 zheevd_(const char *jobz,
49 const char *uplo,
50 const unsigned int *n,
51 std::complex<double> *A,
52 const unsigned int *lda,
53 double *w,
54 std::complex<double> *work,
55 const unsigned int *lwork,
56 double *rwork,
57 const unsigned int *lrwork,
58 int *iwork,
59 const unsigned int *liwork,
60 int *info);
61 }
62
63
64 /**
65 * @brief Contains functions for linear algebra operations on Device
66 *
67 * @author Sambit Das
68 */
70 {
71 /** @brief Apply Chebyshev filter to a given subspace
72 *
73 * @param[in] operatorMatrix An object which has access to the given matrix
74 * @param[in,out] X Given subspace as a dealii array representing multiple
75 * fields as a flattened array. In-place update of the given subspace.
76 * @param[in] numberComponents Number of multiple-fields
77 * @param[in] m Chebyshev polynomial degree
78 * @param[in] a lower bound of unwanted spectrum
79 * @param[in] b upper bound of unwanted spectrum
80 * @param[in] a0 lower bound of wanted spectrum
81 */
82 void
83 chebyshevFilterOverlapComputeCommunication(
84 operatorDFTClass<dftfe::utils::MemorySpace::DEVICE> &operatorMatrix,
85 dftfe::linearAlgebra::MultiVector<dataTypes::number,
87 dftfe::linearAlgebra::MultiVector<dataTypes::number,
89 dftfe::linearAlgebra::MultiVector<dataTypes::number,
91 dftfe::linearAlgebra::MultiVector<dataTypes::number,
93 const dftfe::uInt m,
94 const double a,
95 const double b,
96 const double a0);
97 template <typename T1, typename T2>
98 void
99 reformulatedChebyshevFilterOverlapComputeCommunication(
100 std::shared_ptr<
101 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
102 &BLASWrapperPtr,
103 operatorDFTClass<dftfe::utils::MemorySpace::DEVICE> &operatorMatrix,
104 dftfe::linearAlgebra::MultiVector<T1, dftfe::utils::MemorySpace::DEVICE>
105 &X1,
106 dftfe::linearAlgebra::MultiVector<T1, dftfe::utils::MemorySpace::DEVICE>
107 &Y1,
108 dftfe::linearAlgebra::MultiVector<T1, dftfe::utils::MemorySpace::DEVICE>
109 &X2,
110 dftfe::linearAlgebra::MultiVector<T1, dftfe::utils::MemorySpace::DEVICE>
111 &Y2,
112 dftfe::linearAlgebra::MultiVector<T2, dftfe::utils::MemorySpace::DEVICE>
113 &X1_SP,
114 dftfe::linearAlgebra::MultiVector<T2, dftfe::utils::MemorySpace::DEVICE>
115 &Y1_SP,
116 dftfe::linearAlgebra::MultiVector<T2, dftfe::utils::MemorySpace::DEVICE>
117 &X2_SP,
118 dftfe::linearAlgebra::MultiVector<T2, dftfe::utils::MemorySpace::DEVICE>
119 &Y2_SP,
120 std::vector<double> eigenvalues,
121 const dftfe::uInt m,
122 const double a,
123 const double b,
124 const double a0,
125 const bool approxOverlapMatrix);
126
127 /** @brief Computes Sc=X^{T}*Xc.
128 *
129 *
130 */
131 void
132 fillParallelOverlapMatScalapack(
133 operatorDFTClass<dftfe::utils::MemorySpace::DEVICE> &operatorMatrix,
134 const dataTypes::number *X,
135 distributedDeviceVec<dataTypes::number> &XBlock,
136 distributedDeviceVec<dataTypes::number> &OXBlock,
137 const dftfe::uInt M,
138 const dftfe::uInt N,
139 std::shared_ptr<
140 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
141 &BLASWrapperPtr,
142 const MPI_Comm &mpiCommDomain,
143 utils::DeviceCCLWrapper &devicecclMpiCommDomain,
144 const MPI_Comm &interBandGroupComm,
145 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
146 dftfe::ScaLAPACKMatrix<dataTypes::number> &overlapMatPar,
147 const dftParameters &dftParams);
148
149
150
151 /** @brief Computes Sc=X^{T}*Xc.
152 *
153 *
154 */
155 void
156 fillParallelOverlapMatScalapackAsyncComputeCommun(
157 operatorDFTClass<dftfe::utils::MemorySpace::DEVICE> &operatorMatrix,
158 const dataTypes::number *X,
159 distributedDeviceVec<dataTypes::number> &XBlock,
160 distributedDeviceVec<dataTypes::number> &OXBlock,
161 const dftfe::uInt M,
162 const dftfe::uInt N,
163 std::shared_ptr<
164 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
165 &BLASWrapperPtr,
166 const MPI_Comm &mpiCommDomain,
167 utils::DeviceCCLWrapper &devicecclMpiCommDomain,
168 const MPI_Comm &interBandGroupComm,
169 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
170 dftfe::ScaLAPACKMatrix<dataTypes::number> &overlapMatPar,
171 const dftParameters &dftParams);
172
173
174
175 /** @brief Computes Sc=X^{T}*Xc.
176 *
177 *
178 */
179 void
180 fillParallelOverlapMatMixedPrecScalapackAsyncComputeCommun(
181 operatorDFTClass<dftfe::utils::MemorySpace::DEVICE> &operatorMatrix,
182 const dataTypes::number *X,
183 distributedDeviceVec<dataTypes::number> &XBlock,
184 distributedDeviceVec<dataTypes::number> &OXBlock,
185 const dftfe::uInt M,
186 const dftfe::uInt N,
187 const dftfe::uInt Noc,
188 std::shared_ptr<
189 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
190 &BLASWrapperPtr,
191 const MPI_Comm &mpiCommDomain,
192 utils::DeviceCCLWrapper &devicecclMpiCommDomain,
193 const MPI_Comm &interBandGroupComm,
194 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
195 dftfe::ScaLAPACKMatrix<dataTypes::number> &overlapMatPar,
196 const dftParameters &dftParams);
197
198 /** @brief Computes Sc=X^{T}*Xc.
199 *
200 *
201 */
202 void
203 fillParallelOverlapMatMixedPrecCommunScalapackAsyncComputeCommun(
204 operatorDFTClass<dftfe::utils::MemorySpace::DEVICE> &operatorMatrix,
205 const dataTypes::number *X,
206 distributedDeviceVec<dataTypes::number> &XBlock,
207 distributedDeviceVec<dataTypes::number> &OXBlock,
208 const dftfe::uInt M,
209 const dftfe::uInt N,
210 const dftfe::uInt Noc,
211 std::shared_ptr<
212 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
213 &BLASWrapperPtr,
214 const MPI_Comm &mpiCommDomain,
215 utils::DeviceCCLWrapper &devicecclMpiCommDomain,
216 const MPI_Comm &interBandGroupComm,
217 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
218 dftfe::ScaLAPACKMatrix<dataTypes::number> &overlapMatPar,
219 const dftParameters &dftParams);
220
221 /** @brief Computes Sc=X^{T}*Xc.
222 *
223 *
224 */
225 void
226 fillParallelOverlapMatMixedPrecScalapack(
227 operatorDFTClass<dftfe::utils::MemorySpace::DEVICE> &operatorMatrix,
228 const dataTypes::number *X,
229 distributedDeviceVec<dataTypes::number> &XBlock,
230 distributedDeviceVec<dataTypes::number> &OXBlock,
231 const dftfe::uInt M,
232 const dftfe::uInt N,
233 const dftfe::uInt Noc,
234 std::shared_ptr<
235 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
236 &BLASWrapperPtr,
237 const MPI_Comm &mpiCommDomain,
238 utils::DeviceCCLWrapper &devicecclMpiCommDomain,
239 const MPI_Comm &interBandGroupComm,
240 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
241 dftfe::ScaLAPACKMatrix<dataTypes::number> &overlapMatPar,
242 const dftParameters &dftParams);
243
244
245
246 /** @brief CGS orthogonalization
247 */
248 void
250 operatorDFTClass<dftfe::utils::MemorySpace::DEVICE> &operatorMatrix,
251 elpaScalaManager &elpaScala,
253 distributedDeviceVec<dataTypes::number> &Xb,
254 distributedDeviceVec<dataTypes::number> &HXb,
255 const dftfe::uInt M,
256 const dftfe::uInt N,
257 const MPI_Comm &mpiCommParent,
258 const MPI_Comm &mpiCommDomain,
259 utils::DeviceCCLWrapper &devicecclMpiCommDomain,
260 const MPI_Comm &interBandGroupComm,
261 std::shared_ptr<
262 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
263 &BLASWrapperPtr,
264 const dftParameters &dftParams,
265 const bool useMixedPrecOverall = false);
266
267 void
268 subspaceRotationScalapack(
270 const dftfe::uInt M,
271 const dftfe::uInt N,
272 std::shared_ptr<
273 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
274 &BLASWrapperPtr,
275 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
276 const MPI_Comm &mpiCommDomain,
277 utils::DeviceCCLWrapper &devicecclMpiCommDomain,
278 const MPI_Comm &interBandGroupComm,
279 const dftfe::ScaLAPACKMatrix<dataTypes::number> &rotationMatPar,
280 const dftParameters &dftParams,
281 const bool rotationMatTranspose = false,
282 const bool isRotationMatLowerTria = false);
283
284
285
286 void
287 subspaceRotationCGSMixedPrecScalapack(
289 const dftfe::uInt M,
290 const dftfe::uInt N,
291 std::shared_ptr<
292 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
293 &BLASWrapperPtr,
294 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
295 const MPI_Comm &mpiCommDomain,
296 utils::DeviceCCLWrapper &devicecclMpiCommDomain,
297 const MPI_Comm &interBandGroupComm,
298 const dftfe::ScaLAPACKMatrix<dataTypes::number> &rotationMatPar,
299 const dftParameters &dftParams,
300 const bool rotationMatTranspose = false);
301
302
303 void
304 subspaceRotationRRMixedPrecScalapack(
306 const dftfe::uInt M,
307 const dftfe::uInt N,
308 std::shared_ptr<
309 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
310 &BLASWrapperPtr,
311 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
312 const MPI_Comm &mpiCommDomain,
313 utils::DeviceCCLWrapper &devicecclMpiCommDomain,
314 const MPI_Comm &interBandGroupComm,
315 const dftfe::ScaLAPACKMatrix<dataTypes::number> &rotationMatPar,
316 const dftParameters &dftParams,
317 const bool rotationMatTranspose = false);
318
319
320 void
322 operatorDFTClass<dftfe::utils::MemorySpace::DEVICE> &operatorMatrix,
323 elpaScalaManager &elpaScala,
325 distributedDeviceVec<dataTypes::number> &Xb,
326 distributedDeviceVec<dataTypes::number> &HXb,
327 const dftfe::uInt M,
328 const dftfe::uInt N,
329 const MPI_Comm &mpiCommParent,
330 const MPI_Comm &mpiCommDomain,
331 utils::DeviceCCLWrapper &devicecclMpiCommDomain,
332 const MPI_Comm &interBandGroupComm,
333 std::vector<double> &eigenValues,
334 std::shared_ptr<
335 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
336 &BLASWrapperPtr,
337 const dftParameters &dftParams,
338 const bool useMixedPrecOverall = false);
339
340 void
342 operatorDFTClass<dftfe::utils::MemorySpace::DEVICE> &operatorMatrix,
343 elpaScalaManager &elpaScala,
345 distributedDeviceVec<dataTypes::number> &Xb,
346 distributedDeviceVec<dataTypes::number> &HXb,
347 const dftfe::uInt M,
348 const dftfe::uInt N,
349 const MPI_Comm &mpiCommParent,
350 const MPI_Comm &mpiCommDomain,
351 utils::DeviceCCLWrapper &devicecclMpiCommDomain,
352 const MPI_Comm &interBandGroupComm,
353 std::vector<double> &eigenValues,
354 std::shared_ptr<
355 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
356 &BLASWrapperPtr,
357 const dftParameters &dftParams,
358 const bool useMixedPrecOverall = false);
359
360
361
362 void
364 operatorDFTClass<dftfe::utils::MemorySpace::DEVICE> &operatorMatrix,
366 distributedDeviceVec<dataTypes::number> &Xb,
367 distributedDeviceVec<dataTypes::number> &HXb,
368 const dftfe::uInt M,
369 const dftfe::uInt N,
370 const MPI_Comm &mpiCommParent,
371 const MPI_Comm &mpiCommDomain,
372 utils::DeviceCCLWrapper &devicecclMpiCommDomain,
373 const MPI_Comm &interBandGroupComm,
374 const std::vector<double> &eigenValues,
375 const double fermiEnergy,
376 std::vector<double> &densityMatDerFermiEnergy,
377 dftfe::elpaScalaManager &elpaScala,
378 std::shared_ptr<
379 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
380 &BLASWrapperPtr,
381 const dftParameters &dftParams);
382
383 void
385 operatorDFTClass<dftfe::utils::MemorySpace::DEVICE> &operatorMatrix,
387 distributedDeviceVec<dataTypes::number> &Xb,
388 distributedDeviceVec<dataTypes::number> &HXb,
389 const dftfe::uInt M,
390 const dftfe::uInt N,
391 const std::vector<double> &eigenValues,
392 const MPI_Comm &mpiCommParent,
393 const MPI_Comm &mpiCommDomain,
394 const MPI_Comm &interBandGroupComm,
395 std::shared_ptr<
396 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
397 &BLASWrapperPtr,
398 std::vector<double> &residualNorm,
399 const dftParameters &dftParams,
400 const bool useBandParal = false);
401
402 void
403 XtHX(operatorDFTClass<dftfe::utils::MemorySpace::DEVICE> &operatorMatrix,
404 const dataTypes::number *X,
405 distributedDeviceVec<dataTypes::number> &XBlock,
406 distributedDeviceVec<dataTypes::number> &HXBlock,
407 const dftfe::uInt M,
408 const dftfe::uInt N,
409 std::shared_ptr<
410 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
411 &BLASWrapperPtr,
412 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
413 dftfe::ScaLAPACKMatrix<dataTypes::number> &projHamPar,
414 utils::DeviceCCLWrapper &devicecclMpiCommDomain,
415 const MPI_Comm &mpiCommDomain,
416 const MPI_Comm &interBandGroupComm,
417 const dftParameters &dftParams,
418 const bool onlyHPrimePartForFirstOrderDensityMatResponse = false);
419
420 void
421 XtHXMixedPrecOverlapComputeCommun(
422 operatorDFTClass<dftfe::utils::MemorySpace::DEVICE> &operatorMatrix,
423 const dataTypes::number *X,
424 distributedDeviceVec<dataTypes::number> &XBlock,
425 distributedDeviceVec<dataTypes::number> &HXBlock,
426 const dftfe::uInt M,
427 const dftfe::uInt N,
428 const dftfe::uInt Noc,
429 std::shared_ptr<
430 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
431 &BLASWrapperPtr,
432 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
433 dftfe::ScaLAPACKMatrix<dataTypes::number> &projHamPar,
434 utils::DeviceCCLWrapper &devicecclMpiCommDomain,
435 const MPI_Comm &mpiCommDomain,
436 const MPI_Comm &interBandGroupComm,
437 const dftParameters &dftParams,
438 const bool onlyHPrimePartForFirstOrderDensityMatResponse = false);
439
440 void
441 XtHXOverlapComputeCommun(
442 operatorDFTClass<dftfe::utils::MemorySpace::DEVICE> &operatorMatrix,
443 const dataTypes::number *X,
444 distributedDeviceVec<dataTypes::number> &XBlock,
445 distributedDeviceVec<dataTypes::number> &HXBlock,
446 const dftfe::uInt M,
447 const dftfe::uInt N,
448 std::shared_ptr<
449 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
450 &BLASWrapperPtr,
451 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
452 dftfe::ScaLAPACKMatrix<dataTypes::number> &projHamPar,
453 utils::DeviceCCLWrapper &devicecclMpiCommDomain,
454 const MPI_Comm &mpiCommDomain,
455 const MPI_Comm &interBandGroupComm,
456 const dftParameters &dftParams,
457 const bool onlyHPrimePartForFirstOrderDensityMatResponse = false);
458
459 void
460 XtHXMixedPrecCommunOverlapComputeCommun(
461 operatorDFTClass<dftfe::utils::MemorySpace::DEVICE> &operatorMatrix,
462 const dataTypes::number *X,
463 distributedDeviceVec<dataTypes::number> &XBlock,
464 distributedDeviceVec<dataTypes::number> &HXBlock,
465 const dftfe::uInt M,
466 const dftfe::uInt N,
467 const dftfe::uInt Noc,
468 std::shared_ptr<
469 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
470 &BLASWrapperPtr,
471 const std::shared_ptr<const dftfe::ProcessGrid> &processGrid,
472 dftfe::ScaLAPACKMatrix<dataTypes::number> &projHamPar,
473 utils::DeviceCCLWrapper &devicecclMpiCommDomain,
474 const MPI_Comm &mpiCommDomain,
475 const MPI_Comm &interBandGroupComm,
476 const dftParameters &dftParams,
477 const bool onlyHPrimePartForFirstOrderDensityMatResponse = false);
478
479 } // namespace linearAlgebraOperationsDevice
480} // namespace dftfe
481# endif
482#endif
double number
Definition dftfeDataTypes.h:42
Definition linearAlgebraOperationsDeviceKernels.h:9
void rayleighRitzGEP(operatorDFTClass< dftfe::utils::MemorySpace::HOST > &operatorMatrix, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, elpaScalaManager &elpaScala, T *X, const dftfe::uInt numberComponents, const dftfe::uInt numberDofs, const MPI_Comm &mpiCommParent, const MPI_Comm &interBandGroupComm, const MPI_Comm &mpiCommDomain, std::vector< double > &eigenValues, const bool useMixedPrec, const dftParameters &dftParams)
Compute Rayleigh-Ritz projection (serial version using LAPACK, parallel version using ScaLAPACK)
void computeEigenResidualNorm(operatorDFTClass< dftfe::utils::MemorySpace::HOST > &operatorMatrix, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, T *X, const std::vector< double > &eigenValues, const dftfe::uInt numberComponents, const dftfe::uInt numberDofs, const MPI_Comm &mpiCommParent, const MPI_Comm &mpiCommDomain, const MPI_Comm &interBandGroupComm, std::vector< double > &residualNorm, const dftParameters &dftParams)
Compute residual norm associated with eigenValue problem of the given operator.
void rayleighRitz(operatorDFTClass< dftfe::utils::MemorySpace::HOST > &operatorMatrix, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, elpaScalaManager &elpaScala, T *X, const dftfe::uInt numberComponents, const dftfe::uInt numberDofs, const MPI_Comm &mpiCommParent, const MPI_Comm &interBandGroupComm, const MPI_Comm &mpiCommDomain, std::vector< double > &eigenValues, const dftParameters &dftParams, const bool doCommAfterBandParal=true)
Compute Rayleigh-Ritz projection (serial version using LAPACK, parallel version using ScaLAPACK)
dftfe::uInt pseudoGramSchmidtOrthogonalization(elpaScalaManager &elpaScala, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, T *X, const dftfe::uInt numberComponents, const dftfe::uInt numberDofs, const MPI_Comm &mpiCommParent, const MPI_Comm &interBandGroupComm, const MPI_Comm &mpiCommDomain, const bool useMixedPrec, const dftParameters &dftParams)
Orthogonalize given subspace using Pseudo-Gram-Schmidt orthogonalization (serial version using LAPACK...
void XtHX(operatorDFTClass< dftfe::utils::MemorySpace::HOST > &operatorMatrix, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, const dataTypes::number *X, const dftfe::uInt numberComponents, const dftfe::uInt numberLocalDofs, const MPI_Comm &mpiCommDomain, const MPI_Comm &interBandGroupComm, const dftParameters &dftParams, std::vector< dataTypes::number > &ProjHam)
Compute projection of the operator into a subspace spanned by a given orthogonal basis HProjConj=X^{T...
void densityMatrixEigenBasisFirstOrderResponse(operatorDFTClass< dftfe::utils::MemorySpace::HOST > &operatorMatrix, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, T *X, const dftfe::uInt N, const dftfe::uInt 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)
Compute first order response in density matrix with respect to perturbation in the Hamiltonian....
@ DEVICE
Definition MemorySpaceType.h:36
Definition pseudoPotentialToDftfeConverter.cc:34
std::uint32_t uInt
Definition TypeConfig.h:10