From 91bc06334e33c3a2c74144deba9965e5f96772a4 Mon Sep 17 00:00:00 2001 From: Wolfgang Bangerth Date: Wed, 7 Dec 2016 10:14:19 -0700 Subject: [PATCH] Rename DataPostprocessor functions. 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. --- doc/news/changes.h | 3 +- examples/step-29/step-29.cc | 10 +++--- examples/step-32/step-32.cc | 10 +++--- examples/step-33/step-33.cc | 6 ++-- examples/step-47/step-47.cc | 6 ++-- include/deal.II/numerics/data_postprocessor.h | 33 +++++++++---------- source/numerics/data_out.cc | 8 ++--- source/numerics/data_out_faces.cc | 8 ++--- source/numerics/data_out_rotation.cc | 8 ++--- source/numerics/data_postprocessor.cc | 8 ++--- .../numerics/data_out_postprocessor_01-new.cc | 4 +-- .../data_out_postprocessor_scalar_01-new.cc | 4 +-- .../data_out_postprocessor_scalar_01.cc | 2 +- .../data_out_postprocessor_vector_01-new.cc | 4 +-- .../data_out_postprocessor_vector_01.cc | 2 +- 15 files changed, 54 insertions(+), 62 deletions(-) diff --git a/doc/news/changes.h b/doc/news/changes.h index 401747de27..47a7e2f067 100644 --- a/doc/news/changes.h +++ b/doc/news/changes.h @@ -270,7 +270,8 @@ inconvenience this causes. 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 diff --git a/examples/step-29/step-29.cc b/examples/step-29/step-29.cc index 4732319ecc..1d6acef8ba 100644 --- a/examples/step-29/step-29.cc +++ b/examples/step-29/step-29.cc @@ -273,8 +273,8 @@ namespace Step29 // 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 @@ -298,11 +298,9 @@ namespace Step29 public: ComputeIntensity (); - using DataPostprocessorScalar::compute_derived_quantities_vector; - virtual void - compute_derived_quantities_vector + evaluate_vector_field (const DataPostprocessorInputs::Vector &inputs, std::vector > &computed_quantities) const; }; @@ -343,7 +341,7 @@ namespace Step29 // are provided. template void - ComputeIntensity::compute_derived_quantities_vector + ComputeIntensity::evaluate_vector_field (const DataPostprocessorInputs::Vector &inputs, std::vector > &computed_quantities) const { diff --git a/examples/step-32/step-32.cc b/examples/step-32/step-32.cc index 6d0812ce23..00e7ee0e1f 100644 --- a/examples/step-32/step-32.cc +++ b/examples/step-32/step-32.cc @@ -3136,8 +3136,8 @@ namespace Step32 // 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. @@ -3148,11 +3148,9 @@ namespace Step32 Postprocessor (const unsigned int partition, const double minimal_pressure); - using DataPostprocessor::compute_derived_quantities_vector; - virtual void - compute_derived_quantities_vector + evaluate_vector_field (const DataPostprocessorInputs::Vector &inputs, std::vector > &computed_quantities) const; @@ -3240,7 +3238,7 @@ namespace Step32 template void BoussinesqFlowProblem::Postprocessor:: - compute_derived_quantities_vector + evaluate_vector_field (const DataPostprocessorInputs::Vector &inputs, std::vector > &computed_quantities) const { diff --git a/examples/step-33/step-33.cc b/examples/step-33/step-33.cc index 1ce831f0fa..9aa53c9235 100644 --- a/examples/step-33/step-33.cc +++ b/examples/step-33/step-33.cc @@ -541,11 +541,9 @@ namespace Step33 public: Postprocessor (const bool do_schlieren_plot); - using DataPostprocessor::compute_derived_quantities_vector; - virtual void - compute_derived_quantities_vector + evaluate_vector_field (const DataPostprocessorInputs::Vector &inputs, std::vector > &computed_quantities) const; @@ -588,7 +586,7 @@ namespace Step33 template void EulerEquations::Postprocessor:: - compute_derived_quantities_vector + evaluate_vector_field (const DataPostprocessorInputs::Vector &inputs, std::vector > &computed_quantities) const { diff --git a/examples/step-47/step-47.cc b/examples/step-47/step-47.cc index 226c03f70e..98cd5b4f07 100644 --- a/examples/step-47/step-47.cc +++ b/examples/step-47/step-47.cc @@ -903,11 +903,9 @@ namespace Step47 class Postprocessor : public DataPostprocessor { public: - using DataPostprocessor::compute_derived_quantities_vector; - virtual void - compute_derived_quantities_vector + evaluate_vector_field (const dealii::DataPostprocessorInputs::Vector &inputs, std::vector > &computed_quantities) const; @@ -954,7 +952,7 @@ namespace Step47 template void Postprocessor:: - compute_derived_quantities_vector + evaluate_vector_field (const dealii::DataPostprocessorInputs::Vector &inputs, std::vector > &computed_quantities) const { diff --git a/include/deal.II/numerics/data_postprocessor.h b/include/deal.II/numerics/data_postprocessor.h index de652967af..4348172a75 100644 --- a/include/deal.II/numerics/data_postprocessor.h +++ b/include/deal.II/numerics/data_postprocessor.h @@ -39,7 +39,7 @@ namespace DataPostprocessorInputs { /** * 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 @@ -122,7 +122,7 @@ namespace DataPostprocessorInputs /** * 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 @@ -255,8 +255,8 @@ namespace DataPostprocessorInputs * 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 @@ -279,8 +279,7 @@ namespace DataPostprocessorInputs * 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 @@ -326,12 +325,12 @@ public: */ virtual void - compute_derived_quantities_scalar (const DataPostprocessorInputs::Scalar &input_data, - std::vector > &computed_quantities) const; + evaluate_scalar_field (const DataPostprocessorInputs::Scalar &input_data, + std::vector > &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. * @@ -359,18 +358,18 @@ public: std::vector > &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 &input_data, - std::vector > &computed_quantities) const; + evaluate_vector_field (const DataPostprocessorInputs::Vector &input_data, + std::vector > &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. * @@ -452,8 +451,8 @@ public: * 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. * @@ -526,8 +525,8 @@ private: * 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. diff --git a/source/numerics/data_out.cc b/source/numerics/data_out.cc index 516b096519..eaa4d63f53 100644 --- a/source/numerics/data_out.cc +++ b/source/numerics/data_out.cc @@ -169,8 +169,8 @@ build_one_patch 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 { @@ -192,8 +192,8 @@ build_one_patch 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 - compute_derived_quantities_scalar(data.patch_values_scalar, - data.postprocessed_values[dataset]); + evaluate_scalar_field(data.patch_values_scalar, + data.postprocessed_values[dataset]); } else { @@ -195,8 +195,8 @@ build_one_patch (const FaceDescriptor *cell_and_face, 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 - compute_derived_quantities_scalar(data.patch_values_scalar, - data.postprocessed_values[dataset]); + evaluate_scalar_field(data.patch_values_scalar, + data.postprocessed_values[dataset]); } else { @@ -241,8 +241,8 @@ build_one_patch (const cell_iterator std::vector > 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; diff --git a/source/numerics/data_postprocessor.cc b/source/numerics/data_postprocessor.cc index a3ab1cb092..8d0154db19 100644 --- a/source/numerics/data_postprocessor.cc +++ b/source/numerics/data_postprocessor.cc @@ -30,8 +30,8 @@ DataPostprocessor::~DataPostprocessor() template void DataPostprocessor:: -compute_derived_quantities_scalar (const DataPostprocessorInputs::Scalar &inputs, - std::vector > &computed_quantities) const +evaluate_scalar_field (const DataPostprocessorInputs::Scalar &inputs, + std::vector > &computed_quantities) const { // for backward compatibility, call the old function. // this also requires converting the accidental use @@ -67,8 +67,8 @@ compute_derived_quantities_scalar (const std::vector &/*solution template void DataPostprocessor:: -compute_derived_quantities_vector (const DataPostprocessorInputs::Vector &inputs, - std::vector > &computed_quantities) const +evaluate_vector_field (const DataPostprocessorInputs::Vector &inputs, + std::vector > &computed_quantities) const { // for backward compatibility, call the old function. // this also requires converting the accidental use diff --git a/tests/numerics/data_out_postprocessor_01-new.cc b/tests/numerics/data_out_postprocessor_01-new.cc index 228a73522a..dc5d6f8119 100644 --- a/tests/numerics/data_out_postprocessor_01-new.cc +++ b/tests/numerics/data_out_postprocessor_01-new.cc @@ -150,8 +150,8 @@ public: virtual void - compute_derived_quantities_vector (const DataPostprocessorInputs::Vector &input_data, - std::vector > &computed_quantities) const + evaluate_vector_field (const DataPostprocessorInputs::Vector &input_data, + std::vector > &computed_quantities) const { for (unsigned int q=0; q &input_data, - std::vector > &computed_quantities) const + evaluate_vector_field (const DataPostprocessorInputs::Vector &input_data, + std::vector > &computed_quantities) const { for (unsigned int q=0; q &input_data, - std::vector > &computed_quantities) const + evaluate_vector_field (const DataPostprocessorInputs::Vector &input_data, + std::vector > &computed_quantities) const { for (unsigned int q=0; q