]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix initialization of protected/private solver members.
authorDavid Wells <wellsd2@rpi.edu>
Sat, 27 Jan 2018 18:23:39 +0000 (13:23 -0500)
committerDavid Wells <wellsd2@rpi.edu>
Sat, 27 Jan 2018 22:47:54 +0000 (17:47 -0500)
This was found by cppcheck.

include/deal.II/lac/solver_bicgstab.h
include/deal.II/lac/solver_minres.h
include/deal.II/lac/solver_qmrs.h

index 43093e076d89a6d56b5b14af518b94aed5335925..89089636cb641ef3f0a3f29fbd2c5da911ade528 100644 (file)
@@ -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 <typename VectorType = Vector<double> >
-class SolverBicgstab : public Solver<VectorType>
+class SolverBicgstab : public Solver<VectorType>,
+  protected internal::SolverBicgstabData
 {
 public:
   /**
@@ -137,89 +186,66 @@ public:
          const PreconditionerType &preconditioner);
 
 protected:
-  /**
-   * Computation of the stopping criterion.
-   */
-  template <typename MatrixType>
-  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<VectorType>::Pointer Vr;
+
   /**
    * Auxiliary vector.
    */
   typename VectorMemory<VectorType>::Pointer Vrbar;
+
   /**
    * Auxiliary vector.
    */
   typename VectorMemory<VectorType>::Pointer Vp;
+
   /**
    * Auxiliary vector.
    */
   typename VectorMemory<VectorType>::Pointer Vy;
+
   /**
    * Auxiliary vector.
    */
   typename VectorMemory<VectorType>::Pointer Vz;
+
   /**
    * Auxiliary vector.
    */
   typename VectorMemory<VectorType>::Pointer Vt;
+
   /**
    * Auxiliary vector.
    */
   typename VectorMemory<VectorType>::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 <typename MatrixType>
+  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<VectorType>::IterationResult::IterationResult
 {}
 
 
+
+internal::SolverBicgstabData::SolverBicgstabData ()
+  :
+  alpha(0.),
+  beta(0.),
+  omega(0.),
+  rho(0.),
+  rhobar(0.),
+  step(numbers::invalid_unsigned_int),
+  res(numbers::signaling_nan<double>())
+{}
+
+
+
 template <typename VectorType>
 SolverBicgstab<VectorType>::SolverBicgstab (SolverControl            &cn,
                                             VectorMemory<VectorType> &mem,
                                             const AdditionalData     &data)
   :
   Solver<VectorType>(cn,mem),
+  Vx (nullptr),
+  Vb (nullptr),
   additional_data(data)
 {}
 
@@ -296,13 +338,8 @@ SolverBicgstab<VectorType>::SolverBicgstab (SolverControl        &cn,
                                             const AdditionalData &data)
   :
   Solver<VectorType>(cn),
-  alpha(0.),
-  beta(0.),
-  omega(0.),
-  rho(0.),
-  rhobar(0.),
-  step(numbers::invalid_unsigned_int),
-  res(numbers::signaling_nan<double>()),
+  Vx (nullptr),
+  Vb (nullptr),
   additional_data(data)
 {}
 
index 3c7f79566e3266abed7765a89a05783975ffb147..e1fadba6e9dfb6ce0cccfae38f95627e31633934 100644 (file)
@@ -152,7 +152,8 @@ SolverMinRes<VectorType>::SolverMinRes (SolverControl            &cn,
                                         VectorMemory<VectorType> &mem,
                                         const AdditionalData &)
   :
-  Solver<VectorType>(cn,mem)
+  Solver<VectorType>(cn,mem),
+  res2(numbers::signaling_nan<double>())
 {}
 
 
index cb885e8a8eb562ff66a972c12b3adc2ff0c7d790..dd70e588c792c3803d317f4aa218322dd55c022d 100644 (file)
@@ -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 <class VectorType>
 SolverQMRS<VectorType>::IterationResult::IterationResult (const SolverControl::State state,
                                                           const double last_residual)
   :
   state (state),
   last_residual (last_residual)
-{
-}
+{}
+
+
 
 template <class VectorType>
-SolverQMRS<VectorType>::SolverQMRS (SolverControl &cn,
+SolverQMRS<VectorType>::SolverQMRS (SolverControl            &cn,
                                     VectorMemory<VectorType> &mem,
-                                    const AdditionalData &data)
+                                    const AdditionalData     &data)
   :
   Solver<VectorType> (cn, mem),
-  additional_data (data)
+  additional_data (data),
+  step (0)
 {
 }
 
 template <class VectorType>
-SolverQMRS<VectorType>::SolverQMRS (SolverControl &cn,
+SolverQMRS<VectorType>::SolverQMRS (SolverControl        &cn,
                                     const AdditionalData &data)
   :
   Solver<VectorType> (cn),
-  additional_data (data)
+  additional_data (data),
+  step (0)
 {
 }
 

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.