Class to update triangulation under affine transformation.
More...
#include <meshMovementAffineTransform.h>
|
| meshMovementAffineTransform (const MPI_Comm &mpi_comm_parent, const MPI_Comm &mpi_comm_domain, const dftParameters &dftParams) |
| Constructor.
|
|
std::pair< bool, double > | transform (const dealii::Tensor< 2, 3, double > &deformationGradient) |
| Performs affine transformation of the triangulation.
|
|
std::pair< bool, double > | moveMesh (const std::vector< dealii::Point< 3 > > &controlPointLocations, const std::vector< dealii::Tensor< 1, 3, double > > &controlPointDisplacements, const double controllingParameter, const bool moveSubdivided=false) |
|
| 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 () |
| internal function which computes the nodal increment field in the local processor
|
|
|
dealii::Tensor< 2, 3, double > | d_deformationGradient |
| storage for the deformation gradient to be applied to the triangulation
|
|
Class to update triangulation under affine transformation.
- Author
- Sambit Das
◆ meshMovementAffineTransform()
dftfe::meshMovementAffineTransform::meshMovementAffineTransform |
( |
const MPI_Comm & | mpi_comm_parent, |
|
|
const MPI_Comm & | mpi_comm_domain, |
|
|
const dftParameters & | dftParams ) |
Constructor.
- Parameters
-
mpi_comm_parent | parent mpi communicator |
mpi_comm_domain | mpi communicator domain decomposition |
◆ computeIncrement()
void dftfe::meshMovementAffineTransform::computeIncrement |
( |
| ) |
|
|
private |
internal function which computes the nodal increment field in the local processor
◆ moveMesh()
std::pair< bool, double > dftfe::meshMovementAffineTransform::moveMesh |
( |
const std::vector< dealii::Point< 3 > > & | controlPointLocations, |
|
|
const std::vector< dealii::Tensor< 1, 3, double > > & | controlPointDisplacements, |
|
|
const double | controllingParameter, |
|
|
const bool | moveSubdivided = false ) |
Not implemented, just present to override the pure virtual from base class
◆ transform()
std::pair< bool, double > dftfe::meshMovementAffineTransform::transform |
( |
const dealii::Tensor< 2, 3, double > & | deformationGradient | ) |
|
Performs affine transformation of the triangulation.
- Parameters
-
- Returns
- std::pair<bool,double> mesh quality metrics pair(bool for is negative jacobian, maximum jacobian ratio)
◆ d_deformationGradient
dealii::Tensor<2, 3, double> dftfe::meshMovementAffineTransform::d_deformationGradient |
|
private |
storage for the deformation gradient to be applied to the triangulation
The documentation for this class was generated from the following file: