41 namespace DataOutFacesImplementation
43 template <
int dim,
int spacedim>
45 const unsigned int n_datasets,
46 const unsigned int n_subdivisions,
47 const std::vector<unsigned int> &n_postprocessor_outputs,
53 :
internal::DataOutImplementation::ParallelDataBase<dim, spacedim>(
56 n_postprocessor_outputs,
69 template <
int dim,
int spacedim>
80 patches.push_back(patch);
81 patches.back().patch_index = patches.size() - 1;
88template <
int dim,
int spacedim>
95template <
int dim,
int spacedim>
103 const unsigned int face_number = cell_and_face->second;
115 for (
const unsigned int vertex : cell->face(face_number)->vertex_indices())
117 const Point<dim> vertex_reference_coordinates =
118 cell->reference_cell().template vertex<dim>(
119 cell->reference_cell().face_to_cell_vertices(
120 face_number, vertex, cell->combined_face_orientation(face_number)));
124 cell, vertex_reference_coordinates);
126 patch.
vertices[vertex] = vertex_real_coordinates;
139 const std::vector<Point<dim>> &q_points =
148 for (
unsigned int i = 0; i < dim; ++i)
149 for (
unsigned int q = 0; q < n_q_points; ++q)
150 patch.
data(patch.
data.size(0) - dim + i, q) = q_points[q][i];
153 unsigned int offset = 0;
156 for (
unsigned int dataset = 0; dataset < this->dof_data.size(); ++dataset)
160 const unsigned int n_components =
163 this->dof_data[dataset]->postprocessor;
164 if (postprocessor !=
nullptr)
171 if (n_components == 1)
176 this->dof_data[dataset]->get_function_values(
177 this_fe_patch_values,
182 this->dof_data[dataset]->get_function_gradients(
183 this_fe_patch_values,
188 this->dof_data[dataset]->get_function_hessians(
189 this_fe_patch_values,
203 dh_cell(&cell->get_triangulation(),
206 this->dof_data[dataset]->dof_handler);
208 dh_cell, face_number);
220 this->dof_data[dataset]->get_function_values(
221 this_fe_patch_values,
226 this->dof_data[dataset]->get_function_gradients(
227 this_fe_patch_values,
232 this->dof_data[dataset]->get_function_hessians(
233 this_fe_patch_values,
247 dh_cell(&cell->get_triangulation(),
250 this->dof_data[dataset]->dof_handler);
252 dh_cell, face_number);
259 for (
unsigned int q = 0; q < n_q_points; ++q)
260 for (
unsigned int component = 0;
261 component < this->dof_data[dataset]->n_output_variables;
263 patch.
data(offset + component, q) =
270 if (n_components == 1)
272 this->dof_data[dataset]->get_function_values(
273 this_fe_patch_values,
276 for (
unsigned int q = 0; q < n_q_points; ++q)
277 patch.
data(offset, q) =
283 this->dof_data[dataset]->get_function_values(
284 this_fe_patch_values,
287 for (
unsigned int component = 0; component < n_components;
289 for (
unsigned int q = 0; q < n_q_points; ++q)
290 patch.
data(offset + component, q) =
294 offset += this->dof_data[dataset]->n_output_variables;
298 for (
unsigned int dataset = 0; dataset < this->cell_data.size();
307 "The current function is trying to generate cell-data output "
308 "for a face that does not belong to an active cell. This is "
310 const unsigned int cell_number = std::distance(
314 const double value = this->cell_data[dataset]->get_cell_data_value(
317 for (
unsigned int q = 0; q < n_q_points; ++q)
318 patch.
data(dataset + offset, q) = value;
325template <
int dim,
int spacedim>
331 .template get_default_linear_mapping<dim, spacedim>(),
335 ExcMessage(
"The DataOutFaces class can currently not be "
336 "used on meshes that do not have the same cell type "
342template <
int dim,
int spacedim>
346 const unsigned int n_subdivisions_)
348 const unsigned int n_subdivisions =
349 (n_subdivisions_ != 0) ? n_subdivisions_ : this->default_subdivisions;
351 Assert(n_subdivisions >= 1,
358 this->validate_dataset_names();
360 unsigned int n_datasets = this->cell_data.size();
361 for (
unsigned int i = 0; i < this->dof_data.size(); ++i)
362 n_datasets += this->dof_data[i]->n_output_variables;
368 std::vector<FaceDescriptor> all_faces;
372 face = next_face(face))
373 all_faces.push_back(face);
376 this->patches.clear();
377 this->patches.reserve(all_faces.size());
381 std::vector<unsigned int> n_postprocessor_outputs(this->dof_data.size());
382 for (
unsigned int dataset = 0; dataset < this->dof_data.size(); ++dataset)
383 if (this->dof_data[dataset]->postprocessor)
384 n_postprocessor_outputs[dataset] =
385 this->dof_data[dataset]->n_output_variables;
387 n_postprocessor_outputs[dataset] = 0;
390 for (
unsigned int i = 0; i < this->dof_data.size(); ++i)
391 if (this->dof_data[i]->postprocessor)
393 this->dof_data[i]->postprocessor->get_needed_update_flags();
399 n_postprocessor_outputs,
409 all_faces.data() + all_faces.size(),
414 this->build_one_patch(cell_and_face, data, patch);
417 internal::DataOutFacesImplementation::append_patch_to_list<dim, spacedim>(
418 patch, this->patches);
426template <
int dim,
int spacedim>
431 for (
const auto &cell : this->
triangulation->active_cell_iterators())
432 if (cell->is_locally_owned())
433 for (
const unsigned int f : cell->face_indices())
434 if ((surface_only ==
false) || cell->face(f)->at_boundary())
445template <
int dim,
int spacedim>
454 for (
unsigned int f = face.second + 1; f < face.first->n_faces(); ++f)
455 if (!surface_only || face.first->face(f)->at_boundary())
477 if (active_cell->is_locally_owned())
478 for (
const unsigned int f : face.first->face_indices())
479 if (!surface_only || active_cell->face(f)->at_boundary())
481 face.first = active_cell;
500#include "data_out_faces.inst"
void build_one_patch(const FaceDescriptor *cell_and_face, internal::DataOutFacesImplementation::ParallelData< dim, spacedim > &data, DataOutBase::Patch< patch_dim, patch_spacedim > &patch)
typename DataOut_DoFData< dim, patch_dim, spacedim, patch_spacedim >::cell_iterator cell_iterator
virtual FaceDescriptor next_face(const FaceDescriptor &face)
typename std::pair< cell_iterator, unsigned int > FaceDescriptor
virtual void build_patches(const unsigned int n_subdivisions=0)
DataOutFaces(const bool surface_only=true)
virtual FaceDescriptor first_face()
virtual UpdateFlags get_needed_update_flags() const =0
virtual void evaluate_vector_field(const DataPostprocessorInputs::Vector< dim > &input_data, std::vector< Vector< double > > &computed_quantities) const
virtual void evaluate_scalar_field(const DataPostprocessorInputs::Scalar< dim > &input_data, std::vector< Vector< double > > &computed_quantities) const
const std::vector< Point< spacedim > > & get_quadrature_points() const
const unsigned int n_quadrature_points
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors() const
const FiniteElement< dim, spacedim > & get_fe() const
unsigned int n_components() const
Abstract base class for mapping classes.
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcInvalidNumberOfSubdivisions(int arg1)
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcNoTriangulationSelected()
static ::ExceptionBase & ExcMessage(std::string arg1)
typename ActiveSelector::active_cell_iterator active_cell_iterator
@ update_hessians
Second derivatives of shape functions.
@ update_values
Shape function values.
@ update_normal_vectors
Normal vectors.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
void run(const std::vector< std::vector< Iterator > > &colored_iterators, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8)
void append_patch_to_list(const DataOutBase::Patch< DataOutFaces< dim, spacedim >::patch_dim, DataOutFaces< dim, spacedim >::patch_spacedim > &patch, std::vector< DataOutBase::Patch< DataOutFaces< dim, spacedim >::patch_dim, DataOutFaces< dim, spacedim >::patch_spacedim > > &patches)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
ReferenceCell reference_cell
static const unsigned int space_dim
unsigned int n_subdivisions
std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > vertices
bool points_are_available
ParallelData(const unsigned int n_datasets, const unsigned int n_subdivisions, const std::vector< unsigned int > &n_postprocessor_outputs, const Mapping< dim, spacedim > &mapping, const std::vector< std::shared_ptr<::hp::FECollection< dim, spacedim > > > &finite_elements, const UpdateFlags update_flags)
const FEValuesBase< dim, spacedim > & get_present_fe_values(const unsigned int dataset) const
const unsigned int n_datasets
DataPostprocessorInputs::Scalar< spacedim > patch_values_scalar
std::vector< std::vector<::Vector< double > > > postprocessed_values
const ::hp::MappingCollection< dim, spacedim > mapping_collection
DataPostprocessorInputs::Vector< spacedim > patch_values_system
void reinit_all_fe_values(std::vector< std::shared_ptr< DataEntryBase< dim, spacedim > > > &dof_data, const typename ::Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face=numbers::invalid_unsigned_int)
void resize_system_vectors(const unsigned int n_components)