DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
ExcDFTPlusU.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#ifndef DFTFE_EXE_EXCDFTPLUSU_H
19#define DFTFE_EXE_EXCDFTPLUSU_H
20
21
22
24#include "hubbardClass.h"
25namespace dftfe
26{
27 template <typename ValueType, dftfe::utils::MemorySpace memorySpace>
28 class ExcDFTPlusU : public ExcSSDFunctionalBaseClass<memorySpace>
29 {
30 public:
32 std::shared_ptr<ExcSSDFunctionalBaseClass<memorySpace>> excSSDObjPtr,
33 dftfe::uInt numSpins);
34
36
37 void
40 &src,
42 const dftfe::uInt inputVecSize,
43 const dftfe::uInt kPointIndex,
44 const dftfe::uInt spinIndex) override;
45
46 /*
47 * @brief The apply function that will be called in HXCheby() with single precision.
48 * The distribute() and updateGhostValues() for src
49 * has to be called before this function.
50 * Similarly for dst, accumulateLocallyOwned() should be called in HX()
51 * after this function is called. param[in] src The input vector param[out]
52 * dst The output vector param[in] inputVecSize The size of the input vector
53 * param[in] kPointIndex the k point for which the HX() is called
54 * param[in] spinIndex the spin index for which the HX() is called
55 */
56 void
60 memorySpace> &src,
63 memorySpace> &dst,
64 const dftfe::uInt inputVecSize,
65 const dftfe::uInt kPointIndex,
66 const dftfe::uInt spinIndex) override;
67
68 void
70 const std::shared_ptr<AuxDensityMatrix<memorySpace>> &auxDensityMatrixPtr,
71 const std::vector<double> &kPointWeights) override;
72
73 void
75 const std::shared_ptr<AuxDensityMatrix<memorySpace>> &auxDensityMatrix,
76 const std::vector<double> &kPointWeights) override;
77
78 double
80
81 double
83
84 /**
85 * x and c denotes exchange and correlation respectively
86 */
87 void
89 AuxDensityMatrix<memorySpace> &auxDensityMatrix,
90 const std::pair<dftfe::uInt, dftfe::uInt> &quadIndexRange,
91 std::unordered_map<xcRemainderOutputDataAttributes, std::vector<double>>
92 &xDataOut,
93 std::unordered_map<xcRemainderOutputDataAttributes, std::vector<double>>
94 &cDataout) const override;
95 void
97 const std::vector<xcRemainderOutputDataAttributes> &outputDataAttributes)
98 const override;
99
100 void
102
103 void
105 const MPI_Comm &mpi_comm_parent,
106 const MPI_Comm &mpi_comm_domain,
107 const MPI_Comm &mpi_comm_interPool,
108 const MPI_Comm &mpi_comm_interBandGroup,
109 std::shared_ptr<
111 basisOperationsMemPtr,
112 std::shared_ptr<
114 FEBasisOperations<ValueType, double, dftfe::utils::MemorySpace::HOST>>
115 basisOperationsHostPtr,
117 BLASWrapperMemPtr,
118 std::shared_ptr<
120 BLASWrapperHostPtr,
121 const dftfe::uInt matrixFreeVectorComponent,
122 const dftfe::uInt densityQuadratureId,
123 const dftfe::uInt sparsityPatternQuadratureId,
124 const dftfe::uInt numberWaveFunctions,
125 const dftfe::uInt numSpins,
126 const dftParameters &dftParam,
127 const std::string &scratchFolderName,
128 const bool singlePrecNonLocalOperator,
129 const bool updateNonlocalSparsity,
130 const std::vector<std::vector<double>> &atomLocations,
131 const std::vector<std::vector<double>> &atomLocationsFrac,
132 const std::vector<dftfe::Int> &imageIds,
133 const std::vector<std::vector<double>> &imagePositions,
134 std::vector<double> &kPointCoordinates,
135 const std::vector<double> &kPointWeights,
136 const std::vector<std::vector<double>> &domainBoundaries);
137
138 std::shared_ptr<hubbard<ValueType, memorySpace>> &
140
141 public:
142 std::shared_ptr<ExcSSDFunctionalBaseClass<memorySpace>> d_excSSDObjPtr;
143 std::shared_ptr<hubbard<ValueType, memorySpace>> d_hubbardClassPtr;
144 };
145} // namespace dftfe
146#endif // DFTFE_EXE_EXCDFTPLUSU_H
Definition AuxDensityMatrix.h:40
void computeRhoTauDependentXCData(AuxDensityMatrix< memorySpace > &auxDensityMatrix, const std::pair< dftfe::uInt, dftfe::uInt > &quadIndexRange, std::unordered_map< xcRemainderOutputDataAttributes, std::vector< double > > &xDataOut, std::unordered_map< xcRemainderOutputDataAttributes, std::vector< double > > &cDataout) const override
void initialiseHubbardClass(const MPI_Comm &mpi_comm_parent, const MPI_Comm &mpi_comm_domain, const MPI_Comm &mpi_comm_interPool, const MPI_Comm &mpi_comm_interBandGroup, std::shared_ptr< dftfe::basis::FEBasisOperations< ValueType, double, memorySpace > > basisOperationsMemPtr, std::shared_ptr< dftfe::basis::FEBasisOperations< ValueType, double, dftfe::utils::MemorySpace::HOST > > basisOperationsHostPtr, std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > BLASWrapperMemPtr, std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > BLASWrapperHostPtr, const dftfe::uInt matrixFreeVectorComponent, const dftfe::uInt densityQuadratureId, const dftfe::uInt sparsityPatternQuadratureId, const dftfe::uInt numberWaveFunctions, const dftfe::uInt numSpins, const dftParameters &dftParam, const std::string &scratchFolderName, const bool singlePrecNonLocalOperator, const bool updateNonlocalSparsity, const std::vector< std::vector< double > > &atomLocations, const std::vector< std::vector< double > > &atomLocationsFrac, const std::vector< dftfe::Int > &imageIds, const std::vector< std::vector< double > > &imagePositions, std::vector< double > &kPointCoordinates, const std::vector< double > &kPointWeights, const std::vector< std::vector< double > > &domainBoundaries)
void reinitKPointDependentVariables(dftfe::uInt kPointIndex) override
std::shared_ptr< ExcSSDFunctionalBaseClass< memorySpace > > d_excSSDObjPtr
Definition ExcDFTPlusU.h:142
void applyWaveFunctionDependentFuncDerWrtPsi(const dftfe::linearAlgebra::MultiVector< typename dataTypes::singlePrecType< ValueType >::type, memorySpace > &src, dftfe::linearAlgebra::MultiVector< typename dataTypes::singlePrecType< ValueType >::type, memorySpace > &dst, const dftfe::uInt inputVecSize, const dftfe::uInt kPointIndex, const dftfe::uInt spinIndex) override
void computeWaveFunctionDependentExcEnergy(const std::shared_ptr< AuxDensityMatrix< memorySpace > > &auxDensityMatrix, const std::vector< double > &kPointWeights) override
void updateWaveFunctionDependentFuncDerWrtPsi(const std::shared_ptr< AuxDensityMatrix< memorySpace > > &auxDensityMatrixPtr, const std::vector< double > &kPointWeights) override
double getWaveFunctionDependentExcEnergy() override
std::shared_ptr< hubbard< ValueType, memorySpace > > & getHubbardClass()
void checkInputOutputDataAttributesConsistency(const std::vector< xcRemainderOutputDataAttributes > &outputDataAttributes) const override
double getExpectationOfWaveFunctionDependentExcFuncDerWrtPsi() override
ExcDFTPlusU(std::shared_ptr< ExcSSDFunctionalBaseClass< memorySpace > > excSSDObjPtr, dftfe::uInt numSpins)
void applyWaveFunctionDependentFuncDerWrtPsi(const dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > &src, dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > &dst, const dftfe::uInt inputVecSize, const dftfe::uInt kPointIndex, const dftfe::uInt spinIndex) override
std::shared_ptr< hubbard< ValueType, memorySpace > > d_hubbardClassPtr
Definition ExcDFTPlusU.h:143
ExcSSDFunctionalBaseClass(const ExcFamilyType excFamType, const densityFamilyType densityFamType, const std::vector< DensityDescriptorDataAttributes > &densityDescriptorAttributesList)
Definition ExcSSDFunctionalBaseClass.t.cc:25
Definition FEBasisOperations.h:85
Namespace which declares the input parameters and the functions to parse them from the input paramete...
Definition dftParameters.h:36
Definition BLASWrapper.h:35
An class template to encapsulate a MultiVector. A MultiVector is a collection of vectors belonging t...
Definition MultiVector.h:127
Definition FEBasisOperations.h:30
Definition pseudoPotentialToDftfeConverter.cc:34
xcRemainderOutputDataAttributes
Definition ExcSSDFunctionalBaseClass.h:58
std::uint32_t uInt
Definition TypeConfig.h:10
T type
Definition dftfeDataTypes.h:113