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

Class to move triangulation nodes using Gaussian functions attached to control points. More...

#include <meshMovementGaussian.h>

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

Public Member Functions

 meshMovementGaussianClass (const MPI_Comm &mpi_comm_parent, const MPI_Comm &mpi_comm_domaim, const dftParameters &dftParams)
 Constructor.
 
std::pair< bool, double > moveMesh (const std::vector< dealii::Point< 3 > > &controlPointLocations, const std::vector< dealii::Tensor< 1, 3, double > > &controlPointDisplacements, const std::vector< double > &gaussianWidthParameter, const std::vector< double > &flatTopWidthParameter, const bool moveSubdivided=false)
 Moves the triangulation corresponding to Gaussians attached to control points.
 
std::pair< bool, double > moveMeshTwoStep (const std::vector< dealii::Point< 3 > > &controlPointLocations1, const std::vector< dealii::Point< 3 > > &controlPointLocations2, const std::vector< dealii::Tensor< 1, 3, double > > &controlPointDisplacements1, const std::vector< dealii::Tensor< 1, 3, double > > &controlPointDisplacements2, const std::vector< double > &controllingParameter1, const std::vector< double > &controllingParameter2, const std::vector< double > &flatTopWidthParameter, const bool moveSubdivided=false)
 
void moveMeshTwoLevelElectro ()
 
- Public Member Functions inherited from dftfe::meshMovementClass
 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.
 

Private Member Functions

void computeIncrement (const std::vector< dealii::Point< 3 > > &controlPointLocations, const std::vector< dealii::Tensor< 1, 3, double > > &controlPointDisplacements, const std::vector< double > &gaussianWidthParameter, const std::vector< double > &flatTopWidthParameter)
 internal function which computes the nodal increment field in the local processor
 
void computeIncrementTwoStep (const std::vector< dealii::Point< 3 > > &controlPointLocations1, const std::vector< dealii::Point< 3 > > &controlPointLocations2, const std::vector< dealii::Tensor< 1, 3, double > > &controlPointDisplacements1, const std::vector< dealii::Tensor< 1, 3, double > > &controlPointDisplacements2, const std::vector< double > &gaussianWidthParameter1, const std::vector< double > &gaussianWidthParameter2, const std::vector< double > &flatTopWidthParameter)
 

Additional Inherited Members

- Protected Member Functions inherited from dftfe::meshMovementClass
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 inherited from dftfe::meshMovementClass
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

Class to move triangulation nodes using Gaussian functions attached to control points.

Author
Sambit Das

Constructor & Destructor Documentation

◆ meshMovementGaussianClass()

dftfe::meshMovementGaussianClass::meshMovementGaussianClass ( const MPI_Comm & mpi_comm_parent,
const MPI_Comm & mpi_comm_domaim,
const dftParameters & dftParams )

Constructor.

Parameters
mpi_comm_parentparent mpi communicator
mpi_comm_domainmpi communicator for domain decomposition

Member Function Documentation

◆ computeIncrement()

void dftfe::meshMovementGaussianClass::computeIncrement ( const std::vector< dealii::Point< 3 > > & controlPointLocations,
const std::vector< dealii::Tensor< 1, 3, double > > & controlPointDisplacements,
const std::vector< double > & gaussianWidthParameter,
const std::vector< double > & flatTopWidthParameter )
private

internal function which computes the nodal increment field in the local processor

◆ computeIncrementTwoStep()

void dftfe::meshMovementGaussianClass::computeIncrementTwoStep ( const std::vector< dealii::Point< 3 > > & controlPointLocations1,
const std::vector< dealii::Point< 3 > > & controlPointLocations2,
const std::vector< dealii::Tensor< 1, 3, double > > & controlPointDisplacements1,
const std::vector< dealii::Tensor< 1, 3, double > > & controlPointDisplacements2,
const std::vector< double > & gaussianWidthParameter1,
const std::vector< double > & gaussianWidthParameter2,
const std::vector< double > & flatTopWidthParameter )
private

◆ moveMesh()

std::pair< bool, double > dftfe::meshMovementGaussianClass::moveMesh ( const std::vector< dealii::Point< 3 > > & controlPointLocations,
const std::vector< dealii::Tensor< 1, 3, double > > & controlPointDisplacements,
const std::vector< double > & gaussianWidthParameter,
const std::vector< double > & flatTopWidthParameter,
const bool moveSubdivided = false )

Moves the triangulation corresponding to Gaussians attached to control points.

This functions takes into account the hanging node and periodic constraints when computing the nodal increment field.

Parameters
controlPointLocationsvector of coordinates of control points
controlPointDisplacementsvector of displacements of control points @ controllingParameter constant in the Gaussian function: exp(-(r/controllingParameter)^pow) where pow is controlled via the input file parameter (.prm)
Returns
std::pair<bool,double> mesh quality metrics pair(bool for is negative jacobian, maximum jacobian ratio)

◆ moveMeshTwoLevelElectro()

void dftfe::meshMovementGaussianClass::moveMeshTwoLevelElectro ( )

◆ moveMeshTwoStep()

std::pair< bool, double > dftfe::meshMovementGaussianClass::moveMeshTwoStep ( const std::vector< dealii::Point< 3 > > & controlPointLocations1,
const std::vector< dealii::Point< 3 > > & controlPointLocations2,
const std::vector< dealii::Tensor< 1, 3, double > > & controlPointDisplacements1,
const std::vector< dealii::Tensor< 1, 3, double > > & controlPointDisplacements2,
const std::vector< double > & controllingParameter1,
const std::vector< double > & controllingParameter2,
const std::vector< double > & flatTopWidthParameter,
const bool moveSubdivided = false )

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