PolyDEAL
 
Loading...
Searching...
No Matches
poly_utils.h File Reference
#include <deal.II/base/config.h>
#include <deal.II/base/conditional_ostream.h>
#include <deal.II/base/point.h>
#include <deal.II/base/quadrature.h>
#include <deal.II/base/std_cxx20/iota_view.h>
#include <deal.II/boost_adaptors/bounding_box.h>
#include <deal.II/boost_adaptors/point.h>
#include <deal.II/boost_adaptors/segment.h>
#include <deal.II/distributed/tria.h>
#include <deal.II/dofs/dof_handler.h>
#include <deal.II/fe/fe_dgq.h>
#include <deal.II/fe/fe_values.h>
#include <deal.II/grid/grid_tools.h>
#include <deal.II/lac/dynamic_sparsity_pattern.h>
#include <deal.II/lac/sparse_matrix.h>
#include <deal.II/lac/sparsity_pattern.h>
#include <deal.II/lac/sparsity_tools.h>
#include <deal.II/lac/trilinos_sparse_matrix.h>
#include <deal.II/numerics/vector_tools_common.h>
#include <boost/geometry/algorithms/distance.hpp>
#include <boost/geometry/index/detail/rtree/utilities/print.hpp>
#include <boost/geometry/index/rtree.hpp>
#include <boost/geometry/strategies/strategies.hpp>
#include <deal.II/cgal/point_conversion.h>
#include <memory>

Go to the source code of this file.

Classes

struct  dealii::PolyUtils::Rtree_visitor< Value, Options, Translator, Box, Allocators >
 

Namespaces

namespace  dealii
 
namespace  dealii::PolyUtils
 
namespace  dealii::PolyUtils::internal
 

Functions

types::global_cell_index dealii::PolyUtils::internal::get_index (const std::vector< types::global_cell_index > &v, const types::global_cell_index index)
 
template<int dim, int spacedim>
void dealii::PolyUtils::internal::get_face_connectivity_of_cells (const parallel::fullydistributed::Triangulation< dim, spacedim > &triangulation, DynamicSparsityPattern &cell_connectivity, const std::vector< types::global_cell_index > locally_owned_cells)
 
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)
 
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)
 
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)
 
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)
 
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)
 
template<int dim>
std::tuple< std::vector< double >, std::vector< double >, std::vector< double >, double > dealii::PolyUtils::compute_quality_metrics (const AgglomerationHandler< dim > &ah)
 
template<int dim>
void dealii::PolyUtils::export_polygon_to_csv_file (const AgglomerationHandler< dim > &agglomeration_handler, const std::string &filename)
 
template<typename T >
constexpr T dealii::PolyUtils::constexpr_pow (T num, unsigned int pow)
 
void dealii::PolyUtils::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 dealii::PolyUtils::internal::interpolate_to_fine_grid (const AgglomerationHandler< dim, spacedim > &agglomeration_handler, VectorType &dst, const VectorType &src)
 
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)
 
template<int dim, int spacedim, typename MatrixType >
void dealii::PolyUtils::fill_interpolation_matrix (const AgglomerationHandler< dim, spacedim > &agglomeration_handler, MatrixType &interpolation_matrix)
 
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)
 
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)
 
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)
 
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)
 
template<int dim, typename MatrixType >
void dealii::PolyUtils::assemble_dg_matrix (MatrixType &system_matrix, const FiniteElement< dim > &fe_dg, const AgglomerationHandler< dim > &ah)
 
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)