54 namespace DataOutRotationImplementation
56 template <
int dim,
int spacedim>
58 const unsigned int n_datasets,
59 const unsigned int n_subdivisions,
60 const unsigned int n_patches_per_circle,
61 const std::vector<unsigned int> &n_postprocessor_outputs,
67 :
internal::DataOutImplementation::ParallelDataBase<dim, spacedim>(
70 n_postprocessor_outputs,
75 , n_patches_per_circle(n_patches_per_circle)
84 template <
int dim,
int spacedim>
90 for (
unsigned int i = 0; i < new_patches.size(); ++i)
92 patches.push_back(new_patches[i]);
93 patches.back().patch_index = patches.size() - 1;
101template <
int dim,
typename DoFHandlerType>
106 space_dimension> &data,
120 const unsigned int n_patches_per_circle = data.n_patches_per_circle;
123 const unsigned int n_points = data.n_subdivisions + 1;
129 std::vector<Point<dimension + 1>> angle_directions(n_patches_per_circle + 1);
130 for (
unsigned int i = 0; i <= n_patches_per_circle; ++i)
132 angle_directions[i][dimension - 1] =
134 angle_directions[i][dimension] =
147 const double r1 = (*cell)->vertex(0)(0),
148 r2 = (*cell)->vertex(1)(0);
149 Assert(r1 >= 0, ExcRadialVariableHasNegativeValues(r1));
150 Assert(r2 >= 0, ExcRadialVariableHasNegativeValues(r2));
152 my_patches[
angle].vertices[0] = r1 * angle_directions[
angle];
153 my_patches[
angle].vertices[1] = r2 * angle_directions[
angle];
154 my_patches[
angle].vertices[2] = r1 * angle_directions[
angle + 1];
155 my_patches[
angle].vertices[3] = r2 * angle_directions[
angle + 1];
162 for (
const unsigned int vertex :
168 Assert(v(0) >= 0, ExcRadialVariableHasNegativeValues(v(0)));
171 my_patches[
angle].vertices[vertex] =
172 v(0) * angle_directions[
angle];
173 my_patches[
angle].vertices[vertex][0] = v(1);
178 v(0) * angle_directions[
angle + 1];
193 if (data.n_datasets > 0)
195 unsigned int offset = 0;
197 data.reinit_all_fe_values(this->dof_data, *cell);
199 for (
unsigned int dataset = 0; dataset < this->dof_data.size();
203 data.get_present_fe_values(dataset);
204 const unsigned int n_components =
207 this->dof_data[dataset]->postprocessor;
208 if (postprocessor !=
nullptr)
216 if (n_components == 1)
222 this->dof_data[dataset]->get_function_values(
226 data.patch_values_scalar.solution_values);
228 this->dof_data[dataset]->get_function_gradients(
232 data.patch_values_scalar.solution_gradients);
234 this->dof_data[dataset]->get_function_hessians(
238 data.patch_values_scalar.solution_hessians);
241 data.patch_values_scalar.evaluation_points =
244 const typename DoFHandlerType::active_cell_iterator
245 dh_cell(&(*cell)->get_triangulation(),
248 this->dof_data[dataset]->dof_handler);
249 data.patch_values_scalar
250 .template set_cell<DoFHandlerType>(dh_cell);
253 data.patch_values_scalar,
254 data.postprocessed_values[dataset]);
258 data.resize_system_vectors(n_components);
263 this->dof_data[dataset]->get_function_values(
267 data.patch_values_system.solution_values);
269 this->dof_data[dataset]->get_function_gradients(
273 data.patch_values_system.solution_gradients);
275 this->dof_data[dataset]->get_function_hessians(
279 data.patch_values_system.solution_hessians);
282 data.patch_values_system.evaluation_points =
285 const typename DoFHandlerType::active_cell_iterator
286 dh_cell(&(*cell)->get_triangulation(),
289 this->dof_data[dataset]->dof_handler);
290 data.patch_values_system
291 .template set_cell<DoFHandlerType>(dh_cell);
294 data.patch_values_system,
295 data.postprocessed_values[dataset]);
298 for (
unsigned int component = 0;
299 component < this->dof_data[dataset]->n_output_variables;
305 for (
unsigned int x = 0; x < n_points; ++x)
306 for (
unsigned int y = 0; y < n_points; ++y)
307 my_patches[
angle].data(offset + component,
309 data.postprocessed_values[dataset][x](
314 for (
unsigned int x = 0; x < n_points; ++x)
315 for (
unsigned int y = 0; y < n_points; ++y)
316 for (
unsigned int z = 0; z < n_points; ++z)
317 my_patches[
angle].data(offset + component,
321 data.postprocessed_values[dataset]
331 else if (n_components == 1)
333 this->dof_data[dataset]->get_function_values(
337 data.patch_values_scalar.solution_values);
342 for (
unsigned int x = 0; x < n_points; ++x)
343 for (
unsigned int y = 0; y < n_points; ++y)
344 my_patches[
angle].data(offset, x * n_points + y) =
345 data.patch_values_scalar.solution_values[x];
349 for (
unsigned int x = 0; x < n_points; ++x)
350 for (
unsigned int y = 0; y < n_points; ++y)
351 for (
unsigned int z = 0; z < n_points; ++z)
352 my_patches[
angle].data(offset,
353 x * n_points * n_points +
355 data.patch_values_scalar
356 .solution_values[x * n_points + z];
366 data.resize_system_vectors(n_components);
367 this->dof_data[dataset]->get_function_values(
371 data.patch_values_system.solution_values);
373 for (
unsigned int component = 0; component < n_components;
379 for (
unsigned int x = 0; x < n_points; ++x)
380 for (
unsigned int y = 0; y < n_points; ++y)
381 my_patches[
angle].data(offset + component,
383 data.patch_values_system.solution_values[x](
388 for (
unsigned int x = 0; x < n_points; ++x)
389 for (
unsigned int y = 0; y < n_points; ++y)
390 for (
unsigned int z = 0; z < n_points; ++z)
391 my_patches[
angle].data(offset + component,
395 data.patch_values_system
396 .solution_values[x * n_points + z](
405 offset += this->dof_data[dataset]->n_output_variables;
409 for (
unsigned int dataset = 0; dataset < this->cell_data.size();
415 Assert((*cell)->is_active(),
416 ExcMessage(
"Cell must be active for cell data"));
417 const unsigned int cell_number = std::distance(
420 active_cell_iterator(*cell));
422 this->cell_data[dataset]->get_cell_data_value(
429 for (
unsigned int x = 0; x < n_points; ++x)
430 for (
unsigned int y = 0; y < n_points; ++y)
431 my_patches[
angle].data(dataset + offset,
432 x * n_points + y) = value;
436 for (
unsigned int x = 0; x < n_points; ++x)
437 for (
unsigned int y = 0; y < n_points; ++y)
438 for (
unsigned int z = 0; z < n_points; ++z)
439 my_patches[
angle].data(dataset + offset,
440 x * n_points * n_points +
441 y * n_points + z) = value;
454template <
int dim,
typename DoFHandlerType>
457 const unsigned int n_patches_per_circle,
458 const unsigned int nnnn_subdivisions)
466 const unsigned int n_subdivisions =
467 (nnnn_subdivisions != 0) ? nnnn_subdivisions : this->default_subdivisions;
468 Assert(n_subdivisions >= 1,
472 this->validate_dataset_names();
474 unsigned int n_datasets = this->cell_data.size();
475 for (
unsigned int i = 0; i < this->dof_data.size(); ++i)
476 n_datasets += this->dof_data[i]->n_output_variables;
479 for (
unsigned int i = 0; i < this->dof_data.size(); ++i)
480 if (this->dof_data[i]->postprocessor)
482 this->dof_data[i]->postprocessor->get_needed_update_flags();
487 ExcMessage(
"The update of normal vectors may not be requested for "
488 "evaluation of data on cells via DataPostprocessor."));
493 std::vector<cell_iterator> all_cells;
495 cell = next_cell(cell))
496 all_cells.push_back(cell);
503 this->patches.clear();
504 this->patches.reserve(all_cells.size() * n_patches_per_circle);
507 std::vector<unsigned int> n_postprocessor_outputs(this->dof_data.size());
508 for (
unsigned int dataset = 0; dataset < this->dof_data.size(); ++dataset)
509 if (this->dof_data[dataset]->postprocessor)
510 n_postprocessor_outputs[dataset] =
511 this->dof_data[dataset]->n_output_variables;
513 n_postprocessor_outputs[dataset] = 0;
517 thread_data(n_datasets,
519 n_patches_per_circle,
520 n_postprocessor_outputs,
524 std::vector<DataOutBase::Patch<dimension + 1, space_dimension + 1>>
525 new_patches(n_patches_per_circle);
526 for (
unsigned int i = 0; i < new_patches.size(); ++i)
528 new_patches[i].n_subdivisions = n_subdivisions;
529 new_patches[i].data.reinit(
530 n_datasets, Utilities::fixed_power<dimension + 1>(n_subdivisions + 1));
536 all_cells.data() + all_cells.size(),
539 ParallelData<dimension, space_dimension> &data,
541 &my_patches) { this->build_one_patch(cell, data, my_patches); },
545 internal::DataOutRotationImplementation::
546 append_patch_to_list<dimension, space_dimension>(new_patches,
555template <
int dim,
typename DoFHandlerType>
563template <
int dim,
typename DoFHandlerType>
579#include "data_out_rotation.inst"
void build_one_patch(const cell_iterator *cell, internal::DataOutRotationImplementation::ParallelData< dimension, space_dimension > &data, std::vector< DataOutBase::Patch< dimension+1, space_dimension+1 > > &my_patches)
virtual cell_iterator next_cell(const cell_iterator &cell)
virtual void build_patches(const unsigned int n_patches_per_circle, const unsigned int n_subdivisions=0)
virtual cell_iterator first_cell()
typename DataOut_DoFData< DoFHandlerType, dimension+1 >::cell_iterator cell_iterator
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 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)
static ::ExceptionBase & ExcNoTriangulationSelected()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
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 std::vector< DataOutBase::Patch< dim+1, spacedim+1 > > &new_patches, std::vector< DataOutBase::Patch< dim+1, spacedim+1 > > &patches)
static constexpr double PI
::VectorizedArray< Number, width > cos(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sin(const ::VectorizedArray< Number, width > &)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
ParallelData(const unsigned int n_datasets, const unsigned int n_subdivisions, const unsigned int n_patches_per_circle, 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)