From: Wolfgang Bangerth Date: Wed, 21 Dec 2016 16:28:59 +0000 (-0700) Subject: Provide the cell via DataOutFaces. Add tests. X-Git-Tag: v8.5.0-rc1~308^2~2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=ca0b331be3e675f93c167fb7ce1700df0b5a0595;p=dealii.git Provide the cell via DataOutFaces. Add tests. --- diff --git a/source/numerics/data_out_faces.cc b/source/numerics/data_out_faces.cc index cc999ba3f0..a81fd310bc 100644 --- a/source/numerics/data_out_faces.cc +++ b/source/numerics/data_out_faces.cc @@ -169,6 +169,12 @@ build_one_patch (const FaceDescriptor *cell_and_face, if (update_flags & update_normal_vectors) data.patch_values_scalar.normals = this_fe_patch_values.get_all_normal_vectors(); + const typename DoFHandlerType::active_cell_iterator dh_cell(&cell_and_face->first->get_triangulation(), + cell_and_face->first->level(), + cell_and_face->first->index(), + this->dof_data[dataset]->dof_handler); + data.patch_values_scalar.template set_cell (dh_cell); + postprocessor-> evaluate_scalar_field(data.patch_values_scalar, data.postprocessed_values[dataset]); @@ -194,6 +200,12 @@ build_one_patch (const FaceDescriptor *cell_and_face, if (update_flags & update_normal_vectors) data.patch_values_system.normals = this_fe_patch_values.get_all_normal_vectors(); + const typename DoFHandlerType::active_cell_iterator dh_cell(&cell_and_face->first->get_triangulation(), + cell_and_face->first->level(), + cell_and_face->first->index(), + this->dof_data[dataset]->dof_handler); + data.patch_values_system.template set_cell (dh_cell); + postprocessor-> evaluate_vector_field(data.patch_values_system, data.postprocessed_values[dataset]); diff --git a/tests/numerics/data_out_faces_postprocessor_scalar_02.cc b/tests/numerics/data_out_faces_postprocessor_scalar_02.cc new file mode 100644 index 0000000000..2345b90fce --- /dev/null +++ b/tests/numerics/data_out_faces_postprocessor_scalar_02.cc @@ -0,0 +1,125 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2016 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 at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + + +// tests DataPostprocessor: access the cell we are currently working +// on for a scalar finite element field and DataOutFaces + + +#include "../tests.h" +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include + +#include +#include +#include + +#include + + +std::ofstream logfile("output"); + + +template +class MyPostprocessor : public DataPostprocessorScalar +{ +public: + MyPostprocessor () + : + DataPostprocessorScalar ("data", update_values) + {} + + virtual + void + evaluate_scalar_field (const DataPostprocessorInputs::Scalar &input_data, + std::vector > &computed_quantities) const + { + for (unsigned int q=0; q::cell_iterator + cell = input_data.template get_cell >(); + + Assert (input_data.solution_values[q] + == + double(cell->active_cell_index()), + ExcInternalError()); + + computed_quantities[q][0] = input_data.solution_values[q]; + } + } +}; + + + +template +void test () +{ + Triangulation triangulation; + FE_DGQ fe(1); + DoFHandler dof_handler(triangulation); + + GridGenerator::hyper_cube (triangulation, 0, 1); + triangulation.refine_global (1); + triangulation.begin_active()->set_refine_flag (); + triangulation.execute_coarsening_and_refinement (); + + dof_handler.distribute_dofs (fe); + + // create a vector with piecewise constants, and set the elements of + // the vector so that on each cell the field has values equal to the + // cell's active_cell_index + Vector solution (dof_handler.n_dofs()); + for (typename DoFHandler::active_cell_iterator cell = dof_handler.begin_active(); + cell != dof_handler.end(); ++cell) + { + std::vector local_dof_indices (cell->get_fe().dofs_per_cell); + cell->get_dof_indices (local_dof_indices); + for (unsigned int i=0; iactive_cell_index(); + } + + MyPostprocessor p; + DataOutFaces data_out; + data_out.attach_dof_handler (dof_handler); + data_out.add_data_vector (solution, p); + data_out.build_patches (); + data_out.write_gnuplot (logfile); +} + + +int main () +{ + logfile << std::setprecision(2); + deallog << std::setprecision(2); + + test<2>(); + test<3>(); + + return 0; +} diff --git a/tests/numerics/data_out_faces_postprocessor_scalar_02.output b/tests/numerics/data_out_faces_postprocessor_scalar_02.output new file mode 100644 index 0000000000..68365512fc --- /dev/null +++ b/tests/numerics/data_out_faces_postprocessor_scalar_02.output @@ -0,0 +1,285 @@ +# This file was generated by the deal.II library. + + +# +# For a description of the GNUPLOT format see the GNUPLOT manual. +# +# +1 0 0 +1 0.5 0 + + +0.5 0 0 +1 0 0 + + +0 0.5 1 +0 1 1 + + +0 1 1 +0.5 1 1 + + +1 0.5 2 +1 1 2 + + +0.5 1 2 +1 1 2 + + +0 0 3 +0 0.25 3 + + +0 0 3 +0.25 0 3 + + +0.25 0 4 +0.5 0 4 + + +0 0.25 5 +0 0.5 5 + + +# This file was generated by the deal.II library. + + +# +# For a description of the GNUPLOT format see the GNUPLOT manual. +# +# +1 0 0 0 +1 0.5 0 0 + +1 0 0.5 0 +1 0.5 0.5 0 + + +0.5 0 0 0 +0.5 0 0.5 0 + +1 0 0 0 +1 0 0.5 0 + + +0.5 0 0 0 +1 0 0 0 + +0.5 0.5 0 0 +1 0.5 0 0 + + +0 0.5 0 1 +0 1 0 1 + +0 0.5 0.5 1 +0 1 0.5 1 + + +0 1 0 1 +0 1 0.5 1 + +0.5 1 0 1 +0.5 1 0.5 1 + + +0 0.5 0 1 +0.5 0.5 0 1 + +0 1 0 1 +0.5 1 0 1 + + +1 0.5 0 2 +1 1 0 2 + +1 0.5 0.5 2 +1 1 0.5 2 + + +0.5 1 0 2 +0.5 1 0.5 2 + +1 1 0 2 +1 1 0.5 2 + + +0.5 0.5 0 2 +1 0.5 0 2 + +0.5 1 0 2 +1 1 0 2 + + +0 0 0.5 3 +0 0.5 0.5 3 + +0 0 1 3 +0 0.5 1 3 + + +0 0 0.5 3 +0 0 1 3 + +0.5 0 0.5 3 +0.5 0 1 3 + + +0 0 1 3 +0.5 0 1 3 + +0 0.5 1 3 +0.5 0.5 1 3 + + +1 0 0.5 4 +1 0.5 0.5 4 + +1 0 1 4 +1 0.5 1 4 + + +0.5 0 0.5 4 +0.5 0 1 4 + +1 0 0.5 4 +1 0 1 4 + + +0.5 0 1 4 +1 0 1 4 + +0.5 0.5 1 4 +1 0.5 1 4 + + +0 0.5 0.5 5 +0 1 0.5 5 + +0 0.5 1 5 +0 1 1 5 + + +0 1 0.5 5 +0 1 1 5 + +0.5 1 0.5 5 +0.5 1 1 5 + + +0 0.5 1 5 +0.5 0.5 1 5 + +0 1 1 5 +0.5 1 1 5 + + +1 0.5 0.5 6 +1 1 0.5 6 + +1 0.5 1 6 +1 1 1 6 + + +0.5 1 0.5 6 +0.5 1 1 6 + +1 1 0.5 6 +1 1 1 6 + + +0.5 0.5 1 6 +1 0.5 1 6 + +0.5 1 1 6 +1 1 1 6 + + +0 0 0 7 +0 0.25 0 7 + +0 0 0.25 7 +0 0.25 0.25 7 + + +0 0 0 7 +0 0 0.25 7 + +0.25 0 0 7 +0.25 0 0.25 7 + + +0 0 0 7 +0.25 0 0 7 + +0 0.25 0 7 +0.25 0.25 0 7 + + +0.25 0 0 8 +0.25 0 0.25 8 + +0.5 0 0 8 +0.5 0 0.25 8 + + +0.25 0 0 8 +0.5 0 0 8 + +0.25 0.25 0 8 +0.5 0.25 0 8 + + +0 0.25 0 9 +0 0.5 0 9 + +0 0.25 0.25 9 +0 0.5 0.25 9 + + +0 0.25 0 9 +0.25 0.25 0 9 + +0 0.5 0 9 +0.25 0.5 0 9 + + +0.25 0.25 0 10 +0.5 0.25 0 10 + +0.25 0.5 0 10 +0.5 0.5 0 10 + + +0 0 0.25 11 +0 0.25 0.25 11 + +0 0 0.5 11 +0 0.25 0.5 11 + + +0 0 0.25 11 +0 0 0.5 11 + +0.25 0 0.25 11 +0.25 0 0.5 11 + + +0.25 0 0.25 12 +0.25 0 0.5 12 + +0.5 0 0.25 12 +0.5 0 0.5 12 + + +0 0.25 0.25 13 +0 0.5 0.25 13 + +0 0.25 0.5 13 +0 0.5 0.5 13 + + diff --git a/tests/numerics/data_out_faces_postprocessor_vector_02.cc b/tests/numerics/data_out_faces_postprocessor_vector_02.cc new file mode 100644 index 0000000000..53edf76282 --- /dev/null +++ b/tests/numerics/data_out_faces_postprocessor_vector_02.cc @@ -0,0 +1,128 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2016 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 at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + + +// tests DataPostprocessor: access the cell we are currently working +// on for a vector finite element field and DataOutFaces + + +#include "../tests.h" +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include + +#include +#include +#include + +#include + + +std::ofstream logfile("output"); + + +template +class MyPostprocessor : public DataPostprocessorVector +{ +public: + MyPostprocessor () + : + DataPostprocessorVector ("data", update_values) + {} + + virtual + void + evaluate_vector_field (const DataPostprocessorInputs::Vector &input_data, + std::vector > &computed_quantities) const + { + for (unsigned int q=0; q::cell_iterator + cell = input_data.template get_cell >(); + + Assert (input_data.solution_values[q](0) + == + double(cell->active_cell_index()), + ExcInternalError()); + + computed_quantities[q][0] = input_data.solution_values[q](0); + } + } +}; + + + +template +void test () +{ + Triangulation triangulation; + FESystem fe(FE_DGQ(0),2); + DoFHandler dof_handler(triangulation); + + GridGenerator::hyper_cube (triangulation, 0, 1); + triangulation.refine_global (1); + triangulation.begin_active()->set_refine_flag (); + triangulation.execute_coarsening_and_refinement (); + + dof_handler.distribute_dofs (fe); + + // create a vector with piecewise constants, and set the elements of + // the vector so that on each cell the field has values equal to the + // cell's active_cell_index + Vector solution (dof_handler.n_dofs()); + for (typename DoFHandler::active_cell_iterator cell = dof_handler.begin_active(); + cell != dof_handler.end(); ++cell) + { + std::vector local_dof_indices (cell->get_fe().dofs_per_cell); + cell->get_dof_indices (local_dof_indices); + for (unsigned int i=0; iactive_cell_index(); + } + + MyPostprocessor p; + DataOutFaces data_out; + data_out.attach_dof_handler (dof_handler); + data_out.add_data_vector (solution, p); + data_out.build_patches (); + data_out.write_gnuplot (logfile); +} + + +int main () +{ + logfile << std::setprecision(2); + deallog << std::setprecision(2); + + test<2>(); + test<3>(); + + return 0; +} diff --git a/tests/numerics/data_out_faces_postprocessor_vector_02.output b/tests/numerics/data_out_faces_postprocessor_vector_02.output new file mode 100644 index 0000000000..85f25875da --- /dev/null +++ b/tests/numerics/data_out_faces_postprocessor_vector_02.output @@ -0,0 +1,285 @@ +# This file was generated by the deal.II library. + + +# +# For a description of the GNUPLOT format see the GNUPLOT manual. +# +# +1 0 0 0 +1 0.5 0 0 + + +0.5 0 0 0 +1 0 0 0 + + +0 0.5 1 0 +0 1 1 0 + + +0 1 1 0 +0.5 1 1 0 + + +1 0.5 2 0 +1 1 2 0 + + +0.5 1 2 0 +1 1 2 0 + + +0 0 3 0 +0 0.25 3 0 + + +0 0 3 0 +0.25 0 3 0 + + +0.25 0 4 0 +0.5 0 4 0 + + +0 0.25 5 0 +0 0.5 5 0 + + +# This file was generated by the deal.II library. + + +# +# For a description of the GNUPLOT format see the GNUPLOT manual. +# +# +1 0 0 0 0 0 +1 0.5 0 0 0 0 + +1 0 0.5 0 0 0 +1 0.5 0.5 0 0 0 + + +0.5 0 0 0 0 0 +0.5 0 0.5 0 0 0 + +1 0 0 0 0 0 +1 0 0.5 0 0 0 + + +0.5 0 0 0 0 0 +1 0 0 0 0 0 + +0.5 0.5 0 0 0 0 +1 0.5 0 0 0 0 + + +0 0.5 0 1 0 0 +0 1 0 1 0 0 + +0 0.5 0.5 1 0 0 +0 1 0.5 1 0 0 + + +0 1 0 1 0 0 +0 1 0.5 1 0 0 + +0.5 1 0 1 0 0 +0.5 1 0.5 1 0 0 + + +0 0.5 0 1 0 0 +0.5 0.5 0 1 0 0 + +0 1 0 1 0 0 +0.5 1 0 1 0 0 + + +1 0.5 0 2 0 0 +1 1 0 2 0 0 + +1 0.5 0.5 2 0 0 +1 1 0.5 2 0 0 + + +0.5 1 0 2 0 0 +0.5 1 0.5 2 0 0 + +1 1 0 2 0 0 +1 1 0.5 2 0 0 + + +0.5 0.5 0 2 0 0 +1 0.5 0 2 0 0 + +0.5 1 0 2 0 0 +1 1 0 2 0 0 + + +0 0 0.5 3 0 0 +0 0.5 0.5 3 0 0 + +0 0 1 3 0 0 +0 0.5 1 3 0 0 + + +0 0 0.5 3 0 0 +0 0 1 3 0 0 + +0.5 0 0.5 3 0 0 +0.5 0 1 3 0 0 + + +0 0 1 3 0 0 +0.5 0 1 3 0 0 + +0 0.5 1 3 0 0 +0.5 0.5 1 3 0 0 + + +1 0 0.5 4 0 0 +1 0.5 0.5 4 0 0 + +1 0 1 4 0 0 +1 0.5 1 4 0 0 + + +0.5 0 0.5 4 0 0 +0.5 0 1 4 0 0 + +1 0 0.5 4 0 0 +1 0 1 4 0 0 + + +0.5 0 1 4 0 0 +1 0 1 4 0 0 + +0.5 0.5 1 4 0 0 +1 0.5 1 4 0 0 + + +0 0.5 0.5 5 0 0 +0 1 0.5 5 0 0 + +0 0.5 1 5 0 0 +0 1 1 5 0 0 + + +0 1 0.5 5 0 0 +0 1 1 5 0 0 + +0.5 1 0.5 5 0 0 +0.5 1 1 5 0 0 + + +0 0.5 1 5 0 0 +0.5 0.5 1 5 0 0 + +0 1 1 5 0 0 +0.5 1 1 5 0 0 + + +1 0.5 0.5 6 0 0 +1 1 0.5 6 0 0 + +1 0.5 1 6 0 0 +1 1 1 6 0 0 + + +0.5 1 0.5 6 0 0 +0.5 1 1 6 0 0 + +1 1 0.5 6 0 0 +1 1 1 6 0 0 + + +0.5 0.5 1 6 0 0 +1 0.5 1 6 0 0 + +0.5 1 1 6 0 0 +1 1 1 6 0 0 + + +0 0 0 7 0 0 +0 0.25 0 7 0 0 + +0 0 0.25 7 0 0 +0 0.25 0.25 7 0 0 + + +0 0 0 7 0 0 +0 0 0.25 7 0 0 + +0.25 0 0 7 0 0 +0.25 0 0.25 7 0 0 + + +0 0 0 7 0 0 +0.25 0 0 7 0 0 + +0 0.25 0 7 0 0 +0.25 0.25 0 7 0 0 + + +0.25 0 0 8 0 0 +0.25 0 0.25 8 0 0 + +0.5 0 0 8 0 0 +0.5 0 0.25 8 0 0 + + +0.25 0 0 8 0 0 +0.5 0 0 8 0 0 + +0.25 0.25 0 8 0 0 +0.5 0.25 0 8 0 0 + + +0 0.25 0 9 0 0 +0 0.5 0 9 0 0 + +0 0.25 0.25 9 0 0 +0 0.5 0.25 9 0 0 + + +0 0.25 0 9 0 0 +0.25 0.25 0 9 0 0 + +0 0.5 0 9 0 0 +0.25 0.5 0 9 0 0 + + +0.25 0.25 0 10 0 0 +0.5 0.25 0 10 0 0 + +0.25 0.5 0 10 0 0 +0.5 0.5 0 10 0 0 + + +0 0 0.25 11 0 0 +0 0.25 0.25 11 0 0 + +0 0 0.5 11 0 0 +0 0.25 0.5 11 0 0 + + +0 0 0.25 11 0 0 +0 0 0.5 11 0 0 + +0.25 0 0.25 11 0 0 +0.25 0 0.5 11 0 0 + + +0.25 0 0.25 12 0 0 +0.25 0 0.5 12 0 0 + +0.5 0 0.25 12 0 0 +0.5 0 0.5 12 0 0 + + +0 0.25 0.25 13 0 0 +0 0.5 0.25 13 0 0 + +0 0.25 0.5 13 0 0 +0 0.5 0.5 13 0 0 + + diff --git a/tests/numerics/data_out_postprocessor_scalar_02.cc b/tests/numerics/data_out_postprocessor_scalar_02.cc index 560edbdedc..be28e10de2 100644 --- a/tests/numerics/data_out_postprocessor_scalar_02.cc +++ b/tests/numerics/data_out_postprocessor_scalar_02.cc @@ -63,7 +63,7 @@ public: // get the cell this all belongs to typename DoFHandler::cell_iterator - cell = input_data.template get_cell >(); + cell = input_data.template get_cell >(); Assert (input_data.solution_values[q] == @@ -103,7 +103,7 @@ void test () for (unsigned int i=0; iactive_cell_index(); } - + MyPostprocessor p; DataOut data_out; data_out.attach_dof_handler (dof_handler);