DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
chebyshevOrthogonalizedSubspaceIterationSolverDevice.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
18#if defined(DFTFE_WITH_DEVICE)
19# ifndef chebyshevOrthogonalizedSubspaceIterationSolverDevice_h
20# define chebyshevOrthogonalizedSubspaceIterationSolverDevice_h
21
22
24# include "headers.h"
25# include "operator.h"
26# include "elpaScalaManager.h"
27# include "dftParameters.h"
28# include <BLASWrapper.h>
29
30namespace dftfe
31{
32 /**
33 * @brief Concrete class implementing Chebyshev filtered orthogonalized subspace
34 * iteration solver.
35 * @author Sambit Das, Phani Motamarri
36 */
37 class chebyshevOrthogonalizedSubspaceIterationSolverDevice
38 {
39 public:
40 /**
41 * @brief Constructor.
42 *
43 * @param mpi_comm_parent parent mpi communicator
44 * @param mpi_comm_domain domain decomposition mpi communicator
45 * @param lowerBoundWantedSpectrum Lower Bound of the Wanted Spectrum.
46 * @param lowerBoundUnWantedSpectrum Lower Bound of the UnWanted Spectrum.
47 */
48 chebyshevOrthogonalizedSubspaceIterationSolverDevice(
49 const MPI_Comm & mpi_comm_parent,
50 const MPI_Comm & mpi_comm_domain,
51 double lowerBoundWantedSpectrum,
52 double lowerBoundUnWantedSpectrum,
53 double upperBoundUnWantedSpectrum,
54 const dftParameters &dftParams);
55
56
57
58 /**
59 * @brief Solve a generalized eigen problem.
60 */
61 double
62 solve(operatorDFTClass<dftfe::utils::MemorySpace::DEVICE> &operatorMatrix,
63 std::shared_ptr<dftfe::linearAlgebra::BLASWrapper<
64 dftfe::utils::MemorySpace::DEVICE>> & BLASWrapperPtr,
65 elpaScalaManager & elpaScala,
66 dataTypes::number * eigenVectorsFlattenedDevice,
67 const unsigned int flattenedSize,
68 const unsigned int totalNumberWaveFunctions,
69 std::vector<double> & eigenValues,
70 std::vector<double> & residuals,
71 utils::DeviceCCLWrapper &devicecclMpiCommDomain,
72 const MPI_Comm & interBandGroupComm,
73 const bool isFirstFilteringCall,
74 const bool computeResidual,
75 const bool useMixedPrecOverall = false,
76 const bool isFirstScf = false);
77
78
79
80 /**
81 * @brief Used for LRD preconditioner, also required for XL-BOMD
82 */
83 void
85 operatorDFTClass<dftfe::utils::MemorySpace::DEVICE> &operatorMatrix,
86 std::shared_ptr<
87 dftfe::linearAlgebra::BLASWrapper<dftfe::utils::MemorySpace::DEVICE>>
88 & BLASWrapperPtr,
89 dataTypes::number * eigenVectorsFlattenedDevice,
90 const unsigned int flattenedSize,
91 const unsigned int totalNumberWaveFunctions,
92 const std::vector<double> &eigenValues,
93 const double fermiEnergy,
94 std::vector<double> & densityMatDerFermiEnergy,
95 utils::DeviceCCLWrapper & devicecclMpiCommDomain,
96 const MPI_Comm & interBandGroupComm,
97 dftfe::elpaScalaManager & elpaScala);
98
99
100 /**
101 * @brief reinit spectrum bounds
102 */
103 void
104 reinitSpectrumBounds(double lowerBoundWantedSpectrum,
105 double lowerBoundUnWantedSpectrum,
106 double upperBoundUnWantedSpectrum);
107
108 private:
109 const MPI_Comm d_mpiCommParent;
110 //
111 // stores lower bound of wanted spectrum
112 //
113 double d_lowerBoundWantedSpectrum;
114
115 //
116 // stores lower bound of unwanted spectrum
117 //
118 double d_lowerBoundUnWantedSpectrum;
119
120 //
121 // stores upper bound of unwanted spectrum
122 //
123 double d_upperBoundUnWantedSpectrum;
124
125 const dftParameters &d_dftParams;
126
127 //
128 // variables for printing out and timing
129 //
130 dealii::ConditionalOStream pcout;
131 dealii::TimerOutput computing_timer;
132 };
133} // namespace dftfe
134# endif
135#endif
void densityMatrixEigenBasisFirstOrderResponse(operatorDFTClass< dftfe::utils::MemorySpace::HOST > &operatorMatrix, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, 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)
Compute first order response in density matrix with respect to perturbation in the Hamiltonian....
@ DEVICE
Definition MemorySpaceType.h:36
Definition pseudoPotentialToDftfeConverter.cc:34