24#include <deal.II/numerics/data_out_dof_data.templates.h>
33template <
int dim,
int patch_dim,
int spacedim>
37 : patch_dof_handler(patch_tria)
38 , patch_mapping(&patch_mapping)
43template <
int dim,
int patch_dim,
int spacedim>
47 const unsigned int n_subdivisions)
49 this->mapping = &mapping;
50 this->point_to_local_vector_indices.clear();
53 std::max<unsigned int>(1, n_subdivisions));
54 patch_dof_handler.distribute_dofs(fe);
56 std::vector<std::pair<types::global_dof_index, Point<spacedim>>> points_all;
69 partitioner = std::make_shared<Utilities::MPI::Partitioner>(
70 patch_dof_handler.locally_owned_dofs(), active_dofs, MPI_COMM_WORLD);
72 for (
const auto &cell : patch_dof_handler.active_cell_iterators() |
77 cell->get_dof_indices(dof_indices);
79 for (
unsigned int i = 0; i < dof_indices.size(); ++i)
80 points_all.emplace_back(partitioner->global_to_local(dof_indices[i]),
84 std::sort(points_all.begin(),
86 [](
const auto &a,
const auto &b) { return a.first < b.first; });
87 points_all.erase(std::unique(points_all.begin(),
89 [](
const auto &a,
const auto &b) {
90 return a.first == b.first;
94 std::vector<Point<spacedim>> points;
96 for (
const auto &i : points_all)
98 point_to_local_vector_indices.push_back(i.first);
99 points.push_back(i.second);
107template <
int dim,
int patch_dim,
int spacedim>
111 const unsigned int n_subdivisions,
115 this->build_patches(curved_region);
120template <
int dim,
int patch_dim,
int spacedim>
125 patch_data_out.clear();
127 if (rpe.is_ready() ==
false)
132 "Mapping is not valid anymore! Please register a new mapping via "
133 "update_mapping() or the other build_patches() function."));
134 update_mapping(*this->mapping, patch_dof_handler.get_fe().degree);
137 std::vector<std::shared_ptr<LinearAlgebra::distributed::Vector<double>>>
140 patch_data_out.attach_dof_handler(patch_dof_handler);
142 unsigned int counter = 0;
144 for (
const auto &data : this->dof_data)
146 const auto data_ptr =
dynamic_cast<
147 internal::DataOutImplementation::DataEntry<dim, spacedim, double> *
>(
152 const auto &dh = *data_ptr->dof_handler;
155 for (
const auto &fe : dh.get_fe_collection())
157 fe.n_base_elements() == 1,
159 "This class currently only supports scalar elements and elements "
160 "with a single base element."));
163 for (
unsigned int comp = 0; comp < dh.get_fe_collection().n_components();
166 const auto values = VectorTools::point_values<1>(
169 vectors.emplace_back(
173 for (
unsigned int j = 0; j < values.size(); ++j)
174 vectors.back()->local_element(point_to_local_vector_indices[j]) =
177 vectors.back()->set_ghost_state(
true);
182 patch_data_out.add_data_vector(
184 std::string(
"temp_" + std::to_string(counter)),
186 DataVectorType::type_dof_data);
192 patch_data_out.build_patches(*patch_mapping,
193 patch_dof_handler.get_fe().degree,
199template <
int dim,
int patch_dim,
int spacedim>
200const std::vector<typename DataOutBase::Patch<patch_dim, spacedim>> &
203 return patch_data_out.get_patches();
208#include "data_out_resample.inst"
virtual const std::vector< typename DataOutBase::Patch< patch_dim, spacedim > > & get_patches() const override
DataOutResample(const Triangulation< patch_dim, spacedim > &patch_tria, const Mapping< patch_dim, spacedim > &patch_mapping)
void build_patches(const Mapping< dim, spacedim > &mapping, const unsigned int n_subdivisions=0, const typename DataOut< patch_dim, spacedim >::CurvedCellRegion curved_region=DataOut< patch_dim, spacedim >::CurvedCellRegion::curved_boundary)
void update_mapping(const Mapping< dim, spacedim > &mapping, const unsigned int n_subdivisions=0)
const Point< spacedim > & quadrature_point(const unsigned int q_point) const
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
unsigned int n_dofs_per_cell() const
const std::vector< Point< dim > > & get_unit_support_points() const
Abstract base class for mapping classes.
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcMessage(std::string arg1)
@ update_quadrature_points
Transformed quadrature points.
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation