16 #include <deal.II/base/work_stream.h> 17 #include <deal.II/numerics/data_out.h> 18 #include <deal.II/grid/tria.h> 19 #include <deal.II/dofs/dof_handler.h> 20 #include <deal.II/dofs/dof_accessor.h> 21 #include <deal.II/hp/dof_handler.h> 22 #include <deal.II/grid/tria_iterator.h> 23 #include <deal.II/fe/fe.h> 24 #include <deal.II/fe/fe_dgq.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> 31 DEAL_II_NAMESPACE_OPEN
38 template <
int dim,
int spacedim>
39 ParallelData<dim,spacedim>::
40 ParallelData (
const unsigned int n_datasets,
41 const unsigned int n_subdivisions,
42 const std::vector<unsigned int> &n_postprocessor_outputs,
46 const std::vector<std::vector<unsigned int> > &cell_to_patch_index_map)
48 ParallelDataBase<dim,spacedim> (n_datasets,
50 n_postprocessor_outputs,
55 cell_to_patch_index_map (&cell_to_patch_index_map)
62 template <
int dim,
typename DoFHandlerType>
66 (
const std::pair<cell_iterator, unsigned int> *cell_and_index,
68 const unsigned int n_subdivisions,
80 for (
unsigned int vertex=0; vertex<GeometryInfo<DoFHandlerType::dimension>::vertices_per_cell; ++vertex)
81 if (scratch_data.mapping_collection[0].preserves_vertex_locations())
82 patch.
vertices[vertex] = cell_and_index->first->vertex(vertex);
84 patch.
vertices[vertex] = scratch_data.mapping_collection[0].transform_unit_to_real_cell
85 (cell_and_index->first,
88 if (scratch_data.n_datasets > 0)
91 scratch_data.reinit_all_fe_values(this->dof_data, cell_and_index->first);
94 = scratch_data.get_present_fe_values (0);
105 if (curved_cell_region==curved_inner_cells
107 (curved_cell_region==curved_boundary
109 (cell_and_index->first->at_boundary()
111 (DoFHandlerType::dimension != DoFHandlerType::space_dimension))))
113 Assert(patch.
space_dim==DoFHandlerType::space_dimension, ExcInternalError());
114 const std::vector<Point<DoFHandlerType::space_dimension> > &q_points=fe_patch_values.
get_quadrature_points();
117 patch.
data.
reinit (scratch_data.n_datasets+DoFHandlerType::space_dimension, n_q_points);
122 for (
unsigned int i=0; i<DoFHandlerType::space_dimension; ++i)
123 for (
unsigned int q=0; q<n_q_points; ++q)
124 patch.
data(patch.
data.
size(0)-DoFHandlerType::space_dimension+i,q)=q_points[q][i];
128 patch.
data.
reinit(scratch_data.n_datasets, n_q_points);
134 unsigned int offset=0;
137 for (
unsigned int dataset=0; dataset<this->dof_data.size(); ++dataset)
140 = scratch_data.get_present_fe_values (dataset);
141 const unsigned int n_components =
142 this_fe_patch_values.
get_fe().n_components();
146 if (postprocessor != 0)
151 if (n_components == 1)
156 this->dof_data[dataset]->get_function_values (this_fe_patch_values,
157 scratch_data.patch_values);
159 this->dof_data[dataset]->get_function_gradients (this_fe_patch_values,
160 scratch_data.patch_gradients);
162 this->dof_data[dataset]->get_function_hessians (this_fe_patch_values,
163 scratch_data.patch_hessians);
169 std::vector<Point<DoFHandlerType::space_dimension> > dummy_normals;
171 compute_derived_quantities_scalar(scratch_data.patch_values,
172 scratch_data.patch_gradients,
173 scratch_data.patch_hessians,
175 scratch_data.patch_evaluation_points,
176 scratch_data.postprocessed_values[dataset]);
180 scratch_data.resize_system_vectors (n_components);
185 this->dof_data[dataset]->get_function_values (this_fe_patch_values,
186 scratch_data.patch_values_system);
188 this->dof_data[dataset]->get_function_gradients (this_fe_patch_values,
189 scratch_data.patch_gradients_system);
191 this->dof_data[dataset]->get_function_hessians (this_fe_patch_values,
192 scratch_data.patch_hessians_system);
197 std::vector<Point<DoFHandlerType::space_dimension> > dummy_normals;
200 compute_derived_quantities_vector(scratch_data.patch_values_system,
201 scratch_data.patch_gradients_system,
202 scratch_data.patch_hessians_system,
204 scratch_data.patch_evaluation_points,
205 scratch_data.postprocessed_values[dataset]);
208 for (
unsigned int q=0; q<n_q_points; ++q)
209 for (
unsigned int component=0;
210 component<this->dof_data[dataset]->n_output_variables;
212 patch.
data(offset+component,q)
213 = scratch_data.postprocessed_values[dataset][q](component);
219 if (n_components == 1)
221 this->dof_data[dataset]->get_function_values (this_fe_patch_values,
222 scratch_data.patch_values);
223 for (
unsigned int q=0; q<n_q_points; ++q)
224 patch.
data(offset,q) = scratch_data.patch_values[q];
228 scratch_data.resize_system_vectors(n_components);
229 this->dof_data[dataset]->get_function_values (this_fe_patch_values,
230 scratch_data.patch_values_system);
231 for (
unsigned int component=0; component<n_components;
233 for (
unsigned int q=0; q<n_q_points; ++q)
234 patch.
data(offset+component,q) =
235 scratch_data.patch_values_system[q](component);
238 offset+=this->dof_data[dataset]->n_output_variables;
244 if (this->cell_data.size() != 0)
246 Assert (!cell_and_index->first->has_children(), ExcNotImplemented());
248 for (
unsigned int dataset=0; dataset<this->cell_data.size(); ++dataset)
251 = this->cell_data[dataset]->get_cell_data_value (cell_and_index->second);
252 for (
unsigned int q=0; q<n_q_points; ++q)
253 patch.
data(offset+dataset,q) = value;
259 for (
unsigned int f=0; f<GeometryInfo<DoFHandlerType::dimension>::faces_per_cell; ++f)
271 if (cell_and_index->first->at_boundary(f)
273 (cell_and_index->first->neighbor(f)->level() != cell_and_index->first->level()))
279 const cell_iterator neighbor = cell_and_index->first->neighbor(f);
280 Assert (static_cast<unsigned int>(neighbor->level()) <
281 scratch_data.cell_to_patch_index_map->size(),
283 if ((static_cast<unsigned int>(neighbor->index()) >=
284 (*scratch_data.cell_to_patch_index_map)[neighbor->level()].size())
286 ((*scratch_data.cell_to_patch_index_map)[neighbor->level()][neighbor->index()]
297 = (*scratch_data.cell_to_patch_index_map)[neighbor->level()][neighbor->index()];
300 const unsigned int patch_idx =
301 (*scratch_data.cell_to_patch_index_map)[cell_and_index->first->level()][cell_and_index->first->index()];
303 Assert(patch_idx < patches.size(), ExcInternalError());
309 patches[patch_idx].swap (patch);
314 template <
int dim,
typename DoFHandlerType>
318 n_subdivisions, no_curved_cells);
323 template <
int dim,
typename DoFHandlerType>
326 const unsigned int n_subdivisions_,
330 Assert (dim==DoFHandlerType::dimension, ExcDimensionMismatch(dim, DoFHandlerType::dimension));
332 Assert (this->triangulation != 0,
333 Exceptions::DataOut::ExcNoTriangulationSelected());
335 const unsigned int n_subdivisions = (n_subdivisions_ != 0)
337 : this->default_subdivisions;
338 Assert (n_subdivisions >= 1,
339 Exceptions::DataOut::ExcInvalidNumberOfSubdivisions(n_subdivisions));
354 std::vector<std::vector<unsigned int> > cell_to_patch_index_map;
355 cell_to_patch_index_map.resize (this->triangulation->n_levels());
356 for (
unsigned int l=0; l<this->triangulation->n_levels(); ++l)
359 unsigned int max_index = 0;
360 for (
cell_iterator cell=first_locally_owned_cell(); cell != this->triangulation->end();
361 cell = next_locally_owned_cell(cell))
362 if (static_cast<unsigned int>(cell->level()) == l)
363 max_index = std::max (max_index,
364 static_cast<unsigned int>(cell->index()));
366 cell_to_patch_index_map[l].resize (max_index+1,
368 DoFHandlerType::space_dimension>::no_neighbor);
372 std::vector<std::pair<cell_iterator, unsigned int> > all_cells;
379 active_cell_iterator active_cell = this->triangulation->begin_active();
380 unsigned int active_index = 0;
382 for (; cell != this->triangulation->end();
383 cell = next_locally_owned_cell(cell))
387 while (active_cell!=this->triangulation->end()
389 && active_cell_iterator(cell) != active_cell)
395 Assert (static_cast<unsigned int>(cell->level()) <
396 cell_to_patch_index_map.size(),
398 Assert (static_cast<unsigned int>(cell->index()) <
399 cell_to_patch_index_map[cell->level()].size(),
401 Assert (active_index < this->triangulation->n_active_cells(),
403 cell_to_patch_index_map[cell->level()][cell->index()] = all_cells.size();
405 all_cells.push_back (std::make_pair(cell, active_index));
409 this->patches.clear ();
410 this->patches.resize(all_cells.size());
413 unsigned int n_datasets=this->cell_data.size();
414 for (
unsigned int i=0; i<this->dof_data.size(); ++i)
415 n_datasets += this->dof_data[i]->n_output_variables;
417 std::vector<unsigned int> n_postprocessor_outputs (this->dof_data.size());
418 for (
unsigned int dataset=0; dataset<this->dof_data.size(); ++dataset)
419 if (this->dof_data[dataset]->postprocessor)
420 n_postprocessor_outputs[dataset] = this->dof_data[dataset]->n_output_variables;
422 n_postprocessor_outputs[dataset] = 0;
425 = (n_subdivisions<2 ? no_curved_cells : curved_region);
428 if (curved_cell_region != no_curved_cells)
431 for (
unsigned int i=0; i<this->dof_data.size(); ++i)
432 if (this->dof_data[i]->postprocessor)
433 update_flags |= this->dof_data[i]->postprocessor->get_needed_update_flags();
437 ExcMessage(
"The update of normal vectors may not be requested for evaluation of " 438 "data on cells via DataPostprocessor."));
441 thread_data (n_datasets, n_subdivisions,
442 n_postprocessor_outputs,
444 this->get_finite_elements(),
446 cell_to_patch_index_map);
449 if (all_cells.size() > 0)
451 &all_cells[0]+all_cells.size(),
461 std_cxx11::ref(this->patches)),
463 std_cxx11::function<void (const int &)>(),
478 template <
int dim,
typename DoFHandlerType>
482 return this->triangulation->begin_active ();
487 template <
int dim,
typename DoFHandlerType>
495 active_cell_iterator active_cell = cell;
502 template <
int dim,
typename DoFHandlerType>
511 while ((cell != this->triangulation->end()) &&
512 (cell->has_children() ==
false) &&
513 !cell->is_locally_owned())
514 cell = next_cell(cell);
521 template <
int dim,
typename DoFHandlerType>
527 cell = next_cell(old_cell);
528 while ((cell != this->triangulation->end()) &&
529 (cell->has_children() ==
false) &&
530 !cell->is_locally_owned())
531 cell = next_cell(cell);
537 #include "data_out.inst" 539 DEAL_II_NAMESPACE_CLOSE
static const unsigned int invalid_unsigned_int
void build_one_patch(const std::pair< cell_iterator, unsigned int > *cell_and_index, internal::DataOut::ParallelData< DoFHandlerType::dimension, DoFHandlerType::space_dimension > &scratch_data, const unsigned int n_subdivisions, const CurvedCellRegion curved_cell_region, std::vector< DataOutBase::Patch< DoFHandlerType::dimension, DoFHandlerType::space_dimension > > &patches)
::ExceptionBase & ExcMessage(std::string arg1)
unsigned int neighbors[dim > 0 ? GeometryInfo< dim >::faces_per_cell :1]
virtual cell_iterator next_cell(const cell_iterator &cell)
const FiniteElement< dim, spacedim > & get_fe() const
Transformed quadrature points.
Point< spacedim > vertices[GeometryInfo< dim >::vertices_per_cell]
virtual void build_patches(const unsigned int n_subdivisions=0)
const std::vector< Point< spacedim > > & get_quadrature_points() const
virtual cell_iterator first_locally_owned_cell()
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
#define Assert(cond, exc)
Abstract base class for mapping classes.
virtual cell_iterator first_cell()
unsigned int n_subdivisions
DataOut_DoFData< DoFHandler< dim >, DoFHandler< dim > ::dimension, DoFHandler< dim > ::space_dimension >::cell_iterator cell_iterator
Second derivatives of shape functions.
const unsigned int n_quadrature_points
static const unsigned int space_dim
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 unsigned int n_threads()
unsigned int size(const unsigned int i) const
bool points_are_available
virtual cell_iterator next_locally_owned_cell(const cell_iterator &cell)
virtual UpdateFlags get_needed_update_flags() const =0