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

Overloads dealii's distribute and distribute_local_to_global functions associated with constraints class. Stores the dealii's constraint matrix data into STL vectors for faster memory access costs. More...

#include <constraintMatrixInfo.h>

Public Member Functions

 constraintMatrixInfo ()
 
 ~constraintMatrixInfo ()
 
void initialize (const std::shared_ptr< const dealii::Utilities::MPI::Partitioner > &partitioner, const dealii::AffineConstraints< double > &constraintMatrixData)
 convert a given constraintMatrix to simple arrays (STL) for fast access
 
void distribute (distributedCPUVec< double > &fieldVector) const
 overloaded dealii internal function "distribute" which sets the slave node field values from master nodes
 
template<typename T>
void distribute (distributedCPUVec< T > &fieldVector, const unsigned int blockSize) const
 overloaded dealii internal function distribute for flattened dealii array which sets the slave node field values from master nodes
 
template<typename T>
void distribute (dftfe::linearAlgebra::MultiVector< T, memorySpace > &fieldVector) const
 
template<typename T>
void distribute_slave_to_master (distributedCPUVec< T > &fieldVector, const unsigned int blockSize) const
 transfers the contributions of slave nodes to master nodes using the constraint equation slave nodes are the nodes which are to the right of the constraint equation and master nodes are the nodes which are left of the constraint equation.
 
template<typename T>
void distribute_slave_to_master (dftfe::linearAlgebra::MultiVector< T, memorySpace > &fieldVector) const
 
void initializeScaledConstraints (const distributedCPUVec< double > &invSqrtMassVec)
 Scales the constraints with the inverse diagonal mass matrix so that the scaling of the vector can be done at the cell level.
 
void initializeScaledConstraints (const dftfe::utils::MemoryStorage< double, memorySpace > &invSqrtMassVec)
 
template<typename T>
void set_zero (distributedCPUVec< T > &fieldVector, const unsigned int blockSize) const
 sets field values at constrained nodes to be zero
 
template<typename T>
void set_zero (dftfe::linearAlgebra::MultiVector< T, memorySpace > &fieldVector) const
 
void clear ()
 

Private Attributes

std::vector< dealii::types::global_dof_index > d_rowIdsGlobal
 
std::vector< dealii::types::global_dof_index > d_rowIdsLocal
 
std::vector< dealii::types::global_dof_index > d_columnIdsLocal
 
std::vector< dealii::types::global_dof_index > d_columnIdsGlobal
 
std::vector< double > d_columnValues
 
std::vector< double > d_inhomogenities
 
std::vector< dealii::types::global_dof_index > d_rowSizes
 
std::vector< dealii::types::global_dof_index > d_localIndexMapUnflattenedToFlattened
 

Detailed Description

template<dftfe::utils::MemorySpace memorySpace>
class dftfe::dftUtils::constraintMatrixInfo< memorySpace >

Overloads dealii's distribute and distribute_local_to_global functions associated with constraints class. Stores the dealii's constraint matrix data into STL vectors for faster memory access costs.

Author
Phani Motamarri

Constructor & Destructor Documentation

◆ constraintMatrixInfo()

template<dftfe::utils::MemorySpace memorySpace>
dftfe::dftUtils::constraintMatrixInfo< memorySpace >::constraintMatrixInfo ( )

class constructor

◆ ~constraintMatrixInfo()

class destructor

Member Function Documentation

◆ clear()

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::dftUtils::constraintMatrixInfo< memorySpace >::clear ( )

clear data members

◆ distribute() [1/3]

template<dftfe::utils::MemorySpace memorySpace>
template<typename T>
void dftfe::dftUtils::constraintMatrixInfo< memorySpace >::distribute ( dftfe::linearAlgebra::MultiVector< T, memorySpace > & fieldVector) const

◆ distribute() [2/3]

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::dftUtils::constraintMatrixInfo< memorySpace >::distribute ( distributedCPUVec< double > & fieldVector) const

overloaded dealii internal function "distribute" which sets the slave node field values from master nodes

Parameters
fieldVectorparallel dealii vector

◆ distribute() [3/3]

template<dftfe::utils::MemorySpace memorySpace>
template<typename T>
void dftfe::dftUtils::constraintMatrixInfo< memorySpace >::distribute ( distributedCPUVec< T > & fieldVector,
const unsigned int blockSize ) const

overloaded dealii internal function distribute for flattened dealii array which sets the slave node field values from master nodes

Parameters
blockSizenumber of components for a given node

◆ distribute_slave_to_master() [1/2]

template<dftfe::utils::MemorySpace memorySpace>
template<typename T>
void dftfe::dftUtils::constraintMatrixInfo< memorySpace >::distribute_slave_to_master ( dftfe::linearAlgebra::MultiVector< T, memorySpace > & fieldVector) const

◆ distribute_slave_to_master() [2/2]

template<dftfe::utils::MemorySpace memorySpace>
template<typename T>
void dftfe::dftUtils::constraintMatrixInfo< memorySpace >::distribute_slave_to_master ( distributedCPUVec< T > & fieldVector,
const unsigned int blockSize ) const

transfers the contributions of slave nodes to master nodes using the constraint equation slave nodes are the nodes which are to the right of the constraint equation and master nodes are the nodes which are left of the constraint equation.

Parameters
fieldVectorparallel dealii vector which is the result of matrix-vector product(vmult) withot taking care of constraints
blockSizenumber of components for a given node

◆ initialize()

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::dftUtils::constraintMatrixInfo< memorySpace >::initialize ( const std::shared_ptr< const dealii::Utilities::MPI::Partitioner > & partitioner,
const dealii::AffineConstraints< double > & constraintMatrixData )

convert a given constraintMatrix to simple arrays (STL) for fast access

Parameters
partitionerassociated with the dealii vector
constraintMatrixDatadealii constraint matrix from which the data is extracted

◆ initializeScaledConstraints() [1/2]

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::dftUtils::constraintMatrixInfo< memorySpace >::initializeScaledConstraints ( const dftfe::utils::MemoryStorage< double, memorySpace > & invSqrtMassVec)

◆ initializeScaledConstraints() [2/2]

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::dftUtils::constraintMatrixInfo< memorySpace >::initializeScaledConstraints ( const distributedCPUVec< double > & invSqrtMassVec)

Scales the constraints with the inverse diagonal mass matrix so that the scaling of the vector can be done at the cell level.

Parameters
invSqrtMassVecthe inverse diagonal mass matrix

◆ set_zero() [1/2]

template<dftfe::utils::MemorySpace memorySpace>
template<typename T>
void dftfe::dftUtils::constraintMatrixInfo< memorySpace >::set_zero ( dftfe::linearAlgebra::MultiVector< T, memorySpace > & fieldVector) const

◆ set_zero() [2/2]

template<dftfe::utils::MemorySpace memorySpace>
template<typename T>
void dftfe::dftUtils::constraintMatrixInfo< memorySpace >::set_zero ( distributedCPUVec< T > & fieldVector,
const unsigned int blockSize ) const

sets field values at constrained nodes to be zero

Parameters
fieldVectorparallel dealii vector with fields stored in a flattened format
blockSizenumber of field components for a given node

Member Data Documentation

◆ d_columnIdsGlobal

template<dftfe::utils::MemorySpace memorySpace>
std::vector<dealii::types::global_dof_index> dftfe::dftUtils::constraintMatrixInfo< memorySpace >::d_columnIdsGlobal
private

◆ d_columnIdsLocal

template<dftfe::utils::MemorySpace memorySpace>
std::vector<dealii::types::global_dof_index> dftfe::dftUtils::constraintMatrixInfo< memorySpace >::d_columnIdsLocal
private

◆ d_columnValues

template<dftfe::utils::MemorySpace memorySpace>
std::vector<double> dftfe::dftUtils::constraintMatrixInfo< memorySpace >::d_columnValues
private

◆ d_inhomogenities

template<dftfe::utils::MemorySpace memorySpace>
std::vector<double> dftfe::dftUtils::constraintMatrixInfo< memorySpace >::d_inhomogenities
private

◆ d_localIndexMapUnflattenedToFlattened

template<dftfe::utils::MemorySpace memorySpace>
std::vector<dealii::types::global_dof_index> dftfe::dftUtils::constraintMatrixInfo< memorySpace >::d_localIndexMapUnflattenedToFlattened
private

◆ d_rowIdsGlobal

template<dftfe::utils::MemorySpace memorySpace>
std::vector<dealii::types::global_dof_index> dftfe::dftUtils::constraintMatrixInfo< memorySpace >::d_rowIdsGlobal
private

◆ d_rowIdsLocal

template<dftfe::utils::MemorySpace memorySpace>
std::vector<dealii::types::global_dof_index> dftfe::dftUtils::constraintMatrixInfo< memorySpace >::d_rowIdsLocal
private

◆ d_rowSizes

template<dftfe::utils::MemorySpace memorySpace>
std::vector<dealii::types::global_dof_index> dftfe::dftUtils::constraintMatrixInfo< memorySpace >::d_rowSizes
private

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