Class to move triangulation nodes using Gaussian functions attached to control points.
More...
#include <meshMovementGaussian.h>
|
| 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 () |
|
| 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.
|
|
|
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) |
|
Class to move triangulation nodes using Gaussian functions attached to control points.
- Author
- Sambit Das
◆ meshMovementGaussianClass()
dftfe::meshMovementGaussianClass::meshMovementGaussianClass |
( |
const MPI_Comm & | mpi_comm_parent, |
|
|
const MPI_Comm & | mpi_comm_domaim, |
|
|
const dftParameters & | dftParams ) |
Constructor.
- Parameters
-
mpi_comm_parent | parent mpi communicator |
mpi_comm_domain | mpi communicator for domain decomposition |
◆ 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
-
controlPointLocations | vector of coordinates of control points |
controlPointDisplacements | vector 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: