[P]arallel [Hi]gh-order [Li]brary for [P]DEs
Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
|
Free form deformation class from Sederberg 1986. More...
#include <free_form_deformation.h>
Public Member Functions | |
FreeFormDeformation (const dealii::Point< dim > &_origin, const std::array< dealii::Tensor< 1, dim, double >, dim > _parallepiped_vectors, const std::array< unsigned int, dim > &_ndim_control) | |
Constructor for an oblique parallepiped. | |
FreeFormDeformation (const dealii::Point< dim > &_origin, const std::array< double, dim > &rectangle_lengths, const std::array< unsigned int, dim > &_ndim_control) | |
Constructor for a rectangular FFD box. | |
template<typename real > | |
dealii::Point< dim, real > | new_point_location (const dealii::Point< dim, double > &initial_point, const std::vector< dealii::Point< dim, real >> &control_pts) const |
dealii::Point< dim, double > | new_point_location (const dealii::Point< dim, double > &initial_point) const |
dealii::LinearAlgebra::distributed::Vector< double > | get_surface_displacement (const HighOrderGrid< dim, double > &high_order_grid) const |
void | deform_mesh (HighOrderGrid< dim, double > &high_order_grid) const |
Deform HighOrderGrid using its initial volume_nodes to retrieve the deformed set of volume_nodes. | |
dealii::Point< dim, double > | dXdXp (const dealii::Point< dim, double > &initial_point, const unsigned int ctl_index, const unsigned int ctl_axis) const |
void | get_dXvsdXp (const HighOrderGrid< dim, double > &high_order_grid, const std::vector< std::pair< unsigned int, unsigned int > > &ffd_design_variables_indices_dim, dealii::TrilinosWrappers::SparseMatrix &dXvsdXp) const |
std::vector< dealii::LinearAlgebra::distributed::Vector< double > > | get_dXvsdXp (const HighOrderGrid< dim, double > &high_order_grid, const std::vector< std::pair< unsigned int, unsigned int > > &ffd_design_variables_indices_dim) const |
std::vector< dealii::LinearAlgebra::distributed::Vector< double > > | get_dXvsdXp_FD (const HighOrderGrid< dim, double > &high_order_grid, const std::vector< std::pair< unsigned int, unsigned int > > &ffd_design_variables_indices_dim, const double eps) |
void | get_dXvdXp (const HighOrderGrid< dim, double > &high_order_grid, const std::vector< std::pair< unsigned int, unsigned int > > ffd_design_variables_indices_dim, dealii::TrilinosWrappers::SparseMatrix &dXvdXp) const |
void | get_dXvdXp_FD (HighOrderGrid< dim, double > &high_order_grid, const std::vector< std::pair< unsigned int, unsigned int > > ffd_design_variables_indices_dim, dealii::TrilinosWrappers::SparseMatrix &dXvdXp_FD, const double eps) |
template<typename real > | |
dealii::Tensor< 1, dim, real > | get_displacement (const dealii::Point< dim, double > &initial_point, const std::vector< dealii::Point< dim, real >> &control_pts) const |
template<typename real > | |
std::vector< dealii::Tensor< 1, dim, real > > | get_displacement (const std::vector< dealii::Point< dim, double >> &initial_point, const std::vector< dealii::Point< dim, real >> &control_pts) const |
template<typename real > | |
dealii::Point< dim, real > | evaluate_ffd (const dealii::Point< dim, double > &s_t_u_point, const std::vector< dealii::Point< dim, real >> &control_pts) const |
std::array< unsigned int, dim > | global_to_grid (const unsigned int global_ictl) const |
Given a control points' global index return its ijk coordinate. More... | |
unsigned int | grid_to_global (const std::array< unsigned int, dim > &ijk_index) const |
Given a control points' ijk coordinate return its global indexing. More... | |
void | move_ctl_dx (const unsigned i, const dealii::Tensor< 1, dim, double >) |
Move control point with global index i. | |
void | move_ctl_dx (const std::array< unsigned int, dim > ijk, const dealii::Tensor< 1, dim, double >) |
Move control point with grid index (i,j,k). | |
void | get_design_variables (const std::vector< std::pair< unsigned int, unsigned int > > &ffd_design_variables_indices_dim, dealii::LinearAlgebra::distributed::Vector< double > &vector_to_copy_into) const |
Copies the desired control points from FreeFormDeformation object into vector_to_copy_into. | |
void | set_design_variables (const std::vector< std::pair< unsigned int, unsigned int > > &ffd_design_variables_indices_dim, const dealii::LinearAlgebra::distributed::Vector< double > &vector_to_copy_from) |
Copies the desired control points from vector_to_copy_from into FreeFormDeformation object. | |
void | output_ffd_vtu (const unsigned int cycle) const |
Output a .vtu file of the FFD box to visualize. | |
Public Attributes | |
std::vector< dealii::Point< dim > > | control_pts |
Control points of the FFD box used to deform the geometry. | |
const unsigned int | n_control_pts |
Total number of control points. | |
Protected Member Functions | |
dealii::Point< dim, double > | get_local_coordinates (const dealii::Point< dim, double > p) const |
Returns the local coordinates s-t-u within the FFD box. More... | |
Protected Attributes | |
const dealii::Point< dim > | origin |
Parallepiped origin. | |
const std::array< dealii::Tensor< 1, dim, double >, dim > | parallepiped_vectors |
Parallepiped vectors. More... | |
const std::array< unsigned int, dim > | ndim_control_pts |
Number of control points in each direction. | |
Private Member Functions | |
std::array< dealii::Tensor< 1, dim, double >, dim > | get_rectangular_parallepiped_vectors (const std::array< double, dim > &rectangle_lengths) const |
Returns rectangular vector box given the box lengths. | |
unsigned int | compute_total_ctl_pts () const |
Used by constructor to evaluate total number of control points. | |
void | init_msg () const |
Initial message. | |
Private Attributes | |
dealii::ConditionalOStream | pcout |
Outputs if MPI rank is 0. | |
Free form deformation class from Sederberg 1986.
Definition at line 10 of file free_form_deformation.h.
dealii::Point< dim, double > PHiLiP::FreeFormDeformation< dim >::dXdXp | ( | const dealii::Point< dim, double > & | initial_point, |
const unsigned int | ctl_index, | ||
const unsigned int | ctl_axis | ||
) | const |
Given an initial point in the undeformed initial parallepiped and the index a control point, return the derivative dXdXp of the new point location point_i with respect to that control_point_j.
Definition at line 240 of file free_form_deformation.cpp.
dealii::Point< dim, real > PHiLiP::FreeFormDeformation< dim >::evaluate_ffd | ( | const dealii::Point< dim, double > & | s_t_u_point, |
const std::vector< dealii::Point< dim, real >> & | control_pts | ||
) | const |
Given the s,t,u reference location within the FFD box, return its position in the actual domain.
Definition at line 263 of file free_form_deformation.cpp.
dealii::Tensor< 1, dim, real > PHiLiP::FreeFormDeformation< dim >::get_displacement | ( | const dealii::Point< dim, double > & | initial_point, |
const std::vector< dealii::Point< dim, real >> & | control_pts | ||
) | const |
Given an initial point in the undeformed initial parallepiped, return the its displacement due to the free-form deformation.
Definition at line 214 of file free_form_deformation.cpp.
std::vector< dealii::Tensor< 1, dim, real > > PHiLiP::FreeFormDeformation< dim >::get_displacement | ( | const std::vector< dealii::Point< dim, double >> & | initial_point, |
const std::vector< dealii::Point< dim, real >> & | control_pts | ||
) | const |
Given a vector of initial points in the undeformed initial parallepiped, return the its corresponding displacements due to the free-form deformation.
Definition at line 227 of file free_form_deformation.cpp.
void PHiLiP::FreeFormDeformation< dim >::get_dXvdXp | ( | const HighOrderGrid< dim, double > & | high_order_grid, |
const std::vector< std::pair< unsigned int, unsigned int > > | ffd_design_variables_indices_dim, | ||
dealii::TrilinosWrappers::SparseMatrix & | dXvdXp | ||
) | const |
For the given list of FFD indices and direction, return the finite-differenced derivatives of the HighOrderGrid's initial volume points with respect to the FFD.
Definition at line 554 of file free_form_deformation.cpp.
void PHiLiP::FreeFormDeformation< dim >::get_dXvdXp_FD | ( | HighOrderGrid< dim, double > & | high_order_grid, |
const std::vector< std::pair< unsigned int, unsigned int > > | ffd_design_variables_indices_dim, | ||
dealii::TrilinosWrappers::SparseMatrix & | dXvdXp_FD, | ||
const double | eps | ||
) |
For the given list of FFD indices and direction, return the analytical derivatives of the HighOrderGrid's initial volume points with respect to the FFD.
Definition at line 579 of file free_form_deformation.cpp.
void PHiLiP::FreeFormDeformation< dim >::get_dXvsdXp | ( | const HighOrderGrid< dim, double > & | high_order_grid, |
const std::vector< std::pair< unsigned int, unsigned int > > & | ffd_design_variables_indices_dim, | ||
dealii::TrilinosWrappers::SparseMatrix & | dXvsdXp | ||
) | const |
For the given list of FFD indices and direction, return the analytical derivatives of the HighOrderGrid's initial surface points with respect to the FFD. The result is written into the given dXvsdXp SparseMatrix.
Definition at line 496 of file free_form_deformation.cpp.
std::vector< dealii::LinearAlgebra::distributed::Vector< double > > PHiLiP::FreeFormDeformation< dim >::get_dXvsdXp | ( | const HighOrderGrid< dim, double > & | high_order_grid, |
const std::vector< std::pair< unsigned int, unsigned int > > & | ffd_design_variables_indices_dim | ||
) | const |
For the given list of FFD indices and direction, return the analytical derivatives of the HighOrderGrid's initial surface points with respect to the FFD. Note that the result is returned as a vector of vector because it will likely be used with a MeshMover's apply_dXvdXvs() function.
Definition at line 454 of file free_form_deformation.cpp.
std::vector< dealii::LinearAlgebra::distributed::Vector< double > > PHiLiP::FreeFormDeformation< dim >::get_dXvsdXp_FD | ( | const HighOrderGrid< dim, double > & | high_order_grid, |
const std::vector< std::pair< unsigned int, unsigned int > > & | ffd_design_variables_indices_dim, | ||
const double | eps | ||
) |
For the given list of FFD indices and direction, return the finite-differenced derivatives of the HighOrderGrid's initial surface points with respect to the FFD. Note that the result is returned as a vector of vector because it will likely be used with a MeshMover's apply_dXvdXvs() function.
Definition at line 389 of file free_form_deformation.cpp.
|
protected |
Returns the local coordinates s-t-u within the FFD box.
s,t,u should be in [0,1] for a point inside the box.
Definition at line 57 of file free_form_deformation.cpp.
dealii::LinearAlgebra::distributed::Vector< double > PHiLiP::FreeFormDeformation< dim >::get_surface_displacement | ( | const HighOrderGrid< dim, double > & | high_order_grid | ) | const |
Using the initial surface nodes from the given HighOrderGrid, return the surface displacements based on the free-form deformation displacements.
Definition at line 360 of file free_form_deformation.cpp.
std::array< unsigned int, dim > PHiLiP::FreeFormDeformation< dim >::global_to_grid | ( | const unsigned int | global_ictl | ) | const |
Given a control points' global index return its ijk coordinate.
Opposite of grid_to_global
Definition at line 106 of file free_form_deformation.cpp.
unsigned int PHiLiP::FreeFormDeformation< dim >::grid_to_global | ( | const std::array< unsigned int, dim > & | ijk_index | ) | const |
Given a control points' ijk coordinate return its global indexing.
Opposite of global_to_grid
Definition at line 121 of file free_form_deformation.cpp.
dealii::Point< dim, real > PHiLiP::FreeFormDeformation< dim >::new_point_location | ( | const dealii::Point< dim, double > & | initial_point, |
const std::vector< dealii::Point< dim, real >> & | control_pts | ||
) | const |
Given an initial point in the undeformed initial parallepiped, return the position of the new point location.
Definition at line 194 of file free_form_deformation.cpp.
dealii::Point< dim, double > PHiLiP::FreeFormDeformation< dim >::new_point_location | ( | const dealii::Point< dim, double > & | initial_point | ) | const |
Given an initial point in the undeformed initial parallepiped, return the position of the new point location using the current control point locations.
Definition at line 185 of file free_form_deformation.cpp.
|
protected |
Parallepiped vectors.
Not unit vector since the parallepiped length will be determined by those vectors' magnitude.
Definition at line 160 of file free_form_deformation.h.