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,
78 for (
unsigned int vertex=0; vertex<GeometryInfo<DoFHandlerType::dimension>::vertices_per_cell; ++vertex)
79 if (scratch_data.mapping_collection[0].preserves_vertex_locations())
80 patch.
vertices[vertex] = cell_and_index->first->vertex(vertex);
82 patch.
vertices[vertex] = scratch_data.mapping_collection[0].transform_unit_to_real_cell
83 (cell_and_index->first,
87 scratch_data.reinit_all_fe_values(this->dof_data, cell_and_index->first);
90 = scratch_data.get_present_fe_values (0);
101 if (curved_cell_region==curved_inner_cells
103 (curved_cell_region==curved_boundary
105 (cell_and_index->first->at_boundary()
107 (DoFHandlerType::dimension != DoFHandlerType::space_dimension))))
117 const std::vector<Point<DoFHandlerType::space_dimension> > &q_points
120 patch.
data.
reinit (scratch_data.n_datasets+DoFHandlerType::space_dimension, n_q_points);
121 for (
unsigned int i=0; i<DoFHandlerType::space_dimension; ++i)
122 for (
unsigned int q=0; q<n_q_points; ++q)
123 patch.
data(patch.
data.
size(0)-DoFHandlerType::space_dimension+i,q) = q_points[q][i];
127 patch.
data.
reinit(scratch_data.n_datasets, n_q_points);
132 if (scratch_data.n_datasets > 0)
135 unsigned int offset=0;
138 for (
unsigned int dataset=0; dataset<this->dof_data.size(); ++dataset)
141 = scratch_data.get_present_fe_values (dataset);
143 this_fe_patch_values.
get_fe().n_components();
146 = this->dof_data[dataset]->postprocessor;
148 if (postprocessor != 0)
158 this->dof_data[dataset]->get_function_values (this_fe_patch_values,
159 scratch_data.patch_values_scalar.solution_values);
161 this->dof_data[dataset]->get_function_gradients (this_fe_patch_values,
162 scratch_data.patch_values_scalar.solution_gradients);
164 this->dof_data[dataset]->get_function_hessians (this_fe_patch_values,
165 scratch_data.patch_values_scalar.solution_hessians);
170 const typename DoFHandlerType::active_cell_iterator dh_cell(&cell_and_index->first->get_triangulation(),
171 cell_and_index->first->level(),
172 cell_and_index->first->index(),
173 this->dof_data[dataset]->dof_handler);
174 scratch_data.patch_values_scalar.template set_cell<DoFHandlerType> (dh_cell);
177 evaluate_scalar_field(scratch_data.patch_values_scalar,
178 scratch_data.postprocessed_values[dataset]);
187 this->dof_data[dataset]->get_function_values (this_fe_patch_values,
188 scratch_data.patch_values_system.solution_values);
190 this->dof_data[dataset]->get_function_gradients (this_fe_patch_values,
191 scratch_data.patch_values_system.solution_gradients);
193 this->dof_data[dataset]->get_function_hessians (this_fe_patch_values,
194 scratch_data.patch_values_system.solution_hessians);
199 const typename DoFHandlerType::active_cell_iterator dh_cell(&cell_and_index->first->get_triangulation(),
200 cell_and_index->first->level(),
201 cell_and_index->first->index(),
202 this->dof_data[dataset]->dof_handler);
203 scratch_data.patch_values_system.template set_cell<DoFHandlerType> (dh_cell);
206 evaluate_vector_field(scratch_data.patch_values_system,
207 scratch_data.postprocessed_values[dataset]);
210 for (
unsigned int q=0; q<n_q_points; ++q)
211 for (
unsigned int component=0;
212 component<this->dof_data[dataset]->n_output_variables;
214 patch.
data(offset+component,q)
215 = scratch_data.postprocessed_values[dataset][q](component);
223 this->dof_data[dataset]->get_function_values (this_fe_patch_values,
224 scratch_data.patch_values_scalar.solution_values);
225 for (
unsigned int q=0; q<n_q_points; ++q)
226 patch.
data(offset,q) = scratch_data.patch_values_scalar.solution_values[q];
231 this->dof_data[dataset]->get_function_values (this_fe_patch_values,
232 scratch_data.patch_values_system.solution_values);
235 for (
unsigned int q=0; q<n_q_points; ++q)
236 patch.
data(offset+component,q) =
237 scratch_data.patch_values_system.solution_values[q](component);
240 offset+=this->dof_data[dataset]->n_output_variables;
246 if (this->cell_data.size() != 0)
250 for (
unsigned int dataset=0; dataset<this->cell_data.size(); ++dataset)
253 = this->cell_data[dataset]->get_cell_data_value (cell_and_index->second);
254 for (
unsigned int q=0; q<n_q_points; ++q)
255 patch.
data(offset+dataset,q) = value;
261 for (
unsigned int f=0; f<GeometryInfo<DoFHandlerType::dimension>::faces_per_cell; ++f)
273 if (cell_and_index->first->at_boundary(f)
275 (cell_and_index->first->neighbor(f)->level() != cell_and_index->first->level()))
281 const cell_iterator neighbor = cell_and_index->first->neighbor(f);
282 Assert (static_cast<unsigned int>(neighbor->level()) <
283 scratch_data.cell_to_patch_index_map->size(),
285 if ((static_cast<unsigned int>(neighbor->index()) >=
286 (*scratch_data.cell_to_patch_index_map)[neighbor->level()].size())
288 ((*scratch_data.cell_to_patch_index_map)[neighbor->level()][neighbor->index()]
299 = (*scratch_data.cell_to_patch_index_map)[neighbor->level()][neighbor->index()];
302 const unsigned int patch_idx =
303 (*scratch_data.cell_to_patch_index_map)[cell_and_index->first->level()][cell_and_index->first->index()];
311 this->patches[patch_idx].swap (patch);
316 template <
int dim,
typename DoFHandlerType>
320 n_subdivisions, no_curved_cells);
325 template <
int dim,
typename DoFHandlerType>
328 const unsigned int n_subdivisions_,
334 Assert (this->triangulation != 0,
337 const unsigned int n_subdivisions = (n_subdivisions_ != 0)
339 : this->default_subdivisions;
340 Assert (n_subdivisions >= 1,
343 this->validate_dataset_names();
358 std::vector<std::vector<unsigned int> > cell_to_patch_index_map;
359 cell_to_patch_index_map.resize (this->triangulation->n_levels());
360 for (
unsigned int l=0; l<this->triangulation->n_levels(); ++l)
363 unsigned int max_index = 0;
364 for (
cell_iterator cell=first_locally_owned_cell(); cell != this->triangulation->end();
365 cell = next_locally_owned_cell(cell))
366 if (static_cast<unsigned int>(cell->level()) == l)
367 max_index = std::max (max_index,
368 static_cast<unsigned int>(cell->index()));
370 cell_to_patch_index_map[l].resize (max_index+1,
372 DoFHandlerType::space_dimension>::no_neighbor);
376 std::vector<std::pair<cell_iterator, unsigned int> > all_cells;
383 active_cell_iterator active_cell = this->triangulation->begin_active();
384 unsigned int active_index = 0;
386 for (; cell != this->triangulation->end();
387 cell = next_locally_owned_cell(cell))
391 while (active_cell!=this->triangulation->end()
393 && active_cell_iterator(cell) != active_cell)
399 Assert (static_cast<unsigned int>(cell->level()) <
400 cell_to_patch_index_map.size(),
402 Assert (static_cast<unsigned int>(cell->index()) <
403 cell_to_patch_index_map[cell->level()].size(),
405 Assert (active_index < this->triangulation->n_active_cells(),
407 cell_to_patch_index_map[cell->level()][cell->index()] = all_cells.size();
409 all_cells.push_back (std::make_pair(cell, active_index));
413 this->patches.clear ();
414 this->patches.resize(all_cells.size());
417 unsigned int n_datasets=this->cell_data.size();
418 for (
unsigned int i=0; i<this->dof_data.size(); ++i)
419 n_datasets += this->dof_data[i]->n_output_variables;
421 std::vector<unsigned int> n_postprocessor_outputs (this->dof_data.size());
422 for (
unsigned int dataset=0; dataset<this->dof_data.size(); ++dataset)
423 if (this->dof_data[dataset]->postprocessor)
424 n_postprocessor_outputs[dataset] = this->dof_data[dataset]->n_output_variables;
426 n_postprocessor_outputs[dataset] = 0;
429 = (n_subdivisions<2 ? no_curved_cells : curved_region);
432 if (curved_cell_region != no_curved_cells)
435 for (
unsigned int i=0; i<this->dof_data.size(); ++i)
436 if (this->dof_data[i]->postprocessor)
437 update_flags |= this->dof_data[i]->postprocessor->get_needed_update_flags();
441 ExcMessage(
"The update of normal vectors may not be requested for evaluation of " 442 "data on cells via DataPostprocessor."));
445 thread_data (n_datasets, n_subdivisions,
446 n_postprocessor_outputs,
448 this->get_finite_elements(),
450 cell_to_patch_index_map);
453 if (all_cells.size() > 0)
455 &all_cells[0]+all_cells.size(),
466 std_cxx11::function<void (const int &)>(),
481 template <
int dim,
typename DoFHandlerType>
485 return this->triangulation->begin_active ();
490 template <
int dim,
typename DoFHandlerType>
498 active_cell_iterator active_cell = cell;
505 template <
int dim,
typename DoFHandlerType>
514 while ((cell != this->triangulation->end()) &&
515 (cell->has_children() ==
false) &&
516 !cell->is_locally_owned())
517 cell = next_cell(cell);
524 template <
int dim,
typename DoFHandlerType>
530 cell = next_cell(old_cell);
531 while ((cell != this->triangulation->end()) &&
532 (cell->has_children() ==
false) &&
533 !cell->is_locally_owned())
534 cell = next_cell(cell);
540 #include "data_out.inst" 542 DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcInvalidNumberOfSubdivisions(int arg1)
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)
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
static ::ExceptionBase & ExcMessage(std::string arg1)
virtual cell_iterator first_locally_owned_cell()
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
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 ::ExceptionBase & ExcNotImplemented()
static unsigned int n_threads()
unsigned int size(const unsigned int i) const
bool points_are_available
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
virtual cell_iterator next_locally_owned_cell(const cell_iterator &cell)
static ::ExceptionBase & ExcNoTriangulationSelected()
static ::ExceptionBase & ExcInternalError()
virtual UpdateFlags get_needed_update_flags() const =0