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>
75 patches.push_back(patch);
76 patches.back().patch_index = patches.size() - 1;
83template <
int dim,
typename DoFHandlerType>
92template <
int dim,
typename DoFHandlerType>
106 for (
unsigned int vertex = 0;
107 vertex <
GeometryInfo<dimension - 1>::vertices_per_cell;
111 cell_and_face->first,
114 cell_and_face->second,
116 cell_and_face->first->face_orientation(cell_and_face->second),
117 cell_and_face->first->face_flip(cell_and_face->second),
118 cell_and_face->first->face_rotation(cell_and_face->second))));
123 cell_and_face->first,
124 cell_and_face->second);
132 const std::vector<Point<dimension>> &q_points =
141 for (
unsigned int i = 0; i < dimension; ++i)
142 for (
unsigned int q = 0; q < n_q_points; ++q)
143 patch.
data(patch.
data.size(0) - dimension + i, q) = q_points[q][i];
146 unsigned int offset = 0;
149 for (
unsigned int dataset = 0; dataset < this->dof_data.size(); ++dataset)
153 const unsigned int n_components =
156 this->dof_data[dataset]->postprocessor;
157 if (postprocessor !=
nullptr)
164 if (n_components == 1)
169 this->dof_data[dataset]->get_function_values(
170 this_fe_patch_values,
175 this->dof_data[dataset]->get_function_gradients(
176 this_fe_patch_values,
181 this->dof_data[dataset]->get_function_hessians(
182 this_fe_patch_values,
195 const typename DoFHandlerType::active_cell_iterator dh_cell(
196 &cell_and_face->first->get_triangulation(),
197 cell_and_face->first->level(),
198 cell_and_face->first->index(),
199 this->dof_data[dataset]->dof_handler);
213 this->dof_data[dataset]->get_function_values(
214 this_fe_patch_values,
219 this->dof_data[dataset]->get_function_gradients(
220 this_fe_patch_values,
225 this->dof_data[dataset]->get_function_hessians(
226 this_fe_patch_values,
239 const typename DoFHandlerType::active_cell_iterator dh_cell(
240 &cell_and_face->first->get_triangulation(),
241 cell_and_face->first->level(),
242 cell_and_face->first->index(),
243 this->dof_data[dataset]->dof_handler);
252 for (
unsigned int q = 0; q < n_q_points; ++q)
253 for (
unsigned int component = 0;
254 component < this->dof_data[dataset]->n_output_variables;
256 patch.
data(offset + component, q) =
263 if (n_components == 1)
265 this->dof_data[dataset]->get_function_values(
266 this_fe_patch_values,
269 for (
unsigned int q = 0; q < n_q_points; ++q)
270 patch.
data(offset, q) =
276 this->dof_data[dataset]->get_function_values(
277 this_fe_patch_values,
280 for (
unsigned int component = 0; component < n_components;
282 for (
unsigned int q = 0; q < n_q_points; ++q)
283 patch.
data(offset + component, q) =
287 offset += this->dof_data[dataset]->n_output_variables;
291 for (
unsigned int dataset = 0; dataset < this->cell_data.size();
298 cell_and_face->first->is_active(),
300 "The current function is trying to generate cell-data output "
301 "for a face that does not belong to an active cell. This is "
303 const unsigned int cell_number =
306 active_cell_iterator(cell_and_face->first));
308 const double value = this->cell_data[dataset]->get_cell_data_value(
311 for (
unsigned int q = 0; q < n_q_points; ++q)
312 patch.
data(dataset + offset, q) = value;
319template <
int dim,
typename DoFHandlerType>
322 const unsigned int n_subdivisions_)
329template <
int dim,
typename DoFHandlerType>
333 const unsigned int n_subdivisions_)
338 const unsigned int n_subdivisions =
339 (n_subdivisions_ != 0) ? n_subdivisions_ : this->default_subdivisions;
341 Assert(n_subdivisions >= 1,
348 this->validate_dataset_names();
350 unsigned int n_datasets = this->cell_data.size();
351 for (
unsigned int i = 0; i < this->dof_data.size(); ++i)
352 n_datasets += this->dof_data[i]->n_output_variables;
358 std::vector<FaceDescriptor> all_faces;
362 face = next_face(face))
363 all_faces.push_back(face);
366 this->patches.clear();
367 this->patches.reserve(all_faces.size());
371 std::vector<unsigned int> n_postprocessor_outputs(this->dof_data.size());
372 for (
unsigned int dataset = 0; dataset < this->dof_data.size(); ++dataset)
373 if (this->dof_data[dataset]->postprocessor)
374 n_postprocessor_outputs[dataset] =
375 this->dof_data[dataset]->n_output_variables;
377 n_postprocessor_outputs[dataset] = 0;
380 for (
unsigned int i = 0; i < this->dof_data.size(); ++i)
381 if (this->dof_data[i]->postprocessor)
383 this->dof_data[i]->postprocessor->get_needed_update_flags();
387 thread_data(n_datasets,
389 n_postprocessor_outputs,
395 sample_patch.data.reinit(
396 n_datasets, Utilities::fixed_power<dimension - 1>(n_subdivisions + 1));
401 all_faces.data() + all_faces.size(),
406 this->build_one_patch(cell_and_face, data, patch);
409 internal::DataOutFacesImplementation::
410 append_patch_to_list<dim, space_dimension>(patch, this->patches);
418template <
int dim,
typename DoFHandlerType>
426 if (cell->is_locally_owned())
428 if (!surface_only || cell->face(f)->at_boundary())
439template <
int dim,
typename DoFHandlerType>
448 for (
unsigned int f = face.second + 1;
451 if (!surface_only || face.first->face(f)->at_boundary())
463 active_cell = face.first;
473 if (active_cell->is_locally_owned())
475 if (!surface_only || active_cell->face(f)->at_boundary())
477 face.first = active_cell;
496#include "data_out_faces.inst"
virtual FaceDescriptor first_face()
DataOutFaces(const bool surface_only=true)
virtual void build_patches(const unsigned int n_subdivisions=0)
virtual FaceDescriptor next_face(const FaceDescriptor &face)
void build_one_patch(const FaceDescriptor *cell_and_face, internal::DataOutFacesImplementation::ParallelData< dimension, dimension > &data, DataOutBase::Patch< dimension - 1, space_dimension > &patch)
typename std::pair< cell_iterator, unsigned int > FaceDescriptor
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
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
@ 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.
static ::ExceptionBase & ExcInvalidNumberOfSubdivisions(int arg1)
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
Point< spacedim > vertices[GeometryInfo< dim >::vertices_per_cell]
static const unsigned int space_dim
static ::ExceptionBase & ExcInternalError()
unsigned int n_subdivisions
static ::ExceptionBase & ExcNoTriangulationSelected()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
bool points_are_available
static ::ExceptionBase & ExcMessage(std::string arg1)
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< dim - 1, spacedim > &patch, std::vector< DataOutBase::Patch< dim - 1, spacedim > > &patches)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
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)