<h3>Specific improvements</h3>
<ol>
+<li> Changed: The DataPostprocessor class previously required users of this
+class to overload DataPostprocessor::get_names(),
+DataPostprocessor::get_data_component_interpretation()
+and DataPostprocessor::n_output_variables(). The latter function is redundant
+since its output must equal the length of the arrays returned by the
+first two of these functions. It has therefore been removed.
+<br>
+(Wolfgang Bangerth, 2011/12/15)
<li> Improved: Objects of the type LogStream::Prefix can now be used
as a safe implementation of the push and pop mechanism for log
virtual std::vector<std::string> get_names () const;
virtual UpdateFlags get_needed_update_flags () const;
- virtual unsigned int n_output_variables () const;
};
// The <code>get_names</code>
return update_values;
}
- // To allow the caller to find out
- // how many derived quantities are
- // returned by the postprocessor, the
- // <code>n_output_variables</code>
- // function is used. Since we compute
- // only $|u|$, the correct value to
- // return in our case is just 1:
- template <int dim>
- unsigned int
- ComputeIntensity<dim>::n_output_variables () const
- {
- return 1;
- }
-
// The actual prostprocessing happens
// in the following function. Its
virtual std::vector<std::string> get_names () const;
- virtual unsigned int n_output_variables() const;
-
virtual
std::vector<DataComponentInterpretation::DataComponentInterpretation>
get_data_component_interpretation () const;
}
- template <int dim>
- unsigned int
- BoussinesqFlowProblem<dim>::Postprocessor::n_output_variables() const
- {
- return get_names().size();
- }
-
-
template <int dim>
std::vector<DataComponentInterpretation::DataComponentInterpretation>
BoussinesqFlowProblem<dim>::Postprocessor::
Assert (duh.size() == n_quadrature_points, ExcInternalError());
Assert (computed_quantities.size() == n_quadrature_points, ExcInternalError());
Assert (uh[0].size() == dim+2, ExcInternalError());
- Assert (computed_quantities[0].size()==n_output_variables(),ExcInternalError());
for (unsigned int q=0; q<n_quadrature_points; ++q)
{
virtual UpdateFlags get_needed_update_flags () const;
- virtual unsigned int n_output_variables() const;
-
private:
const bool do_schlieren_plot;
};
}
-
- template <int dim>
- unsigned int
- EulerEquations<dim>::Postprocessor::
- n_output_variables () const
- {
- if (do_schlieren_plot == true)
- return dim+2;
- else
- return dim+1;
- }
-
-
// @sect3{Run time parameter handling}
// Our next job is to define a few
virtual std::vector<std::string> get_names () const;
- virtual unsigned int n_output_variables() const;
-
virtual
std::vector<DataComponentInterpretation::DataComponentInterpretation>
get_data_component_interpretation () const;
}
- template <int dim>
- unsigned int
- Postprocessor<dim>::n_output_variables() const
- {
- return get_names().size();
- }
-
-
template <int dim>
std::vector<DataComponentInterpretation::DataComponentInterpretation>
Postprocessor<dim>::
const unsigned int n_quadrature_points = uh.size();
Assert (computed_quantities.size() == n_quadrature_points, ExcInternalError());
Assert (uh[0].size() == 2, ExcInternalError());
- Assert (computed_quantities[0].size()==n_output_variables(),ExcInternalError());
for (unsigned int q=0; q<n_quadrature_points; ++q)
{
* Return the vector of strings describing
* the names of the computed quantities.
*/
- virtual std::vector<std::string> get_names () const=0;
+ virtual std::vector<std::string> get_names () const = 0;
/**
* This functions returns
* also ask for a update of normals via the
* @p update_normal_vectors flag.
*/
- virtual UpdateFlags get_needed_update_flags () const=0;
-
- /**
- * Number of postprocessed
- * variables. Has to match the
- * number of entries filled by
- * compute_derived_quantities_scalar()
- * or
- * compute_derived_quantities_vector()
- * as well as the size of the
- * vector of names returned by
- * get_names().
- */
- virtual unsigned int n_output_variables() const=0;
-
+ virtual UpdateFlags get_needed_update_flags () const = 0;
};
DEAL_II_NAMESPACE_CLOSE
names(data_postprocessor->get_names()),
data_component_interpretation (data_postprocessor->get_data_component_interpretation()),
postprocessor(data_postprocessor, typeid(*this).name()),
- n_output_variables(data_postprocessor->n_output_variables())
+ n_output_variables(names.size())
{
- // if there is a post processor, then we
- // should have gotten the names from the
- // postprocessor. check that the number of
- // elements in the names vector is
- // correct. otherwise there is nothing for
- // us to check
- Assert(data_postprocessor->n_output_variables()==names.size(),
- ExcDimensionMismatch(data_postprocessor->n_output_variables(),
- names.size()));
- Assert (names.size() == data_component_interpretation.size(),
- ExcDimensionMismatch(data_component_interpretation.size(),
- names.size()));
-
+ Assert (data_postprocessor->get_names().size()
+ ==
+ data_postprocessor->get_data_component_interpretation().size(),
+ ExcDimensionMismatch (data_postprocessor->get_names().size(),
+ data_postprocessor->get_data_component_interpretation().size()));
+
// check that the names use only allowed
// characters
- // check names for invalid characters
for (unsigned int i=0; i<names.size(); ++i)
Assert (names[i].find_first_not_of("abcdefghijklmnopqrstuvwxyz"
"ABCDEFGHIJKLMNOPQRSTUVWXYZ"
// components are independent scalars
return
std::vector<DataComponentInterpretation::DataComponentInterpretation>
- (n_output_variables(),
+ (get_names().size(),
DataComponentInterpretation::component_is_scalar);
}