// output preprocessor defines as is:
if (has_prefix(whole_file, "#"))
- {
- std::cout << get_substring_with_delim(whole_file, "\n") << '\n';
- skip_space(whole_file);
- continue;
- }
+ {
+ std::cout << get_substring_with_delim(whole_file, "\n") << '\n';
+ skip_space(whole_file);
+ continue;
+ }
if (!has_prefix(whole_file, "for"))
{
const unsigned int n_datasets,
const unsigned int n_subdivisions,
const std::vector<unsigned int> &n_postprocessor_outputs,
- const Mapping<dim, spacedim> & mapping,
+ const dealii::hp::MappingCollection<dim, spacedim> &mapping,
const std::vector<
std::shared_ptr<dealii::hp::FECollection<dim, spacedim>>>
& finite_elements,
* even if the mesh is internally stored in its undeformed configuration and
* the deformation is only tracked by an additional vector that holds the
* deformation of each vertex.
- *
- * @todo The @p mapping argument should be replaced by a
- * hp::MappingCollection in case of a hp::DoFHandler.
*/
virtual void
build_patches(const Mapping<DoFHandlerType::dimension,
const unsigned int n_subdivisions = 0,
const CurvedCellRegion curved_region = curved_boundary);
+ /**
+ * Same as above, but for hp::MappingCollection.
+ */
+ virtual void
+ build_patches(
+ const hp::MappingCollection<DoFHandlerType::dimension,
+ DoFHandlerType::space_dimension> &mapping,
+ const unsigned int n_subdivisions = 0,
+ const CurvedCellRegion curved_region = curved_boundary);
+
/**
* A function that allows selecting for which cells output should be
* generated. This function takes two arguments, both `std::function`
const UpdateFlags update_flags,
const bool use_face_values);
+ ParallelDataBase(
+ const unsigned int n_datasets,
+ const unsigned int n_subdivisions,
+ const std::vector<unsigned int> &n_postprocessor_outputs,
+ const dealii::hp::MappingCollection<dim, spacedim> &mapping,
+ const std::vector<
+ std::shared_ptr<dealii::hp::FECollection<dim, spacedim>>>
+ & finite_elements,
+ const UpdateFlags update_flags,
+ const bool use_face_values);
+
ParallelDataBase(const ParallelDataBase &data);
template <typename DoFHandlerType>
& finite_elements,
const UpdateFlags update_flags,
const bool use_face_values)
+ : ParallelDataBase<dim, spacedim>(
+ n_datasets,
+ n_subdivisions,
+ n_postprocessor_outputs,
+ dealii::hp::MappingCollection<dim, spacedim>(mapping),
+ finite_elements,
+ update_flags,
+ use_face_values)
+ {}
+
+
+ template <int dim, int spacedim>
+ ParallelDataBase<dim, spacedim>::ParallelDataBase(
+ const unsigned int n_datasets,
+ const unsigned int n_subdivisions,
+ const std::vector<unsigned int> &n_postprocessor_outputs,
+ const dealii::hp::MappingCollection<dim, spacedim> &mapping,
+ const std::vector<
+ std::shared_ptr<dealii::hp::FECollection<dim, spacedim>>>
+ & finite_elements,
+ const UpdateFlags update_flags,
+ const bool use_face_values)
: n_datasets(n_datasets)
, n_subdivisions(n_subdivisions)
, postprocessed_values(n_postprocessor_outputs.size())
const unsigned int n_datasets,
const unsigned int n_subdivisions,
const std::vector<unsigned int> &n_postprocessor_outputs,
- const Mapping<dim, spacedim> & mapping,
+ const dealii::hp::MappingCollection<dim, spacedim> &mapping,
const std::vector<
std::shared_ptr<dealii::hp::FECollection<dim, spacedim>>>
& finite_elements,
patch;
patch.n_subdivisions = n_subdivisions;
+ // initialize FEValues
+ scratch_data.reinit_all_fe_values(this->dof_data, cell_and_index->first);
+
+ const FEValuesBase<DoFHandlerType::dimension, DoFHandlerType::space_dimension>
+ &fe_patch_values = scratch_data.get_present_fe_values(0);
+
// set the vertices of the patch. if the mapping does not preserve locations
// (e.g. MappingQEulerian), we need to compute the offset of the vertex for
// the graphical output. Otherwise, we can just use the vertex info.
for (const unsigned int vertex :
GeometryInfo<DoFHandlerType::dimension>::vertex_indices())
- if (scratch_data.mapping_collection[0].preserves_vertex_locations())
+ if (fe_patch_values.get_mapping().preserves_vertex_locations())
patch.vertices[vertex] = cell_and_index->first->vertex(vertex);
else
patch.vertices[vertex] =
- scratch_data.mapping_collection[0].transform_unit_to_real_cell(
+ fe_patch_values.get_mapping().transform_unit_to_real_cell(
cell_and_index->first,
GeometryInfo<DoFHandlerType::dimension>::unit_cell_vertex(vertex));
- // initialize FEValues
- scratch_data.reinit_all_fe_values(this->dof_data, cell_and_index->first);
-
- const FEValuesBase<DoFHandlerType::dimension, DoFHandlerType::space_dimension>
- &fe_patch_values = scratch_data.get_present_fe_values(0);
-
const unsigned int n_q_points = fe_patch_values.n_quadrature_points;
// First fill the geometric information for the patch: Where are the
& mapping,
const unsigned int n_subdivisions_,
const CurvedCellRegion curved_region)
+{
+ hp::MappingCollection<DoFHandlerType::dimension,
+ DoFHandlerType::space_dimension>
+ mapping_collection(mapping);
+
+ build_patches(mapping_collection, n_subdivisions_, curved_region);
+}
+
+
+
+template <int dim, typename DoFHandlerType>
+void
+DataOut<dim, DoFHandlerType>::build_patches(
+ const hp::MappingCollection<DoFHandlerType::dimension,
+ DoFHandlerType::space_dimension> &mapping,
+ const unsigned int n_subdivisions_,
+ const CurvedCellRegion curved_region)
{
// Check consistency of redundant template parameter
Assert(dim == DoFHandlerType::dimension,
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+
+// Test DataOut::build_patch() that takes MappingCollection.
+
+
+#include <deal.II/dofs/dof_accessor.h>
+#include <deal.II/dofs/dof_tools.h>
+
+#include <deal.II/fe/fe_nothing.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_system.h>
+#include <deal.II/fe/mapping_fe_field.h>
+#include <deal.II/fe/mapping_q.h>
+#include <deal.II/fe/mapping_q_cache.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_refinement.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/tria_iterator.h>
+
+#include <deal.II/hp/dof_handler.h>
+#include <deal.II/hp/fe_collection.h>
+#include <deal.II/hp/fe_values.h>
+
+#include <deal.II/numerics/data_out.h>
+#include <deal.II/numerics/vector_tools.h>
+
+#include "../tests.h"
+
+
+
+template <int dim, int spacedim = dim>
+void
+test()
+{
+ // create quarter hyper ball with 3 cells
+ Triangulation<dim> triangulation;
+ GridGenerator::quarter_hyper_ball(triangulation);
+
+ // create first mapping class, that does NOT preserve position of vertices
+ const FE_Q<dim, spacedim> feq(1);
+ const FESystem<dim, spacedim> fesystem(feq, spacedim);
+ DoFHandler<dim, spacedim> dhq(triangulation);
+ dhq.distribute_dofs(fesystem);
+ const ComponentMask mask(spacedim, true);
+ Vector<double> eulerq(dhq.n_dofs());
+ VectorTools::get_position_vector(dhq, eulerq, mask);
+ eulerq.add(-0.01);
+ MappingFEField<dim, spacedim> mapping_1(dhq, eulerq, mask);
+
+ // create first mapping class, that does preserve position of vertices
+ MappingQGeneric<dim, spacedim> mapping_2(1);
+
+ // create mapping collection
+ hp::FECollection<dim> fe_collection(FE_Q<dim>(1), FE_Q<dim>(1));
+ hp::MappingCollection<dim> mapping_collection(mapping_1, mapping_2);
+
+ // create dof-handler and assign cells to different fes/manifolds
+ hp::DoFHandler<dim> dof_handler(triangulation);
+
+ for (const auto &cell : dof_handler.active_cell_iterators())
+ cell->set_active_fe_index(std::min(cell->active_cell_index(), 1u));
+
+ dof_handler.distribute_dofs(fe_collection);
+
+ // output vector
+ DataOut<dim> data_out;
+ data_out.attach_dof_handler(dof_handler);
+
+ Vector<double> solution(dof_handler.n_dofs());
+
+ data_out.add_data_vector(solution, "solution");
+
+ data_out.build_patches(mapping_collection, 3);
+
+#if false
+ std::ofstream output("test.vtk");
+ data_out.write_vtk(output);
+#else
+ data_out.write_vtk(deallog.get_file_stream());
+#endif
+}
+
+
+
+int
+main()
+{
+ initlog();
+ deallog.get_file_stream().precision(2);
+
+ test<2>();
+}
--- /dev/null
+
+# vtk DataFile Version 3.0
+#This file was generated
+ASCII
+DATASET UNSTRUCTURED_GRID
+
+POINTS 48 double
+-0.010 -0.010 0
+0.18 -0.010 0
+0.36 -0.010 0
+0.55 -0.010 0
+-0.010 0.18 0
+0.16 0.16 0
+0.33 0.15 0
+0.50 0.13 0
+-0.010 0.36 0
+0.15 0.33 0
+0.30 0.30 0
+0.46 0.28 0
+-0.010 0.55 0
+0.13 0.50 0
+0.28 0.46 0
+0.42 0.42 0
+1.0 0.0 0
+0.90 0.24 0
+0.80 0.47 0
+0.71 0.71 0
+0.85 0.0 0
+0.77 0.20 0
+0.69 0.41 0
+0.61 0.61 0
+0.70 0.0 0
+0.64 0.17 0
+0.58 0.35 0
+0.52 0.52 0
+0.56 0.0 0
+0.51 0.14 0
+0.47 0.29 0
+0.43 0.43 0
+0.0 1.0 0
+0.0 0.85 0
+0.0 0.70 0
+0.0 0.56 0
+0.24 0.90 0
+0.20 0.77 0
+0.17 0.64 0
+0.14 0.51 0
+0.47 0.80 0
+0.41 0.69 0
+0.35 0.58 0
+0.29 0.47 0
+0.71 0.71 0
+0.61 0.61 0
+0.52 0.52 0
+0.43 0.43 0
+
+CELLS 27 135
+4 0 1 5 4
+4 1 2 6 5
+4 2 3 7 6
+4 4 5 9 8
+4 5 6 10 9
+4 6 7 11 10
+4 8 9 13 12
+4 9 10 14 13
+4 10 11 15 14
+4 16 17 21 20
+4 17 18 22 21
+4 18 19 23 22
+4 20 21 25 24
+4 21 22 26 25
+4 22 23 27 26
+4 24 25 29 28
+4 25 26 30 29
+4 26 27 31 30
+4 32 33 37 36
+4 33 34 38 37
+4 34 35 39 38
+4 36 37 41 40
+4 37 38 42 41
+4 38 39 43 42
+4 40 41 45 44
+4 41 42 46 45
+4 42 43 47 46
+
+CELL_TYPES 27
+ 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9
+POINT_DATA 48
+SCALARS solution double 1
+LOOKUP_TABLE default
+0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0