In #3333, I added virtual functions to DataPostprocessor with the same name as the existing
functions. This leads to issues where we get a warning in every derived class that only
overloads one of these functions, because that hides the other function. This is, well,
suboptimal.
This patch is therefore a redo of my earlier attempt in which I continue to deprecate
the old functions, but the new functions have a different name. I think they also
have a better name (for a discussion of the naming, see
https://github.com/geodynamics/aspect/issues/1284 ). The different names avoid the
problem of getting the warning and should therefore lead to less discontent. They
also avoid the need to try to work around the warnings using 'using' declarations, like
in #3528.
Rather than adding each possible argument anyone may want to use
individually to the list of the postprocessor function arguments, the
existing functions have been deprecated in favor of a new set of
- functions that
+ functions DataPostprocessor::evaluate_scalar_field() and
+ DataPostprocessor::evaluate_vector_field() that
take a reference to a structure that contains these individual pieces
of information. We can extend the members of these structures without
backward compatibility issues because the functions still get a
// version of this function that in addition to the data vector has an
// additional argument of type DataPostprocessor. What happens when this
// function is used for output is that at each point where output data is to
- // be generated, the DataPostprocessor::compute_derived_quantities_scalar or
- // DataPostprocessor::compute_derived_quantities_vector function of the
+ // be generated, the DataPostprocessor::evaluate_scalar_field() or
+ // DataPostprocessor::evaluate_vector_field() function of the
// specified DataPostprocessor object is invoked to compute the output
// quantities from the values, the gradients and the second derivatives of
// the finite element function represented by the data vector (in the case
public:
ComputeIntensity ();
- using DataPostprocessorScalar<dim>::compute_derived_quantities_vector;
-
virtual
void
- compute_derived_quantities_vector
+ evaluate_vector_field
(const DataPostprocessorInputs::Vector<dim> &inputs,
std::vector<Vector<double> > &computed_quantities) const;
};
// are provided.
template <int dim>
void
- ComputeIntensity<dim>::compute_derived_quantities_vector
+ ComputeIntensity<dim>::evaluate_vector_field
(const DataPostprocessorInputs::Vector<dim> &inputs,
std::vector<Vector<double> > &computed_quantities) const
{
// inherits from the class DataPostprocessor, which can be attached to
// DataOut. This allows us to output derived quantities from the solution,
// like the friction heating included in this example. It overloads the
- // virtual function DataPostprocessor::compute_derived_quantities_vector,
- // which is then internally called from DataOut::build_patches. We have to
+ // virtual function DataPostprocessor::evaluate_vector_field(),
+ // which is then internally called from DataOut::build_patches(). We have to
// give it values of the numerical solution, its derivatives, normals to the
// cell, the actual evaluation points and any additional quantities. This
// follows the same procedure as discussed in step-29 and other programs.
Postprocessor (const unsigned int partition,
const double minimal_pressure);
- using DataPostprocessor<dim>::compute_derived_quantities_vector;
-
virtual
void
- compute_derived_quantities_vector
+ evaluate_vector_field
(const DataPostprocessorInputs::Vector<dim> &inputs,
std::vector<Vector<double> > &computed_quantities) const;
template <int dim>
void
BoussinesqFlowProblem<dim>::Postprocessor::
- compute_derived_quantities_vector
+ evaluate_vector_field
(const DataPostprocessorInputs::Vector<dim> &inputs,
std::vector<Vector<double> > &computed_quantities) const
{
public:
Postprocessor (const bool do_schlieren_plot);
- using DataPostprocessor<dim>::compute_derived_quantities_vector;
-
virtual
void
- compute_derived_quantities_vector
+ evaluate_vector_field
(const DataPostprocessorInputs::Vector<dim> &inputs,
std::vector<Vector<double> > &computed_quantities) const;
template <int dim>
void
EulerEquations<dim>::Postprocessor::
- compute_derived_quantities_vector
+ evaluate_vector_field
(const DataPostprocessorInputs::Vector<dim> &inputs,
std::vector<Vector<double> > &computed_quantities) const
{
class Postprocessor : public DataPostprocessor<dim>
{
public:
- using DataPostprocessor<dim>::compute_derived_quantities_vector;
-
virtual
void
- compute_derived_quantities_vector
+ evaluate_vector_field
(const dealii::DataPostprocessorInputs::Vector<dim> &inputs,
std::vector<Vector<double> > &computed_quantities) const;
template <int dim>
void
Postprocessor<dim>::
- compute_derived_quantities_vector
+ evaluate_vector_field
(const dealii::DataPostprocessorInputs::Vector<dim> &inputs,
std::vector<Vector<double> > &computed_quantities) const
{
{
/**
* A structure that is used to pass information to
- * DataPostprocessor::compute_derived_quantities_scalar(). It contains
+ * DataPostprocessor::evaluate_scalar_field(). It contains
* the values and (if requested) derivatives of a scalar solution
* variable at the evaluation points on a cell or face. If appropriate,
* it also contains the normal vectors to the geometry on which output
/**
* A structure that is used to pass information to
- * DataPostprocessor::compute_derived_quantities_vector(). It contains
+ * DataPostprocessor::evaluate_vector_field(). It contains
* the values and (if requested) derivatives of a vector-valued solution
* variable at the evaluation points on a cell or face. If appropriate,
* it also contains the normal vectors to the geometry on which output
* function get_needed_update_flags(). It is your responsibility to use only
* those values which were updated in the calculation of derived quantities.
* The DataOut object will provide references to the requested data in the
- * call to compute_derived_quantities_scalar() or
- * compute_derived_quantities_vector() (DataOut decides which of the two
+ * call to evaluate_scalar_field() or
+ * evaluate_vector_field() (DataOut decides which of the two
* functions to call depending on whether the finite element in use has only a
* single, or multiple vector components; note that this is only determined by
* the number of components in the finite element in use, and not by whether
* hand, in step-29 we implement a postprocessor that only computes the
* magnitude of a complex number given by a two-component finite element. It
* seems silly to have to implement four virtual functions for this
- * (compute_derived_quantities_scalar() or
- * compute_derived_quantities_vector(), get_names(), get_update_flags() and
+ * (evaluate_scalar_field() or evaluate_vector_field(), get_names(), get_update_flags() and
* get_data_component_interpretation()).
*
* To this end there are two classes DataPostprocessorScalar and
*/
virtual
void
- compute_derived_quantities_scalar (const DataPostprocessorInputs::Scalar<dim> &input_data,
- std::vector<Vector<double> > &computed_quantities) const;
+ evaluate_scalar_field (const DataPostprocessorInputs::Scalar<dim> &input_data,
+ std::vector<Vector<double> > &computed_quantities) const;
/**
* @deprecated This function is deprecated. It has been superseded by
- * function of same name above that receives a superset of the
+ * the evaluate_scalar_field() function that receives a superset of the
* information provided to the current function through the members
* of the structure it receives as the first argument.
*
std::vector<Vector<double> > &computed_quantities) const DEAL_II_DEPRECATED;
/**
- * Same as the compute_derived_quantities_scalar() function, but this
+ * Same as the evaluate_scalar_field() function, but this
* function is called when the original data vector represents vector data,
* i.e. the finite element in use has multiple vector components.
*/
virtual
void
- compute_derived_quantities_vector (const DataPostprocessorInputs::Vector<dim> &input_data,
- std::vector<Vector<double> > &computed_quantities) const;
+ evaluate_vector_field (const DataPostprocessorInputs::Vector<dim> &input_data,
+ std::vector<Vector<double> > &computed_quantities) const;
/**
* @deprecated This function is deprecated. It has been superseded by
- * function of same name above that receives a superset of the
+ * the evaluate_vector_field() function that receives a superset of the
* information provided to the current function through the members
* of the structure it receives as the first argument.
*
* functions by hand.
*
* All derived classes have to do is implement a constructor and overload
- * either DataPostprocessor::compute_derived_quantities_scalar() or
- * DataPostprocessor::compute_derived_quantities_vector().
+ * either DataPostprocessor::evaluate_scalar_field() or
+ * DataPostprocessor::evaluate_vector_field().
*
* An example of how this class can be used can be found in step-29.
*
* functions by hand.
*
* All derived classes have to do is implement a constructor and overload
- * either DataPostprocessor::compute_derived_quantities_scalar() or
- * DataPostprocessor::compute_derived_quantities_vector().
+ * either DataPostprocessor::evaluate_scalar_field() or
+ * DataPostprocessor::evaluate_vector_field().
*
* An example of how the closely related class DataPostprocessorScalar is used
* can be found in step-29.
postprocessor->
- compute_derived_quantities_scalar(scratch_data.patch_values_scalar,
- scratch_data.postprocessed_values[dataset]);
+ evaluate_scalar_field(scratch_data.patch_values_scalar,
+ scratch_data.postprocessed_values[dataset]);
}
else
{
scratch_data.patch_values_system.evaluation_points = this_fe_patch_values.get_quadrature_points();
postprocessor->
- compute_derived_quantities_vector(scratch_data.patch_values_system,
- scratch_data.postprocessed_values[dataset]);
+ evaluate_vector_field(scratch_data.patch_values_system,
+ scratch_data.postprocessed_values[dataset]);
}
for (unsigned int q=0; q<n_q_points; ++q)
data.patch_values_scalar.normals = this_fe_patch_values.get_all_normal_vectors();
postprocessor->
- compute_derived_quantities_scalar(data.patch_values_scalar,
- data.postprocessed_values[dataset]);
+ evaluate_scalar_field(data.patch_values_scalar,
+ data.postprocessed_values[dataset]);
}
else
{
data.patch_values_system.normals = this_fe_patch_values.get_all_normal_vectors();
postprocessor->
- compute_derived_quantities_vector(data.patch_values_system,
- data.postprocessed_values[dataset]);
+ evaluate_vector_field(data.patch_values_system,
+ data.postprocessed_values[dataset]);
}
for (unsigned int q=0; q<n_q_points; ++q)
data.patch_values_scalar.evaluation_points = fe_patch_values.get_quadrature_points();
postprocessor->
- compute_derived_quantities_scalar(data.patch_values_scalar,
- data.postprocessed_values[dataset]);
+ evaluate_scalar_field(data.patch_values_scalar,
+ data.postprocessed_values[dataset]);
}
else
{
std::vector<Point<space_dimension> > dummy_normals;
postprocessor->
- compute_derived_quantities_vector(data.patch_values_system,
- data.postprocessed_values[dataset]);
+ evaluate_vector_field(data.patch_values_system,
+ data.postprocessed_values[dataset]);
}
for (unsigned int component=0;
template <int dim>
void
DataPostprocessor<dim>::
-compute_derived_quantities_scalar (const DataPostprocessorInputs::Scalar<dim> &inputs,
- std::vector<Vector<double> > &computed_quantities) const
+evaluate_scalar_field (const DataPostprocessorInputs::Scalar<dim> &inputs,
+ std::vector<Vector<double> > &computed_quantities) const
{
// for backward compatibility, call the old function.
// this also requires converting the accidental use
template <int dim>
void
DataPostprocessor<dim>::
-compute_derived_quantities_vector (const DataPostprocessorInputs::Vector<dim> &inputs,
- std::vector<Vector<double> > &computed_quantities) const
+evaluate_vector_field (const DataPostprocessorInputs::Vector<dim> &inputs,
+ std::vector<Vector<double> > &computed_quantities) const
{
// for backward compatibility, call the old function.
// this also requires converting the accidental use
virtual
void
- compute_derived_quantities_vector (const DataPostprocessorInputs::Vector<dim> &input_data,
- std::vector<Vector<double> > &computed_quantities) const
+ evaluate_vector_field (const DataPostprocessorInputs::Vector<dim> &input_data,
+ std::vector<Vector<double> > &computed_quantities) const
{
for (unsigned int q=0; q<input_data.solution_values.size(); ++q)
{
virtual
void
- compute_derived_quantities_vector (const DataPostprocessorInputs::Vector<dim> &input_data,
- std::vector<Vector<double> > &computed_quantities) const
+ evaluate_vector_field (const DataPostprocessorInputs::Vector<dim> &input_data,
+ std::vector<Vector<double> > &computed_quantities) const
{
for (unsigned int q=0; q<input_data.solution_values.size(); ++q)
{
// ---------------------------------------------------------------------
//
-// Copyright (C) 2000 - 2015 by the deal.II authors
+// Copyright (C) 2000 - 2016 by the deal.II authors
//
// This file is part of the deal.II library.
//
virtual
void
- compute_derived_quantities_vector (const DataPostprocessorInputs::Vector<dim> &input_data,
- std::vector<Vector<double> > &computed_quantities) const
+ evaluate_vector_field (const DataPostprocessorInputs::Vector<dim> &input_data,
+ std::vector<Vector<double> > &computed_quantities) const
{
for (unsigned int q=0; q<input_data.solution_values.size(); ++q)
{
// ---------------------------------------------------------------------
//
-// Copyright (C) 2000 - 2015 by the deal.II authors
+// Copyright (C) 2000 - 2016 by the deal.II authors
//
// This file is part of the deal.II library.
//