From: hartmann Date: Mon, 5 Jul 1999 11:19:10 +0000 (+0000) Subject: add the weight function argument also in the scalar #integrate_difference# function X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=928ddafda6a550f1e9289936ec7460a4c5b16d44;p=dealii-svn.git add the weight function argument also in the scalar #integrate_difference# function git-svn-id: https://svn.dealii.org/trunk@1533 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/deal.II/include/numerics/vectors.h b/deal.II/deal.II/include/numerics/vectors.h index 95fda10a0e..2e4a8f1e37 100644 --- a/deal.II/deal.II/include/numerics/vectors.h +++ b/deal.II/deal.II/include/numerics/vectors.h @@ -434,20 +434,6 @@ class VectorTools // * accuracy of the #double# data type is * used. * - * See the general documentation of this - * class for more information. - */ - static void integrate_difference (const DoFHandler &dof, - const Vector &fe_function, - const Function &exact_solution, - Vector &difference, - const Quadrature &q, - const NormType &norm); - - /** - * Compute the error for the solution of a system. - * See the other #integrate_difference#. - * * The additional argument #weight# allows * to evaluate weighted norms. This is useful * for weighting the error of different parts @@ -460,6 +446,21 @@ class VectorTools // * solutions. * By default, no weighting function is given, * i.e. weight=1 in the whole domain. + * + * See the general documentation of this + * class for more information. + */ + static void integrate_difference (const DoFHandler &dof, + const Vector &fe_function, + const Function &exact_solution, + Vector &difference, + const Quadrature &q, + const NormType &norm, + const Function *weight=0); + + /** + * Compute the error for the solution of a system. + * See the other #integrate_difference#. */ static void integrate_difference (const DoFHandler &dof, const Vector &fe_function, diff --git a/deal.II/deal.II/source/numerics/vectors.cc b/deal.II/deal.II/source/numerics/vectors.cc index 2dc3a163e0..5bcdaeb20e 100644 --- a/deal.II/deal.II/source/numerics/vectors.cc +++ b/deal.II/deal.II/source/numerics/vectors.cc @@ -693,7 +693,9 @@ void VectorTools::integrate_difference (const DoFHandler &dof, const Function &exact_solution, Vector &difference, const Quadrature &q, - const NormType &norm) { + const NormType &norm, + const Function *weight=0) +{ const FiniteElement &fe = dof.get_fe(); difference.reinit (dof.get_tria().n_active_cells()); @@ -782,6 +784,17 @@ void VectorTools::integrate_difference (const DoFHandler &dof, Assert (false, ExcNotImplemented()); }; + // now weight the values + // at the quadrature points due + // to the weighting function + if (weight) + { + vector w(n_q_points); + weight->value_list(fe_values.get_quadrature_points(),w); + for (unsigned int q=0; q::integrate_difference (const DoFHandler &dof, for (unsigned int i=0; i w(n_q_points); + weight->value_list(fe_values.get_quadrature_points(),w); + for (unsigned int q=0; q