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

Categorizes atoms into bins for efficient solution of nuclear electrostatic self-potential. template parameter FEOrderElectro is the finite element polynomial order. 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< 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< 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 unsigned int offset, const unsigned int matrixFreeQuadratureIdAX, const dealii::AffineConstraints< double > &hangingPeriodicConstraintMatrix, const std::vector< std::vector< double > > &imagePositions, const std::vector< 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< int > > &bQuadAtomIdsAllAtoms, std::map< dealii::CellId, std::vector< int > > &bQuadAtomIdsAllAtomsImages, std::map< dealii::CellId, std::vector< unsigned int > > &bCellNonTrivialAtomIds, std::vector< std::map< dealii::CellId, std::vector< unsigned int > > > &bCellNonTrivialAtomIdsBins, std::map< dealii::CellId, std::vector< unsigned int > > &bCellNonTrivialAtomImageIds, std::vector< std::map< dealii::CellId, std::vector< unsigned int > > > &bCellNonTrivialAtomImageIdsBins, const std::vector< double > &smearingWidths, std::vector< double > &smearedChargeScaling, const unsigned int 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 unsigned int mfBaseDofHandlerIndex, const unsigned int matrixFreeQuadratureIdAX, const unsigned int offset, const dealii::AffineConstraints< double > &hangingPeriodicConstraintMatrix, const std::vector< std::vector< double > > &imagePositions, const std::vector< int > &imageIds, const std::vector< double > &imageCharges, const std::vector< double > &smearingWidths, const unsigned int 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< int, std::set< int > > & getAtomIdsBins () const
 get const reference map of binIds and atomIds
 
const std::map< int, std::set< int > > & getAtomImageIdsBins () const
 get const reference map of binIds and atomIds
 
const std::vector< std::map< dealii::types::global_dof_index, int > > & getBoundaryFlagsBins () const
 
const std::vector< std::map< dealii::types::global_dof_index, int > > & getBoundaryFlagsBinsOnlyChargeId () const
 
const std::vector< std::map< dealii::types::global_dof_index, 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< unsigned int, unsigned int > & 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< int, std::set< int > > d_bins
 map of binIds and atomIds
 
std::map< int, std::set< int > > d_binsImages
 map of binIds and atomIds and imageIds
 
std::vector< std::map< dealii::types::global_dof_index, int > > d_boundaryFlag
 map of global dof index and vself solve boundary flag (chargeId or
 
std::vector< std::map< dealii::types::global_dof_index, 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, int > > d_closestAtomBin
 map of global dof index and vself field initial value in each bin
 
std::map< unsigned int, unsigned int > d_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 unsigned int n_mpi_processes
 
const unsigned int this_mpi_process
 
dealii::ConditionalOStream pcout
 

Detailed Description

template<unsigned int FEOrder, unsigned int FEOrderElectro>
class dftfe::vselfBinsManager< FEOrder, FEOrderElectro >

Categorizes atoms into bins for efficient solution of nuclear electrostatic self-potential. template parameter FEOrderElectro is the finite element polynomial order.

Author
Sambit Das, Phani Motamarri

Constructor & Destructor Documentation

◆ vselfBinsManager()

template<unsigned int FEOrder, unsigned int FEOrderElectro>
dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::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()

template<unsigned int FEOrder, unsigned int FEOrderElectro>
void dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::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< 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()

template<unsigned int FEOrder, unsigned int FEOrderElectro>
void dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::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()

template<unsigned int FEOrder, unsigned int FEOrderElectro>
const std::map< unsigned int, unsigned int > & dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::getAtomIdBinIdMapLocalAllImages ( ) const

get const reference to d_atomIdBinIdMapLocalAllImages

◆ getAtomIdsBins()

template<unsigned int FEOrder, unsigned int FEOrderElectro>
const std::map< int, std::set< int > > & dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::getAtomIdsBins ( ) const

get const reference map of binIds and atomIds

◆ getAtomImageIdsBins()

template<unsigned int FEOrder, unsigned int FEOrderElectro>
const std::map< int, std::set< int > > & dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::getAtomImageIdsBins ( ) const

get const reference map of binIds and atomIds

◆ getBoundaryFlagsBins()

template<unsigned int FEOrder, unsigned int FEOrderElectro>
const std::vector< std::map< dealii::types::global_dof_index, int > > & dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::getBoundaryFlagsBins ( ) const

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

◆ getBoundaryFlagsBinsOnlyChargeId()

template<unsigned int FEOrder, unsigned int FEOrderElectro>
const std::vector< std::map< dealii::types::global_dof_index, int > > & dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::getBoundaryFlagsBinsOnlyChargeId ( ) const

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

◆ getClosestAtomIdsBins()

template<unsigned int FEOrder, unsigned int FEOrderElectro>
const std::vector< std::map< dealii::types::global_dof_index, int > > & dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::getClosestAtomIdsBins ( ) const

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

◆ getClosestAtomLocationsBins()

template<unsigned int FEOrder, unsigned int FEOrderElectro>
const std::vector< std::map< dealii::types::global_dof_index, dealii::Point< 3 > > > & dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::getClosestAtomLocationsBins ( ) const

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

◆ getPerturbedVselfFieldBins()

template<unsigned int FEOrder, unsigned int FEOrderElectro>
const std::vector< distributedCPUVec< double > > & dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::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()

template<unsigned int FEOrder, unsigned int FEOrderElectro>
double dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::getStoredAdaptiveBallRadius ( ) const

get stored adaptive ball radius

◆ getVselfFieldBins()

template<unsigned int FEOrder, unsigned int FEOrderElectro>
const std::vector< distributedCPUVec< double > > & dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::getVselfFieldBins ( ) const

get const reference to solved vself fields

◆ getVselfFieldDerRBins()

template<unsigned int FEOrder, unsigned int FEOrderElectro>
const std::vector< distributedCPUVec< double > > & dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::getVselfFieldDerRBins ( ) const

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

◆ locateAtomsInBins()

template<unsigned int FEOrder, unsigned int FEOrderElectro>
void dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::locateAtomsInBins ( const dealii::DoFHandler< 3 > & dofHandler)
private

locate underlying fem nodes for atoms in bins.

◆ solveVselfInBins()

template<unsigned int FEOrder, unsigned int FEOrderElectro>
void dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::solveVselfInBins ( const std::shared_ptr< dftfe::basis::FEBasisOperations< double, double, dftfe::utils::MemorySpace::HOST > > & basisOperationsPtr,
const unsigned int offset,
const unsigned int matrixFreeQuadratureIdAX,
const dealii::AffineConstraints< double > & hangingPeriodicConstraintMatrix,
const std::vector< std::vector< double > > & imagePositions,
const std::vector< 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< int > > & bQuadAtomIdsAllAtoms,
std::map< dealii::CellId, std::vector< int > > & bQuadAtomIdsAllAtomsImages,
std::map< dealii::CellId, std::vector< unsigned int > > & bCellNonTrivialAtomIds,
std::vector< std::map< dealii::CellId, std::vector< unsigned int > > > & bCellNonTrivialAtomIdsBins,
std::map< dealii::CellId, std::vector< unsigned int > > & bCellNonTrivialAtomImageIds,
std::vector< std::map< dealii::CellId, std::vector< unsigned int > > > & bCellNonTrivialAtomImageIdsBins,
const std::vector< double > & smearingWidths,
std::vector< double > & smearedChargeScaling,
const unsigned int 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()

template<unsigned int FEOrder, unsigned int FEOrderElectro>
void dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::solveVselfInBinsPerturbedDomain ( const std::shared_ptr< dftfe::basis::FEBasisOperations< double, double, dftfe::utils::MemorySpace::HOST > > & basisOperationsPtr,
const unsigned int mfBaseDofHandlerIndex,
const unsigned int matrixFreeQuadratureIdAX,
const unsigned int offset,
const dealii::AffineConstraints< double > & hangingPeriodicConstraintMatrix,
const std::vector< std::vector< double > > & imagePositions,
const std::vector< int > & imageIds,
const std::vector< double > & imageCharges,
const std::vector< double > & smearingWidths,
const unsigned int 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()

template<unsigned int FEOrder, unsigned int FEOrderElectro>
void dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::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< 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

template<unsigned int FEOrder, unsigned int FEOrderElectro>
std::map<unsigned int, unsigned int> dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::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

template<unsigned int FEOrder, unsigned int FEOrderElectro>
std::vector<std::vector<double> > dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::d_atomLocations
private

storage for input atomLocations argument in createAtomBins function

◆ d_atomsInBin

template<unsigned int FEOrder, unsigned int FEOrderElectro>
std::vector<std::map<dealii::types::global_dof_index, double> > dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::d_atomsInBin
private

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

◆ d_bins

template<unsigned int FEOrder, unsigned int FEOrderElectro>
std::map<int, std::set<int> > dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::d_bins
private

map of binIds and atomIds

◆ d_binsImages

template<unsigned int FEOrder, unsigned int FEOrderElectro>
std::map<int, std::set<int> > dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::d_binsImages
private

map of binIds and atomIds and imageIds

◆ d_boundaryFlag

template<unsigned int FEOrder, unsigned int FEOrderElectro>
std::vector<std::map<dealii::types::global_dof_index, int> > dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::d_boundaryFlag
private

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

◆ d_boundaryFlagOnlyChargeId

template<unsigned int FEOrder, unsigned int FEOrderElectro>
std::vector<std::map<dealii::types::global_dof_index, int> > dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::d_boundaryFlagOnlyChargeId
private

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

◆ d_closestAtomBin

template<unsigned int FEOrder, unsigned int FEOrderElectro>
std::vector<std::map<dealii::types::global_dof_index, int> > dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::d_closestAtomBin
private

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

◆ d_constraintsOnlyHangingInfo

template<unsigned int FEOrder, unsigned int FEOrderElectro>
dftUtils::constraintMatrixInfo<dftfe::utils::MemorySpace::HOST> dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::d_constraintsOnlyHangingInfo
private

storage for optimized constraints handling object

◆ d_dftParams

template<unsigned int FEOrder, unsigned int FEOrderElectro>
const dftParameters& dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::d_dftParams
private

◆ d_dofClosestChargeLocationMap

template<unsigned int FEOrder, unsigned int FEOrderElectro>
std::vector<std::map<dealii::types::global_dof_index, dealii::Point<3> > > dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::d_dofClosestChargeLocationMap
private

map of global dof index to location of closest charge

◆ d_mpiCommParent

template<unsigned int FEOrder, unsigned int FEOrderElectro>
const MPI_Comm dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::d_mpiCommParent
private

◆ d_mpiInterCommKpts

template<unsigned int FEOrder, unsigned int FEOrderElectro>
const MPI_Comm dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::d_mpiInterCommKpts
private

◆ d_storedAdaptiveBallRadius

template<unsigned int FEOrder, unsigned int FEOrderElectro>
double dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::d_storedAdaptiveBallRadius
private

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

◆ d_vselfBinConstraintMatrices

template<unsigned int FEOrder, unsigned int FEOrderElectro>
std::vector<dealii::AffineConstraints<double> > dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::d_vselfBinConstraintMatrices
private

vector of constraint matrices for vself bins

◆ d_vselfBinField

template<unsigned int FEOrder, unsigned int FEOrderElectro>
std::vector<std::map<dealii::types::global_dof_index, double> > dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::d_vselfBinField
private

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

◆ d_vselfFieldBins

template<unsigned int FEOrder, unsigned int FEOrderElectro>
std::vector<distributedCPUVec<double> > dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::d_vselfFieldBins
private

solved vself solution field for each bin

◆ d_vselfFieldDerRBins

template<unsigned int FEOrder, unsigned int FEOrderElectro>
std::vector<distributedCPUVec<double> > dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::d_vselfFieldDerRBins
private

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

◆ d_vselfFieldPerturbedBins

template<unsigned int FEOrder, unsigned int FEOrderElectro>
std::vector<distributedCPUVec<double> > dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::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

template<unsigned int FEOrder, unsigned int FEOrderElectro>
const MPI_Comm dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::mpi_communicator
private

◆ n_mpi_processes

template<unsigned int FEOrder, unsigned int FEOrderElectro>
const unsigned int dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::n_mpi_processes
private

◆ pcout

template<unsigned int FEOrder, unsigned int FEOrderElectro>
dealii::ConditionalOStream dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::pcout
private

◆ this_mpi_process

template<unsigned int FEOrder, unsigned int FEOrderElectro>
const unsigned int dftfe::vselfBinsManager< FEOrder, FEOrderElectro >::this_mpi_process
private

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