Skip to content

Commit

Permalink
Merge pull request #13027 from bangerth/face_number
Browse files Browse the repository at this point in the history
Allow accessing face numbers in data postprocessors
  • Loading branch information
drwells committed Dec 6, 2021
2 parents e1eec08 + 34d1275 commit 13a4730
Show file tree
Hide file tree
Showing 7 changed files with 354 additions and 5 deletions.
6 changes: 6 additions & 0 deletions doc/news/changes/minor/20211203Bangerth
Original file line number Diff line number Diff line change
@@ -0,0 +1,6 @@
New: The DataPostprocessorInputs::CommonInputs class now also
stores the face number that is currently being worked on when
using a DataOutFaces object. This allows accessing face-related
data like the boundary id from a DataPostprocessor.
<br>
(Wolfgang Bangerth, 2021/12/03)
100 changes: 100 additions & 0 deletions include/deal.II/numerics/data_postprocessor.h
Original file line number Diff line number Diff line change
Expand Up @@ -66,6 +66,7 @@ namespace DataPostprocessorInputs
* then the @p normal_vectors member variable does not contain anything
* useful.
*
*
* <h4>Cell access</h4>
*
* DataPostprocessor is typically called from classes such as DataOut
Expand Down Expand Up @@ -144,10 +145,58 @@ namespace DataPostprocessorInputs
* }
* };
* @endcode
*
*
* <h4>Face access</h4>
*
* When a DataPostprocessor object is used for output via the DataOutFaces
* class, it is sometimes necessary to also know which face is currently
* being worked on. Like accessing the cell as shown above, postprocessors
* can then query the face number via the CommonInputs::get_face_number()
* function. An example postprocessor that ignores its input and
* only puts the @ref GlossBoundaryIndicator "Boundary indicator"
* into the output file would then look as follows:
* @code
* template <int dim>
* class BoundaryIds : public DataPostprocessorScalar<dim>
* {
* public:
* BoundaryIds()
* : DataPostprocessorScalar<dim>("boundary_id", update_quadrature_points)
* {}
*
*
* virtual void
* evaluate_scalar_field(
* const DataPostprocessorInputs::Scalar<dim> &inputs,
* std::vector<Vector<double>> &computed_quantities) const override
* {
* AssertDimension(computed_quantities.size(),
* inputs.solution_values.size());
*
* // Get the cell and face we are currently dealing with:
* const typename DoFHandler<dim>::active_cell_iterator cell =
* inputs.template get_cell<dim>();
* const unsigned int face = inputs.get_face_number();
*
* // Then fill the output fields with the boundary_id of the face
* for (auto &output : computed_quantities)
* {
* AssertDimension(output.size(), 1);
* output(0) = cell->face(face)->boundary_id();
* }
* }
* };
* @endcode
*/
template <int spacedim>
struct CommonInputs
{
/**
* Constructor.
*/
CommonInputs();

/**
* An array of vectors normal to the faces of cells, evaluated at the points
* at which we are generating graphical output. This array is only used by
Expand Down Expand Up @@ -192,6 +241,22 @@ namespace DataPostprocessorInputs
void
set_cell(const typename DoFHandler<dim, spacedim>::cell_iterator &cell);

/**
* Set the cell and face number that is currently being used in evaluating
* the data for which the DataPostprocessor object is being called. Given
* that a face is required, this function is meant to be called by a class
* such as DataOutFaces.
*
* This function is not usually called from user space, but is instead
* called by DataOutFaces and similar classes when creating the object that
* is then passed to DataPostprocessor.
*/
template <int dim>
void
set_cell_and_face(
const typename DoFHandler<dim, spacedim>::cell_iterator &cell,
const unsigned int face_number);

/**
* Set the cell that is currently being used in evaluating the data
* for which the DataPostprocessor object is being called.
Expand Down Expand Up @@ -228,6 +293,18 @@ namespace DataPostprocessorInputs
DEAL_II_DEPRECATED typename DoFHandlerType::cell_iterator
get_cell() const;

/**
* Query the face number on which we currently produce graphical output.
* See the documentation of the current class for an example on how
* to use the related get_cell() function that is meant to query the cell
* currently being worked on.
*
* This function is intended for use when producing graphical output on
* faces, for example through the DataOutFaces class.
*/
unsigned int
get_face_number() const;

private:
/**
* The place where set_cell() stores the cell. Since the actual data
Expand All @@ -237,6 +314,12 @@ namespace DataPostprocessorInputs
* get_cell().
*/
boost::any cell;

/**
* The place where set_cell_and_face() stores the number of the face
* being worked on.
*/
unsigned int face_number;
};

/**
Expand Down Expand Up @@ -1224,6 +1307,23 @@ namespace DataPostprocessorInputs
// if we had nothing stored before, or if we had stored a different
// data type, just let boost::any replace things
cell = new_cell;

// Also reset the face number, just to make sure nobody
// accidentally uses an outdated value.
face_number = numbers::invalid_unsigned_int;
}



template <int spacedim>
template <int dim>
void
CommonInputs<spacedim>::set_cell_and_face(
const typename DoFHandler<dim, spacedim>::cell_iterator &new_cell,
const unsigned int new_face_number)
{
set_cell<dim>(new_cell);
face_number = new_face_number;
}


Expand Down
6 changes: 4 additions & 2 deletions source/numerics/data_out_faces.cc
Original file line number Diff line number Diff line change
Expand Up @@ -204,7 +204,8 @@ DataOutFaces<dim, spacedim>::build_one_patch(
cell->level(),
cell->index(),
this->dof_data[dataset]->dof_handler);
data.patch_values_scalar.template set_cell<dim>(dh_cell);
data.patch_values_scalar.template set_cell_and_face<dim>(
dh_cell, face_number);

postprocessor->evaluate_scalar_field(
data.patch_values_scalar,
Expand Down Expand Up @@ -247,7 +248,8 @@ DataOutFaces<dim, spacedim>::build_one_patch(
cell->level(),
cell->index(),
this->dof_data[dataset]->dof_handler);
data.patch_values_system.template set_cell<dim>(dh_cell);
data.patch_values_system.template set_cell_and_face<dim>(
dh_cell, face_number);

postprocessor->evaluate_vector_field(
data.patch_values_system,
Expand Down
31 changes: 28 additions & 3 deletions source/numerics/data_postprocessor.cc
Original file line number Diff line number Diff line change
Expand Up @@ -19,6 +19,29 @@ DEAL_II_NAMESPACE_OPEN



namespace DataPostprocessorInputs
{
template <int spacedim>
CommonInputs<spacedim>::CommonInputs()
: face_number(numbers::invalid_unsigned_int)
{}



template <int spacedim>
unsigned int
CommonInputs<spacedim>::get_face_number() const
{
Assert(
face_number != numbers::invalid_unsigned_int,
ExcMessage(
"This function can only be called if set_cell_and_face() has "
"previously been called. Typically, this would be by using DataOutFaces "
"or a related class."));
return face_number;
}
} // namespace DataPostprocessorInputs

// -------------------------- DataPostprocessor ---------------------------

template <int dim>
Expand Down Expand Up @@ -54,7 +77,7 @@ DataPostprocessor<dim>::get_data_component_interpretation() const
}


// -------------------------- DataPostprocessorScalar -------------------------
// -------------------- DataPostprocessorScalar -------------------------

template <int dim>
DataPostprocessorScalar<dim>::DataPostprocessorScalar(
Expand Down Expand Up @@ -93,7 +116,8 @@ DataPostprocessorScalar<dim>::get_needed_update_flags() const



// -------------------------- DataPostprocessorVector -------------------------
// -------------------------- DataPostprocessorVector
// -------------------------

template <int dim>
DataPostprocessorVector<dim>::DataPostprocessorVector(
Expand Down Expand Up @@ -132,7 +156,8 @@ DataPostprocessorVector<dim>::get_needed_update_flags() const



// -------------------------- DataPostprocessorTensor -------------------------
// -------------------------- DataPostprocessorTensor
// -------------------------

template <int dim>
DataPostprocessorTensor<dim>::DataPostprocessorTensor(
Expand Down
4 changes: 4 additions & 0 deletions source/numerics/data_postprocessor.inst.in
Original file line number Diff line number Diff line change
Expand Up @@ -16,6 +16,10 @@

for (deal_II_dimension : DIMENSIONS)
{
namespace DataPostprocessorInputs
\{
template struct CommonInputs<deal_II_dimension>;
\}
template class DataPostprocessor<deal_II_dimension>;
template class DataPostprocessorScalar<deal_II_dimension>;
template class DataPostprocessorVector<deal_II_dimension>;
Expand Down
118 changes: 118 additions & 0 deletions tests/data_out/data_out_faces_postprocessor_boundary_id.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,118 @@
// ---------------------------------------------------------------------
//
// Copyright (C) 2016 - 2018 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.
//
// ---------------------------------------------------------------------


// Tests DataPostprocessor: access the cell and face we are currently
// working on for a scalar finite element field and DataOutFaces. This
// is then used to output the boundary id of a face.


#include <deal.II/base/function.h>

#include <deal.II/dofs/dof_accessor.h>
#include <deal.II/dofs/dof_handler.h>

#include <deal.II/fe/fe_dgq.h>
#include <deal.II/fe/fe_system.h>

#include <deal.II/grid/grid_generator.h>
#include <deal.II/grid/tria.h>
#include <deal.II/grid/tria_accessor.h>
#include <deal.II/grid/tria_iterator.h>

#include <deal.II/lac/vector.h>

#include <deal.II/numerics/data_out_faces.h>
#include <deal.II/numerics/data_postprocessor.h>
#include <deal.II/numerics/matrix_tools.h>
#include <deal.II/numerics/vector_tools.h>

#include "../tests.h"



std::ofstream logfile("output");


template <int dim>
class BoundaryIds : public DataPostprocessorScalar<dim>
{
public:
BoundaryIds()
: DataPostprocessorScalar<dim>("boundary_id", update_quadrature_points)
{}


virtual void
evaluate_scalar_field(
const DataPostprocessorInputs::Scalar<dim> &inputs,
std::vector<Vector<double>> &computed_quantities) const override
{
AssertDimension(computed_quantities.size(), inputs.solution_values.size());

// Get the cell and face we are currently dealing with:
const typename DoFHandler<dim>::active_cell_iterator cell =
inputs.template get_cell<dim>();
const unsigned int face = inputs.get_face_number();

// Then fill the output fields with the boundary_id of the face
for (auto &output : computed_quantities)
{
AssertDimension(output.size(), 1);
output(0) = cell->face(face)->boundary_id();
}
}
};



template <int dim>
void
test()
{
Triangulation<dim> triangulation;
FE_DGQ<dim> fe(1);
DoFHandler<dim> dof_handler(triangulation);

GridGenerator::hyper_cube(triangulation, 0, 1, true);
triangulation.refine_global(1);
triangulation.begin_active()->set_refine_flag();
triangulation.execute_coarsening_and_refinement();

dof_handler.distribute_dofs(fe);

// Create a dummy vector. We will ignore its contents.
Vector<double> solution(dof_handler.n_dofs());

BoundaryIds<dim> p;
DataOutFaces<dim> 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<2>();

return 0;
}

0 comments on commit 13a4730

Please sign in to comment.