DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
|
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 |
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.
dftfe::dftUtils::constraintMatrixInfo< memorySpace >::constraintMatrixInfo | ( | ) |
class constructor
dftfe::dftUtils::constraintMatrixInfo< memorySpace >::~constraintMatrixInfo | ( | ) |
class destructor
void dftfe::dftUtils::constraintMatrixInfo< memorySpace >::clear | ( | ) |
clear data members
void dftfe::dftUtils::constraintMatrixInfo< memorySpace >::distribute | ( | dftfe::linearAlgebra::MultiVector< T, memorySpace > & | fieldVector | ) | const |
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
fieldVector | parallel dealii vector |
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
blockSize | number of components for a given node |
void dftfe::dftUtils::constraintMatrixInfo< memorySpace >::distribute_slave_to_master | ( | dftfe::linearAlgebra::MultiVector< T, memorySpace > & | fieldVector | ) | const |
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.
fieldVector | parallel dealii vector which is the result of matrix-vector product(vmult) withot taking care of constraints |
blockSize | number of components for a given node |
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
partitioner | associated with the dealii vector |
constraintMatrixData | dealii constraint matrix from which the data is extracted |
void dftfe::dftUtils::constraintMatrixInfo< memorySpace >::initializeScaledConstraints | ( | const dftfe::utils::MemoryStorage< double, memorySpace > & | invSqrtMassVec | ) |
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.
invSqrtMassVec | the inverse diagonal mass matrix |
void dftfe::dftUtils::constraintMatrixInfo< memorySpace >::set_zero | ( | dftfe::linearAlgebra::MultiVector< T, memorySpace > & | fieldVector | ) | const |
void dftfe::dftUtils::constraintMatrixInfo< memorySpace >::set_zero | ( | distributedCPUVec< T > & | fieldVector, |
const unsigned int | blockSize ) const |
sets field values at constrained nodes to be zero
fieldVector | parallel dealii vector with fields stored in a flattened format |
blockSize | number of field components for a given node |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |