From: David Wells Date: Sat, 27 Jan 2018 18:23:39 +0000 (-0500) Subject: Fix initialization of protected/private solver members. X-Git-Tag: v9.0.0-rc1~505^2~3 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=95ac5c1f4a3b42a37b14e66cdfc001dd5e66edc4;p=dealii.git Fix initialization of protected/private solver members. This was found by cppcheck. --- diff --git a/include/deal.II/lac/solver_bicgstab.h b/include/deal.II/lac/solver_bicgstab.h index 43093e076d..89089636cb 100644 --- a/include/deal.II/lac/solver_bicgstab.h +++ b/include/deal.II/lac/solver_bicgstab.h @@ -30,6 +30,54 @@ DEAL_II_NAMESPACE_OPEN /*!@addtogroup Solvers */ /*@{*/ +namespace internal +{ + /** + * Class containing the non-parameter non-template values used by the + * SolverBicgstab class. + */ + class SolverBicgstabData + { + protected: + /** + * Auxiliary value. + */ + double alpha; + /** + * Auxiliary value. + */ + double beta; + /** + * Auxiliary value. + */ + double omega; + /** + * Auxiliary value. + */ + double rho; + /** + * Auxiliary value. + */ + double rhobar; + + /** + * Current iteration step. + */ + unsigned int step; + + /** + * Residual. + */ + double res; + + /** + * Default constructor. This is protected so that only SolverBicgstab can + * create instances. + */ + SolverBicgstabData(); + }; +} + /** * Bicgstab algorithm by van der Vorst. * @@ -70,7 +118,8 @@ DEAL_II_NAMESPACE_OPEN * */ template > -class SolverBicgstab : public Solver +class SolverBicgstab : public Solver, + protected internal::SolverBicgstabData { public: /** @@ -137,89 +186,66 @@ public: const PreconditionerType &preconditioner); protected: - /** - * Computation of the stopping criterion. - */ - template - double criterion (const MatrixType &A, const VectorType &x, const VectorType &b); - - /** - * Interface for derived class. This function gets the current iteration - * vector, the residual and the update vector in each step. It can be used - * for graphical output of the convergence history. - */ - virtual void print_vectors(const unsigned int step, - const VectorType &x, - const VectorType &r, - const VectorType &d) const; - /** * Auxiliary vector. */ VectorType *Vx; + /** * Auxiliary vector. */ typename VectorMemory::Pointer Vr; + /** * Auxiliary vector. */ typename VectorMemory::Pointer Vrbar; + /** * Auxiliary vector. */ typename VectorMemory::Pointer Vp; + /** * Auxiliary vector. */ typename VectorMemory::Pointer Vy; + /** * Auxiliary vector. */ typename VectorMemory::Pointer Vz; + /** * Auxiliary vector. */ typename VectorMemory::Pointer Vt; + /** * Auxiliary vector. */ typename VectorMemory::Pointer Vv; + /** * Right hand side vector. */ const VectorType *Vb; /** - * Auxiliary value. - */ - double alpha; - /** - * Auxiliary value. - */ - double beta; - /** - * Auxiliary value. - */ - double omega; - /** - * Auxiliary value. - */ - double rho; - /** - * Auxiliary value. - */ - double rhobar; - - /** - * Current iteration step. + * Computation of the stopping criterion. */ - unsigned int step; + template + double criterion (const MatrixType &A, const VectorType &x, const VectorType &b); /** - * Residual. + * Interface for derived class. This function gets the current iteration + * vector, the residual and the update vector in each step. It can be used + * for graphical output of the convergence history. */ - double res; + virtual void print_vectors(const unsigned int step, + const VectorType &x, + const VectorType &r, + const VectorType &d) const; /** * Additional parameters. @@ -280,12 +306,28 @@ SolverBicgstab::IterationResult::IterationResult {} + +internal::SolverBicgstabData::SolverBicgstabData () + : + alpha(0.), + beta(0.), + omega(0.), + rho(0.), + rhobar(0.), + step(numbers::invalid_unsigned_int), + res(numbers::signaling_nan()) +{} + + + template SolverBicgstab::SolverBicgstab (SolverControl &cn, VectorMemory &mem, const AdditionalData &data) : Solver(cn,mem), + Vx (nullptr), + Vb (nullptr), additional_data(data) {} @@ -296,13 +338,8 @@ SolverBicgstab::SolverBicgstab (SolverControl &cn, const AdditionalData &data) : Solver(cn), - alpha(0.), - beta(0.), - omega(0.), - rho(0.), - rhobar(0.), - step(numbers::invalid_unsigned_int), - res(numbers::signaling_nan()), + Vx (nullptr), + Vb (nullptr), additional_data(data) {} diff --git a/include/deal.II/lac/solver_minres.h b/include/deal.II/lac/solver_minres.h index 3c7f79566e..e1fadba6e9 100644 --- a/include/deal.II/lac/solver_minres.h +++ b/include/deal.II/lac/solver_minres.h @@ -152,7 +152,8 @@ SolverMinRes::SolverMinRes (SolverControl &cn, VectorMemory &mem, const AdditionalData &) : - Solver(cn,mem) + Solver(cn,mem), + res2(numbers::signaling_nan()) {} diff --git a/include/deal.II/lac/solver_qmrs.h b/include/deal.II/lac/solver_qmrs.h index cb885e8a8e..dd70e588c7 100644 --- a/include/deal.II/lac/solver_qmrs.h +++ b/include/deal.II/lac/solver_qmrs.h @@ -197,7 +197,6 @@ protected: AdditionalData additional_data; private: - /** * A structure returned by the iterate() function representing what it found * is happening during the iteration. @@ -238,31 +237,35 @@ private: #ifndef DOXYGEN + template SolverQMRS::IterationResult::IterationResult (const SolverControl::State state, const double last_residual) : state (state), last_residual (last_residual) -{ -} +{} + + template -SolverQMRS::SolverQMRS (SolverControl &cn, +SolverQMRS::SolverQMRS (SolverControl &cn, VectorMemory &mem, - const AdditionalData &data) + const AdditionalData &data) : Solver (cn, mem), - additional_data (data) + additional_data (data), + step (0) { } template -SolverQMRS::SolverQMRS (SolverControl &cn, +SolverQMRS::SolverQMRS (SolverControl &cn, const AdditionalData &data) : Solver (cn), - additional_data (data) + additional_data (data), + step (0) { }