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

This class generates and stores adaptive finite element meshes for the real-space dft problem. More...

#include <triangulationManager.h>

Public Member Functions

 triangulationManager (const MPI_Comm &mpi_comm_parent, const MPI_Comm &mpi_comm_domain, const MPI_Comm &interpoolcomm, const MPI_Comm &interBandGroupComm, const unsigned int FEOrder, const dftParameters &dftParams)
 Constructor.
 
 ~triangulationManager ()
 
void generateSerialUnmovedAndParallelMovedUnmovedMesh (const std::vector< std::vector< double > > &atomLocations, const std::vector< std::vector< double > > &imageAtomLocations, const std::vector< int > &imageIds, const std::vector< double > &nearestAtomDistances, const std::vector< std::vector< double > > &domainBoundingVectors, const bool generateSerialTria)
 generates parallel moved and unmoved meshes, and serial unmoved mesh.
 
void generateCoarseMeshesForRestart (const std::vector< std::vector< double > > &atomLocations, const std::vector< std::vector< double > > &imageAtomLocations, const std::vector< int > &imageIds, const std::vector< double > &nearestAtomDistances, const std::vector< std::vector< double > > &domainBoundingVectors, const bool generateSerialTria)
 generates the coarse meshes for restart.
 
void generateAutomaticMeshApriori (const dealii::DoFHandler< 3 > &dofHandler, dealii::parallel::distributed::Triangulation< 3 > &parallelTriangulation, const std::vector< distributedCPUVec< double > > &eigenVectorsArrayIn, const unsigned int FEOrder)
 returns generates A-posteriori refined mesh
 
dealii::parallel::distributed::Triangulation< 3 > & getSerialMeshUnmoved ()
 returns constant reference to serial unmoved triangulation
 
dealii::parallel::distributed::Triangulation< 3 > & getParallelMeshMoved ()
 returns reference to parallel moved triangulation
 
dealii::parallel::distributed::Triangulation< 3 > & getParallelMeshUnmoved ()
 returns constant reference to parallel unmoved triangulation
 
void resetMesh (dealii::parallel::distributed::Triangulation< 3 > &parallelTriangulationA, dealii::parallel::distributed::Triangulation< 3 > &parallelTriangulationB)
 resets the vertices of meshB moved to vertices of meshA.
 
void generateResetMeshes (const std::vector< std::vector< double > > &domainBoundingVectors, const bool generateSerialTria)
 generates reset meshes to the last mesh generated by auto mesh approach. This is used in Gaussian update of meshes during electrostatics
 
void saveTriangulationsSolutionVectors (std::string path, const unsigned int feOrder, const unsigned int nComponents, const std::vector< const distributedCPUVec< double > * > &solutionVectors, const MPI_Comm &interpoolComm, const MPI_Comm &interBandGroupComm)
 serialize the triangulations and the associated solution vectors
 
void loadTriangulationsSolutionVectors (std::string path, const unsigned int feOrder, const unsigned int nComponents, std::vector< distributedCPUVec< double > * > &solutionVectors)
 de-serialize the triangulations and the associated solution vectors
 
void generateMesh (dealii::parallel::distributed::Triangulation< 3 > &parallelTriangulation, dealii::parallel::distributed::Triangulation< 3 > &serialTriangulation, std::vector< std::vector< bool > > &parallelTriaCurrentRefinement, std::vector< std::vector< bool > > &serialTriaCurrentRefinement, const bool generateSerialTria=false, const bool enableManualRepartitioning=false)
 internal function which generates a parallel and serial mesh using a adaptive refinement strategy.
 

Private Member Functions

void generateCoarseMesh (dealii::parallel::distributed::Triangulation< 3 > &parallelTriangulation)
 internal function which generates a coarse mesh which is required for the load function call in restarts.
 
bool refinementAlgorithmA (dealii::parallel::distributed::Triangulation< 3 > &parallelTriangulation, std::vector< unsigned int > &locallyOwnedCellsRefineFlags, std::map< dealii::CellId, unsigned int > &cellIdToCellRefineFlagMapLocal, const bool smoothenCellsOnPeriodicBoundary=false, const double smootheningFactor=2.0)
 internal function which sets refinement flags based on a custom created algorithm
 
bool consistentPeriodicBoundaryRefinement (dealii::parallel::distributed::Triangulation< 3 > &parallelTriangulation, std::vector< unsigned int > &locallyOwnedCellsRefineFlags, std::map< dealii::CellId, unsigned int > &cellIdToCellRefineFlagMapLocal)
 internal function which sets refinement flags to have consistent refinement across periodic boundary
 
bool checkPeriodicSurfaceRefinementConsistency (dealii::parallel::distributed::Triangulation< 3 > &parallelTriangulation)
 check that triangulation has consistent refinement across periodic boundary including for ghost cells
 
bool checkConstraintsConsistency (dealii::parallel::distributed::Triangulation< 3 > &parallelTriangulation)
 check that FEOrder=1 dofHandler using the triangulation has parallel consistent combined hanging node and periodic constraints
 
void refineSerialMesh (const std::map< dealii::CellId, unsigned int > &cellIdToCellRefineFlagMapLocal, const MPI_Comm &mpi_comm, dealii::parallel::distributed::Triangulation< 3 > &serialTriangulation, const dealii::parallel::distributed::Triangulation< 3 > &parallelTriangulation, std::vector< bool > &serialTriaCurrentRefinement)
 internal function which refines the serial mesh based on refinement flags from parallel mesh. This ensures that we get the same mesh in serial and parallel.
 
void saveSupportTriangulations (std::string path)
 internal function to serialize support triangulations. No solution data is attached to them
 
void loadSupportTriangulations (std::string path)
 internal function to de-serialize support triangulations. No solution data is read from them
 

Private Attributes

dealii::parallel::distributed::Triangulation< 3 > d_parallelTriangulationUnmoved
 
dealii::parallel::distributed::Triangulation< 3 > d_parallelTriangulationMoved
 
dealii::parallel::distributed::Triangulation< 3 > d_serialTriangulationUnmoved
 
std::vector< std::vector< bool > > d_parallelTriaCurrentRefinement
 
std::vector< std::vector< bool > > d_serialTriaCurrentRefinement
 
std::vector< std::vector< double > > d_atomPositions
 
std::vector< std::vector< double > > d_imageAtomPositions
 
std::vector< std::vector< double > > d_meshSizes
 
std::vector< int > d_imageIds
 
std::vector< double > d_nearestAtomDistances
 
std::vector< std::vector< double > > d_domainBoundingVectors
 
const unsigned int d_max_refinement_steps = 40
 
const unsigned int d_FEOrder
 
const dftParametersd_dftParams
 
const MPI_Comm d_mpiCommParent
 
const MPI_Comm mpi_communicator
 
const MPI_Comm interpoolcomm
 
const MPI_Comm interBandGroupComm
 
const unsigned int this_mpi_process
 
const unsigned int n_mpi_processes
 
dealii::ConditionalOStream pcout
 
dealii::TimerOutput computing_timer
 

Detailed Description

This class generates and stores adaptive finite element meshes for the real-space dft problem.

The class uses an adpative mesh generation strategy to generate finite element mesh for given domain based on five input parameters: BASE MESH SIZE, ATOM BALL RADIUS, MESH SIZE ATOM BALL, MESH SIZE NEAR ATOM and MAX REFINEMENT STEPS (Refer to utils/dftParameters.cc for their corresponding internal variable names). Additionaly, this class also applies periodicity to mesh. The class stores two types of meshes: moved and unmoved. They are essentially the same meshes, except that we move the nodes of the moved mesh (in the meshMovement class) such that the atoms lie on the nodes. However, once the mesh is moved, dealii has issues using that mesh for further refinement, which is why we also carry an unmoved triangulation.

Author
Phani Motamarri, Sambit Das, Krishnendu Ghosh

Constructor & Destructor Documentation

◆ triangulationManager()

dftfe::triangulationManager::triangulationManager ( const MPI_Comm & mpi_comm_parent,
const MPI_Comm & mpi_comm_domain,
const MPI_Comm & interpoolcomm,
const MPI_Comm & interBandGroupComm,
const unsigned int FEOrder,
const dftParameters & dftParams )

Constructor.

Parameters
mpi_comm_parentparent mpi communicator
mpi_comm_domaindomain decomposition mpi communicator
interpool_commmpi interpool communicator over k points
interBandGroupCommmpi interpool communicator over band groups

◆ ~triangulationManager()

dftfe::triangulationManager::~triangulationManager ( )

triangulationManager destructor

Member Function Documentation

◆ checkConstraintsConsistency()

bool dftfe::triangulationManager::checkConstraintsConsistency ( dealii::parallel::distributed::Triangulation< 3 > & parallelTriangulation)
private

check that FEOrder=1 dofHandler using the triangulation has parallel consistent combined hanging node and periodic constraints

◆ checkPeriodicSurfaceRefinementConsistency()

bool dftfe::triangulationManager::checkPeriodicSurfaceRefinementConsistency ( dealii::parallel::distributed::Triangulation< 3 > & parallelTriangulation)
private

check that triangulation has consistent refinement across periodic boundary including for ghost cells

◆ consistentPeriodicBoundaryRefinement()

bool dftfe::triangulationManager::consistentPeriodicBoundaryRefinement ( dealii::parallel::distributed::Triangulation< 3 > & parallelTriangulation,
std::vector< unsigned int > & locallyOwnedCellsRefineFlags,
std::map< dealii::CellId, unsigned int > & cellIdToCellRefineFlagMapLocal )
private

internal function which sets refinement flags to have consistent refinement across periodic boundary

Returns
bool boolean flag is any local cell has refinement flag set

◆ generateAutomaticMeshApriori()

void dftfe::triangulationManager::generateAutomaticMeshApriori ( const dealii::DoFHandler< 3 > & dofHandler,
dealii::parallel::distributed::Triangulation< 3 > & parallelTriangulation,
const std::vector< distributedCPUVec< double > > & eigenVectorsArrayIn,
const unsigned int FEOrder )

returns generates A-posteriori refined mesh

Parameters
dofHandlercorresponds to starting mesh which has to refined
parallelTriangulationcorresponds to starting triangulation
eigenVectorsArrayInsolution vectors used to compute errors in each cell required for refinement
FEOrderfinite-element interpolating polynomial

◆ generateCoarseMesh()

void dftfe::triangulationManager::generateCoarseMesh ( dealii::parallel::distributed::Triangulation< 3 > & parallelTriangulation)
private

internal function which generates a coarse mesh which is required for the load function call in restarts.

◆ generateCoarseMeshesForRestart()

void dftfe::triangulationManager::generateCoarseMeshesForRestart ( const std::vector< std::vector< double > > & atomLocations,
const std::vector< std::vector< double > > & imageAtomLocations,
const std::vector< int > & imageIds,
const std::vector< double > & nearestAtomDistances,
const std::vector< std::vector< double > > & domainBoundingVectors,
const bool generateSerialTria )

generates the coarse meshes for restart.

Parameters
atomLocationsvector containing cartesian coordinates at atoms with respect to origin (center of domain).
imageAtomLocationsvector containing cartesian coordinates of image atoms with respect to origin.
domainBoundingVectorsvector of domain bounding vectors (refer to description of input parameters.
generateSerialMeshbool to toggle to generation of serial tria

◆ generateMesh()

void dftfe::triangulationManager::generateMesh ( dealii::parallel::distributed::Triangulation< 3 > & parallelTriangulation,
dealii::parallel::distributed::Triangulation< 3 > & serialTriangulation,
std::vector< std::vector< bool > > & parallelTriaCurrentRefinement,
std::vector< std::vector< bool > > & serialTriaCurrentRefinement,
const bool generateSerialTria = false,
const bool enableManualRepartitioning = false )

internal function which generates a parallel and serial mesh using a adaptive refinement strategy.

◆ generateResetMeshes()

void dftfe::triangulationManager::generateResetMeshes ( const std::vector< std::vector< double > > & domainBoundingVectors,
const bool generateSerialTria )

generates reset meshes to the last mesh generated by auto mesh approach. This is used in Gaussian update of meshes during electrostatics

◆ generateSerialUnmovedAndParallelMovedUnmovedMesh()

void dftfe::triangulationManager::generateSerialUnmovedAndParallelMovedUnmovedMesh ( const std::vector< std::vector< double > > & atomLocations,
const std::vector< std::vector< double > > & imageAtomLocations,
const std::vector< int > & imageIds,
const std::vector< double > & nearestAtomDistances,
const std::vector< std::vector< double > > & domainBoundingVectors,
const bool generateSerialTria )

generates parallel moved and unmoved meshes, and serial unmoved mesh.

Parameters
atomLocationsvector containing cartesian coordinates at atoms with respect to origin (center of domain).
imageAtomLocationsvector containing cartesian coordinates of image atoms with respect to origin.
domainBoundingVectorsvector of domain bounding vectors (refer to description of input parameters.
generateSerialMeshbool to toggle to generation of serial tria

◆ getParallelMeshMoved()

dealii::parallel::distributed::Triangulation< 3 > & dftfe::triangulationManager::getParallelMeshMoved ( )

returns reference to parallel moved triangulation

◆ getParallelMeshUnmoved()

dealii::parallel::distributed::Triangulation< 3 > & dftfe::triangulationManager::getParallelMeshUnmoved ( )

returns constant reference to parallel unmoved triangulation

◆ getSerialMeshUnmoved()

dealii::parallel::distributed::Triangulation< 3 > & dftfe::triangulationManager::getSerialMeshUnmoved ( )

returns constant reference to serial unmoved triangulation

◆ loadSupportTriangulations()

void dftfe::triangulationManager::loadSupportTriangulations ( std::string path)
private

internal function to de-serialize support triangulations. No solution data is read from them

◆ loadTriangulationsSolutionVectors()

void dftfe::triangulationManager::loadTriangulationsSolutionVectors ( std::string path,
const unsigned int feOrder,
const unsigned int nComponents,
std::vector< distributedCPUVec< double > * > & solutionVectors )

de-serialize the triangulations and the associated solution vectors

Parameters
[input]feOrderfinite element polynomial order of the dofHandler on which solution vectors to be de-serialized are based upon
[input]nComponentsnumber of components of the dofHandler on which solution vectors to be de-serialized are based upon
[output]solutionVectorsvector of parallel distributed de-serialized solution vectors. The vector length must match the input vector length used in the call to saveTriangulationSolutionVectors

◆ refinementAlgorithmA()

bool dftfe::triangulationManager::refinementAlgorithmA ( dealii::parallel::distributed::Triangulation< 3 > & parallelTriangulation,
std::vector< unsigned int > & locallyOwnedCellsRefineFlags,
std::map< dealii::CellId, unsigned int > & cellIdToCellRefineFlagMapLocal,
const bool smoothenCellsOnPeriodicBoundary = false,
const double smootheningFactor = 2.0 )
private

internal function which sets refinement flags based on a custom created algorithm

Returns
bool boolean flag is any local cell has refinement flag set

◆ refineSerialMesh()

void dftfe::triangulationManager::refineSerialMesh ( const std::map< dealii::CellId, unsigned int > & cellIdToCellRefineFlagMapLocal,
const MPI_Comm & mpi_comm,
dealii::parallel::distributed::Triangulation< 3 > & serialTriangulation,
const dealii::parallel::distributed::Triangulation< 3 > & parallelTriangulation,
std::vector< bool > & serialTriaCurrentRefinement )
private

internal function which refines the serial mesh based on refinement flags from parallel mesh. This ensures that we get the same mesh in serial and parallel.

◆ resetMesh()

void dftfe::triangulationManager::resetMesh ( dealii::parallel::distributed::Triangulation< 3 > & parallelTriangulationA,
dealii::parallel::distributed::Triangulation< 3 > & parallelTriangulationB )

resets the vertices of meshB moved to vertices of meshA.

◆ saveSupportTriangulations()

void dftfe::triangulationManager::saveSupportTriangulations ( std::string path)
private

internal function to serialize support triangulations. No solution data is attached to them

◆ saveTriangulationsSolutionVectors()

void dftfe::triangulationManager::saveTriangulationsSolutionVectors ( std::string path,
const unsigned int feOrder,
const unsigned int nComponents,
const std::vector< const distributedCPUVec< double > * > & solutionVectors,
const MPI_Comm & interpoolComm,
const MPI_Comm & interBandGroupComm )

serialize the triangulations and the associated solution vectors

Parameters
[input]feOrderfinite element polynomial order of the dofHandler on which solution vectors are based upon
[input]nComponentsnumber of components of the dofHandler on which solution vectors are based upon
[input]solutionVectorsvector of parallel distributed solution vectors to be serialized
[input]interpoolCommThis communicator is used to ensure serialization happens only in k point pool
[input]interBandGroupCommThis communicator to ensure serialization happens only in band group

Member Data Documentation

◆ computing_timer

dealii::TimerOutput dftfe::triangulationManager::computing_timer
private

◆ d_atomPositions

std::vector<std::vector<double> > dftfe::triangulationManager::d_atomPositions
private

◆ d_dftParams

const dftParameters& dftfe::triangulationManager::d_dftParams
private

◆ d_domainBoundingVectors

std::vector<std::vector<double> > dftfe::triangulationManager::d_domainBoundingVectors
private

◆ d_FEOrder

const unsigned int dftfe::triangulationManager::d_FEOrder
private

FEOrder to be used for checking parallel consistency of periodic+hanging node constraints

◆ d_imageAtomPositions

std::vector<std::vector<double> > dftfe::triangulationManager::d_imageAtomPositions
private

◆ d_imageIds

std::vector<int> dftfe::triangulationManager::d_imageIds
private

◆ d_max_refinement_steps

const unsigned int dftfe::triangulationManager::d_max_refinement_steps = 40
private

◆ d_meshSizes

std::vector<std::vector<double> > dftfe::triangulationManager::d_meshSizes
private

◆ d_mpiCommParent

const MPI_Comm dftfe::triangulationManager::d_mpiCommParent
private

◆ d_nearestAtomDistances

std::vector<double> dftfe::triangulationManager::d_nearestAtomDistances
private

◆ d_parallelTriaCurrentRefinement

std::vector<std::vector<bool> > dftfe::triangulationManager::d_parallelTriaCurrentRefinement
private

◆ d_parallelTriangulationMoved

dealii::parallel::distributed::Triangulation<3> dftfe::triangulationManager::d_parallelTriangulationMoved
private

◆ d_parallelTriangulationUnmoved

dealii::parallel::distributed::Triangulation<3> dftfe::triangulationManager::d_parallelTriangulationUnmoved
private

◆ d_serialTriaCurrentRefinement

std::vector<std::vector<bool> > dftfe::triangulationManager::d_serialTriaCurrentRefinement
private

◆ d_serialTriangulationUnmoved

dealii::parallel::distributed::Triangulation<3> dftfe::triangulationManager::d_serialTriangulationUnmoved
private

◆ interBandGroupComm

const MPI_Comm dftfe::triangulationManager::interBandGroupComm
private

◆ interpoolcomm

const MPI_Comm dftfe::triangulationManager::interpoolcomm
private

◆ mpi_communicator

const MPI_Comm dftfe::triangulationManager::mpi_communicator
private

◆ n_mpi_processes

const unsigned int dftfe::triangulationManager::n_mpi_processes
private

◆ pcout

dealii::ConditionalOStream dftfe::triangulationManager::pcout
private

◆ this_mpi_process

const unsigned int dftfe::triangulationManager::this_mpi_process
private

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