]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove deprecated constructors of MGSmoother* classes that take VectorMemory objects...
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Tue, 30 Dec 2014 02:36:42 +0000 (20:36 -0600)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Sun, 4 Jan 2015 21:13:42 +0000 (15:13 -0600)
doc/news/changes.h
examples/step-16/step-16.cc
include/deal.II/multigrid/mg_smoother.h
tests/bits/step-16.cc
tests/bits/step-16b.cc

index ea5e2533e6b0c4d2a2dfadad4766ccc7adb3b555..bd92a5a8918bbfc7361205f1e84508b5be6034ef 100644 (file)
@@ -50,9 +50,9 @@ inconvenience this causes.
     data postprocessor classes implemented in your program that overload these
     functions, you will have to change it in a way that they overload the
     functions of same name but with the evaluation point argument instead.
-  <br>
-  (Wolfgang Bangerth, 2015/01/04)
-  </li>
+  - The constructors of classes MGSmoother, MGSmootherRelaxation and
+    MGSmootherPrecondition that take a VectorMemory object.
+</ol>
 
   <li> Removed: The config.h file no longer exports HAVE_* definitions.
   Those are either entirely removed (for the blas/lapack symbols) or
index a3204e4fe800ead6717285c9d71e5233c4ca81c4..48fc783fd83af3bc26fda10e0815eff009d43252 100644 (file)
@@ -504,11 +504,6 @@ namespace Step16
     // iteration. To this end, we define an appropriate <code>typedef</code>
     // and then setup a smoother object.
     //
-    // Since this smoother needs temporary vectors to store intermediate
-    // results, we need to provide a VectorMemory object. Since these vectors
-    // will be reused over and over, the GrowingVectorMemory is more time
-    // efficient than the PrimitiveVectorMemory class in the current case.
-    //
     // The last step is to initialize the smoother object with our level
     // matrices and to set some smoothing parameters.  The
     // <code>initialize()</code> function can optionally take additional
index 1622faa3b788fd6db10ea08af195329cff1f3990..24fdd64c6c653be20b83af17c0bed4cd14c6a6f9 100644 (file)
@@ -47,26 +47,13 @@ class MGSmoother : public MGSmootherBase<VECTOR>
 {
 public:
   /**
-   * Constructor. Sets smoothing parameters and creates a private
-   * GrowingVectorMemory object to be used to retrieve vectors.
+   * Constructor.
    */
   MGSmoother(const unsigned int steps = 1,
              const bool variable = false,
              const bool symmetric = false,
              const bool transpose = false);
 
-  /**
-   * @deprecated Since GrowingVectorMemory now uses a joint memory pool, it is
-   * recommended to use the constructor without the memory object.
-   *
-   * Constructor. Sets memory and smoothing parameters.
-   */
-  MGSmoother(VectorMemory<VECTOR> &mem,
-             const unsigned int steps = 1,
-             const bool variable = false,
-             const bool symmetric = false,
-             const bool transpose = false) DEAL_II_DEPRECATED;
-
   /**
    * Modify the number of smoothing steps on finest level.
    */
@@ -201,7 +188,7 @@ namespace mg
   {
   public:
     /**
-     * Constructor. Sets memory and smoothing parameters.
+     * Constructor. Sets smoothing parameters.
      */
     SmootherRelaxation(const unsigned int steps = 1,
                        const bool variable = false,
@@ -311,17 +298,6 @@ public:
                        const bool symmetric = false,
                        const bool transpose = false);
 
-  /**
-   * Constructor. Sets memory and smoothing parameters.
-   *
-   * @deprecated Use the constructor without the vector memory object
-   */
-  MGSmootherRelaxation(VectorMemory<VECTOR> &mem,
-                       const unsigned int steps = 1,
-                       const bool variable = false,
-                       const bool symmetric = false,
-                       const bool transpose = false) DEAL_II_DEPRECATED;
-
   /**
    * Initialize for matrices. This function stores pointers to the level
    * matrices and initializes the smoothing operator with the same smoother
@@ -451,17 +427,6 @@ public:
                          const bool symmetric = false,
                          const bool transpose = false);
 
-  /**
-   * Constructor. Sets memory and smoothing parameters.
-   *
-   * @deprecated Use the constructor without the vector memory object
-   */
-  MGSmootherPrecondition(VectorMemory<VECTOR> &mem,
-                         const unsigned int steps = 1,
-                         const bool variable = false,
-                         const bool symmetric = false,
-                         const bool transpose = false) DEAL_II_DEPRECATED;
-
   /**
    * Initialize for matrices. This function stores pointers to the level
    * matrices and initializes the smoothing operator with the same smoother
@@ -584,24 +549,6 @@ MGSmoother<VECTOR>::MGSmoother(
 {}
 
 
-template <class VECTOR>
-inline
-MGSmoother<VECTOR>::MGSmoother(
-  VectorMemory<VECTOR> &mem,
-  const unsigned int steps,
-  const bool variable,
-  const bool symmetric,
-  const bool transpose)
-  :
-  steps(steps),
-  variable(variable),
-  symmetric(symmetric),
-  transpose(transpose),
-  debug(0),
-  mem(&mem)
-{}
-
-
 template <class VECTOR>
 inline void
 MGSmoother<VECTOR>::set_steps (const unsigned int s)
@@ -745,20 +692,6 @@ namespace mg
 
 //----------------------------------------------------------------------//
 
-template <class MATRIX, class RELAX, class VECTOR>
-inline
-MGSmootherRelaxation<MATRIX, RELAX, VECTOR>::MGSmootherRelaxation(
-  VectorMemory<VECTOR> &mem,
-  const unsigned int steps,
-  const bool variable,
-  const bool symmetric,
-  const bool transpose)
-  :
-  MGSmoother<VECTOR>(mem, steps, variable, symmetric, transpose)
-{}
-
-
-
 template <class MATRIX, class RELAX, class VECTOR>
 inline
 MGSmootherRelaxation<MATRIX, RELAX, VECTOR>::MGSmootherRelaxation(
@@ -927,20 +860,6 @@ memory_consumption () const
 
 //----------------------------------------------------------------------//
 
-template <class MATRIX, class PRECONDITIONER, class VECTOR>
-inline
-MGSmootherPrecondition<MATRIX, PRECONDITIONER, VECTOR>::MGSmootherPrecondition(
-  VectorMemory<VECTOR> &mem,
-  const unsigned int steps,
-  const bool variable,
-  const bool symmetric,
-  const bool transpose)
-  :
-  MGSmoother<VECTOR>(mem, steps, variable, symmetric, transpose)
-{}
-
-
-
 template <class MATRIX, class PRECONDITIONER, class VECTOR>
 inline
 MGSmootherPrecondition<MATRIX, PRECONDITIONER, VECTOR>::MGSmootherPrecondition(
index c7b3439b460431803ec83eda01c6d26508308282..273ac93cd843bcc4a725154cf29b8ca29f3c2d97 100644 (file)
@@ -242,8 +242,6 @@ void LaplaceProblem<dim>::assemble_multigrid ()
 template <int dim>
 void LaplaceProblem<dim>::solve ()
 {
-  GrowingVectorMemory<>   vector_memory;
-
   MGTransferPrebuilt<Vector<double> > mg_transfer;
   mg_transfer.build_matrices(mg_dof_handler);
 
@@ -254,7 +252,7 @@ void LaplaceProblem<dim>::solve ()
 
   typedef PreconditionSOR<SparseMatrix<float> > RELAXATION;
   MGSmootherRelaxation<SparseMatrix<float>, RELAXATION, Vector<double> >
-  mg_smoother(vector_memory);
+  mg_smoother;
 
   RELAXATION::AdditionalData smoother_data;
   mg_smoother.initialize(mg_matrices, smoother_data);
index 410f89bdcecbd06875caeb3f02a8630636a72f23..25ba59ddd7515bbf4bd92ce37f247617558f0871 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2005 - 2013 by the deal.II authors
+// Copyright (C) 2005 - 2014 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -245,8 +245,6 @@ void LaplaceProblem<dim>::assemble_multigrid ()
 template <int dim>
 void LaplaceProblem<dim>::solve ()
 {
-  GrowingVectorMemory<>   vector_memory;
-
   MGTransferPrebuilt<Vector<double> > mg_transfer;
   mg_transfer.build_matrices(mg_dof_handler);
 
@@ -257,7 +255,7 @@ void LaplaceProblem<dim>::solve ()
 
   typedef PreconditionSOR<SparseMatrix<float> > RELAXATION;
   MGSmootherRelaxation<SparseMatrix<float>, RELAXATION, Vector<double> >
-  mg_smoother(vector_memory);
+  mg_smoother;
 
   RELAXATION::AdditionalData smoother_data;
   mg_smoother.initialize(mg_matrices, smoother_data);

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.