DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace > Class Template Reference

computes configurational forces in KSDFT More...

#include <force.h>

Public Member Functions

 forceClass (dftClass< FEOrder, FEOrderElectro, memorySpace > *_dftPtr, const MPI_Comm &mpi_comm_parent, const MPI_Comm &mpi_comm_domain, const dftParameters &dftParams)
 Constructor.
 
void initUnmoved (const dealii::Triangulation< 3, 3 > &triangulation, const dealii::Triangulation< 3, 3 > &serialTriangulation, const std::vector< std::vector< double > > &domainBoundingVectors, const bool isElectrostaticsMesh)
 initializes data structures inside forceClass assuming unmoved triangulation.
 
void initMoved (std::vector< const dealii::DoFHandler< 3 > * > &dofHandlerVectorMatrixFree, std::vector< const dealii::AffineConstraints< double > * > &constraintsVectorMatrixFree, const bool isElectrostaticsMesh)
 initializes data structures inside forceClass which depend on the moved mesh.
 
void initPseudoData ()
 initializes and precomputes pseudopotential related data structuers required for configurational force and stress computation.
 
void computeAtomsForces (const dealii::MatrixFree< 3, double > &matrixFreeData, const dispersionCorrection &dispersionCorr, const unsigned int eigenDofHandlerIndex, const unsigned int smearedChargeQuadratureId, const unsigned int lpspQuadratureIdElectro, const dealii::MatrixFree< 3, double > &matrixFreeDataElectro, const unsigned int phiTotDofHandlerIndexElectro, const distributedCPUVec< double > &phiTotRhoOutElectro, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &rhoOutValues, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &gradRhoOutValues, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoTotalOutValuesLpsp, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &gradRhoTotalOutValuesLpsp, const std::map< dealii::CellId, std::vector< double > > &rhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &gradRhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &hessianRhoCoreValues, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &gradRhoCoreAtoms, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &hessianRhoCoreAtoms, const std::map< dealii::CellId, std::vector< double > > &pseudoVLocElectro, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &pseudoVLocAtomsElectro, const dealii::AffineConstraints< double > &hangingPlusPBCConstraintsElectro, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManagerElectro)
 computes the configurational force on all atoms corresponding to a Gaussian generator, which represents perturbation of the underlying space.
 
std::vector< double > & getAtomsForces ()
 returns a copy of the configurational force on all global atoms.
 
void printAtomsForces ()
 prints the currently stored configurational forces on atoms and the Gaussian generator constant used to compute them.
 
void computeStress (const dealii::MatrixFree< 3, double > &matrixFreeData, const dispersionCorrection &dispersionCorr, const unsigned int eigenDofHandlerIndex, const unsigned int smearedChargeQuadratureId, const unsigned int lpspQuadratureIdElectro, const dealii::MatrixFree< 3, double > &matrixFreeDataElectro, const unsigned int phiTotDofHandlerIndexElectro, const distributedCPUVec< double > &phiTotRhoOutElectro, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &rhoOutValues, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &gradRhoOutValues, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoTotalOutValuesLpsp, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &gradRhoTotalOutValuesLpsp, const std::map< dealii::CellId, std::vector< double > > &pseudoVLocElectro, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &pseudoVLocAtomsElectro, const std::map< dealii::CellId, std::vector< double > > &rhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &gradRhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &hessianRhoCoreValues, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &gradRhoCoreAtoms, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &hessianRhoCoreAtoms, const dealii::AffineConstraints< double > &hangingPlusPBCConstraintsElectro, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManagerElectro)
 Update force generator Gaussian constant.
 
void printStress ()
 prints the currently stored configurational stress tensor.
 
dealii::Tensor< 2, 3, double > & getStress ()
 returns a copy of the current stress tensor value.
 

Private Member Functions

void locateAtomCoreNodesForce (const dealii::DoFHandler< 3 > &dofHandlerForce, const dealii::IndexSet &locally_owned_dofsForce, std::map< std::pair< unsigned int, unsigned int >, unsigned int > &atomsForceDofs)
 get the value of Gaussian generator parameter (d_gaussianConstant). Gaussian generator: Gamma(r)= exp(-d_gaussianConstant*r^2).
 
void createBinObjectsForce (const dealii::DoFHandler< 3 > &dofHandler, const dealii::DoFHandler< 3 > &dofHandlerForce, const dealii::AffineConstraints< double > &hangingPlusPBCConstraints, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManager, std::vector< std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > > &cellsVselfBallsDofHandler, std::vector< std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > > &cellsVselfBallsDofHandlerForce, std::vector< std::map< dealii::CellId, unsigned int > > &cellsVselfBallsClosestAtomIdDofHandler, std::map< unsigned int, unsigned int > &AtomIdBinIdLocalDofHandler, std::vector< std::map< dealii::DoFHandler< 3 >::active_cell_iterator, std::vector< unsigned int > > > &cellFacesVselfBallSurfacesDofHandler, std::vector< std::map< dealii::DoFHandler< 3 >::active_cell_iterator, std::vector< unsigned int > > > &cellFacesVselfBallSurfacesDofHandlerForce)
 
void configForceLinFEInit (const dealii::MatrixFree< 3, double > &matrixFreeData, const dealii::MatrixFree< 3, double > &matrixFreeDataElectro)
 
void configForceLinFEFinalize ()
 
void computeConfigurationalForceEEshelbyTensorFPSPFnlLinFE (const dealii::MatrixFree< 3, double > &matrixFreeData, const unsigned int eigenDofHandlerIndex, const unsigned int smearedChargeQuadratureId, const unsigned int lpspQuadratureIdElectro, const dealii::MatrixFree< 3, double > &matrixFreeDataElectro, const unsigned int phiTotDofHandlerIndexElectro, const distributedCPUVec< double > &phiTotRhoOutElectro, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &rhoOutValues, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &gradRhoOutValues, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoTotalOutValuesLpsp, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &gradRhoTotalOutValuesLpsp, const std::map< dealii::CellId, std::vector< double > > &rhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &gradRhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &hessianRhoCoreValues, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &gradRhoCoreAtoms, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &hessianRhoCoreAtoms, const std::map< dealii::CellId, std::vector< double > > &pseudoVLocElectro, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &pseudoVLocAtomsElectro, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManagerElectro)
 
void computeConfigurationalForceEEshelbyEElectroPhiTot (const dealii::MatrixFree< 3, double > &matrixFreeDataElectro, const unsigned int phiTotDofHandlerIndexElectro, const unsigned int smearedChargeQuadratureId, const unsigned int lpspQuadratureIdElectro, const distributedCPUVec< double > &phiTotRhoOutElectro, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoTotalOutValues, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoTotalOutValuesLpsp, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &gradRhoTotalOutValues, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &gradRhoTotalOutValuesLpsp, const std::map< dealii::CellId, std::vector< double > > &pseudoVLocElectro, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &pseudoVLocAtomsElectro, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManagerElectro)
 
void computeConfigurationalForcePhiExtLinFE ()
 
void computeConfigurationalForceEselfLinFE (const dealii::DoFHandler< 3 > &dofHandlerElectro, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManagerElectro, const dealii::MatrixFree< 3, double > &matrixFreeDataElectro, const unsigned int smearedChargeQuadratureId)
 
void computeConfigurationalForceEselfNoSurfaceLinFE ()
 
void computeConfigurationalForceTotalLinFE (const dealii::MatrixFree< 3, double > &matrixFreeData, const unsigned int eigenDofHandlerIndex, const unsigned int smearedChargeQuadratureId, const unsigned int lpspQuadratureIdElectro, const dealii::MatrixFree< 3, double > &matrixFreeDataElectro, const unsigned int phiTotDofHandlerIndexElectro, const distributedCPUVec< double > &phiTotRhoOutElectro, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &rhoOutValues, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &gradRhoOutValues, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoTotalOutValuesLpsp, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &gradRhoTotalOutValuesLpsp, const std::map< dealii::CellId, std::vector< double > > &rhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &gradRhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &hessianRhoCoreValues, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &gradRhoCoreAtoms, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &hessianRhoCoreAtoms, const std::map< dealii::CellId, std::vector< double > > &pseudoVLocElectro, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &pseudoVLocAtomsElectro, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManagerElectro)
 
void FPSPLocalGammaAtomsElementalContribution (std::map< unsigned int, std::vector< double > > &forceContributionFPSPLocalGammaAtoms, dealii::FEValues< 3 > &feValues, dealii::FEFaceValues< 3 > &feFaceValues, dealii::FEEvaluation< 3, 1, C_num1DQuadLPSP< FEOrder >() *C_numCopies1DQuadLPSP(), 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const unsigned int phiTotDofHandlerIndexElectro, const unsigned int cell, const dealii::AlignedVector< dealii::VectorizedArray< double > > &rhoQuads, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &gradRhoQuads, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &pseudoVLocAtoms, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManager, const std::vector< std::map< dealii::CellId, unsigned int > > &cellsVselfBallsClosestAtomIdDofHandler)
 
void FPhiTotSmearedChargesGammaAtomsElementalContribution (std::map< unsigned int, std::vector< double > > &forceContributionSmearedChargesGammaAtoms, dealii::FEEvaluation< 3, -1, 1, 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const unsigned int cell, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &gradPhiTotQuads, const std::vector< unsigned int > &nonTrivialAtomIdsMacroCell, const std::map< dealii::CellId, std::vector< int > > &bQuadAtomIdsAllAtoms, const dealii::AlignedVector< dealii::VectorizedArray< double > > &smearedbQuads)
 
void FVselfSmearedChargesGammaAtomsElementalContribution (std::map< unsigned int, std::vector< double > > &forceContributionSmearedChargesGammaAtoms, dealii::FEEvaluation< 3, -1, 1, 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const unsigned int cell, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &gradVselfBinQuads, const std::vector< unsigned int > &nonTrivialAtomIdsMacroCell, const std::map< dealii::CellId, std::vector< int > > &bQuadAtomIdsAllAtoms, const dealii::AlignedVector< dealii::VectorizedArray< double > > &smearedbQuads)
 
void FNonlinearCoreCorrectionGammaAtomsElementalContribution (std::map< unsigned int, std::vector< double > > &forceContributionFNonlinearCoreCorrectionGammaAtoms, dealii::FEEvaluation< 3, 1, C_num1DQuad< C_rhoNodalPolyOrder< FEOrder, FEOrderElectro >()>(), 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const unsigned int cell, const dealii::AlignedVector< dealii::VectorizedArray< double > > &vxcQuads, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &gradRhoCoreAtoms)
 
void FNonlinearCoreCorrectionGammaAtomsElementalContribution (std::map< unsigned int, std::vector< double > > &forceContributionFNonlinearCoreCorrectionGammaAtoms, dealii::FEEvaluation< 3, 1, C_num1DQuad< C_rhoNodalPolyOrder< FEOrder, FEOrderElectro >()>(), 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const unsigned int cell, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &derExcGradRho, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &hessianRhoCoreAtoms)
 
void FNonlinearCoreCorrectionGammaAtomsElementalContributionSpinPolarized (std::map< unsigned int, std::vector< double > > &forceContributionFNonlinearCoreCorrectionGammaAtoms, dealii::FEEvaluation< 3, 1, C_num1DQuad< C_rhoNodalPolyOrder< FEOrder, FEOrderElectro >()>(), 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const unsigned int cell, const dealii::AlignedVector< dealii::VectorizedArray< double > > &vxcQuadsSpin0, const dealii::AlignedVector< dealii::VectorizedArray< double > > &vxcQuadsSpin1, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &derExcGradRhoSpin0, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &derExcGradRhoSpin1, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &gradRhoCoreAtoms, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &hessianRhoCoreAtoms, const bool isXCGGA=false)
 
void distributeForceContributionFPSPLocalGammaAtoms (const std::map< unsigned int, std::vector< double > > &forceContributionFPSPLocalGammaAtoms, const std::map< std::pair< unsigned int, unsigned int >, unsigned int > &atomsForceDofs, const dealii::AffineConstraints< double > &constraintsNoneForce, distributedCPUVec< double > &configForceVectorLinFE)
 
void accumulateForceContributionGammaAtomsFloating (const std::map< unsigned int, std::vector< double > > &forceContributionLocalGammaAtoms, std::vector< double > &accumForcesVector)
 
void FnlGammaAtomsElementalContribution (std::map< unsigned int, std::vector< double > > &forceContributionFnlGammaAtoms, const dealii::MatrixFree< 3, double > &matrixFreeData, dealii::FEEvaluation< 3, 1, C_num1DQuadNLPSP< FEOrder >() *C_numCopies1DQuadNLPSP(), 3 > &forceEvalNLP, const std::shared_ptr< AtomicCenteredNonLocalOperator< dataTypes::number, memorySpace > > nonLocalOp, unsigned int numNonLocalAtomsCurrentProcess, const std::vector< int > &globalChargeIdNonLocalAtoms, const std::vector< unsigned int > &numberPseudoWaveFunctionsPerAtom, const unsigned int cell, const std::map< dealii::CellId, unsigned int > &cellIdToCellNumberMap, const std::vector< dataTypes::number > &zetaDeltaVQuadsFlattened, const std::vector< dataTypes::number > &projectorKetTimesPsiTimesVTimesPartOccContractionGradPsiQuadsFlattened)
 
void FnlGammaxElementalContribution (dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &FVectQuads, const dealii::MatrixFree< 3, double > &matrixFreeData, const unsigned int numQuadPoints, const std::shared_ptr< AtomicCenteredNonLocalOperator< dataTypes::number, memorySpace > > nonLocalOp, const unsigned int numNonLocalAtomsCurrentProcess, const std::vector< int > &globalChargeIdNonLocalAtoms, const std::vector< unsigned int > &numberPseudoWaveFunctionsPerAtom, const unsigned int cell, const std::map< dealii::CellId, unsigned int > &cellIdToCellNumberMap, const std::vector< dataTypes::number > &zetaDeltaVQuadsFlattened, const std::vector< dataTypes::number > &projectorKetTimesPsiTimesVTimesPartOccContractionGradPsiQuadsFlattened)
 
void distributeForceContributionFnlGammaAtoms (const std::map< unsigned int, std::vector< double > > &forceContributionFnlGammaAtoms)
 
void stressEnlElementalContribution (dealii::Tensor< 2, 3, double > &stressContribution, const dealii::MatrixFree< 3, double > &matrixFreeData, const unsigned int numQuadPoints, const std::vector< double > &jxwQuadsSubCells, const unsigned int cell, const unsigned int numNonLocalAtomsCurrentProcess, const std::shared_ptr< AtomicCenteredNonLocalOperator< dataTypes::number, memorySpace > > nonLocalOp, const std::vector< unsigned int > &numberPseudoWaveFunctionsPerAtom, const std::map< dealii::CellId, unsigned int > &cellIdToCellNumberMap, const std::vector< dataTypes::number > &zetalmDeltaVlProductDistImageAtoms, const std::vector< dataTypes::number > &projectorKetTimesPsiTimesVTimesPartOccContractionGradPsiQuadsFlattened, const bool isSpinPolarized)
 
void computeAtomsForcesGaussianGenerator (bool allowGaussianOverlapOnAtoms=false)
 
void computeFloatingAtomsForces ()
 
void computeStressEself (const dealii::DoFHandler< 3 > &dofHandlerElectro, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManagerElectro, const dealii::MatrixFree< 3, double > &matrixFreeDataElectro, const unsigned int smearedChargeQuadratureId)
 
void computeStressEEshelbyEPSPEnlEk (const dealii::MatrixFree< 3, double > &matrixFreeData, const unsigned int eigenDofHandlerIndex, const unsigned int smearedChargeQuadratureId, const unsigned int lpspQuadratureIdElectro, const dealii::MatrixFree< 3, double > &matrixFreeDataElectro, const unsigned int phiTotDofHandlerIndexElectro, const distributedCPUVec< double > &phiTotRhoOutElectro, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &rhoOutValues, const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &gradRhoOutValues, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoTotalOutValuesLpsp, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &gradRhoTotalOutValuesLpsp, const std::map< dealii::CellId, std::vector< double > > &pseudoVLocElectro, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &pseudoVLocAtomsElectro, const std::map< dealii::CellId, std::vector< double > > &rhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &gradRhoCoreValues, const std::map< dealii::CellId, std::vector< double > > &hessianRhoCoreValues, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &gradRhoCoreAtoms, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &hessianRhoCoreAtoms, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManagerElectro)
 
void computeStressEEshelbyEElectroPhiTot (const dealii::MatrixFree< 3, double > &matrixFreeDataElectro, const unsigned int phiTotDofHandlerIndexElectro, const unsigned int smearedChargeQuadratureId, const unsigned int lpspQuadratureIdElectro, const distributedCPUVec< double > &phiTotRhoOutElectro, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoTotalOutValues, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &rhoTotalOutValuesLpsp, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &gradRhoTotalOutValuesElectro, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &gradRhoTotalOutValuesElectroLpsp, const std::map< dealii::CellId, std::vector< double > > &pseudoVLocElectro, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &pseudoVLocAtomsElectro, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManagerElectro)
 
void addEPSPStressContribution (dealii::FEValues< 3 > &feValues, dealii::FEFaceValues< 3 > &feFaceValues, dealii::FEEvaluation< 3, 1, C_num1DQuadLPSP< FEOrder >() *C_numCopies1DQuadLPSP(), 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const unsigned int phiTotDofHandlerIndexElectro, const unsigned int cell, const dealii::AlignedVector< dealii::VectorizedArray< double > > &rhoQuads, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &gradRhoQuads, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &pseudoVLocAtoms, const vselfBinsManager< FEOrder, FEOrderElectro > &vselfBinsManager, const std::vector< std::map< dealii::CellId, unsigned int > > &cellsVselfBallsClosestAtomIdDofHandler)
 
void addENonlinearCoreCorrectionStressContribution (dealii::FEEvaluation< 3, 1, C_num1DQuad< C_rhoNodalPolyOrder< FEOrder, FEOrderElectro >()>(), 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const unsigned int cell, const dealii::AlignedVector< dealii::VectorizedArray< double > > &vxcQuads, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &derExcGradRho, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &gradRhoCoreAtoms, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &hessianRhoCoreAtoms)
 
void addENonlinearCoreCorrectionStressContributionSpinPolarized (dealii::FEEvaluation< 3, 1, C_num1DQuad< C_rhoNodalPolyOrder< FEOrder, FEOrderElectro >()>(), 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const unsigned int cell, const dealii::AlignedVector< dealii::VectorizedArray< double > > &vxcQuadsSpin0, const dealii::AlignedVector< dealii::VectorizedArray< double > > &vxcQuadsSpin1, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &derExcGradRhoSpin0, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &derExcGradRhoSpin1, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &gradRhoCoreAtoms, const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > &hessianRhoCoreAtoms, const bool isXCGGA=false)
 
void addEPhiTotSmearedStressContribution (dealii::FEEvaluation< 3, -1, 1, 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const unsigned int cell, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &gradPhiTotQuads, const std::vector< unsigned int > &nonTrivialAtomImageIdsMacroCell, const std::map< dealii::CellId, std::vector< int > > &bQuadAtomIdsAllAtomsImages, const dealii::AlignedVector< dealii::VectorizedArray< double > > &smearedbQuads)
 
void addEVselfSmearedStressContribution (dealii::FEEvaluation< 3, -1, 1, 3 > &forceEval, const dealii::MatrixFree< 3, double > &matrixFreeData, const unsigned int cell, const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > &gradVselfQuads, const std::vector< unsigned int > &nonTrivialAtomImageIdsMacroCell, const std::map< dealii::CellId, std::vector< int > > &bQuadAtomIdsAllAtomsImages, const dealii::AlignedVector< dealii::VectorizedArray< double > > &smearedbQuads)
 
void computeElementalNonLocalPseudoOVDataForce ()
 

Private Attributes

distributedCPUVec< double > d_configForceVectorLinFE
 
distributedCPUVec< double > d_configForceVectorLinFEElectro
 
std::vector< double > d_forceAtomsFloating
 
std::vector< double > d_globalAtomsForces
 Storage for configurational force on all global atoms.
 
dealii::Tensor< 2, 3, double > d_stress
 Storage for configurational stress tensor.
 
dealii::Tensor< 2, 3, double > d_stressKPoints
 
const bool d_allowGaussianOverlapOnAtoms = false
 
dftClass< FEOrder, FEOrderElectro, memorySpace > * dftPtr
 pointer to dft class
 
dealii::FESystem< 3 > FEForce
 
dealii::DoFHandler< 3 > d_dofHandlerForce
 
dealii::DoFHandler< 3 > d_dofHandlerForceElectro
 
unsigned int d_forceDofHandlerIndex
 
unsigned int d_forceDofHandlerIndexElectro
 
dealii::IndexSet d_locally_owned_dofsForce
 
dealii::IndexSet d_locally_owned_dofsForceElectro
 
dealii::IndexSet d_locally_relevant_dofsForce
 
dealii::IndexSet d_locally_relevant_dofsForceElectro
 
dealii::AffineConstraints< double > d_constraintsNoneForce
 
dealii::AffineConstraints< double > d_constraintsNoneForceElectro
 
std::map< std::pair< unsigned int, unsigned int >, unsigned int > d_atomsForceDofs
 
std::map< std::pair< unsigned int, unsigned int >, unsigned int > d_atomsForceDofsElectro
 
std::vector< std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > > d_cellsVselfBallsDofHandler
 
std::vector< std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > > d_cellsVselfBallsDofHandlerForce
 
std::vector< std::map< dealii::CellId, unsigned int > > d_cellsVselfBallsClosestAtomIdDofHandler
 
std::map< unsigned int, unsigned int > d_AtomIdBinIdLocalDofHandler
 
std::vector< std::map< dealii::DoFHandler< 3 >::active_cell_iterator, std::vector< unsigned int > > > d_cellFacesVselfBallSurfacesDofHandler
 
std::vector< std::map< dealii::DoFHandler< 3 >::active_cell_iterator, std::vector< unsigned int > > > d_cellFacesVselfBallSurfacesDofHandlerForce
 
std::vector< std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > > d_cellsVselfBallsDofHandlerElectro
 
std::vector< std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > > d_cellsVselfBallsDofHandlerForceElectro
 
std::vector< std::map< dealii::CellId, unsigned int > > d_cellsVselfBallsClosestAtomIdDofHandlerElectro
 
std::map< unsigned int, unsigned int > d_AtomIdBinIdLocalDofHandlerElectro
 
std::vector< std::map< dealii::DoFHandler< 3 >::active_cell_iterator, std::vector< unsigned int > > > d_cellFacesVselfBallSurfacesDofHandlerElectro
 
std::vector< std::map< dealii::DoFHandler< 3 >::active_cell_iterator, std::vector< unsigned int > > > d_cellFacesVselfBallSurfacesDofHandlerForceElectro
 
std::map< dealii::CellId, dealii::DoFHandler< 3 >::active_cell_iterator > d_cellIdToActiveCellIteratorMapDofHandlerRhoNodalElectro
 
std::vector< distributedCPUVec< double > > d_gaussianWeightsVecAtoms
 
const MPI_Comm d_mpiCommParent
 domain decomposition mpi_communicator
 
const MPI_Comm mpi_communicator
 domain decomposition mpi_communicator
 
const dftParametersd_dftParams
 
const unsigned int n_mpi_processes
 number of mpi processes in the current pool
 
const unsigned int this_mpi_process
 current mpi process id in the current pool
 
dealii::ConditionalOStream pcout
 

Friends

class dftClass< FEOrder, FEOrderElectro, memorySpace >
 

Detailed Description

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
class dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >

computes configurational forces in KSDFT

This class computes and stores the configurational forces corresponding to geometry optimization. It uses the formulation in the paper by Motamarri et.al. (https://link.aps.org/doi/10.1103/PhysRevB.97.165132) which provides an unified approach to atomic forces corresponding to internal atomic relaxation and cell stress corresponding to cell relaxation.

Author
Sambit Das

Constructor & Destructor Documentation

◆ forceClass()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::forceClass ( dftClass< FEOrder, FEOrderElectro, memorySpace > * _dftPtr,
const MPI_Comm & mpi_comm_parent,
const MPI_Comm & mpi_comm_domain,
const dftParameters & dftParams )

Constructor.

Parameters
_dftPtrpointer to dftClass
mpi_comm_parentparent mpi_communicator
mpi_comm_domaindomain decomposition mpi_communicator

Member Function Documentation

◆ accumulateForceContributionGammaAtomsFloating()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::accumulateForceContributionGammaAtomsFloating ( const std::map< unsigned int, std::vector< double > > & forceContributionLocalGammaAtoms,
std::vector< double > & accumForcesVector )
private

◆ addENonlinearCoreCorrectionStressContribution()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::addENonlinearCoreCorrectionStressContribution ( dealii::FEEvaluation< 3, 1, C_num1DQuad< C_rhoNodalPolyOrder< FEOrder, FEOrderElectro >()>(), 3 > & forceEval,
const dealii::MatrixFree< 3, double > & matrixFreeData,
const unsigned int cell,
const dealii::AlignedVector< dealii::VectorizedArray< double > > & vxcQuads,
const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > & derExcGradRho,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & gradRhoCoreAtoms,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & hessianRhoCoreAtoms )
private

◆ addENonlinearCoreCorrectionStressContributionSpinPolarized()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::addENonlinearCoreCorrectionStressContributionSpinPolarized ( dealii::FEEvaluation< 3, 1, C_num1DQuad< C_rhoNodalPolyOrder< FEOrder, FEOrderElectro >()>(), 3 > & forceEval,
const dealii::MatrixFree< 3, double > & matrixFreeData,
const unsigned int cell,
const dealii::AlignedVector< dealii::VectorizedArray< double > > & vxcQuadsSpin0,
const dealii::AlignedVector< dealii::VectorizedArray< double > > & vxcQuadsSpin1,
const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > & derExcGradRhoSpin0,
const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > & derExcGradRhoSpin1,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & gradRhoCoreAtoms,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & hessianRhoCoreAtoms,
const bool isXCGGA = false )
private

◆ addEPhiTotSmearedStressContribution()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::addEPhiTotSmearedStressContribution ( dealii::FEEvaluation< 3, -1, 1, 3 > & forceEval,
const dealii::MatrixFree< 3, double > & matrixFreeData,
const unsigned int cell,
const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > & gradPhiTotQuads,
const std::vector< unsigned int > & nonTrivialAtomImageIdsMacroCell,
const std::map< dealii::CellId, std::vector< int > > & bQuadAtomIdsAllAtomsImages,
const dealii::AlignedVector< dealii::VectorizedArray< double > > & smearedbQuads )
private

◆ addEPSPStressContribution()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::addEPSPStressContribution ( dealii::FEValues< 3 > & feValues,
dealii::FEFaceValues< 3 > & feFaceValues,
dealii::FEEvaluation< 3, 1, C_num1DQuadLPSP< FEOrder >() *C_numCopies1DQuadLPSP(), 3 > & forceEval,
const dealii::MatrixFree< 3, double > & matrixFreeData,
const unsigned int phiTotDofHandlerIndexElectro,
const unsigned int cell,
const dealii::AlignedVector< dealii::VectorizedArray< double > > & rhoQuads,
const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > & gradRhoQuads,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & pseudoVLocAtoms,
const vselfBinsManager< FEOrder, FEOrderElectro > & vselfBinsManager,
const std::vector< std::map< dealii::CellId, unsigned int > > & cellsVselfBallsClosestAtomIdDofHandler )
private

◆ addEVselfSmearedStressContribution()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::addEVselfSmearedStressContribution ( dealii::FEEvaluation< 3, -1, 1, 3 > & forceEval,
const dealii::MatrixFree< 3, double > & matrixFreeData,
const unsigned int cell,
const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > & gradVselfQuads,
const std::vector< unsigned int > & nonTrivialAtomImageIdsMacroCell,
const std::map< dealii::CellId, std::vector< int > > & bQuadAtomIdsAllAtomsImages,
const dealii::AlignedVector< dealii::VectorizedArray< double > > & smearedbQuads )
private

◆ computeAtomsForces()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::computeAtomsForces ( const dealii::MatrixFree< 3, double > & matrixFreeData,
const dispersionCorrection & dispersionCorr,
const unsigned int eigenDofHandlerIndex,
const unsigned int smearedChargeQuadratureId,
const unsigned int lpspQuadratureIdElectro,
const dealii::MatrixFree< 3, double > & matrixFreeDataElectro,
const unsigned int phiTotDofHandlerIndexElectro,
const distributedCPUVec< double > & phiTotRhoOutElectro,
const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > & rhoOutValues,
const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > & gradRhoOutValues,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & rhoTotalOutValuesLpsp,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & gradRhoTotalOutValuesLpsp,
const std::map< dealii::CellId, std::vector< double > > & rhoCoreValues,
const std::map< dealii::CellId, std::vector< double > > & gradRhoCoreValues,
const std::map< dealii::CellId, std::vector< double > > & hessianRhoCoreValues,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & gradRhoCoreAtoms,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & hessianRhoCoreAtoms,
const std::map< dealii::CellId, std::vector< double > > & pseudoVLocElectro,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & pseudoVLocAtomsElectro,
const dealii::AffineConstraints< double > & hangingPlusPBCConstraintsElectro,
const vselfBinsManager< FEOrder, FEOrderElectro > & vselfBinsManagerElectro )

computes the configurational force on all atoms corresponding to a Gaussian generator, which represents perturbation of the underlying space.

The Gaussian generator is taken to be exp(-d_gaussianConstant*r^2), r being the distance from the atom. Currently d_gaussianConstant is hardcoded to be 4.0. To get the computed atomic forces use getAtomsForces

Returns
void.

◆ computeAtomsForcesGaussianGenerator()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::computeAtomsForcesGaussianGenerator ( bool allowGaussianOverlapOnAtoms = false)
private

◆ computeConfigurationalForceEEshelbyEElectroPhiTot()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::computeConfigurationalForceEEshelbyEElectroPhiTot ( const dealii::MatrixFree< 3, double > & matrixFreeDataElectro,
const unsigned int phiTotDofHandlerIndexElectro,
const unsigned int smearedChargeQuadratureId,
const unsigned int lpspQuadratureIdElectro,
const distributedCPUVec< double > & phiTotRhoOutElectro,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & rhoTotalOutValues,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & rhoTotalOutValuesLpsp,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & gradRhoTotalOutValues,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & gradRhoTotalOutValuesLpsp,
const std::map< dealii::CellId, std::vector< double > > & pseudoVLocElectro,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & pseudoVLocAtomsElectro,
const vselfBinsManager< FEOrder, FEOrderElectro > & vselfBinsManagerElectro )
private

◆ computeConfigurationalForceEEshelbyTensorFPSPFnlLinFE()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::computeConfigurationalForceEEshelbyTensorFPSPFnlLinFE ( const dealii::MatrixFree< 3, double > & matrixFreeData,
const unsigned int eigenDofHandlerIndex,
const unsigned int smearedChargeQuadratureId,
const unsigned int lpspQuadratureIdElectro,
const dealii::MatrixFree< 3, double > & matrixFreeDataElectro,
const unsigned int phiTotDofHandlerIndexElectro,
const distributedCPUVec< double > & phiTotRhoOutElectro,
const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > & rhoOutValues,
const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > & gradRhoOutValues,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & rhoTotalOutValuesLpsp,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & gradRhoTotalOutValuesLpsp,
const std::map< dealii::CellId, std::vector< double > > & rhoCoreValues,
const std::map< dealii::CellId, std::vector< double > > & gradRhoCoreValues,
const std::map< dealii::CellId, std::vector< double > > & hessianRhoCoreValues,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & gradRhoCoreAtoms,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & hessianRhoCoreAtoms,
const std::map< dealii::CellId, std::vector< double > > & pseudoVLocElectro,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & pseudoVLocAtomsElectro,
const vselfBinsManager< FEOrder, FEOrderElectro > & vselfBinsManagerElectro )
private

◆ computeConfigurationalForceEselfLinFE()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::computeConfigurationalForceEselfLinFE ( const dealii::DoFHandler< 3 > & dofHandlerElectro,
const vselfBinsManager< FEOrder, FEOrderElectro > & vselfBinsManagerElectro,
const dealii::MatrixFree< 3, double > & matrixFreeDataElectro,
const unsigned int smearedChargeQuadratureId )
private

◆ computeConfigurationalForceEselfNoSurfaceLinFE()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::computeConfigurationalForceEselfNoSurfaceLinFE ( )
private

◆ computeConfigurationalForcePhiExtLinFE()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::computeConfigurationalForcePhiExtLinFE ( )
private

◆ computeConfigurationalForceTotalLinFE()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::computeConfigurationalForceTotalLinFE ( const dealii::MatrixFree< 3, double > & matrixFreeData,
const unsigned int eigenDofHandlerIndex,
const unsigned int smearedChargeQuadratureId,
const unsigned int lpspQuadratureIdElectro,
const dealii::MatrixFree< 3, double > & matrixFreeDataElectro,
const unsigned int phiTotDofHandlerIndexElectro,
const distributedCPUVec< double > & phiTotRhoOutElectro,
const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > & rhoOutValues,
const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > & gradRhoOutValues,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & rhoTotalOutValuesLpsp,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & gradRhoTotalOutValuesLpsp,
const std::map< dealii::CellId, std::vector< double > > & rhoCoreValues,
const std::map< dealii::CellId, std::vector< double > > & gradRhoCoreValues,
const std::map< dealii::CellId, std::vector< double > > & hessianRhoCoreValues,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & gradRhoCoreAtoms,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & hessianRhoCoreAtoms,
const std::map< dealii::CellId, std::vector< double > > & pseudoVLocElectro,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & pseudoVLocAtomsElectro,
const vselfBinsManager< FEOrder, FEOrderElectro > & vselfBinsManagerElectro )
private

◆ computeElementalNonLocalPseudoOVDataForce()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::computeElementalNonLocalPseudoOVDataForce ( )
private

◆ computeFloatingAtomsForces()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::computeFloatingAtomsForces ( )
private

◆ computeStress()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::computeStress ( const dealii::MatrixFree< 3, double > & matrixFreeData,
const dispersionCorrection & dispersionCorr,
const unsigned int eigenDofHandlerIndex,
const unsigned int smearedChargeQuadratureId,
const unsigned int lpspQuadratureIdElectro,
const dealii::MatrixFree< 3, double > & matrixFreeDataElectro,
const unsigned int phiTotDofHandlerIndexElectro,
const distributedCPUVec< double > & phiTotRhoOutElectro,
const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > & rhoOutValues,
const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > & gradRhoOutValues,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & rhoTotalOutValuesLpsp,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & gradRhoTotalOutValuesLpsp,
const std::map< dealii::CellId, std::vector< double > > & pseudoVLocElectro,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & pseudoVLocAtomsElectro,
const std::map< dealii::CellId, std::vector< double > > & rhoCoreValues,
const std::map< dealii::CellId, std::vector< double > > & gradRhoCoreValues,
const std::map< dealii::CellId, std::vector< double > > & hessianRhoCoreValues,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & gradRhoCoreAtoms,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & hessianRhoCoreAtoms,
const dealii::AffineConstraints< double > & hangingPlusPBCConstraintsElectro,
const vselfBinsManager< FEOrder, FEOrderElectro > & vselfBinsManagerElectro )

Update force generator Gaussian constant.

Returns
void.

computes the configurational stress on the domain corresponding to affine deformation of the periodic cell.

This function cannot be called for fully non-periodic calculations.

Returns
void.

◆ computeStressEEshelbyEElectroPhiTot()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::computeStressEEshelbyEElectroPhiTot ( const dealii::MatrixFree< 3, double > & matrixFreeDataElectro,
const unsigned int phiTotDofHandlerIndexElectro,
const unsigned int smearedChargeQuadratureId,
const unsigned int lpspQuadratureIdElectro,
const distributedCPUVec< double > & phiTotRhoOutElectro,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & rhoTotalOutValues,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & rhoTotalOutValuesLpsp,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & gradRhoTotalOutValuesElectro,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & gradRhoTotalOutValuesElectroLpsp,
const std::map< dealii::CellId, std::vector< double > > & pseudoVLocElectro,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & pseudoVLocAtomsElectro,
const vselfBinsManager< FEOrder, FEOrderElectro > & vselfBinsManagerElectro )
private

◆ computeStressEEshelbyEPSPEnlEk()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::computeStressEEshelbyEPSPEnlEk ( const dealii::MatrixFree< 3, double > & matrixFreeData,
const unsigned int eigenDofHandlerIndex,
const unsigned int smearedChargeQuadratureId,
const unsigned int lpspQuadratureIdElectro,
const dealii::MatrixFree< 3, double > & matrixFreeDataElectro,
const unsigned int phiTotDofHandlerIndexElectro,
const distributedCPUVec< double > & phiTotRhoOutElectro,
const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > & rhoOutValues,
const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > & gradRhoOutValues,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & rhoTotalOutValuesLpsp,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & gradRhoTotalOutValuesLpsp,
const std::map< dealii::CellId, std::vector< double > > & pseudoVLocElectro,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & pseudoVLocAtomsElectro,
const std::map< dealii::CellId, std::vector< double > > & rhoCoreValues,
const std::map< dealii::CellId, std::vector< double > > & gradRhoCoreValues,
const std::map< dealii::CellId, std::vector< double > > & hessianRhoCoreValues,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & gradRhoCoreAtoms,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & hessianRhoCoreAtoms,
const vselfBinsManager< FEOrder, FEOrderElectro > & vselfBinsManagerElectro )
private

◆ computeStressEself()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::computeStressEself ( const dealii::DoFHandler< 3 > & dofHandlerElectro,
const vselfBinsManager< FEOrder, FEOrderElectro > & vselfBinsManagerElectro,
const dealii::MatrixFree< 3, double > & matrixFreeDataElectro,
const unsigned int smearedChargeQuadratureId )
private

◆ configForceLinFEFinalize()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::configForceLinFEFinalize ( )
private

◆ configForceLinFEInit()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::configForceLinFEInit ( const dealii::MatrixFree< 3, double > & matrixFreeData,
const dealii::MatrixFree< 3, double > & matrixFreeDataElectro )
private

◆ createBinObjectsForce()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::createBinObjectsForce ( const dealii::DoFHandler< 3 > & dofHandler,
const dealii::DoFHandler< 3 > & dofHandlerForce,
const dealii::AffineConstraints< double > & hangingPlusPBCConstraints,
const vselfBinsManager< FEOrder, FEOrderElectro > & vselfBinsManager,
std::vector< std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > > & cellsVselfBallsDofHandler,
std::vector< std::vector< dealii::DoFHandler< 3 >::active_cell_iterator > > & cellsVselfBallsDofHandlerForce,
std::vector< std::map< dealii::CellId, unsigned int > > & cellsVselfBallsClosestAtomIdDofHandler,
std::map< unsigned int, unsigned int > & AtomIdBinIdLocalDofHandler,
std::vector< std::map< dealii::DoFHandler< 3 >::active_cell_iterator, std::vector< unsigned int > > > & cellFacesVselfBallSurfacesDofHandler,
std::vector< std::map< dealii::DoFHandler< 3 >::active_cell_iterator, std::vector< unsigned int > > > & cellFacesVselfBallSurfacesDofHandlerForce )
private

◆ distributeForceContributionFnlGammaAtoms()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::distributeForceContributionFnlGammaAtoms ( const std::map< unsigned int, std::vector< double > > & forceContributionFnlGammaAtoms)
private

◆ distributeForceContributionFPSPLocalGammaAtoms()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::distributeForceContributionFPSPLocalGammaAtoms ( const std::map< unsigned int, std::vector< double > > & forceContributionFPSPLocalGammaAtoms,
const std::map< std::pair< unsigned int, unsigned int >, unsigned int > & atomsForceDofs,
const dealii::AffineConstraints< double > & constraintsNoneForce,
distributedCPUVec< double > & configForceVectorLinFE )
private

◆ FnlGammaAtomsElementalContribution()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::FnlGammaAtomsElementalContribution ( std::map< unsigned int, std::vector< double > > & forceContributionFnlGammaAtoms,
const dealii::MatrixFree< 3, double > & matrixFreeData,
dealii::FEEvaluation< 3, 1, C_num1DQuadNLPSP< FEOrder >() *C_numCopies1DQuadNLPSP(), 3 > & forceEvalNLP,
const std::shared_ptr< AtomicCenteredNonLocalOperator< dataTypes::number, memorySpace > > nonLocalOp,
unsigned int numNonLocalAtomsCurrentProcess,
const std::vector< int > & globalChargeIdNonLocalAtoms,
const std::vector< unsigned int > & numberPseudoWaveFunctionsPerAtom,
const unsigned int cell,
const std::map< dealii::CellId, unsigned int > & cellIdToCellNumberMap,
const std::vector< dataTypes::number > & zetaDeltaVQuadsFlattened,
const std::vector< dataTypes::number > & projectorKetTimesPsiTimesVTimesPartOccContractionGradPsiQuadsFlattened )
private

◆ FnlGammaxElementalContribution()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::FnlGammaxElementalContribution ( dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > & FVectQuads,
const dealii::MatrixFree< 3, double > & matrixFreeData,
const unsigned int numQuadPoints,
const std::shared_ptr< AtomicCenteredNonLocalOperator< dataTypes::number, memorySpace > > nonLocalOp,
const unsigned int numNonLocalAtomsCurrentProcess,
const std::vector< int > & globalChargeIdNonLocalAtoms,
const std::vector< unsigned int > & numberPseudoWaveFunctionsPerAtom,
const unsigned int cell,
const std::map< dealii::CellId, unsigned int > & cellIdToCellNumberMap,
const std::vector< dataTypes::number > & zetaDeltaVQuadsFlattened,
const std::vector< dataTypes::number > & projectorKetTimesPsiTimesVTimesPartOccContractionGradPsiQuadsFlattened )
private

◆ FNonlinearCoreCorrectionGammaAtomsElementalContribution() [1/2]

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::FNonlinearCoreCorrectionGammaAtomsElementalContribution ( std::map< unsigned int, std::vector< double > > & forceContributionFNonlinearCoreCorrectionGammaAtoms,
dealii::FEEvaluation< 3, 1, C_num1DQuad< C_rhoNodalPolyOrder< FEOrder, FEOrderElectro >()>(), 3 > & forceEval,
const dealii::MatrixFree< 3, double > & matrixFreeData,
const unsigned int cell,
const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > & derExcGradRho,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & hessianRhoCoreAtoms )
private

◆ FNonlinearCoreCorrectionGammaAtomsElementalContribution() [2/2]

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::FNonlinearCoreCorrectionGammaAtomsElementalContribution ( std::map< unsigned int, std::vector< double > > & forceContributionFNonlinearCoreCorrectionGammaAtoms,
dealii::FEEvaluation< 3, 1, C_num1DQuad< C_rhoNodalPolyOrder< FEOrder, FEOrderElectro >()>(), 3 > & forceEval,
const dealii::MatrixFree< 3, double > & matrixFreeData,
const unsigned int cell,
const dealii::AlignedVector< dealii::VectorizedArray< double > > & vxcQuads,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & gradRhoCoreAtoms )
private

◆ FNonlinearCoreCorrectionGammaAtomsElementalContributionSpinPolarized()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::FNonlinearCoreCorrectionGammaAtomsElementalContributionSpinPolarized ( std::map< unsigned int, std::vector< double > > & forceContributionFNonlinearCoreCorrectionGammaAtoms,
dealii::FEEvaluation< 3, 1, C_num1DQuad< C_rhoNodalPolyOrder< FEOrder, FEOrderElectro >()>(), 3 > & forceEval,
const dealii::MatrixFree< 3, double > & matrixFreeData,
const unsigned int cell,
const dealii::AlignedVector< dealii::VectorizedArray< double > > & vxcQuadsSpin0,
const dealii::AlignedVector< dealii::VectorizedArray< double > > & vxcQuadsSpin1,
const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > & derExcGradRhoSpin0,
const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > & derExcGradRhoSpin1,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & gradRhoCoreAtoms,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & hessianRhoCoreAtoms,
const bool isXCGGA = false )
private

◆ FPhiTotSmearedChargesGammaAtomsElementalContribution()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::FPhiTotSmearedChargesGammaAtomsElementalContribution ( std::map< unsigned int, std::vector< double > > & forceContributionSmearedChargesGammaAtoms,
dealii::FEEvaluation< 3, -1, 1, 3 > & forceEval,
const dealii::MatrixFree< 3, double > & matrixFreeData,
const unsigned int cell,
const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > & gradPhiTotQuads,
const std::vector< unsigned int > & nonTrivialAtomIdsMacroCell,
const std::map< dealii::CellId, std::vector< int > > & bQuadAtomIdsAllAtoms,
const dealii::AlignedVector< dealii::VectorizedArray< double > > & smearedbQuads )
private

◆ FPSPLocalGammaAtomsElementalContribution()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::FPSPLocalGammaAtomsElementalContribution ( std::map< unsigned int, std::vector< double > > & forceContributionFPSPLocalGammaAtoms,
dealii::FEValues< 3 > & feValues,
dealii::FEFaceValues< 3 > & feFaceValues,
dealii::FEEvaluation< 3, 1, C_num1DQuadLPSP< FEOrder >() *C_numCopies1DQuadLPSP(), 3 > & forceEval,
const dealii::MatrixFree< 3, double > & matrixFreeData,
const unsigned int phiTotDofHandlerIndexElectro,
const unsigned int cell,
const dealii::AlignedVector< dealii::VectorizedArray< double > > & rhoQuads,
const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > & gradRhoQuads,
const std::map< unsigned int, std::map< dealii::CellId, std::vector< double > > > & pseudoVLocAtoms,
const vselfBinsManager< FEOrder, FEOrderElectro > & vselfBinsManager,
const std::vector< std::map< dealii::CellId, unsigned int > > & cellsVselfBallsClosestAtomIdDofHandler )
private

◆ FVselfSmearedChargesGammaAtomsElementalContribution()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::FVselfSmearedChargesGammaAtomsElementalContribution ( std::map< unsigned int, std::vector< double > > & forceContributionSmearedChargesGammaAtoms,
dealii::FEEvaluation< 3, -1, 1, 3 > & forceEval,
const dealii::MatrixFree< 3, double > & matrixFreeData,
const unsigned int cell,
const dealii::AlignedVector< dealii::Tensor< 1, 3, dealii::VectorizedArray< double > > > & gradVselfBinQuads,
const std::vector< unsigned int > & nonTrivialAtomIdsMacroCell,
const std::map< dealii::CellId, std::vector< int > > & bQuadAtomIdsAllAtoms,
const dealii::AlignedVector< dealii::VectorizedArray< double > > & smearedbQuads )
private

◆ getAtomsForces()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
std::vector< double > & dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::getAtomsForces ( )

returns a copy of the configurational force on all global atoms.

computeAtomsForces must be called prior to this function call.

Returns
std::vector<double> flattened array of the configurational force on all atoms, the three force components on each atom being the leading dimension. Units- Hartree/Bohr

◆ getStress()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
dealii::Tensor< 2, 3, double > & dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::getStress ( )

returns a copy of the current stress tensor value.

computeStress must be call prior to this function call.

Returns
dealii::Tensor<2,3,double> second order stress Tensor in Hartree/Bohr^3

◆ initMoved()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::initMoved ( std::vector< const dealii::DoFHandler< 3 > * > & dofHandlerVectorMatrixFree,
std::vector< const dealii::AffineConstraints< double > * > & constraintsVectorMatrixFree,
const bool isElectrostaticsMesh )

initializes data structures inside forceClass which depend on the moved mesh.

initMoved is the second step (first step call initUnmoved) of the initialization/reinitialization of force class when starting from a new mesh, and the first step when recomputing forces on a perturbed mesh. initMoved assumes that the triangulation whose reference was passed to the forceClass object in the initUnmoved call now has its nodes moved such that all atomic positions lie on nodes.

Returns
void.

◆ initPseudoData()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::initPseudoData ( )

initializes and precomputes pseudopotential related data structuers required for configurational force and stress computation.

This function is only activated for pseudopotential calculations and is currently called when initializing/reinitializing the dftClass object. This function initializes and precomputes the pseudopotential datastructures for local and non-local parts. Separate internal function calls are made for KB and ONCV projectors.

Returns
void.

◆ initUnmoved()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::initUnmoved ( const dealii::Triangulation< 3, 3 > & triangulation,
const dealii::Triangulation< 3, 3 > & serialTriangulation,
const std::vector< std::vector< double > > & domainBoundingVectors,
const bool isElectrostaticsMesh )

initializes data structures inside forceClass assuming unmoved triangulation.

initUnmoved is the first step of the initialization/reinitialization of force class when starting from a new unmoved triangulation. It creates the dofHandler with linear finite elements and three components corresponding to the three force components. It also creates the corresponding constraint matrices which is why an unmoved triangulation is necessary. Finally this function also initializes the gaussianMovePar data member.

Parameters
triangulationreference to unmoved triangulation where the mesh nodes have not been manually moved.
isElectrostaticsMeshboolean parameter specifying whether this triangulatio is to be used for for the electrostatics part of the configurational force.
Returns
void.

◆ locateAtomCoreNodesForce()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::locateAtomCoreNodesForce ( const dealii::DoFHandler< 3 > & dofHandlerForce,
const dealii::IndexSet & locally_owned_dofsForce,
std::map< std::pair< unsigned int, unsigned int >, unsigned int > & atomsForceDofs )
private

get the value of Gaussian generator parameter (d_gaussianConstant). Gaussian generator: Gamma(r)= exp(-d_gaussianConstant*r^2).

Locates and stores the global dof indices of d_dofHandlerForce whose cooridinates match with the atomic positions.

Returns
void.

◆ printAtomsForces()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::printAtomsForces ( )

prints the currently stored configurational forces on atoms and the Gaussian generator constant used to compute them.

Returns
void.

◆ printStress()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::printStress ( )

prints the currently stored configurational stress tensor.

Returns
void.

◆ stressEnlElementalContribution()

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
void dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::stressEnlElementalContribution ( dealii::Tensor< 2, 3, double > & stressContribution,
const dealii::MatrixFree< 3, double > & matrixFreeData,
const unsigned int numQuadPoints,
const std::vector< double > & jxwQuadsSubCells,
const unsigned int cell,
const unsigned int numNonLocalAtomsCurrentProcess,
const std::shared_ptr< AtomicCenteredNonLocalOperator< dataTypes::number, memorySpace > > nonLocalOp,
const std::vector< unsigned int > & numberPseudoWaveFunctionsPerAtom,
const std::map< dealii::CellId, unsigned int > & cellIdToCellNumberMap,
const std::vector< dataTypes::number > & zetalmDeltaVlProductDistImageAtoms,
const std::vector< dataTypes::number > & projectorKetTimesPsiTimesVTimesPartOccContractionGradPsiQuadsFlattened,
const bool isSpinPolarized )
private

Friends And Related Symbol Documentation

◆ dftClass< FEOrder, FEOrderElectro, memorySpace >

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
friend class dftClass< FEOrder, FEOrderElectro, memorySpace >
friend

Member Data Documentation

◆ d_allowGaussianOverlapOnAtoms

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
const bool dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_allowGaussianOverlapOnAtoms = false
private

◆ d_AtomIdBinIdLocalDofHandler

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
std::map<unsigned int, unsigned int> dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_AtomIdBinIdLocalDofHandler
private

Internal data: stores the map of atom Id (only in the local processor) to the vself bin Id.

◆ d_AtomIdBinIdLocalDofHandlerElectro

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
std::map<unsigned int, unsigned int> dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_AtomIdBinIdLocalDofHandlerElectro
private

Internal data: stores the map of atom Id (only in the local processor) to the vself bin Id.

◆ d_atomsForceDofs

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
std::map<std::pair<unsigned int, unsigned int>, unsigned int> dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_atomsForceDofs
private

Internal data: map < <atomId,force component>, globaldof in d_dofHandlerForce>

◆ d_atomsForceDofsElectro

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
std::map<std::pair<unsigned int, unsigned int>, unsigned int> dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_atomsForceDofsElectro
private

Internal data: map < <atomId,force component>, globaldof in d_dofHandlerForceElectro>

◆ d_cellFacesVselfBallSurfacesDofHandler

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
std::vector<std::map<dealii::DoFHandler<3>::active_cell_iterator, std::vector<unsigned int> > > dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_cellFacesVselfBallSurfacesDofHandler
private

◆ d_cellFacesVselfBallSurfacesDofHandlerElectro

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
std::vector<std::map<dealii::DoFHandler<3>::active_cell_iterator, std::vector<unsigned int> > > dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_cellFacesVselfBallSurfacesDofHandlerElectro
private

◆ d_cellFacesVselfBallSurfacesDofHandlerForce

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
std::vector<std::map<dealii::DoFHandler<3>::active_cell_iterator, std::vector<unsigned int> > > dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_cellFacesVselfBallSurfacesDofHandlerForce
private

◆ d_cellFacesVselfBallSurfacesDofHandlerForceElectro

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
std::vector<std::map<dealii::DoFHandler<3>::active_cell_iterator, std::vector<unsigned int> > > dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_cellFacesVselfBallSurfacesDofHandlerForceElectro
private

◆ d_cellIdToActiveCellIteratorMapDofHandlerRhoNodalElectro

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
std::map<dealii::CellId, dealii::DoFHandler<3>::active_cell_iterator> dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_cellIdToActiveCellIteratorMapDofHandlerRhoNodalElectro
private

◆ d_cellsVselfBallsClosestAtomIdDofHandler

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
std::vector<std::map<dealii::CellId, unsigned int> > dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_cellsVselfBallsClosestAtomIdDofHandler
private

Internal data: stores map of vself ball cell Id to the closest atom Id of that cell. Outer vector over vself bins.

◆ d_cellsVselfBallsClosestAtomIdDofHandlerElectro

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
std::vector<std::map<dealii::CellId, unsigned int> > dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_cellsVselfBallsClosestAtomIdDofHandlerElectro
private

Internal data: stores map of vself ball cell Id to the closest atom Id of that cell. Outer vector over vself bins.

◆ d_cellsVselfBallsDofHandler

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
std::vector<std::vector<dealii::DoFHandler<3>::active_cell_iterator> > dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_cellsVselfBallsDofHandler
private

Internal data: stores cell iterators of all cells in dftPtr->d_dofHandler which are part of the vself ball. Outer vector is over vself bins.

◆ d_cellsVselfBallsDofHandlerElectro

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
std::vector<std::vector<dealii::DoFHandler<3>::active_cell_iterator> > dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_cellsVselfBallsDofHandlerElectro
private

Internal data: stores cell iterators of all cells in dftPtr->d_dofHandler which are part of the vself ball. Outer vector is over vself bins.

◆ d_cellsVselfBallsDofHandlerForce

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
std::vector<std::vector<dealii::DoFHandler<3>::active_cell_iterator> > dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_cellsVselfBallsDofHandlerForce
private

Internal data: stores cell iterators of all cells in d_dofHandlerForce which are part of the vself ball. Outer vector is over vself bins.

◆ d_cellsVselfBallsDofHandlerForceElectro

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
std::vector<std::vector<dealii::DoFHandler<3>::active_cell_iterator> > dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_cellsVselfBallsDofHandlerForceElectro
private

Internal data: stores cell iterators of all cells in d_dofHandlerForce which are part of the vself ball. Outer vector is over vself bins.

◆ d_configForceVectorLinFE

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
distributedCPUVec<double> dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_configForceVectorLinFE
private

Parallel distributed vector field which stores the configurational force for each fem node corresponding to linear shape function generator (see equations 52-53 in (https://link.aps.org/doi/10.1103/PhysRevB.97.165132)). This vector doesn't contain contribution from terms which have sums over k points.

◆ d_configForceVectorLinFEElectro

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
distributedCPUVec<double> dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_configForceVectorLinFEElectro
private

Parallel distributed vector field which stores the configurational force for each fem node corresponding to linear shape function generator (see equations 52-53 in (https://link.aps.org/doi/10.1103/PhysRevB.97.165132)). This vector only containts contribution from the electrostatic part.

◆ d_constraintsNoneForce

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
dealii::AffineConstraints<double> dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_constraintsNoneForce
private

Constraint matrix for hanging node and periodic constaints on d_dofHandlerForce.

◆ d_constraintsNoneForceElectro

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
dealii::AffineConstraints<double> dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_constraintsNoneForceElectro
private

Constraint matrix for hanging node and periodic constaints on d_dofHandlerForceElectro.

◆ d_dftParams

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
const dftParameters& dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_dftParams
private

◆ d_dofHandlerForce

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
dealii::DoFHandler<3> dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_dofHandlerForce
private

◆ d_dofHandlerForceElectro

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
dealii::DoFHandler<3> dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_dofHandlerForceElectro
private

◆ d_forceAtomsFloating

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
std::vector<double> dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_forceAtomsFloating
private

◆ d_forceDofHandlerIndex

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_forceDofHandlerIndex
private

Index of the d_dofHandlerForce in the MatrixFree object stored in dftClass. This is required to correctly use FEEvaluation class.

◆ d_forceDofHandlerIndexElectro

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_forceDofHandlerIndexElectro
private

Index of the d_dofHandlerForceElectro in the MatrixFree object stored in dftClass. This is required to correctly use FEEvaluation class.

◆ d_gaussianWeightsVecAtoms

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
std::vector<distributedCPUVec<double> > dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_gaussianWeightsVecAtoms
private

◆ d_globalAtomsForces

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
std::vector<double> dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_globalAtomsForces
private

Storage for configurational force on all global atoms.

Gaussian generator constant. Gaussian generator: Gamma(r)= exp(-d_gaussianConstant*r^2) FIXME: Until the hanging nodes surface integral issue is fixed use a value >=4.0

◆ d_locally_owned_dofsForce

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
dealii::IndexSet dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_locally_owned_dofsForce
private

dealii::IndexSet of locally owned dofs of in d_dofHandlerForce the current processor

◆ d_locally_owned_dofsForceElectro

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
dealii::IndexSet dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_locally_owned_dofsForceElectro
private

dealii::IndexSet of locally owned dofs of in d_dofHandlerForceElectro the current processor

◆ d_locally_relevant_dofsForce

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
dealii::IndexSet dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_locally_relevant_dofsForce
private

dealii::IndexSet of locally relevant dofs of in d_dofHandlerForce the current processor

◆ d_locally_relevant_dofsForceElectro

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
dealii::IndexSet dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_locally_relevant_dofsForceElectro
private

dealii::IndexSet of locally relevant dofs of in d_dofHandlerForceElectro the current processor

◆ d_mpiCommParent

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
const MPI_Comm dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_mpiCommParent
private

domain decomposition mpi_communicator

◆ d_stress

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
dealii::Tensor<2, 3, double> dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_stress
private

Storage for configurational stress tensor.

◆ d_stressKPoints

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
dealii::Tensor<2, 3, double> dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::d_stressKPoints
private

◆ dftPtr

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
dftClass<FEOrder, FEOrderElectro, memorySpace>* dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::dftPtr
private

pointer to dft class

◆ FEForce

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
dealii::FESystem<3> dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::FEForce
private

Finite element object for configurational force computation. Linear finite elements with three force field components are used.

◆ mpi_communicator

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
const MPI_Comm dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::mpi_communicator
private

domain decomposition mpi_communicator

◆ n_mpi_processes

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
const unsigned int dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::n_mpi_processes
private

number of mpi processes in the current pool

◆ pcout

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
dealii::ConditionalOStream dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::pcout
private

conditional stream object to enable printing only on root processor across all pools

◆ this_mpi_process

template<unsigned int FEOrder, unsigned int FEOrderElectro, dftfe::utils::MemorySpace memorySpace>
const unsigned int dftfe::forceClass< FEOrder, FEOrderElectro, memorySpace >::this_mpi_process
private

current mpi process id in the current pool


The documentation for this class was generated from the following file: