16 #include <deal.II/base/quadrature_lib.h> 17 #include <deal.II/base/work_stream.h> 18 #include <deal.II/lac/vector.h> 19 #include <deal.II/lac/block_vector.h> 20 #include <deal.II/grid/tria.h> 21 #include <deal.II/dofs/dof_handler.h> 22 #include <deal.II/dofs/dof_accessor.h> 23 #include <deal.II/grid/tria_iterator.h> 24 #include <deal.II/fe/fe.h> 25 #include <deal.II/fe/fe_values.h> 26 #include <deal.II/hp/fe_values.h> 27 #include <deal.II/fe/mapping_q1.h> 28 #include <deal.II/numerics/data_out_rotation.h> 32 DEAL_II_NAMESPACE_OPEN
50 template <
int dim,
int spacedim>
51 ParallelData<dim,spacedim>::
52 ParallelData (
const unsigned int n_datasets,
53 const unsigned int n_subdivisions,
54 const unsigned int n_patches_per_circle,
55 const std::vector<unsigned int> &n_postprocessor_outputs,
61 ParallelDataBase<dim,spacedim> (n_datasets,
63 n_postprocessor_outputs,
68 n_patches_per_circle (n_patches_per_circle)
77 template <
int dim,
int spacedim>
82 for (
unsigned int i=0; i<new_patches.size(); ++i)
84 patches.push_back (new_patches[i]);
85 patches.back().patch_index = patches.size()-1;
93 template <
int dim,
typename DoFHandlerType>
108 Assert ((*cell)->is_locally_owned(),
111 const unsigned int n_patches_per_circle = data.n_patches_per_circle;
114 const unsigned int n_points = data.n_subdivisions+1;
120 std::vector<Point<dimension+1> > angle_directions (n_patches_per_circle+1);
121 for (
unsigned int i=0; i<=n_patches_per_circle; ++i)
123 angle_directions[i][dimension-1] = std::cos(2*
numbers::PI *
124 i/n_patches_per_circle);
125 angle_directions[i][dimension] = std::sin(2*
numbers::PI *
126 i/n_patches_per_circle);
129 for (
unsigned int angle=0; angle<n_patches_per_circle; ++angle)
138 const double r1 = (*cell)->vertex(0)(0),
139 r2 = (*cell)->vertex(1)(0);
140 Assert (r1 >= 0, ExcRadialVariableHasNegativeValues(r1));
141 Assert (r2 >= 0, ExcRadialVariableHasNegativeValues(r2));
143 my_patches[angle].vertices[0] = r1*angle_directions[angle];
144 my_patches[angle].vertices[1] = r2*angle_directions[angle];
145 my_patches[angle].vertices[2] = r1*angle_directions[angle+1];
146 my_patches[angle].vertices[3] = r2*angle_directions[angle+1];
153 for (
unsigned int vertex=0;
154 vertex<GeometryInfo<dimension>::vertices_per_cell;
161 Assert (v(0) >= 0, ExcRadialVariableHasNegativeValues(v(0)));
164 my_patches[angle].vertices[vertex] = v(0) * angle_directions[angle];
165 my_patches[angle].vertices[vertex][0] = v(1);
168 = v(0) * angle_directions[angle+1];
181 if (data.n_datasets > 0)
183 unsigned int offset=0;
185 data.reinit_all_fe_values(this->dof_data, *cell);
187 for (
unsigned int dataset=0; dataset<this->dof_data.size(); ++dataset)
190 = data.get_present_fe_values(dataset);
192 = fe_patch_values.
get_fe().n_components();
194 if (postprocessor != 0)
207 this->dof_data[dataset]->get_function_values (fe_patch_values,
208 data.patch_values_scalar.solution_values);
210 this->dof_data[dataset]->get_function_gradients (fe_patch_values,
211 data.patch_values_scalar.solution_gradients);
213 this->dof_data[dataset]->get_function_hessians (fe_patch_values,
214 data.patch_values_scalar.solution_hessians);
219 const typename DoFHandlerType::active_cell_iterator dh_cell(&(*cell)->get_triangulation(),
222 this->dof_data[dataset]->dof_handler);
223 data.patch_values_scalar.template set_cell<DoFHandlerType> (dh_cell);
226 evaluate_scalar_field(data.patch_values_scalar,
227 data.postprocessed_values[dataset]);
236 this->dof_data[dataset]->get_function_values (fe_patch_values,
237 data.patch_values_system.solution_values);
239 this->dof_data[dataset]->get_function_gradients (fe_patch_values,
240 data.patch_values_system.solution_gradients);
242 this->dof_data[dataset]->get_function_hessians (fe_patch_values,
243 data.patch_values_system.solution_hessians);
248 const typename DoFHandlerType::active_cell_iterator dh_cell(&(*cell)->get_triangulation(),
251 this->dof_data[dataset]->dof_handler);
252 data.patch_values_system.template set_cell<DoFHandlerType> (dh_cell);
255 evaluate_vector_field(data.patch_values_system,
256 data.postprocessed_values[dataset]);
259 for (
unsigned int component=0;
260 component<this->dof_data[dataset]->n_output_variables;
266 for (
unsigned int x=0; x<n_points; ++x)
267 for (
unsigned int y=0; y<n_points; ++y)
268 my_patches[angle].data(offset+component,
270 = data.postprocessed_values[dataset][x](component);
274 for (
unsigned int x=0; x<n_points; ++x)
275 for (
unsigned int y=0; y<n_points; ++y)
276 for (
unsigned int z=0; z<n_points; ++z)
277 my_patches[angle].data(offset+component,
278 x*n_points*n_points +
281 = data.postprocessed_values[dataset][x*n_points+z](component);
291 this->dof_data[dataset]->get_function_values (fe_patch_values,
292 data.patch_values_scalar.solution_values);
297 for (
unsigned int x=0; x<n_points; ++x)
298 for (
unsigned int y=0; y<n_points; ++y)
299 my_patches[angle].data(offset,
301 = data.patch_values_scalar.solution_values[x];
305 for (
unsigned int x=0; x<n_points; ++x)
306 for (
unsigned int y=0; y<n_points; ++y)
307 for (
unsigned int z=0; z<n_points; ++z)
308 my_patches[angle].data(offset,
309 x*n_points*n_points +
312 = data.patch_values_scalar.solution_values[x*n_points+z];
323 this->dof_data[dataset]->get_function_values (fe_patch_values,
324 data.patch_values_system.solution_values);
332 for (
unsigned int x=0; x<n_points; ++x)
333 for (
unsigned int y=0; y<n_points; ++y)
334 my_patches[angle].data(offset+component,
336 = data.patch_values_system.solution_values[x](component);
340 for (
unsigned int x=0; x<n_points; ++x)
341 for (
unsigned int y=0; y<n_points; ++y)
342 for (
unsigned int z=0; z<n_points; ++z)
343 my_patches[angle].data(offset+component,
344 x*n_points*n_points +
347 = data.patch_values_system.solution_values[x*n_points+z](component);
355 offset+=this->dof_data[dataset]->n_output_variables;
359 for (
unsigned int dataset=0; dataset<this->cell_data.size(); ++dataset)
364 Assert ((*cell)->active(),
365 ExcMessage(
"Cell must be active for cell data"));
366 const unsigned int cell_number
367 = std::distance (this->triangulation->begin_active(),
370 = this->cell_data[dataset]->get_cell_data_value (cell_number);
374 for (
unsigned int x=0; x<n_points; ++x)
375 for (
unsigned int y=0; y<n_points; ++y)
376 my_patches[angle].data(dataset+offset,
383 for (
unsigned int x=0; x<n_points; ++x)
384 for (
unsigned int y=0; y<n_points; ++y)
385 for (
unsigned int z=0; z<n_points; ++z)
386 my_patches[angle].data(dataset+offset,
387 x*n_points*n_points +
403 template <
int dim,
typename DoFHandlerType>
405 const unsigned int nnnn_subdivisions)
410 Assert (this->triangulation != 0,
413 const unsigned int n_subdivisions = (nnnn_subdivisions != 0)
415 : this->default_subdivisions;
416 Assert (n_subdivisions >= 1,
419 this->validate_dataset_names();
421 unsigned int n_datasets=this->cell_data.size();
422 for (
unsigned int i=0; i<this->dof_data.size(); ++i)
423 n_datasets+= this->dof_data[i]->n_output_variables;
426 for (
unsigned int i=0; i<this->dof_data.size(); ++i)
427 if (this->dof_data[i]->postprocessor)
428 update_flags |= this->dof_data[i]->postprocessor->get_needed_update_flags();
433 ExcMessage(
"The update of normal vectors may not be requested for " 434 "evaluation of data on cells via DataPostprocessor."));
439 std::vector<cell_iterator> all_cells;
440 for (
cell_iterator cell=first_cell(); cell != this->triangulation->end();
441 cell = next_cell(cell))
442 all_cells.push_back (cell);
449 this->patches.clear();
450 this->patches.reserve (all_cells.size() * n_patches_per_circle);
453 std::vector<unsigned int> n_postprocessor_outputs (this->dof_data.size());
454 for (
unsigned int dataset=0; dataset<this->dof_data.size(); ++dataset)
455 if (this->dof_data[dataset]->postprocessor)
456 n_postprocessor_outputs[dataset] = this->dof_data[dataset]->n_output_variables;
458 n_postprocessor_outputs[dataset] = 0;
461 thread_data (n_datasets,
462 n_subdivisions, n_patches_per_circle,
463 n_postprocessor_outputs,
465 this->get_finite_elements(),
467 std::vector<DataOutBase::Patch<dimension+1,space_dimension+1> >
468 new_patches (n_patches_per_circle);
469 for (
unsigned int i=0; i<new_patches.size(); ++i)
471 new_patches[i].n_subdivisions = n_subdivisions;
472 new_patches[i].data.reinit (n_datasets,
473 Utilities::fixed_power<dimension+1>(n_subdivisions+1));
478 &all_cells[0]+all_cells.size(),
480 this, std_cxx11::_1, std_cxx11::_2, std_cxx11::_3),
482 ::append_patch_to_list<dim,space_dimension>,
483 std_cxx11::_1, std_cxx11::ref(this->patches)),
490 template <
int dim,
typename DoFHandlerType>
494 return this->triangulation->begin_active ();
498 template <
int dim,
typename DoFHandlerType>
513 #include "data_out_rotation.inst" 516 DEAL_II_NAMESPACE_CLOSE
TriaActiveIterator< CellAccessor< dim, spacedim > > active_cell_iterator
static ::ExceptionBase & ExcInvalidNumberOfSubdivisions(int arg1)
virtual cell_iterator next_cell(const cell_iterator &cell)
const FiniteElement< dim, spacedim > & get_fe() const
Transformed quadrature points.
const std::vector< Point< spacedim > > & get_quadrature_points() const
void build_one_patch(const cell_iterator *cell, internal::DataOutRotation::ParallelData< dimension, space_dimension > &data, std::vector< DataOutBase::Patch< dimension+1, space_dimension+1 > > &my_patches)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
Second derivatives of shape functions.
virtual cell_iterator first_cell()
virtual void build_patches(const unsigned int n_patches_per_circle, const unsigned int n_subdivisions=0)
DataOut_DoFData< DoFHandlerType, dimension+1 >::cell_iterator cell_iterator
Shape function gradients.
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)
static ::ExceptionBase & ExcNotImplemented()
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
static ::ExceptionBase & ExcNoTriangulationSelected()
virtual UpdateFlags get_needed_update_flags() const =0