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

Base class to move triangulation vertices. More...

#include <meshMovement.h>

Inheritance diagram for dftfe::meshMovementClass:
dftfe::meshMovementAffineTransform dftfe::meshMovementGaussianClass

Public Member Functions

 meshMovementClass (const MPI_Comm &mpi_comm_parent, const MPI_Comm &mpi_comm_domain, const dftParameters &dftParams)
 Constructor.
 
virtual ~meshMovementClass ()
 
void init (dealii::Triangulation< 3, 3 > &triangulation, dealii::Triangulation< 3, 3 > &serialTriangulation, const std::vector< std::vector< double > > &domainBoundingVectors)
 Initializes the required data-structures for a given triangulation.
 
void initMoved (const std::vector< std::vector< double > > &domainBoundingVectors)
 Re-initializes the required data-structures for a given triangulation.
 
void findClosestVerticesToDestinationPoints (const std::vector< dealii::Point< 3 > > &destinationPoints, std::vector< dealii::Point< 3 > > &closestTriaVertexToDestPointsLocation, std::vector< dealii::Tensor< 1, 3, double > > &dispClosestTriaVerticesToDestPoints)
 Finds the closest triangulation vertices to a given vector of position coordinates.
 

Protected Member Functions

void initIncrementField ()
 Initializes the parallel layout of d_incrementalDisplacementParallel.
 
void finalizeIncrementField ()
 
void updateTriangulationVertices ()
 Function which updates the locally relevant triangulation vertices.
 
void moveSubdividedMesh ()
 Function which moves subdivided mesh.
 
std::pair< bool, double > movedMeshCheck ()
 

Protected Attributes

distributedCPUVec< double > d_incrementalDisplacement
 vector of displacements of the triangulation vertices
 
bool d_isParallelMesh
 
dealii::FESystem< 3 > FEMoveMesh
 
dealii::DoFHandler< 3 > d_dofHandlerMoveMesh
 
dealii::parallel::distributed::Triangulation< 3 > * d_triaPtr
 
dealii::Triangulation< 3, 3 > * d_triaPtrSerial
 
dealii::IndexSet d_locally_owned_dofs
 
dealii::IndexSet d_locally_relevant_dofs
 
dealii::AffineConstraints< double > d_constraintsMoveMesh
 
std::vector< dealii::GridTools::PeriodicFacePair< typename dealii::DoFHandler< 3 >::cell_iterator > > d_periodicity_vector
 
std::vector< std::vector< double > > d_domainBoundingVectors
 
const dftParametersd_dftParams
 
MPI_Comm d_mpiCommParent
 
MPI_Comm mpi_communicator
 
const unsigned int this_mpi_process
 
dealii::ConditionalOStream pcout
 

Detailed Description

Base class to move triangulation vertices.

Author
Sambit Das

Constructor & Destructor Documentation

◆ meshMovementClass()

dftfe::meshMovementClass::meshMovementClass ( const MPI_Comm & mpi_comm_parent,
const MPI_Comm & mpi_comm_domain,
const dftParameters & dftParams )

Constructor.

Parameters
[in]mpi_comm_parentparent mpi communicator
[in]mpi_comm_domainmpi communicator for domain decomposition

◆ ~meshMovementClass()

virtual dftfe::meshMovementClass::~meshMovementClass ( )
inlinevirtual

Member Function Documentation

◆ finalizeIncrementField()

void dftfe::meshMovementClass::finalizeIncrementField ( )
protected

Takes care of communicating the movement of triangulation vertices on processor boundaries, and also takes care of hanging nodes and periodic constraints

◆ findClosestVerticesToDestinationPoints()

void dftfe::meshMovementClass::findClosestVerticesToDestinationPoints ( const std::vector< dealii::Point< 3 > > & destinationPoints,
std::vector< dealii::Point< 3 > > & closestTriaVertexToDestPointsLocation,
std::vector< dealii::Tensor< 1, 3, double > > & dispClosestTriaVerticesToDestPoints )

Finds the closest triangulation vertices to a given vector of position coordinates.

Parameters
[in]destinationPointsvector of points in cartesian coordinates (origin at center of the domain) to which closest triangulation vertices are desired.
[out]closestTriaVertexToDestPointsLocationvector of positions of the closest triangulation v vertices.
[out]dispClosestTriaVerticesToDestPointsvector of displacements of the destinationPoints from the closest triangulation vertices.

◆ init()

void dftfe::meshMovementClass::init ( dealii::Triangulation< 3, 3 > & triangulation,
dealii::Triangulation< 3, 3 > & serialTriangulation,
const std::vector< std::vector< double > > & domainBoundingVectors )

Initializes the required data-structures for a given triangulation.

Parameters
[in]triangulationtriangulation object whose nodes are to be moved
[in]serialtriangulation to create constraints from serial dofHandler (temporary fix)
[in]domainBoundingVectorsdomain vectors of the domain corresponding to the triangulation object.

◆ initIncrementField()

void dftfe::meshMovementClass::initIncrementField ( )
protected

Initializes the parallel layout of d_incrementalDisplacementParallel.

◆ initMoved()

void dftfe::meshMovementClass::initMoved ( const std::vector< std::vector< double > > & domainBoundingVectors)

Re-initializes the required data-structures for a given triangulation.

Parameters
[in]domainBoundingVectorscurrent domain vectors of the domain corresponding to the triangulation object.

◆ movedMeshCheck()

std::pair< bool, double > dftfe::meshMovementClass::movedMeshCheck ( )
protected

Performs periodic matching sanity check and returns the pair<if negative jacobian, maximum inverse jacobian magnitude>

◆ moveSubdividedMesh()

void dftfe::meshMovementClass::moveSubdividedMesh ( )
protected

Function which moves subdivided mesh.

◆ updateTriangulationVertices()

void dftfe::meshMovementClass::updateTriangulationVertices ( )
protected

Function which updates the locally relevant triangulation vertices.

Member Data Documentation

◆ d_constraintsMoveMesh

dealii::AffineConstraints<double> dftfe::meshMovementClass::d_constraintsMoveMesh
protected

◆ d_dftParams

const dftParameters& dftfe::meshMovementClass::d_dftParams
protected

◆ d_dofHandlerMoveMesh

dealii::DoFHandler<3> dftfe::meshMovementClass::d_dofHandlerMoveMesh
protected

◆ d_domainBoundingVectors

std::vector<std::vector<double> > dftfe::meshMovementClass::d_domainBoundingVectors
protected

◆ d_incrementalDisplacement

distributedCPUVec<double> dftfe::meshMovementClass::d_incrementalDisplacement
protected

vector of displacements of the triangulation vertices

◆ d_isParallelMesh

bool dftfe::meshMovementClass::d_isParallelMesh
protected

◆ d_locally_owned_dofs

dealii::IndexSet dftfe::meshMovementClass::d_locally_owned_dofs
protected

◆ d_locally_relevant_dofs

dealii::IndexSet dftfe::meshMovementClass::d_locally_relevant_dofs
protected

◆ d_mpiCommParent

MPI_Comm dftfe::meshMovementClass::d_mpiCommParent
protected

◆ d_periodicity_vector

std::vector<dealii::GridTools::PeriodicFacePair< typename dealii::DoFHandler<3>::cell_iterator> > dftfe::meshMovementClass::d_periodicity_vector
protected

◆ d_triaPtr

dealii::parallel::distributed::Triangulation<3>* dftfe::meshMovementClass::d_triaPtr
protected

◆ d_triaPtrSerial

dealii::Triangulation<3, 3>* dftfe::meshMovementClass::d_triaPtrSerial
protected

◆ FEMoveMesh

dealii::FESystem<3> dftfe::meshMovementClass::FEMoveMesh
protected

◆ mpi_communicator

MPI_Comm dftfe::meshMovementClass::mpi_communicator
protected

◆ pcout

dealii::ConditionalOStream dftfe::meshMovementClass::pcout
protected

◆ this_mpi_process

const unsigned int dftfe::meshMovementClass::this_mpi_process
protected

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