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/numerics/data_out_faces.h> 21 #include <deal.II/grid/tria.h> 22 #include <deal.II/dofs/dof_handler.h> 23 #include <deal.II/dofs/dof_accessor.h> 24 #include <deal.II/grid/tria_iterator.h> 25 #include <deal.II/fe/fe.h> 26 #include <deal.II/fe/fe_values.h> 27 #include <deal.II/hp/fe_values.h> 28 #include <deal.II/fe/mapping_q1.h> 30 DEAL_II_NAMESPACE_OPEN
35 namespace DataOutFacesImplementation
37 template <
int dim,
int spacedim>
38 ParallelData<dim,spacedim>::
39 ParallelData (
const unsigned int n_datasets,
40 const unsigned int n_subdivisions,
41 const std::vector<unsigned int> &n_postprocessor_outputs,
47 ParallelDataBase<dim,spacedim> (n_datasets,
49 n_postprocessor_outputs,
62 template <
int dim,
int spacedim>
67 patches.push_back (patch);
68 patches.back().patch_index = patches.size()-1;
75 template <
int dim,
typename DoFHandlerType>
80 Assert (dim == DoFHandlerType::dimension,
86 template <
int dim,
typename DoFHandlerType>
93 Assert (cell_and_face->first->is_locally_owned(),
100 for (
unsigned int vertex=0; vertex<
GeometryInfo<dimension-1>::vertices_per_cell; ++vertex)
101 patch.
vertices[vertex] = data.mapping_collection[0].transform_unit_to_real_cell
102 (cell_and_face->first,
105 (cell_and_face->second,
107 cell_and_face->first->face_orientation(cell_and_face->second),
108 cell_and_face->first->face_flip(cell_and_face->second),
109 cell_and_face->first->face_rotation(cell_and_face->second))));
111 if (data.n_datasets > 0)
113 data.reinit_all_fe_values(this->dof_data, cell_and_face->first,
114 cell_and_face->second);
116 = data.get_present_fe_values (0);
131 for (
unsigned int i=0; i<dimension; ++i)
132 for (
unsigned int q=0; q<n_q_points; ++q)
133 patch.
data(patch.
data.
size(0)-dimension+i,q)=q_points[q][i];
136 unsigned int offset=0;
139 for (
unsigned int dataset=0; dataset<this->dof_data.size(); ++dataset)
142 = data.get_present_fe_values (dataset);
144 = this_fe_patch_values.
get_fe().n_components();
146 if (postprocessor !=
nullptr)
157 this->dof_data[dataset]->get_function_values (this_fe_patch_values,
158 internal::DataOutImplementation::ComponentExtractor::real_part,
159 data.patch_values_scalar.solution_values);
161 this->dof_data[dataset]->get_function_gradients (this_fe_patch_values,
162 internal::DataOutImplementation::ComponentExtractor::real_part,
163 data.patch_values_scalar.solution_gradients);
165 this->dof_data[dataset]->get_function_hessians (this_fe_patch_values,
166 internal::DataOutImplementation::ComponentExtractor::real_part,
167 data.patch_values_scalar.solution_hessians);
175 const typename DoFHandlerType::active_cell_iterator dh_cell(&cell_and_face->first->get_triangulation(),
176 cell_and_face->first->level(),
177 cell_and_face->first->index(),
178 this->dof_data[dataset]->dof_handler);
179 data.patch_values_scalar.template set_cell<DoFHandlerType> (dh_cell);
182 evaluate_scalar_field(data.patch_values_scalar,
183 data.postprocessed_values[dataset]);
191 this->dof_data[dataset]->get_function_values (this_fe_patch_values,
192 internal::DataOutImplementation::ComponentExtractor::real_part,
193 data.patch_values_system.solution_values);
195 this->dof_data[dataset]->get_function_gradients (this_fe_patch_values,
196 internal::DataOutImplementation::ComponentExtractor::real_part,
197 data.patch_values_system.solution_gradients);
199 this->dof_data[dataset]->get_function_hessians (this_fe_patch_values,
200 internal::DataOutImplementation::ComponentExtractor::real_part,
201 data.patch_values_system.solution_hessians);
209 const typename DoFHandlerType::active_cell_iterator dh_cell(&cell_and_face->first->get_triangulation(),
210 cell_and_face->first->level(),
211 cell_and_face->first->index(),
212 this->dof_data[dataset]->dof_handler);
213 data.patch_values_system.template set_cell<DoFHandlerType> (dh_cell);
216 evaluate_vector_field(data.patch_values_system,
217 data.postprocessed_values[dataset]);
220 for (
unsigned int q=0; q<n_q_points; ++q)
221 for (
unsigned int component=0;
222 component<this->dof_data[dataset]->n_output_variables; ++component)
223 patch.
data(offset+component,q)
224 = data.postprocessed_values[dataset][q](component);
232 this->dof_data[dataset]->get_function_values (this_fe_patch_values,
233 internal::DataOutImplementation::ComponentExtractor::real_part,
234 data.patch_values_scalar.solution_values);
235 for (
unsigned int q=0; q<n_q_points; ++q)
236 patch.
data(offset,q) = data.patch_values_scalar.solution_values[q];
241 this->dof_data[dataset]->get_function_values (this_fe_patch_values,
242 internal::DataOutImplementation::ComponentExtractor::real_part,
243 data.patch_values_system.solution_values);
246 for (
unsigned int q=0; q<n_q_points; ++q)
247 patch.
data(offset+component,q) =
248 data.patch_values_system.solution_values[q](component);
251 offset += this->dof_data[dataset]->n_output_variables;
255 for (
unsigned int dataset=0; dataset<this->cell_data.size(); ++dataset)
260 Assert (cell_and_face->first->active(),
261 ExcMessage(
"The current function is trying to generate cell-data output " 262 "for a face that does not belong to an active cell. This is " 264 const unsigned int cell_number
265 = std::distance (this->triangulation->begin_active(),
269 = this->cell_data[dataset]->get_cell_data_value (cell_number,
270 internal::DataOutImplementation::ComponentExtractor::real_part);
271 for (
unsigned int q=0; q<n_q_points; ++q)
272 patch.
data(dataset+offset,q) = value;
280 template <
int dim,
typename DoFHandlerType>
288 template <
int dim,
typename DoFHandlerType>
290 const unsigned int n_subdivisions_)
295 const unsigned int n_subdivisions = (n_subdivisions_ != 0)
297 : this->default_subdivisions;
299 Assert (n_subdivisions >= 1,
302 Assert (this->triangulation !=
nullptr,
305 this->validate_dataset_names();
307 unsigned int n_datasets = this->cell_data.size();
308 for (
unsigned int i=0; i<this->dof_data.size(); ++i)
309 n_datasets += this->dof_data[i]->n_output_variables;
315 std::vector<FaceDescriptor> all_faces;
317 ((face.first != this->triangulation->end())
320 face = next_face(face))
321 all_faces.push_back (face);
324 this->patches.clear ();
325 this->patches.reserve (all_faces.size());
329 std::vector<unsigned int> n_postprocessor_outputs (this->dof_data.size());
330 for (
unsigned int dataset=0; dataset<this->dof_data.size(); ++dataset)
331 if (this->dof_data[dataset]->postprocessor)
332 n_postprocessor_outputs[dataset] = this->dof_data[dataset]->n_output_variables;
334 n_postprocessor_outputs[dataset] = 0;
337 for (
unsigned int i=0; i<this->dof_data.size(); ++i)
338 if (this->dof_data[i]->postprocessor)
339 update_flags |= this->dof_data[i]->postprocessor->get_needed_update_flags();
343 thread_data (n_datasets,
345 n_postprocessor_outputs,
351 sample_patch.data.reinit (n_datasets,
352 Utilities::fixed_power<dimension-1>(n_subdivisions+1));
356 &all_faces[0]+all_faces.size(),
358 this, std::placeholders::_1, std::placeholders::_2, std::placeholders::_3),
360 append_patch_to_list<dim,space_dimension>,
361 std::placeholders::_1, std::ref(this->patches)),
368 template <
int dim,
typename DoFHandlerType>
374 for (; cell != this->triangulation->end(); ++cell)
375 if (cell->is_locally_owned())
376 for (
unsigned int f=0; f<GeometryInfo<dimension>::faces_per_cell; ++f)
377 if (!surface_only || cell->face(f)->at_boundary())
388 template <
int dim,
typename DoFHandlerType>
398 if (!surface_only || face.first->face(f)->at_boundary())
415 while (active_cell != this->triangulation->end())
419 if (active_cell->is_locally_owned())
420 for (
unsigned int f=0; f<GeometryInfo<dimension>::faces_per_cell; ++f)
421 if (!surface_only || active_cell->face(f)->at_boundary())
423 face.first = active_cell;
434 face.first = this->triangulation->end();
442 #include "data_out_faces.inst" 444 DEAL_II_NAMESPACE_CLOSE
TriaActiveIterator< CellAccessor< dim, spacedim > > active_cell_iterator
static ::ExceptionBase & ExcNoTriangulationSelected()
void build_one_patch(const FaceDescriptor *cell_and_face, internal::DataOutFacesImplementation::ParallelData< dimension, dimension > &data, DataOutBase::Patch< dimension-1, space_dimension > &patch)
const FiniteElement< dim, spacedim > & get_fe() const
Transformed quadrature points.
static ::ExceptionBase & ExcInvalidNumberOfSubdivisions(int arg1)
Point< spacedim > vertices[GeometryInfo< dim >::vertices_per_cell]
virtual FaceDescriptor first_face()
const std::vector< Point< spacedim > > & get_quadrature_points() const
static ::ExceptionBase & ExcMessage(std::string arg1)
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 FaceDescriptor next_face(const FaceDescriptor &face)
unsigned int n_subdivisions
Second derivatives of shape functions.
const std::vector< Tensor< 1, spacedim > > & get_all_normal_vectors() const
DataOutFaces(const bool surface_only=true)
const unsigned int n_quadrature_points
std::pair< cell_iterator, unsigned int > FaceDescriptor
size_type size(const unsigned int i) const
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()
bool points_are_available
virtual void build_patches(const unsigned int n_subdivisions=0)
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
static ::ExceptionBase & ExcInternalError()
virtual UpdateFlags get_needed_update_flags() const =0