DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
|
Contains repeatedly used functions in the KSDFT calculations. More...
Classes | |
class | CompositeData |
class | constraintMatrixInfo |
Overloads dealii's distribute and distribute_local_to_global functions associated with constraints class. Stores the dealii's constraint matrix data into STL vectors for faster memory access costs. More... | |
class | MPIWriteOnFile |
class | NodalData |
class | Pool |
class | QuadDataCompositeWrite |
Functions | |
template<typename ValueType> | |
void | distributeDevice (const unsigned int contiguousBlockSize, ValueType *xVec, const unsigned int *constraintLocalRowIdsUnflattened, const unsigned int numConstraints, const unsigned int *constraintRowSizes, const unsigned int *constraintRowSizesAccumulated, const unsigned int *constraintLocalColumnIdsAllRowsUnflattened, const double *constraintColumnValuesAllRowsUnflattened, const double *inhomogenities) |
template<typename ValueType> | |
void | distributeSlaveToMasterAtomicAddDevice (const unsigned int contiguousBlockSize, ValueType *xVec, const unsigned int *constraintLocalRowIdsUnflattened, const unsigned int numConstraints, const unsigned int *constraintRowSizes, const unsigned int *constraintRowSizesAccumulated, const unsigned int *constraintLocalColumnIdsAllRowsUnflattened, const double *constraintColumnValuesAllRowsUnflattened) |
template<typename ValueType> | |
void | setzeroDevice (const unsigned int contiguousBlockSize, ValueType *xVec, const unsigned int *constraintLocalRowIdsUnflattened, const unsigned int numConstraints) |
void | scaleConstraintsDevice (const double *xVec, const unsigned int *constraintLocalRowIdsUnflattened, const unsigned int numConstraints, const unsigned int *constraintRowSizes, const unsigned int *constraintRowSizesAccumulated, const unsigned int *constraintLocalColumnIdsAllRowsUnflattened, double *constraintColumnValuesAllRowsUnflattened) |
void | dgesv_ (int *N, int *NRHS, double *A, int *LDA, int *IPIV, double *B, int *LDB, int *INFO) |
double | smearedCharge (double r, double rc) |
double | smearedChargeDr (double r, double rc) |
double | smearedPot (double r, double rc) |
double | smearedPotDr (double r, double rc) |
std::vector< double > | getFractionalCoordinates (const std::vector< double > &latticeVectorsFlattened, const std::vector< double > &coordWithRespectToCellCorner) |
double | getCompositeGeneratorVal (const double rc, const double r, const double a0, const double power) |
Calculates value of composite generator. | |
dealii::BoundingBox< 3 > | createBoundingBoxForSphere (const dealii::Point< 3 > ¢er, const double sphereRadius) |
Create bounding box around a sphere. | |
double | getPartialOccupancy (const double eigenValue, const double fermiEnergy, const double kb, const double T) |
Calculates partial occupancy of the atomic orbital using Fermi-Dirac smearing. | |
double | getPartialOccupancyDer (const double eigenValue, const double fermiEnergy, const double kb, const double T) |
Calculates the derivative of the partial occupancy of the atomic orbital with respect to (x=eigenvalue-fermiEnergy) using Fermi-Dirac smearing. | |
void | cross_product (const std::vector< double > &a, const std::vector< double > &b, std::vector< double > &crossProductVector) |
Calculates cross product of two vectors. | |
void | transformDomainBoundingVectors (std::vector< std::vector< double > > &domainBoundingVectors, const dealii::Tensor< 2, 3, double > &deformationGradient) |
Applies an affine transformation to the domain bounding vectors. | |
void | writeDataVTUParallelLowestPoolId (const dealii::DoFHandler< 3 > &dofHandler, const dealii::DataOut< 3 > &dataOut, const MPI_Comm &mpiCommParent, const MPI_Comm &mpiCommDomain, const MPI_Comm &interpoolcomm, const MPI_Comm &interBandGroupComm, const std::string &folderName, const std::string &fileName) |
Writes to vtu file only from the lowest pool id. | |
void | createBandParallelizationIndices (const MPI_Comm &interBandGroupComm, const unsigned int numBands, std::vector< unsigned int > &bandGroupLowHighPlusOneIndices) |
Create index vector which is used for band parallelization. | |
void | createKpointParallelizationIndices (const MPI_Comm &interKptPoolComm, const int numberIndices, std::vector< int > &kptGroupLowHighPlusOneIndices) |
void | printCurrentMemoryUsage (const MPI_Comm &mpiComm, const std::string message) |
Wrapper to print current memory usage (prints only the maximum across mpiComm) using PetscMemoryGetCurrentUsage. | |
DeclExceptionMsg (ExcNotImplementedYet, "This functionality is not implemented yet or not needed to be implemented.") | |
Exception handler for not implemented functionality. | |
DeclExceptionMsg (ExcInternalError, "DFT-FE internal error.") | |
Exception handler for DFT-FE internal error. | |
void | readFile (const unsigned int numColumns, std::vector< std::vector< double > > &data, const std::string &fileName) |
Read from file containing only double data in columns. | |
void | readFile (std::vector< std::vector< double > > &data, const std::string &fileName) |
Read from file containing only double data in columns. | |
int | readPsiFile (const unsigned int numColumns, std::vector< std::vector< double > > &data, const std::string &fileName) |
Read from file containing only double data in columns. | |
void | writeDataIntoFile (const std::vector< std::vector< double > > &data, const std::string &fileName, const MPI_Comm &mpi_comm_parent) |
Write data into file containing only double data in rows and columns. | |
void | writeDataIntoFile (const std::vector< std::vector< double > > &data, const std::string &fileName) |
Write data into file containing only double data in rows and columns. | |
void | readRelaxationFlagsFile (const unsigned int numColumns, std::vector< std::vector< int > > &data, std::vector< std::vector< double > > &forceData, const std::string &fileName) |
Read from file containing only integer data in columns. | |
void | moveFile (const std::string &old_name, const std::string &new_name) |
Move/rename checkpoint file. | |
void | copyFile (const std::string &pathold, const std::string &pathnew) |
copy file. | |
void | verifyCheckpointFileExists (const std::string &filename) |
Verify if checkpoint file exists. | |
Contains repeatedly used functions in the KSDFT calculations.
void dftfe::dftUtils::copyFile | ( | const std::string & | pathold, |
const std::string & | pathnew ) |
copy file.
void dftfe::dftUtils::createBandParallelizationIndices | ( | const MPI_Comm & | interBandGroupComm, |
const unsigned int | numBands, | ||
std::vector< unsigned int > & | bandGroupLowHighPlusOneIndices ) |
Create index vector which is used for band parallelization.
@[in]param interBandGroupComm mpi communicator across band groups @[in]param numBands @[out]param bandGroupLowHighPlusOneIndices
dealii::BoundingBox< 3 > dftfe::dftUtils::createBoundingBoxForSphere | ( | const dealii::Point< 3 > & | center, |
const double | sphereRadius ) |
Create bounding box around a sphere.
sphere | center |
sphere | radius |
void dftfe::dftUtils::createKpointParallelizationIndices | ( | const MPI_Comm & | interKptPoolComm, |
const int | numberIndices, | ||
std::vector< int > & | kptGroupLowHighPlusOneIndices ) |
void dftfe::dftUtils::cross_product | ( | const std::vector< double > & | a, |
const std::vector< double > & | b, | ||
std::vector< double > & | crossProductVector ) |
Calculates cross product of two vectors.
a | first vector |
b | second vector |
crossProductVector | cross product of a and b |
dftfe::dftUtils::DeclExceptionMsg | ( | ExcInternalError | , |
"DFT-FE internal error." | ) |
Exception handler for DFT-FE internal error.
dftfe::dftUtils::DeclExceptionMsg | ( | ExcNotImplementedYet | , |
"This functionality is not implemented yet or not needed to be implemented." | ) |
Exception handler for not implemented functionality.
void dftfe::dftUtils::dgesv_ | ( | int * | N, |
int * | NRHS, | ||
double * | A, | ||
int * | LDA, | ||
int * | IPIV, | ||
double * | B, | ||
int * | LDB, | ||
int * | INFO ) |
void dftfe::dftUtils::distributeDevice | ( | const unsigned int | contiguousBlockSize, |
ValueType * | xVec, | ||
const unsigned int * | constraintLocalRowIdsUnflattened, | ||
const unsigned int | numConstraints, | ||
const unsigned int * | constraintRowSizes, | ||
const unsigned int * | constraintRowSizesAccumulated, | ||
const unsigned int * | constraintLocalColumnIdsAllRowsUnflattened, | ||
const double * | constraintColumnValuesAllRowsUnflattened, | ||
const double * | inhomogenities ) |
void dftfe::dftUtils::distributeSlaveToMasterAtomicAddDevice | ( | const unsigned int | contiguousBlockSize, |
ValueType * | xVec, | ||
const unsigned int * | constraintLocalRowIdsUnflattened, | ||
const unsigned int | numConstraints, | ||
const unsigned int * | constraintRowSizes, | ||
const unsigned int * | constraintRowSizesAccumulated, | ||
const unsigned int * | constraintLocalColumnIdsAllRowsUnflattened, | ||
const double * | constraintColumnValuesAllRowsUnflattened ) |
double dftfe::dftUtils::getCompositeGeneratorVal | ( | const double | rc, |
const double | r, | ||
const double | a0, | ||
const double | power ) |
Calculates value of composite generator.
|
inline |
double dftfe::dftUtils::getPartialOccupancy | ( | const double | eigenValue, |
const double | fermiEnergy, | ||
const double | kb, | ||
const double | T ) |
Calculates partial occupancy of the atomic orbital using Fermi-Dirac smearing.
eigenValue | |
fermiEnergy | |
kb | Boltzmann constant |
T | smearing temperature |
double dftfe::dftUtils::getPartialOccupancyDer | ( | const double | eigenValue, |
const double | fermiEnergy, | ||
const double | kb, | ||
const double | T ) |
Calculates the derivative of the partial occupancy of the atomic orbital with respect to (x=eigenvalue-fermiEnergy) using Fermi-Dirac smearing.
eigenValue | |
fermiEnergy | |
kb | Boltzmann constant |
T | smearing temperature |
void dftfe::dftUtils::moveFile | ( | const std::string & | old_name, |
const std::string & | new_name ) |
Move/rename checkpoint file.
void dftfe::dftUtils::printCurrentMemoryUsage | ( | const MPI_Comm & | mpiComm, |
const std::string | message ) |
Wrapper to print current memory usage (prints only the maximum across mpiComm) using PetscMemoryGetCurrentUsage.
@[in]param mpiComm mpi communicator across which the memory printing will be synchronized @[in]param message message to be printed alongwith the memory usage
void dftfe::dftUtils::readFile | ( | const unsigned int | numColumns, |
std::vector< std::vector< double > > & | data, | ||
const std::string & | fileName ) |
Read from file containing only double data in columns.
[in] | numColumns | number of data columsn in the file to be read |
[out] | data | output double data in [rows][columns] format |
[in] | fileName |
void dftfe::dftUtils::readFile | ( | std::vector< std::vector< double > > & | data, |
const std::string & | fileName ) |
Read from file containing only double data in columns.
[out] | data | output double data in [rows][columns] format |
[in] | fileName |
int dftfe::dftUtils::readPsiFile | ( | const unsigned int | numColumns, |
std::vector< std::vector< double > > & | data, | ||
const std::string & | fileName ) |
Read from file containing only double data in columns.
void dftfe::dftUtils::readRelaxationFlagsFile | ( | const unsigned int | numColumns, |
std::vector< std::vector< int > > & | data, | ||
std::vector< std::vector< double > > & | forceData, | ||
const std::string & | fileName ) |
Read from file containing only integer data in columns.
void dftfe::dftUtils::scaleConstraintsDevice | ( | const double * | xVec, |
const unsigned int * | constraintLocalRowIdsUnflattened, | ||
const unsigned int | numConstraints, | ||
const unsigned int * | constraintRowSizes, | ||
const unsigned int * | constraintRowSizesAccumulated, | ||
const unsigned int * | constraintLocalColumnIdsAllRowsUnflattened, | ||
double * | constraintColumnValuesAllRowsUnflattened ) |
void dftfe::dftUtils::setzeroDevice | ( | const unsigned int | contiguousBlockSize, |
ValueType * | xVec, | ||
const unsigned int * | constraintLocalRowIdsUnflattened, | ||
const unsigned int | numConstraints ) |
|
inline |
|
inline |
|
inline |
|
inline |
void dftfe::dftUtils::transformDomainBoundingVectors | ( | std::vector< std::vector< double > > & | domainBoundingVectors, |
const dealii::Tensor< 2, 3, double > & | deformationGradient ) |
Applies an affine transformation to the domain bounding vectors.
d_domainBoundingVectors | the bounding vectors of the domain given as a 2d array |
deformationGradient |
void dftfe::dftUtils::verifyCheckpointFileExists | ( | const std::string & | filename | ) |
Verify if checkpoint file exists.
void dftfe::dftUtils::writeDataIntoFile | ( | const std::vector< std::vector< double > > & | data, |
const std::string & | fileName ) |
Write data into file containing only double data in rows and columns.
[in] | data | input double data in [rows][columns] format |
[in] | fileName |
void dftfe::dftUtils::writeDataIntoFile | ( | const std::vector< std::vector< double > > & | data, |
const std::string & | fileName, | ||
const MPI_Comm & | mpi_comm_parent ) |
Write data into file containing only double data in rows and columns.
[in] | data | input double data in [rows][columns] format |
[in] | fileName | |
[in] | mpi_comm_parent | parent communicator |
void dftfe::dftUtils::writeDataVTUParallelLowestPoolId | ( | const dealii::DoFHandler< 3 > & | dofHandler, |
const dealii::DataOut< 3 > & | dataOut, | ||
const MPI_Comm & | mpiCommParent, | ||
const MPI_Comm & | mpiCommDomain, | ||
const MPI_Comm & | interpoolcomm, | ||
const MPI_Comm & | interBandGroupComm, | ||
const std::string & | folderName, | ||
const std::string & | fileName ) |
Writes to vtu file only from the lowest pool id.
dataOut | DataOut class object |
mpiCommParent | parent mpi communicator |
mpiCommDomain | mpi communicator of domain decomposition inside each pool |
interpoolcomm | mpi communicator across k point pools |
interBandGroupComm | mpi communicator across band groups |
fileName |