DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
dftfe::vselfBinsManager Class Reference

Categorizes atoms into bins for efficient solution of nuclear electrostatic self-potential. More...

#include <vselfBinsManager.h>

Public Member Functions

 vselfBinsManager (const MPI_Comm &mpi_comm_parent, const MPI_Comm &mpi_comm_domain, const MPI_Comm &mpi_intercomm_kpts, const dftParameters &dftParams)
 Constructor.
 
void createAtomBins (std::vector< const dealii::AffineConstraints< double > * > &constraintsVector, const dealii::AffineConstraints< double > &onlyHangingNodeConstraints, const dealii::DoFHandler< 3 > &dofHandler, const dealii::AffineConstraints< double > &constraintMatrix, const std::vector< std::vector< double > > &atomLocations, const std::vector< std::vector< double > > &imagePositions, const std::vector< dftfe::Int > &imageIds, const std::vector< double > &imageCharges, const double radiusAtomBall)
 Categorize atoms into bins based on self-potential ball radius around each atom such that no two atoms in each bin has overlapping balls.
 
void updateBinsBc (std::vector< const dealii::AffineConstraints< double > * > &constraintsVector, const dealii::AffineConstraints< double > &onlyHangingNodeConstraints, const dealii::DoFHandler< 3 > &dofHandler, const dealii::AffineConstraints< double > &constraintMatrix, const std::vector< std::vector< double > > &atomLocations, const std::vector< std::vector< double > > &imagePositions, const std::vector< dftfe::Int > &imageIds, const std::vector< double > &imageCharges, const bool vselfPerturbationUpdateForStress=false)
 Categorize atoms into bins based on self-potential ball radius around each atom such that no two atoms in each bin has overlapping balls.
 
void solveVselfInBins (const std::shared_ptr< dftfe::basis::FEBasisOperations< double, double, dftfe::utils::MemorySpace::HOST > > &basisOperationsPtr, const dftfe::uInt offset, const dftfe::uInt matrixFreeQuadratureIdAX, const dealii::AffineConstraints< double > &hangingPeriodicConstraintMatrix, const std::vector< std::vector< double > > &imagePositions, const std::vector< dftfe::Int > &imageIds, const std::vector< double > &imageCharges, std::vector< std::vector< double > > &localVselfs, std::map< dealii::CellId, std::vector< double > > &bQuadValuesAllAtoms, std::map< dealii::CellId, std::vector< dftfe::Int > > &bQuadAtomIdsAllAtoms, std::map< dealii::CellId, std::vector< dftfe::Int > > &bQuadAtomIdsAllAtomsImages, std::map< dealii::CellId, std::vector< dftfe::uInt > > &bCellNonTrivialAtomIds, std::vector< std::map< dealii::CellId, std::vector< dftfe::uInt > > > &bCellNonTrivialAtomIdsBins, std::map< dealii::CellId, std::vector< dftfe::uInt > > &bCellNonTrivialAtomImageIds, std::vector< std::map< dealii::CellId, std::vector< dftfe::uInt > > > &bCellNonTrivialAtomImageIdsBins, const std::vector< double > &smearingWidths, std::vector< double > &smearedChargeScaling, const dftfe::uInt smearedChargeQuadratureId, const bool useSmearedCharges=false, const bool isVselfPerturbationSolve=false)
 Solve nuclear electrostatic self-potential of atoms in each bin one-by-one.
 
void solveVselfInBinsPerturbedDomain (const std::shared_ptr< dftfe::basis::FEBasisOperations< double, double, dftfe::utils::MemorySpace::HOST > > &basisOperationsPtr, const dftfe::uInt mfBaseDofHandlerIndex, const dftfe::uInt matrixFreeQuadratureIdAX, const dftfe::uInt offset, const dealii::AffineConstraints< double > &hangingPeriodicConstraintMatrix, const std::vector< std::vector< double > > &imagePositions, const std::vector< dftfe::Int > &imageIds, const std::vector< double > &imageCharges, const std::vector< double > &smearingWidths, const dftfe::uInt smearedChargeQuadratureId, const bool useSmearedCharges=false)
 Solve nuclear electrostatic self-potential of atoms in each bin in a perturbed domain (used for cell stress calculation)
 
const std::map< dftfe::Int, std::set< dftfe::Int > > & getAtomIdsBins () const
 get const reference map of binIds and atomIds
 
const std::map< dftfe::Int, std::set< dftfe::Int > > & getAtomImageIdsBins () const
 get const reference map of binIds and atomIds
 
const std::vector< std::map< dealii::types::global_dof_index, dftfe::Int > > & getBoundaryFlagsBins () const
 
const std::vector< std::map< dealii::types::global_dof_index, dftfe::Int > > & getBoundaryFlagsBinsOnlyChargeId () const
 
const std::vector< std::map< dealii::types::global_dof_index, dftfe::Int > > & getClosestAtomIdsBins () const
 
const std::vector< std::map< dealii::types::global_dof_index, dealii::Point< 3 > > > & getClosestAtomLocationsBins () const
 
const std::vector< distributedCPUVec< double > > & getVselfFieldBins () const
 get const reference to solved vself fields
 
const std::vector< distributedCPUVec< double > > & getVselfFieldDerRBins () const
 get const reference to del{vself}/del{R_i} fields
 
const std::vector< distributedCPUVec< double > > & getPerturbedVselfFieldBins () const
 
const std::map< dftfe::uInt, dftfe::uInt > & getAtomIdBinIdMapLocalAllImages () const
 get const reference to d_atomIdBinIdMapLocalAllImages
 
double getStoredAdaptiveBallRadius () const
 get stored adaptive ball radius
 

Private Member Functions

void locateAtomsInBins (const dealii::DoFHandler< 3 > &dofHandler)
 locate underlying fem nodes for atoms in bins.
 
void createAtomBinsSanityCheck (const dealii::DoFHandler< 3 > &dofHandler, const dealii::AffineConstraints< double > &onlyHangingNodeConstraints)
 sanity check for Dirichlet boundary conditions on the vself balls in each bin
 

Private Attributes

std::vector< std::vector< double > > d_atomLocations
 storage for input atomLocations argument in createAtomBins function
 
dftUtils::constraintMatrixInfo< dftfe::utils::MemorySpace::HOSTd_constraintsOnlyHangingInfo
 storage for optimized constraints handling object
 
std::vector< dealii::AffineConstraints< double > > d_vselfBinConstraintMatrices
 vector of constraint matrices for vself bins
 
std::map< dftfe::Int, std::set< dftfe::Int > > d_bins
 map of binIds and atomIds
 
std::map< dftfe::Int, std::set< dftfe::Int > > d_binsImages
 map of binIds and atomIds and imageIds
 
std::vector< std::map< dealii::types::global_dof_index, dftfe::Int > > d_boundaryFlag
 map of global dof index and vself solve boundary flag (chargeId or
 
std::vector< std::map< dealii::types::global_dof_index, dftfe::Int > > d_boundaryFlagOnlyChargeId
 
std::vector< std::map< dealii::types::global_dof_index, dealii::Point< 3 > > > d_dofClosestChargeLocationMap
 map of global dof index to location of closest charge
 
std::vector< std::map< dealii::types::global_dof_index, double > > d_vselfBinField
 map of global dof index and vself field initial value in each bin
 
std::vector< std::map< dealii::types::global_dof_index, dftfe::Int > > d_closestAtomBin
 map of global dof index and vself field initial value in each bin
 
std::map< dftfe::uInt, dftfe::uIntd_atomIdBinIdMapLocalAllImages
 
std::vector< distributedCPUVec< double > > d_vselfFieldBins
 solved vself solution field for each bin
 
std::vector< distributedCPUVec< double > > d_vselfFieldDerRBins
 solved del{vself}/del{R_i} solution field for each bin
 
std::vector< distributedCPUVec< double > > d_vselfFieldPerturbedBins
 
std::vector< std::map< dealii::types::global_dof_index, double > > d_atomsInBin
 
double d_storedAdaptiveBallRadius
 
const dftParametersd_dftParams
 
const MPI_Comm d_mpiCommParent
 
const MPI_Comm mpi_communicator
 
const MPI_Comm d_mpiInterCommKpts
 
const dftfe::uInt n_mpi_processes
 
const dftfe::uInt this_mpi_process
 
dealii::ConditionalOStream pcout
 

Detailed Description

Categorizes atoms into bins for efficient solution of nuclear electrostatic self-potential.

Author
Sambit Das, Phani Motamarri

Constructor & Destructor Documentation

◆ vselfBinsManager()

dftfe::vselfBinsManager::vselfBinsManager ( const MPI_Comm & mpi_comm_parent,
const MPI_Comm & mpi_comm_domain,
const MPI_Comm & mpi_intercomm_kpts,
const dftParameters & dftParams )

Constructor.

Parameters
mpi_comm_parentparent mpi communicator
mpi_comm_domaindomain decomposition mpi communicator

Member Function Documentation

◆ createAtomBins()

void dftfe::vselfBinsManager::createAtomBins ( std::vector< const dealii::AffineConstraints< double > * > & constraintsVector,
const dealii::AffineConstraints< double > & onlyHangingNodeConstraints,
const dealii::DoFHandler< 3 > & dofHandler,
const dealii::AffineConstraints< double > & constraintMatrix,
const std::vector< std::vector< double > > & atomLocations,
const std::vector< std::vector< double > > & imagePositions,
const std::vector< dftfe::Int > & imageIds,
const std::vector< double > & imageCharges,
const double radiusAtomBall )

Categorize atoms into bins based on self-potential ball radius around each atom such that no two atoms in each bin has overlapping balls.

Parameters
[out]constraintsVectorconstraintsVector to which the vself bins solve constraint matrices will be pushed back
[in]dofHandlerDofHandler object
[in]constraintMatrixdealii::AffineConstraints<double> which was used for the total electrostatics solve
[in]atomLocationsglobal atom locations and charge values data
[in]imagePositionsimage atoms positions data
[in]imageIdsimage atoms Ids data
[in]imageChargesimage atoms charge values data
[in]radiusAtomBallself-potential ball radius

◆ createAtomBinsSanityCheck()

void dftfe::vselfBinsManager::createAtomBinsSanityCheck ( const dealii::DoFHandler< 3 > & dofHandler,
const dealii::AffineConstraints< double > & onlyHangingNodeConstraints )
private

sanity check for Dirichlet boundary conditions on the vself balls in each bin

◆ getAtomIdBinIdMapLocalAllImages()

const std::map< dftfe::uInt, dftfe::uInt > & dftfe::vselfBinsManager::getAtomIdBinIdMapLocalAllImages ( ) const

get const reference to d_atomIdBinIdMapLocalAllImages

◆ getAtomIdsBins()

const std::map< dftfe::Int, std::set< dftfe::Int > > & dftfe::vselfBinsManager::getAtomIdsBins ( ) const

get const reference map of binIds and atomIds

◆ getAtomImageIdsBins()

const std::map< dftfe::Int, std::set< dftfe::Int > > & dftfe::vselfBinsManager::getAtomImageIdsBins ( ) const

get const reference map of binIds and atomIds

◆ getBoundaryFlagsBins()

const std::vector< std::map< dealii::types::global_dof_index, dftfe::Int > > & dftfe::vselfBinsManager::getBoundaryFlagsBins ( ) const

get const reference to map of global dof index and vself solve boundary flag in each bin

◆ getBoundaryFlagsBinsOnlyChargeId()

const std::vector< std::map< dealii::types::global_dof_index, dftfe::Int > > & dftfe::vselfBinsManager::getBoundaryFlagsBinsOnlyChargeId ( ) const

get const reference to map of global dof index and vself solve boundary flag in each bin

◆ getClosestAtomIdsBins()

const std::vector< std::map< dealii::types::global_dof_index, dftfe::Int > > & dftfe::vselfBinsManager::getClosestAtomIdsBins ( ) const

get const reference to map of global dof index and vself field initial value in each bin

◆ getClosestAtomLocationsBins()

const std::vector< std::map< dealii::types::global_dof_index, dealii::Point< 3 > > > & dftfe::vselfBinsManager::getClosestAtomLocationsBins ( ) const

get const reference to map of global dof index and vself field initial value in each bin

◆ getPerturbedVselfFieldBins()

const std::vector< distributedCPUVec< double > > & dftfe::vselfBinsManager::getPerturbedVselfFieldBins ( ) const

perturbation of vself solution field to be used to evaluate the Gateaux derivative of vself field with respect to affine strain components using central finite difference

◆ getStoredAdaptiveBallRadius()

double dftfe::vselfBinsManager::getStoredAdaptiveBallRadius ( ) const

get stored adaptive ball radius

◆ getVselfFieldBins()

const std::vector< distributedCPUVec< double > > & dftfe::vselfBinsManager::getVselfFieldBins ( ) const

get const reference to solved vself fields

◆ getVselfFieldDerRBins()

const std::vector< distributedCPUVec< double > > & dftfe::vselfBinsManager::getVselfFieldDerRBins ( ) const

get const reference to del{vself}/del{R_i} fields

◆ locateAtomsInBins()

void dftfe::vselfBinsManager::locateAtomsInBins ( const dealii::DoFHandler< 3 > & dofHandler)
private

locate underlying fem nodes for atoms in bins.

◆ solveVselfInBins()

void dftfe::vselfBinsManager::solveVselfInBins ( const std::shared_ptr< dftfe::basis::FEBasisOperations< double, double, dftfe::utils::MemorySpace::HOST > > & basisOperationsPtr,
const dftfe::uInt offset,
const dftfe::uInt matrixFreeQuadratureIdAX,
const dealii::AffineConstraints< double > & hangingPeriodicConstraintMatrix,
const std::vector< std::vector< double > > & imagePositions,
const std::vector< dftfe::Int > & imageIds,
const std::vector< double > & imageCharges,
std::vector< std::vector< double > > & localVselfs,
std::map< dealii::CellId, std::vector< double > > & bQuadValuesAllAtoms,
std::map< dealii::CellId, std::vector< dftfe::Int > > & bQuadAtomIdsAllAtoms,
std::map< dealii::CellId, std::vector< dftfe::Int > > & bQuadAtomIdsAllAtomsImages,
std::map< dealii::CellId, std::vector< dftfe::uInt > > & bCellNonTrivialAtomIds,
std::vector< std::map< dealii::CellId, std::vector< dftfe::uInt > > > & bCellNonTrivialAtomIdsBins,
std::map< dealii::CellId, std::vector< dftfe::uInt > > & bCellNonTrivialAtomImageIds,
std::vector< std::map< dealii::CellId, std::vector< dftfe::uInt > > > & bCellNonTrivialAtomImageIdsBins,
const std::vector< double > & smearingWidths,
std::vector< double > & smearedChargeScaling,
const dftfe::uInt smearedChargeQuadratureId,
const bool useSmearedCharges = false,
const bool isVselfPerturbationSolve = false )

Solve nuclear electrostatic self-potential of atoms in each bin one-by-one.

Parameters
[in]matrix_free_dataMatrixFree object
[in]offsetMatrixFree object starting offset for vself bins solve
[out]phiExtsum of the self-potential fields of all atoms and image atoms
[in]phiExtConstraintMatrixconstraintMatrix corresponding to phiExt
[in]imagePositionsimage atoms positions data
[in]imageIdsimage atoms Ids data
[in]imageChargesimage atoms charge values data *
[out]localVselfspeak self-potential values of atoms in the local processor

◆ solveVselfInBinsPerturbedDomain()

void dftfe::vselfBinsManager::solveVselfInBinsPerturbedDomain ( const std::shared_ptr< dftfe::basis::FEBasisOperations< double, double, dftfe::utils::MemorySpace::HOST > > & basisOperationsPtr,
const dftfe::uInt mfBaseDofHandlerIndex,
const dftfe::uInt matrixFreeQuadratureIdAX,
const dftfe::uInt offset,
const dealii::AffineConstraints< double > & hangingPeriodicConstraintMatrix,
const std::vector< std::vector< double > > & imagePositions,
const std::vector< dftfe::Int > & imageIds,
const std::vector< double > & imageCharges,
const std::vector< double > & smearingWidths,
const dftfe::uInt smearedChargeQuadratureId,
const bool useSmearedCharges = false )

Solve nuclear electrostatic self-potential of atoms in each bin in a perturbed domain (used for cell stress calculation)

◆ updateBinsBc()

void dftfe::vselfBinsManager::updateBinsBc ( std::vector< const dealii::AffineConstraints< double > * > & constraintsVector,
const dealii::AffineConstraints< double > & onlyHangingNodeConstraints,
const dealii::DoFHandler< 3 > & dofHandler,
const dealii::AffineConstraints< double > & constraintMatrix,
const std::vector< std::vector< double > > & atomLocations,
const std::vector< std::vector< double > > & imagePositions,
const std::vector< dftfe::Int > & imageIds,
const std::vector< double > & imageCharges,
const bool vselfPerturbationUpdateForStress = false )

Categorize atoms into bins based on self-potential ball radius around each atom such that no two atoms in each bin has overlapping balls.

Parameters
[out]constraintsVectorconstraintsVector to which the vself bins solve constraint matrices will be pushed back
[in]dofHandlerDofHandler object
[in]constraintMatrixdealii::AffineConstraints<double> which was used for the total electrostatics solve
[in]atomLocationsglobal atom locations and charge values data
[in]imagePositionsimage atoms positions data
[in]imageIdsimage atoms Ids data
[in]imageChargesimage atoms charge values data

Member Data Documentation

◆ d_atomIdBinIdMapLocalAllImages

std::map<dftfe::uInt, dftfe::uInt> dftfe::vselfBinsManager::d_atomIdBinIdMapLocalAllImages
private

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

◆ d_atomLocations

std::vector<std::vector<double> > dftfe::vselfBinsManager::d_atomLocations
private

storage for input atomLocations argument in createAtomBins function

◆ d_atomsInBin

std::vector<std::map<dealii::types::global_dof_index, double> > dftfe::vselfBinsManager::d_atomsInBin
private

Map of locally relevant global dof index and the atomic charge in each bin

◆ d_bins

std::map<dftfe::Int, std::set<dftfe::Int> > dftfe::vselfBinsManager::d_bins
private

map of binIds and atomIds

◆ d_binsImages

std::map<dftfe::Int, std::set<dftfe::Int> > dftfe::vselfBinsManager::d_binsImages
private

map of binIds and atomIds and imageIds

◆ d_boundaryFlag

std::vector<std::map<dealii::types::global_dof_index, dftfe::Int> > dftfe::vselfBinsManager::d_boundaryFlag
private

map of global dof index and vself solve boundary flag (chargeId or

◆ d_boundaryFlagOnlyChargeId

std::vector<std::map<dealii::types::global_dof_index, dftfe::Int> > dftfe::vselfBinsManager::d_boundaryFlagOnlyChargeId
private

map of global dof index and vself solve boundary flag (only chargeId)in each bin

◆ d_closestAtomBin

std::vector<std::map<dealii::types::global_dof_index, dftfe::Int> > dftfe::vselfBinsManager::d_closestAtomBin
private

map of global dof index and vself field initial value in each bin

◆ d_constraintsOnlyHangingInfo

dftUtils::constraintMatrixInfo<dftfe::utils::MemorySpace::HOST> dftfe::vselfBinsManager::d_constraintsOnlyHangingInfo
private

storage for optimized constraints handling object

◆ d_dftParams

const dftParameters& dftfe::vselfBinsManager::d_dftParams
private

◆ d_dofClosestChargeLocationMap

std::vector<std::map<dealii::types::global_dof_index, dealii::Point<3> > > dftfe::vselfBinsManager::d_dofClosestChargeLocationMap
private

map of global dof index to location of closest charge

◆ d_mpiCommParent

const MPI_Comm dftfe::vselfBinsManager::d_mpiCommParent
private

◆ d_mpiInterCommKpts

const MPI_Comm dftfe::vselfBinsManager::d_mpiInterCommKpts
private

◆ d_storedAdaptiveBallRadius

double dftfe::vselfBinsManager::d_storedAdaptiveBallRadius
private

Vself ball radius. This is stored after the first call to createAtomBins and reused for subsequent calls

◆ d_vselfBinConstraintMatrices

std::vector<dealii::AffineConstraints<double> > dftfe::vselfBinsManager::d_vselfBinConstraintMatrices
private

vector of constraint matrices for vself bins

◆ d_vselfBinField

std::vector<std::map<dealii::types::global_dof_index, double> > dftfe::vselfBinsManager::d_vselfBinField
private

map of global dof index and vself field initial value in each bin

◆ d_vselfFieldBins

std::vector<distributedCPUVec<double> > dftfe::vselfBinsManager::d_vselfFieldBins
private

solved vself solution field for each bin

◆ d_vselfFieldDerRBins

std::vector<distributedCPUVec<double> > dftfe::vselfBinsManager::d_vselfFieldDerRBins
private

solved del{vself}/del{R_i} solution field for each bin

◆ d_vselfFieldPerturbedBins

std::vector<distributedCPUVec<double> > dftfe::vselfBinsManager::d_vselfFieldPerturbedBins
private

perturbation of vself solution field to be used to evaluate the Gateaux derivative of vself field with respect to affine strain components using central finite difference

◆ mpi_communicator

const MPI_Comm dftfe::vselfBinsManager::mpi_communicator
private

◆ n_mpi_processes

const dftfe::uInt dftfe::vselfBinsManager::n_mpi_processes
private

◆ pcout

dealii::ConditionalOStream dftfe::vselfBinsManager::pcout
private

◆ this_mpi_process

const dftfe::uInt dftfe::vselfBinsManager::this_mpi_process
private

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