PolyDEAL
 
Loading...
Searching...
No Matches
dealii::PolyUtils Namespace Reference

Namespaces

namespace  internal
 

Classes

struct  Rtree_visitor
 

Functions

template<typename Rtree >
std::pair< std::vector< std::vector< unsigned int > >, std::vector< std::vector< typename Triangulation< boost::geometry::dimension< typename Rtree::indexable_type >::value >::active_cell_iterator > > > extract_children_of_level (const Rtree &tree, const unsigned int level)
 
template<int dim, typename Number = double>
Number compute_h_orthogonal (const unsigned int face_index, const std::vector< typename Triangulation< dim >::active_face_iterator > &polygon_boundary, const Tensor< 1, dim > &deal_normal)
 
template<int dim, int spacedim = dim>
void collect_cells_for_agglomeration (const Triangulation< dim, spacedim > &tria, const std::vector< types::global_cell_index > &cell_idxs, std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > &cells_to_be_agglomerated)
 
template<int dim, int spacedim>
void partition_locally_owned_regions (const unsigned int n_partitions, Triangulation< dim, spacedim > &triangulation, const SparsityTools::Partitioner partitioner)
 
template<int dim, int spacedim>
void partition_locally_owned_regions (AgglomerationHandler< dim > &agglomeration_handler, const unsigned int n_partitions, Triangulation< dim, spacedim > &triangulation, const SparsityTools::Partitioner partitioner)
 
template<int dim>
std::tuple< std::vector< double >, std::vector< double >, std::vector< double >, double > compute_quality_metrics (const AgglomerationHandler< dim > &ah)
 
template<int dim>
void export_polygon_to_csv_file (const AgglomerationHandler< dim > &agglomeration_handler, const std::string &filename)
 
template<typename T >
constexpr T constexpr_pow (T num, unsigned int pow)
 
void write_to_matrix_market_format (const std::string &filename, const std::string &matrix_name, const TrilinosWrappers::SparseMatrix &matrix)
 
template<int dim, int spacedim, typename VectorType >
void interpolate_to_fine_grid (const AgglomerationHandler< dim, spacedim > &agglomeration_handler, VectorType &dst, const VectorType &src, const bool on_the_fly=true)
 
template<int dim, int spacedim, typename MatrixType >
void fill_interpolation_matrix (const AgglomerationHandler< dim, spacedim > &agglomeration_handler, MatrixType &interpolation_matrix)
 
template<int dim, typename Number , typename VectorType >
void compute_global_error (const AgglomerationHandler< dim > &agglomeration_handler, const VectorType &solution, const Function< dim, Number > &exact_solution, const std::vector< VectorTools::NormType > &norms, std::vector< double > &global_errors)
 
template<int dim>
unsigned int construct_agglomerated_levels (const Triangulation< dim > &tria, std::vector< std::unique_ptr< AgglomerationHandler< dim > > > &agglomeration_handlers, const FE_DGQ< dim > &fe_dg, const Mapping< dim > &mapping, const unsigned int starting_tree_level)
 
template<int dim>
void assemble_local_jumps_and_averages (FullMatrix< double > &M11, FullMatrix< double > &M12, FullMatrix< double > &M21, FullMatrix< double > &M22, const FEValuesBase< dim > &fe_faces0, const FEValuesBase< dim > &fe_faces1, const double penalty_constant, const double h_f)
 
template<int dim>
void assemble_local_jumps_and_averages_ghost (FullMatrix< double > &M11, FullMatrix< double > &M12, FullMatrix< double > &M21, FullMatrix< double > &M22, const FEValuesBase< dim > &fe_faces0, const std::vector< std::vector< double > > &recv_values, const std::vector< std::vector< Tensor< 1, dim > > > &recv_gradients, const std::vector< double > &recv_jxws, const double penalty_constant, const double h_f)
 
template<int dim, typename MatrixType >
void assemble_dg_matrix (MatrixType &system_matrix, const FiniteElement< dim > &fe_dg, const AgglomerationHandler< dim > &ah)
 
template<int dim, typename MatrixType , typename VectorType >
void assemble_dg_matrix_on_standard_mesh (MatrixType &system_matrix, VectorType &system_rhs, const Mapping< dim > &mapping, const FiniteElement< dim > &fe_dg, const DoFHandler< dim > &dof_handler)
 

Function Documentation

◆ assemble_dg_matrix()

◆ assemble_dg_matrix_on_standard_mesh()

template<int dim, typename MatrixType , typename VectorType >
void dealii::PolyUtils::assemble_dg_matrix_on_standard_mesh ( MatrixType & system_matrix,
VectorType & system_rhs,
const Mapping< dim > & mapping,
const FiniteElement< dim > & fe_dg,
const DoFHandler< dim > & dof_handler )

Compute SIPDG matrix as well as rhs vector.

Note
Hardcoded for $f=1$ and simplex elements. TODO: Pass Function object for boundary conditions and forcing term.

Definition at line 2012 of file poly_utils.h.

References DoFHandler< int dim, int spacedim >::active_cell_iterators(), VectorOperation::add, Triangulation< int dim, int spacedim >::all_reference_cells_are_simplex(), Assert, AffineConstraints< typename number >::close(), FiniteElement< int dim, int spacedim >::degree, AffineConstraints< typename number >::distribute_local_to_global(), SparsityTools::distribute_sparsity_pattern(), DoFTools::extract_locally_relevant_dofs(), DoFHandler< int dim, int spacedim >::get_communicator(), DoFHandler< int dim, int spacedim >::get_mpi_communicator(), FEFaceValues< int dim, int spacedim >::get_normal_vectors(), DoFHandler< int dim, int spacedim >::get_triangulation(), FEFaceValues< int dim, int spacedim >::JxW(), FEValues< int dim, int spacedim >::JxW(), DoFHandler< int dim, int spacedim >::locally_owned_dofs(), DoFTools::make_flux_sparsity_pattern(), FiniteElement< int dim, int spacedim >::n_dofs_per_cell(), FEFaceValues< int dim, int spacedim >::quadrature_point_indices(), FEValues< int dim, int spacedim >::quadrature_point_indices(), FEFaceValues< int dim, int spacedim >::reinit(), FEValues< int dim, int spacedim >::reinit(), FEFaceValues< int dim, int spacedim >::shape_grad(), FEValues< int dim, int spacedim >::shape_grad(), FEFaceValues< int dim, int spacedim >::shape_value(), FEValues< int dim, int spacedim >::shape_value(), update_gradients, update_JxW_values, update_normal_vectors, update_quadrature_points, and update_values.

◆ assemble_local_jumps_and_averages()

template<int dim>
void dealii::PolyUtils::assemble_local_jumps_and_averages ( FullMatrix< double > & M11,
FullMatrix< double > & M12,
FullMatrix< double > & M21,
FullMatrix< double > & M22,
const FEValuesBase< dim > & fe_faces0,
const FEValuesBase< dim > & fe_faces1,
const double penalty_constant,
const double h_f )

◆ assemble_local_jumps_and_averages_ghost()

template<int dim>
void dealii::PolyUtils::assemble_local_jumps_and_averages_ghost ( FullMatrix< double > & M11,
FullMatrix< double > & M12,
FullMatrix< double > & M21,
FullMatrix< double > & M22,
const FEValuesBase< dim > & fe_faces0,
const std::vector< std::vector< double > > & recv_values,
const std::vector< std::vector< Tensor< 1, dim > > > & recv_gradients,
const std::vector< double > & recv_jxws,
const double penalty_constant,
const double h_f )

◆ collect_cells_for_agglomeration()

template<int dim, int spacedim = dim>
void dealii::PolyUtils::collect_cells_for_agglomeration ( const Triangulation< dim, spacedim > & tria,
const std::vector< types::global_cell_index > & cell_idxs,
std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > & cells_to_be_agglomerated )

Agglomerate cells together based on their global index. This function is not efficient and should be used for testing purposes only.

Definition at line 523 of file poly_utils.h.

References Triangulation< int dim, int spacedim >::active_cell_iterators(), and Assert.

Referenced by Poisson< dim >::setup_agglomeration().

◆ compute_global_error()

template<int dim, typename Number , typename VectorType >
void dealii::PolyUtils::compute_global_error ( const AgglomerationHandler< dim > & agglomeration_handler,
const VectorType & solution,
const Function< dim, Number > & exact_solution,
const std::vector< VectorTools::NormType > & norms,
std::vector< double > & global_errors )

Similar to VectorTools::compute_global_error(), but customized for polytopic elements. Aside from the solution vector and a reference function, this function takes in addition a vector norms with types VectorTools::NormType to be computed and later stored in the last argument global_errors. In case of a parallel vector, the local errors are collected over each processor and later a classical reduction operation is performed.

Definition at line 1455 of file poly_utils.h.

References Assert, b(), d(), Triangulation< int dim, int spacedim >::get_mpi_communicator(), FEValues< int dim, int spacedim >::get_quadrature_points(), AgglomerationHandler< dim, spacedim >::get_triangulation(), Function< int dim, typename RangeNumberType >::gradient_list(), AgglomerationHandler< dim, spacedim >::n_dofs_per_cell(), AgglomerationHandler< dim, spacedim >::polytope_iterators(), std::pow(), Utilities::MPI::reduce(), AgglomerationHandler< dim, spacedim >::reinit(), std::sqrt(), and Function< int dim, typename RangeNumberType >::value_list().

◆ compute_h_orthogonal()

template<int dim, typename Number = double>
Number dealii::PolyUtils::compute_h_orthogonal ( const unsigned int face_index,
const std::vector< typename Triangulation< dim >::active_face_iterator > & polygon_boundary,
const Tensor< 1, dim > & deal_normal )

Definition at line 401 of file poly_utils.h.

References Assert, center, and std::sqrt().

◆ compute_quality_metrics()

◆ constexpr_pow()

template<typename T >
T dealii::PolyUtils::constexpr_pow ( T num,
unsigned int pow )
inlineconstexpr

Definition at line 896 of file poly_utils.h.

References constexpr_pow(), and pow().

Referenced by constexpr_pow(), and construct_agglomerated_levels().

◆ construct_agglomerated_levels()

template<int dim>
unsigned int dealii::PolyUtils::construct_agglomerated_levels ( const Triangulation< dim > & tria,
std::vector< std::unique_ptr< AgglomerationHandler< dim > > > & agglomeration_handlers,
const FE_DGQ< dim > & fe_dg,
const Mapping< dim > & mapping,
const unsigned int starting_tree_level )

Utility function that builds the multilevel hierarchy from the tree level starting_level. This function fills the vector of AgglomerationHandlers objects by distributing degrees of freedom on each level of the hierarchy. It returns the total number of levels in the hierarchy.

Definition at line 1569 of file poly_utils.h.

References Assert, comm, constexpr_pow(), FE_DGQ< int dim, int spacedim >::degree, dealii::CellsAgglomerator< dim, RtreeType >::extract_agglomerates(), Mapping< int dim, int spacedim >::get_bounding_box(), Utilities::MPI::sum(), Utilities::MPI::this_mpi_process(), update_gradients, update_JxW_values, update_quadrature_points, and update_values.

◆ export_polygon_to_csv_file()

template<int dim>
void dealii::PolyUtils::export_polygon_to_csv_file ( const AgglomerationHandler< dim > & agglomeration_handler,
const std::string & filename )

Export each polygon in a csv file as a collection of segments.

Definition at line 863 of file poly_utils.h.

References AgglomerationHandler< dim, spacedim >::polytope_iterators().

◆ extract_children_of_level()

template<typename Rtree >
std::pair< std::vector< std::vector< unsigned int > >, std::vector< std::vector< typename Triangulation< boost::geometry::dimension< typename Rtree::indexable_type >::value >::active_cell_iterator > > > dealii::PolyUtils::extract_children_of_level ( const Rtree & tree,
const unsigned int level )
inline

Definition at line 359 of file poly_utils.h.

References AssertDimension, level, and std::min().

◆ fill_interpolation_matrix()

template<int dim, int spacedim, typename MatrixType >
void dealii::PolyUtils::fill_interpolation_matrix ( const AgglomerationHandler< dim, spacedim > & agglomeration_handler,
MatrixType & interpolation_matrix )

Construct the interpolation matrix from the DG space defined the polytopic elements defined in agglomeration_handler to the DG space defined on the DoFHandler associated to standard shapes. The interpolation matrix is assumed to be default-constructed and is filled inside this function.

Definition at line 1277 of file poly_utils.h.

References DoFHandler< int dim, int spacedim >::active_cell_iterators(), VectorOperation::add, DynamicSparsityPattern::add_entries(), AgglomerationHandler< dim, spacedim >::agglo_dh, Assert, AssertThrow, AgglomerationHandler< dim, spacedim >::cell_to_polytope_index(), AffineConstraints< typename number >::close(), DoFHandler< int dim, int spacedim >::distribute_dofs(), AffineConstraints< typename number >::distribute_local_to_global(), SparsityTools::distribute_sparsity_pattern(), FEValues< int dim, int spacedim >::dof_indices(), FiniteElement< int dim, int spacedim >::dofs_per_cell, DoFTools::extract_locally_relevant_dofs(), Triangulation< int dim, int spacedim >::get_communicator(), AgglomerationHandler< dim, spacedim >::get_fe(), AgglomerationHandler< dim, spacedim >::get_local_bboxes(), AgglomerationHandler< dim, spacedim >::get_mapping(), FEValues< int dim, int spacedim >::get_quadrature_points(), AgglomerationHandler< dim, spacedim >::get_slaves_of_idx(), AgglomerationHandler< dim, spacedim >::get_triangulation(), FiniteElement< int dim, int spacedim >::get_unit_support_points(), AgglomerationHandler< dim, spacedim >::is_master_cell(), DoFHandler< int dim, int spacedim >::locally_owned_dofs(), DoFHandler< int dim, int spacedim >::n_dofs(), AgglomerationHandler< dim, spacedim >::output_dh, BoundingBox< int spacedim, typename Number >::real_to_unit(), DoFHandler< int dim, int spacedim >::reinit(), FEValues< int dim, int spacedim >::reinit(), FiniteElement< int dim, int spacedim >::shape_value(), and update_quadrature_points.

◆ interpolate_to_fine_grid()

template<int dim, int spacedim, typename VectorType >
void dealii::PolyUtils::interpolate_to_fine_grid ( const AgglomerationHandler< dim, spacedim > & agglomeration_handler,
VectorType & dst,
const VectorType & src,
const bool on_the_fly = true )

Given a vector src, typically the solution stemming after the agglomerate problem has been solved, this function interpolates src onto the finer grid and stores the result in vector dst. The last argument on_the_fly does not build any interpolation matrix and allows computing the entries in dst in a matrix-free fashion.

Note
Supported parallel types are TrilinosWrappers::SparseMatrix and TrilinosWrappers::MPI::Vector.

Definition at line 1147 of file poly_utils.h.

References Triangulation< int dim, int spacedim >::all_reference_cells_are_hyper_cube(), Triangulation< int dim, int spacedim >::all_reference_cells_are_simplex(), Assert, AssertThrow, FiniteElement< int dim, int spacedim >::degree, DoFHandler< int dim, int spacedim >::distribute_dofs(), AgglomerationHandler< dim, spacedim >::get_fe(), AgglomerationHandler< dim, spacedim >::get_local_bboxes(), AgglomerationHandler< dim, spacedim >::get_mapping(), FEValues< int dim, int spacedim >::get_quadrature_points(), AgglomerationHandler< dim, spacedim >::get_triangulation(), dealii::PolyUtils::internal::interpolate_to_fine_grid(), DoFHandler< int dim, int spacedim >::locally_owned_dofs(), DoFHandler< int dim, int spacedim >::n_dofs(), AgglomerationHandler< dim, spacedim >::n_dofs_per_cell(), AgglomerationHandler< dim, spacedim >::output_dh, AgglomerationHandler< dim, spacedim >::polytope_iterators(), BoundingBox< int spacedim, typename Number >::real_to_unit(), DoFHandler< int dim, int spacedim >::reinit(), FEValues< int dim, int spacedim >::reinit(), FiniteElement< int dim, int spacedim >::shape_value(), and update_quadrature_points.

◆ partition_locally_owned_regions() [1/2]

template<int dim, int spacedim>
void dealii::PolyUtils::partition_locally_owned_regions ( AgglomerationHandler< dim > & agglomeration_handler,
const unsigned int n_partitions,
Triangulation< dim, spacedim > & triangulation,
const SparsityTools::Partitioner partitioner )

Partition with METIS the locally owned regions of the given triangulation and insert agglomerates in the polytopic grid.

Note
The given triangulation must be a parallel::fullydistributed::Triangulation. This is required as the partitions generated by p4est, the partitioner for parallell::distributed::Triangulation, can generate discontinuous partitions which are not supported by the METIS partitioner.

Definition at line 631 of file poly_utils.h.

References Assert, AssertDimension, SparsityPattern::copy_from(), AgglomerationHandler< dim, spacedim >::define_agglomerate(), dealii::PolyUtils::internal::get_face_connectivity_of_cells(), dealii::PolyUtils::internal::get_index(), AgglomerationHandler< dim, spacedim >::n_agglomerates(), SparsityTools::partition(), and triangulation.

◆ partition_locally_owned_regions() [2/2]

template<int dim, int spacedim>
void dealii::PolyUtils::partition_locally_owned_regions ( const unsigned int n_partitions,
Triangulation< dim, spacedim > & triangulation,
const SparsityTools::Partitioner partitioner )

Partition with METIS the locally owned regions of the given triangulation.

Note
The given triangulation must be a parallel::fullydistributed::Triangulation. This is required as the partitions generated by p4est, the partitioner for parallell::distributed::Triangulation, can generate discontinuous partitions which are not supported by the METIS partitioner.

Definition at line 555 of file poly_utils.h.

References Assert, AssertDimension, SparsityPattern::copy_from(), dealii::PolyUtils::internal::get_face_connectivity_of_cells(), dealii::PolyUtils::internal::get_index(), SparsityTools::partition(), and triangulation.

◆ write_to_matrix_market_format()

void dealii::PolyUtils::write_to_matrix_market_format ( const std::string & filename,
const std::string & matrix_name,
const TrilinosWrappers::SparseMatrix & matrix )

Definition at line 906 of file poly_utils.h.

References AssertThrow, and ExcTrilinosError().