[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType > Class Template Reference

#include <high_order_grid.h>

Collaboration diagram for PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >:

Public Member Functions

 HighOrderGrid (const unsigned int max_degree, const std::shared_ptr< MeshType > triangulation_input, const bool check_valid_metric_Jacobian_input=true, const bool renumber_dof_handler_Cuthill_Mckee_input=true, const bool output_high_order_grid=true)
 Principal constructor that will call delegated constructor.
 
void reinit ()
 Reinitialize high_order_grid after a change in triangulation.
 
void update_mapping_fe_field ()
 Update the MappingFEField. More...
 
void ensure_conforming_mesh ()
 Ensures that hanging nodes are updated for a conforming mesh.
 
void initialize_with_triangulation_manifold (const bool output_mesh=true)
 Sets the volume_nodes to the interpolated position of the Manifold associated to the triangulation.
 
void allocate ()
 Needed to allocate the correct number of volume_nodes when initializing and after the mesh is refined.
 
dealii::MappingFEField< dim, dim, VectorType, DoFHandlerType > get_MappingFEField ()
 Return a MappingFEField that corresponds to the current node locations.
 
void reset_initial_nodes ()
 Sets the initial_volume_nodes and initial_surface_nodes to the current volume_nodes and surface_nodes.
 
void update_map_nodes_surf_to_vol ()
 Sets the initial_volume_nodes and initial_surface_nodes to the current volume_nodes and surface_nodes.
 
void update_surface_nodes ()
 Update list of surface nodes (all_locally_relevant_surface_nodes). More...
 
VectorType transform_surface_nodes (std::function< dealii::Point< dim >(dealii::Point< dim >)> transformation) const
 
void deform_mesh (std::vector< real > local_surface_displacements)
 RBF mesh deformation - To be done.
 
void test_jacobian ()
 Test metric Jacobian.
 
std::vector< real > evaluate_jacobian_at_points (const VectorType &solution, const typename DoFHandlerType::cell_iterator &cell, const std::vector< dealii::Point< dim >> &points) const
 Evaluates Jacobian of a cell given some points using a global solution vector.
 
template<typename real2 >
void evaluate_jacobian_at_points (const std::vector< real2 > &dofs, const dealii::FESystem< dim > &fe, const std::vector< dealii::Point< dim >> &points, std::vector< real2 > &jacobian_determinants) const
 Evaluates Jacobian given some DoF, associated FE, and some points. More...
 
template<typename real2 >
real2 evaluate_jacobian_at_point (const std::vector< real2 > &dofs, const dealii::FESystem< dim > &fe, const dealii::Point< dim > &point) const
 Evaluates Jacobian given some DoF, associated FE, and some point. More...
 
bool check_valid_cell (const typename DoFHandlerType::cell_iterator &cell) const
 Evaluate exact Jacobian determinant polynomial and uses Bernstein polynomials to determine positivity.
 
bool fix_invalid_cell (const typename DoFHandlerType::cell_iterator &cell)
 Evaluate exact Jacobian determinant polynomial and uses Bernstein polynomials to determine positivity.
 
void evaluate_lagrange_to_bernstein_operator (const unsigned int order)
 Evaluates the operator to obtain Bernstein coefficients from a set of Lagrange coefficients. More...
 
void output_results_vtk (const unsigned int cycle) const
 Output mesh with metric informations.
 
void refine_global ()
 Globally refines the high-order mesh. More...
 
void prepare_for_coarsening_and_refinement ()
 Prepares the solution transfer such that the curved refined grid is on top of the curved coarse grid. More...
 
void execute_coarsening_and_refinement (const bool output_mesh=false)
 Executes the solution transfer such that the curved refined grid is on top of the curved coarse grid. More...
 

Public Attributes

const unsigned int max_degree
 Maximum degree of the geometry polynomial representing the grid.
 
const std::shared_ptr< MeshType > triangulation
 Mesh.
 
const bool check_valid_metric_Jacobian
 Flag to check validity of Jacobian.
 
const bool renumber_dof_handler_Cuthill_Mckee
 Flag to renumber dof_handler_grid with Cuthill Mckee.
 
dealii::DoFHandler< dim > dof_handler_grid
 Degrees of freedom handler for the high-order grid.
 
VectorType volume_nodes
 Current nodal coefficients of the high-order grid. More...
 
VectorType surface_nodes
 
VectorType initial_volume_nodes
 Initial nodal coefficients of the high-order grid. More...
 
VectorType initial_surface_nodes
 Distributed ghosted vector of initial surface nodes. More...
 
dealii::LinearAlgebra::distributed::Vector< int > surface_to_volume_indices
 
dealii::TrilinosWrappers::SparseMatrix map_nodes_surf_to_vol
 
dealii::IndexSet locally_owned_surface_nodes_indexset
 
dealii::IndexSet ghost_surface_nodes_indexset
 
std::vector< unsigned int > n_locally_owned_surface_nodes_per_mpi
 
std::vector< real > locally_relevant_surface_nodes
 List of surface nodes. More...
 
std::vector< dealii::types::global_dof_index > locally_relevant_surface_nodes_indices
 List of surface node indices.
 
std::vector< dealii::types::global_dof_index > locally_relevant_surface_nodes_boundary_id
 List of surface node boundary IDs, corresponding to locally_relevant_surface_nodes_indices.
 
std::vector< real > locally_owned_surface_nodes
 List of surface nodes. More...
 
std::vector< real > all_surface_nodes
 List of all surface nodes. More...
 
std::vector< dealii::types::global_dof_index > all_surface_indices
 
std::vector< dealii::types::global_dof_index > locally_owned_surface_nodes_indices
 List of surface node indices.
 
std::vector< dealii::Point< dim > > locally_relevant_surface_points
 Locally relevant surface points. More...
 
std::vector< dealii::Point< dim > > initial_locally_relevant_surface_points
 Initial locally relevant surface points. More...
 
std::map< dealii::types::global_dof_index, std::pair< unsigned int, unsigned int > > global_index_to_point_and_axis
 
std::map< std::pair< unsigned int, unsigned int >, dealii::types::global_dof_index > point_and_axis_to_global_index
 
dealii::FullMatrix< double > lagrange_to_bernstein_operator
 Used to transform coefficients from a Lagrange basis to a Bernstein basis.
 
const dealii::FE_Q< dim > fe_q
 Use Lagrange polynomial to represent the spatial location.
 
const dealii::FESystem< dim > fe_system
 Using system of polynomials to represent the x, y, and z directions.
 
const dealii::FE_DGQ< 1 > oneD_fe_q
 Use oneD Lagrange polynomial to represent the spatial location. More...
 
const dealii::FESystem< 1 > oneD_fe_system
 One-dimensional fe system for each direction.
 
const dealii::QGaussLobatto< 1 > oneD_grid_nodes
 One Dimensional grid nodes in reference space.
 
const dealii::QGaussLobatto< dim > dim_grid_nodes
 Dim-Dimensional grid nodes in reference space.
 
std::shared_ptr< dealii::MappingFEField< dim, dim, VectorType, DoFHandlerType > > mapping_fe_field
 MappingFEField that will provide the polynomial-based grid. More...
 
std::shared_ptr< dealii::MappingFEField< dim, dim, VectorType, DoFHandlerType > > initial_mapping_fe_field
 MappingFEField that will provide the polynomial-based grid for the initial volume_nodes. More...
 
dealii::IndexSet locally_owned_dofs_grid
 Locally own degrees of freedom for the grid.
 
dealii::IndexSet ghost_dofs_grid
 Locally relevant ghost degrees of freedom for the grid.
 
dealii::IndexSet locally_relevant_dofs_grid
 Union of locally owned degrees of freedom and relevant ghost degrees of freedom for the grid.
 

Static Public Attributes

static unsigned int nth_refinement =0
 Used to name the various files outputted.
 

Protected Member Functions

void update_surface_indices ()
 Update list of surface indices (locally_relevant_surface_nodes_indices and locally_relevant_surface_nodes_boundary_id)
 
template<typename real2 >
real2 determinant (const std::array< dealii::Tensor< 1, dim, real2 >, dim > jacobian) const
 Evaluate the determinant of a matrix given in the format of a std::array<dealii::Tensor<1,dim,real2>,dim>. More...
 
void get_position_vector (const DoFHandlerType &dh, VectorType &vector, const dealii::ComponentMask &mask)
 A stripped down copy of dealii::VectorTools::get_position_vector()
 
void get_projected_position_vector (const DoFHandlerType &dh, VectorType &vector, const dealii::ComponentMask &mask)
 

Protected Attributes

int n_mpi
 Number of MPI processes.
 
int mpi_rank
 
VectorType old_volume_nodes
 Used for the SolutionTransfer when performing grid adaptation.
 
SolutionTransfer solution_transfer
 
MPI_Comm mpi_communicator
 MPI communicator.
 
dealii::ConditionalOStream pcout
 Parallel std::cout that only outputs on mpi_rank==0.
 

Private Types

using SolutionTransfer = typename MeshTypeHelper< MeshType >::template SolutionTransfer< dim, VectorType, DoFHandlerType >
 

Detailed Description

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
class PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >

This HighOrderGrid class basically contains all the different part necessary to generate a dealii::MappingFEField that corresponds to the current Triangulation and attached Manifold. Once the high order grid is generated, the mesh can be deformed by assigning different values to the volume_nodes vector. The dof_handler_grid is used to access and loop through those volume_nodes. This will especially be useful when performing shape optimization, where the surface and volume volume_nodes will need to be displaced. Note that there are a lot of pre-processor statements, and that is because the SolutionTransfer class and the Vector class act quite differently between serial and parallel implementation. Hopefully, deal.II will change this one day such that we have one interface for both.

Definition at line 130 of file high_order_grid.h.

Member Function Documentation

◆ determinant()

template<int dim, typename real , typename MeshType , typename VectorType , typename DoFHandlerType >
template<typename real2 >
real2 PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::determinant ( const std::array< dealii::Tensor< 1, dim, real2 >, dim >  jacobian) const
inlineprotected

Evaluate the determinant of a matrix given in the format of a std::array<dealii::Tensor<1,dim,real2>,dim>.

The indices of the array represent the matrix rows, and the indices of the Tensor represents its columns.

Definition at line 439 of file high_order_grid.cpp.

◆ evaluate_jacobian_at_point()

template<int dim, typename real , typename MeshType , typename VectorType , typename DoFHandlerType >
template<typename real2 >
real2 PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::evaluate_jacobian_at_point ( const std::vector< real2 > &  dofs,
const dealii::FESystem< dim > &  fe,
const dealii::Point< dim > &  point 
) const

Evaluates Jacobian given some DoF, associated FE, and some point.

Can be used in conjunction with AD since the type is templated.

Definition at line 504 of file high_order_grid.cpp.

◆ evaluate_jacobian_at_points()

template<int dim, typename real , typename MeshType , typename VectorType , typename DoFHandlerType >
template<typename real2 >
void PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::evaluate_jacobian_at_points ( const std::vector< real2 > &  dofs,
const dealii::FESystem< dim > &  fe,
const std::vector< dealii::Point< dim >> &  points,
std::vector< real2 > &  jacobian_determinants 
) const

Evaluates Jacobian given some DoF, associated FE, and some points.

Calls evaluate_jacobian_at_point() on the vector of points. Can be used in conjunction with AD since the type is templated.

Definition at line 524 of file high_order_grid.cpp.

◆ evaluate_lagrange_to_bernstein_operator()

template<int dim, typename real , typename MeshType , typename VectorType , typename DoFHandlerType >
void PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::evaluate_lagrange_to_bernstein_operator ( const unsigned int  order)

Evaluates the operator to obtain Bernstein coefficients from a set of Lagrange coefficients.

This is used in the evaluation of the Jacobian positivity by checking the convex hull of the resulting Bezier curve.

Definition at line 555 of file high_order_grid.cpp.

◆ execute_coarsening_and_refinement()

template<int dim, typename real , typename MeshType , typename VectorType , typename DoFHandlerType >
void PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::execute_coarsening_and_refinement ( const bool  output_mesh = false)

Executes the solution transfer such that the curved refined grid is on top of the curved coarse grid.

This function needs to be after this->prepare_for_coarsening_and_refinement() and dealii::Triangulation::execute_coarsening_and_refinement() or dealii::Triangulation::refine_global().

Definition at line 1034 of file high_order_grid.cpp.

◆ get_projected_position_vector()

template<int dim, typename real , typename MeshType , typename VectorType , typename DoFHandlerType >
void PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::get_projected_position_vector ( const DoFHandlerType &  dh,
VectorType &  vector,
const dealii::ComponentMask &  mask 
)
protected

Another version of dealii::VectorTools::get_position_vector() that projects the points instead of interpolating them.

Definition at line 274 of file high_order_grid.cpp.

◆ prepare_for_coarsening_and_refinement()

template<int dim, typename real , typename MeshType , typename VectorType , typename DoFHandlerType >
void PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::prepare_for_coarsening_and_refinement ( )

Prepares the solution transfer such that the curved refined grid is on top of the curved coarse grid.

This function needs to be called before dealii::Triangulation::execute_coarsening_and_refinement() or dealii::Triangulation::refine_global() and this->execute_coarsening_and_refinement().

Definition at line 1025 of file high_order_grid.cpp.

◆ refine_global()

template<int dim, typename real , typename MeshType , typename VectorType , typename DoFHandlerType >
void PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::refine_global ( )

Globally refines the high-order mesh.

The metric Jacobian is given by the gradient of the physical location with respect to the reference locations

Definition at line 1016 of file high_order_grid.cpp.

◆ transform_surface_nodes()

template<int dim, typename real , typename MeshType , typename VectorType , typename DoFHandlerType >
VectorType PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::transform_surface_nodes ( std::function< dealii::Point< dim >(dealii::Point< dim >)>  transformation) const

Transforms the surface_nodes vector using a std::function tranformation.

Definition at line 1358 of file high_order_grid.cpp.

◆ update_mapping_fe_field()

template<int dim, typename real , typename MeshType , typename VectorType , typename DoFHandlerType >
void PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::update_mapping_fe_field ( )

Update the MappingFEField.

Note that this rarely needs to be called since MappingFEField stores a pointer to the DoFHandler and to the node Vector.

Definition at line 186 of file high_order_grid.cpp.

◆ update_surface_nodes()

template<int dim, typename real , typename MeshType , typename VectorType , typename DoFHandlerType >
void PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::update_surface_nodes ( )

Update list of surface nodes (all_locally_relevant_surface_nodes).

Given the Point index within the locally_relevant_surface_points and its component, this will return the index of the surface_nodes' index. This allows us to obtain a deformation vector matching the locally_relevant_surface_points, and then copy its results to the surface_nodes.

Definition at line 1088 of file high_order_grid.cpp.

Member Data Documentation

◆ all_surface_indices

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
std::vector<dealii::types::global_dof_index> PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::all_surface_indices

List of global surface node indices including those across different processors. Ordering corresponds to all_surface_nodes. TODO: Might want to make this a std::pair.

Definition at line 260 of file high_order_grid.h.

◆ all_surface_nodes

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
std::vector<real> PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::all_surface_nodes

List of all surface nodes.

Same as locally_owned_surface_nodes except that it stores a global vector of all the surface nodes that will be needed to evaluate the A matrix in the RBF deformation dxv = A * coeff = A * (Minv*dxs)

Definition at line 255 of file high_order_grid.h.

◆ ghost_surface_nodes_indexset

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
dealii::IndexSet PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::ghost_surface_nodes_indexset

Ghost surface nodes dealii::IndexSet

Definition at line 223 of file high_order_grid.h.

◆ global_index_to_point_and_axis

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
std::map<dealii::types::global_dof_index, std::pair<unsigned int, unsigned int> > PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::global_index_to_point_and_axis

Given a global DoF index, this will return the Point index within the locally_relevant_surface_points and its component. This is the inverse map of point_and_axis_to_global_index.

Definition at line 288 of file high_order_grid.h.

◆ initial_locally_relevant_surface_points

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
std::vector<dealii::Point<dim> > PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::initial_locally_relevant_surface_points

Initial locally relevant surface points.

In the optimization framework, we often displace the mesh from the initial mesh.

Definition at line 282 of file high_order_grid.h.

◆ initial_mapping_fe_field

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
std::shared_ptr<dealii::MappingFEField<dim,dim,VectorType,DoFHandlerType> > PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::initial_mapping_fe_field

MappingFEField that will provide the polynomial-based grid for the initial volume_nodes.

Will likely be used for deformations based on initial grids. It is a shared smart pointer because the constructor requires the dof_handler_grid to be properly initialized. See discussion in the following thread.

Definition at line 416 of file high_order_grid.h.

◆ initial_surface_nodes

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
VectorType PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::initial_surface_nodes

Distributed ghosted vector of initial surface nodes.

Used for deformation schemes that always deform the mesh from the initial mesh.

Definition at line 198 of file high_order_grid.h.

◆ initial_volume_nodes

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
VectorType PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::initial_volume_nodes

Initial nodal coefficients of the high-order grid.

Used for deformation schemes that always deform the mesh from the initial mesh.

Definition at line 194 of file high_order_grid.h.

◆ locally_owned_surface_nodes

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
std::vector<real> PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::locally_owned_surface_nodes

List of surface nodes.

Note that this contains all <dim> directions. By convention, the DoF representing the z-direction follows the DoF representing the y-direction, which follows the one representing the x-direction such that the integer division "idof_index / dim" gives the coordinates related to the same point.

Definition at line 249 of file high_order_grid.h.

◆ locally_owned_surface_nodes_indexset

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
dealii::IndexSet PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::locally_owned_surface_nodes_indexset

Locally owned surface nodes dealii::IndexSet

Definition at line 220 of file high_order_grid.h.

◆ locally_relevant_surface_nodes

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
std::vector<real> PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::locally_relevant_surface_nodes

List of surface nodes.

Note that this contains all <dim> directions. By convention, the DoF representing the z-direction follows the DoF representing the y-direction, which follows the one representing the x-direction such that the integer division "idof_index / dim" gives the coordinates related to the same point.

Definition at line 235 of file high_order_grid.h.

◆ locally_relevant_surface_points

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
std::vector<dealii::Point<dim> > PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::locally_relevant_surface_points

Locally relevant surface points.

Might be useful if the transformation operates on a Point. TODO: However, might want to remove this and use the surface_nodes vector for consistency and create a function to retrieve the Point.

Definition at line 277 of file high_order_grid.h.

◆ map_nodes_surf_to_vol

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
dealii::TrilinosWrappers::SparseMatrix PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::map_nodes_surf_to_vol

Maps a vector of surface_nodes to a vector of volume_nodes. This is equivalent to using surface_to_volume_indices to assign the surface_nodes values into a vector of size volume_nodes.size(), except that it uses a matrix-vector multiplication to easily obtain the mapping.

Definition at line 213 of file high_order_grid.h.

◆ mapping_fe_field

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
std::shared_ptr<dealii::MappingFEField<dim,dim,VectorType,DoFHandlerType> > PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::mapping_fe_field

MappingFEField that will provide the polynomial-based grid.

It is a shared smart pointer because the constructor requires the dof_handler_grid to be properly initialized. See discussion in the following thread.

Definition at line 408 of file high_order_grid.h.

◆ mpi_rank

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
int PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::mpi_rank
protected

This processor's MPI rank.

Definition at line 425 of file high_order_grid.h.

◆ n_locally_owned_surface_nodes_per_mpi

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
std::vector<unsigned int> PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::n_locally_owned_surface_nodes_per_mpi

Number of locally_relevant_surface_nodes per process

Definition at line 226 of file high_order_grid.h.

◆ oneD_fe_q

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
const dealii::FE_DGQ<1> PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::oneD_fe_q

Use oneD Lagrange polynomial to represent the spatial location.

We use FE_DGQ instead of FE_Q for this because DGQ orders in a tensor-product way, with x running fastest, then y, and z the slowest; whereas FE_Q order by the vertices first, then proceeds to do the edges, then volume. We need it in a tensor product indexing in order to use sum-factorization on the mapping support nodes.
Lastly, since it is collocated on the 1D Gauss-Legendre-Lobatto quadrature set, the nodal location recovers the FE_Q nodeal set, but indexed in a way we can use sum-factorization on.

Definition at line 394 of file high_order_grid.h.

◆ point_and_axis_to_global_index

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
std::map<std::pair<unsigned int, unsigned int>, dealii::types::global_dof_index> PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::point_and_axis_to_global_index

Given the Point index within the locally_relevant_surface_points and its component, this will return the global DoF index. This is the inverse map of global_index_to_point_and_axis.

Definition at line 293 of file high_order_grid.h.

◆ solution_transfer

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
SolutionTransfer PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::solution_transfer
protected

Transfers the coarse curved curve onto the fine curved grid. Used in prepare_for_coarsening_and_refinement() and execute_coarsening_and_refinement()

Definition at line 435 of file high_order_grid.h.

◆ surface_nodes

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
VectorType PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::surface_nodes

Distributed ghosted vector of surface nodes.

Definition at line 189 of file high_order_grid.h.

◆ surface_to_volume_indices

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
dealii::LinearAlgebra::distributed::Vector<int> PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::surface_to_volume_indices

Distributed ghosted vector of surface indices. Ordering matches the surface_nodes.

Definition at line 206 of file high_order_grid.h.

◆ volume_nodes

template<int dim = PHILIP_DIM, typename real = double, typename MeshType = dealii::parallel::distributed::Triangulation<dim>, typename VectorType = typename MeshTypeHelper<MeshType>::VectorType, typename DoFHandlerType = typename MeshTypeHelper<MeshType>::DoFHandlerType>
VectorType PHiLiP::HighOrderGrid< dim, real, MeshType, VectorType, DoFHandlerType >::volume_nodes

Current nodal coefficients of the high-order grid.

Note that this contains all <dim> directions. By convention, the DoF representing the z-direction follows the DoF representing the y-direction, which follows the one representing the x-direction such that the integer division "idof_index / dim" gives the coordinates related to the same point.

Definition at line 184 of file high_order_grid.h.


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