DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
|
Contains generic utils functions related to custom partitioned flattened dealii vector. More...
Functions | |
void | createParallelConstraintMatrixFromSerial (const dealii::Triangulation< 3, 3 > &serTria, const dealii::DoFHandler< 3 > &dofHandlerPar, const MPI_Comm &mpi_comm_parent, const MPI_Comm &mpi_comm_domain, const std::vector< std::vector< double > > &domainBoundingVectors, dealii::AffineConstraints< double > &periodicHangingConstraints, dealii::AffineConstraints< double > &onlyHangingConstraints, const int verbosity, const bool periodicX, const bool periodicY, const bool periodicZ) |
Create constraint matrix using serial mesh. Temporary fix for a bug (Issue #7053) in deal.ii until it is resolved. | |
template<typename T> | |
void | createDealiiVector (const std::shared_ptr< const dealii::Utilities::MPI::Partitioner > &partitioner, const unsigned int blockSize, distributedCPUVec< T > &flattenedArray) |
Creates a custom partitioned flattened dealii vector. stores multiple components asociated with a node sequentially. | |
void | computeCellLocalIndexSetMap (const std::shared_ptr< const utils::mpi::MPIPatternP2P< dftfe::utils::MemorySpace::HOST > > &partitioner, const dealii::MatrixFree< 3, double > &matrix_free_data, const unsigned int mfDofHandlerIndex, const unsigned int blockSize, std::vector< std::vector< dealii::types::global_dof_index > > &flattenedArrayMacroCellLocalProcIndexId, std::vector< std::vector< dealii::types::global_dof_index > > &flattenedArrayCellLocalProcIndexId) |
Creates a cell local index set map for flattened array. | |
void | computeCellLocalIndexSetMap (const std::shared_ptr< const utils::mpi::MPIPatternP2P< dftfe::utils::MemorySpace::HOST > > &partitioner, const dealii::MatrixFree< 3, double > &matrix_free_data, const unsigned int mfDofHandlerIndex, const unsigned int blockSize, std::vector< dealii::types::global_dof_index > &flattenedArrayMacroCellLocalProcIndexId, std::vector< unsigned int > &normalCellIdToMacroCellIdMap, std::vector< unsigned int > ¯oCellIdToNormalCellIdMap, std::vector< dealii::types::global_dof_index > &flattenedArrayCellLocalProcIndexId) |
Creates a cell local index set map for flattened array. | |
void | computeCellLocalIndexSetMap (const std::shared_ptr< const utils::mpi::MPIPatternP2P< dftfe::utils::MemorySpace::HOST > > &partitioner, const dealii::MatrixFree< 3, double > &matrix_free_data, const unsigned int mfDofHandlerIndex, const unsigned int blockSize, std::vector< dealii::types::global_dof_index > &flattenedArrayCellLocalProcIndexId) |
Creates a cell local index set map for flattened array. | |
void | copyFlattenedSTLVecToSingleCompVec (const double *flattenedArray, const unsigned int totalNumberComponents, const unsigned int localVectorSize, const std::pair< unsigned int, unsigned int > componentIndexRange, std::vector< distributedCPUVec< double > > &componentVectors) |
Copies a single field component from a flattenedArray STL vector containing multiple component fields to a single field parallel distributed vector. | |
void | copyFlattenedDealiiVecToSingleCompVec (const distributedCPUVec< double > &flattenedArray, const unsigned int totalNumberComponents, const std::pair< unsigned int, unsigned int > componentIndexRange, std::vector< distributedCPUVec< double > > &componentVectors, const bool isFlattenedDealiiGhostValuesUpdated=false) |
Copies a single field component from a flattenedArray parallel distributed vector containing multiple component fields to a single field parallel distributed vector. | |
void | copySingleCompVecToFlattenedDealiiVec (distributedCPUVec< double > &flattenedArray, const unsigned int totalNumberComponents, const std::pair< unsigned int, unsigned int > componentIndexRange, const std::vector< distributedCPUVec< double > > &componentVectors) |
Copies to a flattenedArray parallel distributed vector containing multiple component fields from a single field parallel distributed vector. | |
void | copySingleCompVecToFlattenedSTLVec (std::vector< double > &flattenedArray, const unsigned int totalNumberComponents, const std::pair< unsigned int, unsigned int > componentIndexRange, const std::vector< distributedCPUVec< double > > &componentVectors) |
Copies to a flattenedArray stl vector containing multiple component fields from a single field parallel distributed vector. | |
std::pair< dealii::Point< 3 >, dealii::Point< 3 > > | createBoundingBoxTriaLocallyOwned (const dealii::DoFHandler< 3 > &dofHandler) |
void | classifyInteriorSurfaceNodesInCell (const dealii::MatrixFree< 3, double > &matrix_free_data, const unsigned int mfDofHandlerIndex, std::vector< unsigned int > &nodesPerCellClassificationMap) |
void | classifyInteriorSurfaceNodesInGlobalArray (const dealii::MatrixFree< 3, double > &matrix_free_data, const unsigned int mfDofHandlerIndex, const dealii::AffineConstraints< double > &constraintMatrix, std::vector< unsigned int > &nodesPerCellClassificationMap, std::vector< unsigned int > &globalArrayClassificationMap) |
Contains generic utils functions related to custom partitioned flattened dealii vector.
void dftfe::vectorTools::classifyInteriorSurfaceNodesInCell | ( | const dealii::MatrixFree< 3, double > & | matrix_free_data, |
const unsigned int | mfDofHandlerIndex, | ||
std::vector< unsigned int > & | nodesPerCellClassificationMap ) |
void dftfe::vectorTools::classifyInteriorSurfaceNodesInGlobalArray | ( | const dealii::MatrixFree< 3, double > & | matrix_free_data, |
const unsigned int | mfDofHandlerIndex, | ||
const dealii::AffineConstraints< double > & | constraintMatrix, | ||
std::vector< unsigned int > & | nodesPerCellClassificationMap, | ||
std::vector< unsigned int > & | globalArrayClassificationMap ) |
void dftfe::vectorTools::computeCellLocalIndexSetMap | ( | const std::shared_ptr< const utils::mpi::MPIPatternP2P< dftfe::utils::MemorySpace::HOST > > & | partitioner, |
const dealii::MatrixFree< 3, double > & | matrix_free_data, | ||
const unsigned int | mfDofHandlerIndex, | ||
const unsigned int | blockSize, | ||
std::vector< dealii::types::global_dof_index > & | flattenedArrayCellLocalProcIndexId ) |
Creates a cell local index set map for flattened array.
partitioner | associated with the flattened array |
matrix_free_data | object pointer associated with the matrix free data structure |
blockSize | number of components associated with each node |
void dftfe::vectorTools::computeCellLocalIndexSetMap | ( | const std::shared_ptr< const utils::mpi::MPIPatternP2P< dftfe::utils::MemorySpace::HOST > > & | partitioner, |
const dealii::MatrixFree< 3, double > & | matrix_free_data, | ||
const unsigned int | mfDofHandlerIndex, | ||
const unsigned int | blockSize, | ||
std::vector< dealii::types::global_dof_index > & | flattenedArrayMacroCellLocalProcIndexId, | ||
std::vector< unsigned int > & | normalCellIdToMacroCellIdMap, | ||
std::vector< unsigned int > & | macroCellIdToNormalCellIdMap, | ||
std::vector< dealii::types::global_dof_index > & | flattenedArrayCellLocalProcIndexId ) |
Creates a cell local index set map for flattened array.
partitioner | associated with the flattened array |
matrix_free_data | object pointer associated with the matrix free data structure |
blockSize | number of components associated with each node |
void dftfe::vectorTools::computeCellLocalIndexSetMap | ( | const std::shared_ptr< const utils::mpi::MPIPatternP2P< dftfe::utils::MemorySpace::HOST > > & | partitioner, |
const dealii::MatrixFree< 3, double > & | matrix_free_data, | ||
const unsigned int | mfDofHandlerIndex, | ||
const unsigned int | blockSize, | ||
std::vector< std::vector< dealii::types::global_dof_index > > & | flattenedArrayMacroCellLocalProcIndexId, | ||
std::vector< std::vector< dealii::types::global_dof_index > > & | flattenedArrayCellLocalProcIndexId ) |
Creates a cell local index set map for flattened array.
partitioner | associated with the flattened array |
matrix_free_data | object pointer associated with the matrix free data structure |
blockSize | number of components associated with each node |
void dftfe::vectorTools::copyFlattenedDealiiVecToSingleCompVec | ( | const distributedCPUVec< double > & | flattenedArray, |
const unsigned int | totalNumberComponents, | ||
const std::pair< unsigned int, unsigned int > | componentIndexRange, | ||
std::vector< distributedCPUVec< double > > & | componentVectors, | ||
const bool | isFlattenedDealiiGhostValuesUpdated = false ) |
Copies a single field component from a flattenedArray parallel distributed vector containing multiple component fields to a single field parallel distributed vector.
[in] | flattenedArray | flattened parallel distributed vector with multiple component fields |
[in] | totalNumberComponents | total number of component fiels in flattenedArray |
[in] | componentIndexRange | desired range field components [componentIndexRange.first,componentIndexRange.second) |
[out] | componentVectors | vector of parallel distributed vectors with fields corresponding to componentIndexRange. componentVectors is expected to be of the size componentIndexRange.second-componentIndexRange.first. Further, each entry of componentVectors is assumed to be already initialized with the same single component partitioner used in the creation of the flattenedArray partitioner. |
[in] | isFlattenedDealiiGhostValuesUpdated | default is false. Use true for optimization if update ghost values has already been called in the flattened dealii vec. |
void dftfe::vectorTools::copyFlattenedSTLVecToSingleCompVec | ( | const double * | flattenedArray, |
const unsigned int | totalNumberComponents, | ||
const unsigned int | localVectorSize, | ||
const std::pair< unsigned int, unsigned int > | componentIndexRange, | ||
std::vector< distributedCPUVec< double > > & | componentVectors ) |
Copies a single field component from a flattenedArray STL vector containing multiple component fields to a single field parallel distributed vector.
[in] | flattenedArray | flattened STL vector with multiple component fields |
[in] | totalNumberComponents | total number of component fiels in flattenedArray |
[in] | componentIndexRange | desired range field components [componentIndexRange.first,componentIndexRange.second) |
[out] | componentVectors | vector of parallel distributed vectors with fields corresponding to componentIndexRange. componentVectors is expected to be of the size componentIndexRange.second-componentIndexRange.first. Further, each entry of componentVectors is assumed to be already initialized with the same single component partitioner used in the creation of the flattenedArray partitioner. |
void dftfe::vectorTools::copySingleCompVecToFlattenedDealiiVec | ( | distributedCPUVec< double > & | flattenedArray, |
const unsigned int | totalNumberComponents, | ||
const std::pair< unsigned int, unsigned int > | componentIndexRange, | ||
const std::vector< distributedCPUVec< double > > & | componentVectors ) |
Copies to a flattenedArray parallel distributed vector containing multiple component fields from a single field parallel distributed vector.
[out] | flattenedArray | flattened parallel distributed vector with multiple component fields |
[in] | totalNumberComponents | total number of component fiels in flattenedArray |
[in] | componentIndexRange | desired range field components [componentIndexRange.first,componentIndexRange.second) |
[in] | componentVectors | vector of parallel distributed vectors with fields corresponding to componentIndexRange. componentVectors is expected to be of the size componentIndexRange.second-componentIndexRange.first. Further, each entry of componentVectors is assumed to be already initialized with the same single component partitioner used in the creation of the flattenedArray partitioner. |
void dftfe::vectorTools::copySingleCompVecToFlattenedSTLVec | ( | std::vector< double > & | flattenedArray, |
const unsigned int | totalNumberComponents, | ||
const std::pair< unsigned int, unsigned int > | componentIndexRange, | ||
const std::vector< distributedCPUVec< double > > & | componentVectors ) |
Copies to a flattenedArray stl vector containing multiple component fields from a single field parallel distributed vector.
[out] | flattenedArray | flattened stl vector with multiple component fields |
[in] | totalNumberComponents | total number of component fiels in flattenedArray |
[in] | componentIndexRange | desired range field components [componentIndexRange.first,componentIndexRange.second) |
[in] | componentVectors | vector of parallel distributed vectors with fields corresponding to componentIndexRange. componentVectors is expected to be of the size componentIndexRange.second-componentIndexRange.first. Further, each entry of componentVectors is assumed to be already initialized with the same single component partitioner used in the creation of the flattenedArray partitioner. |
std::pair< dealii::Point< 3 >, dealii::Point< 3 > > dftfe::vectorTools::createBoundingBoxTriaLocallyOwned | ( | const dealii::DoFHandler< 3 > & | dofHandler | ) |
void dftfe::vectorTools::createDealiiVector | ( | const std::shared_ptr< const dealii::Utilities::MPI::Partitioner > & | partitioner, |
const unsigned int | blockSize, | ||
distributedCPUVec< T > & | flattenedArray ) |
Creates a custom partitioned flattened dealii vector. stores multiple components asociated with a node sequentially.
partitioner | associated with single component vector |
blockSize | number of components associated with each node |
void dftfe::vectorTools::createParallelConstraintMatrixFromSerial | ( | const dealii::Triangulation< 3, 3 > & | serTria, |
const dealii::DoFHandler< 3 > & | dofHandlerPar, | ||
const MPI_Comm & | mpi_comm_parent, | ||
const MPI_Comm & | mpi_comm_domain, | ||
const std::vector< std::vector< double > > & | domainBoundingVectors, | ||
dealii::AffineConstraints< double > & | periodicHangingConstraints, | ||
dealii::AffineConstraints< double > & | onlyHangingConstraints, | ||
const int | verbosity, | ||
const bool | periodicX, | ||
const bool | periodicY, | ||
const bool | periodicZ ) |
Create constraint matrix using serial mesh. Temporary fix for a bug (Issue #7053) in deal.ii until it is resolved.
[in] | serial | Triangulation which must be exactly same as the parallel triangulation associated with dofHandlerPar |
[in] | parallel | DofHandler |
[out] | periodic | hanging constraints. |
[out] | only | hanging constraints |