]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix default definitions in header files
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Fri, 15 Sep 2017 07:56:57 +0000 (09:56 +0200)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Fri, 15 Sep 2017 08:23:10 +0000 (10:23 +0200)
12 files changed:
include/deal.II/base/parallel.h
include/deal.II/base/point.h
include/deal.II/base/table.h
include/deal.II/grid/tria_boundary.h
include/deal.II/lac/petsc_parallel_block_vector.h
include/deal.II/lac/pointer_matrix.h
include/deal.II/lac/precondition_block.h
include/deal.II/lac/precondition_block.templates.h
include/deal.II/lac/solver_cg.h
include/deal.II/matrix_free/operators.h
source/base/function_parser.cc
source/grid/tria_boundary.cc

index 74a2c3b18bc481850e5f7f4a0f82a5eeff1d8b82..4a683dd446e6d690c602f41afd0023ade7c63fad 100644 (file)
@@ -477,7 +477,7 @@ namespace parallel
      * Destructor. Made virtual to ensure that derived classes also have
      * virtual destructors.
      */
-    virtual ~ParallelForInteger ();
+    virtual ~ParallelForInteger () = default;
 
     /**
      * This function runs the for loop over the given range
@@ -843,11 +843,6 @@ namespace parallel
 #endif
 
 
-  inline
-  ParallelForInteger::~ParallelForInteger () = default;
-
-
-
   inline
   void
   ParallelForInteger::apply_parallel (const std::size_t begin,
index 40b387abcd750d39ece0e96a1a440d35a050cbc8..103582dbde5899815475e6557e3a4e791ff66b44 100644 (file)
@@ -253,7 +253,8 @@ public:
 // and we can't use 'Point<dim,Number>::Point () = default' here.
 template <int dim, typename Number>
 inline
-Point<dim,Number>::Point () = default;
+Point<dim,Number>::Point ()
+{}
 
 
 
index 80228463892134a3b169de61759beab85c3e5219..cfebb63c632cffe249f00197b4fbd80b93684479 100644 (file)
@@ -452,7 +452,7 @@ public:
   /**
    * Destructor. Free allocated memory.
    */
-  ~TableBase ();
+  ~TableBase () = default;
 
   /**
    * Assignment operator. Copy all elements of <tt>src</tt> into the matrix.
@@ -1797,12 +1797,6 @@ namespace internal
 
 
 
-template <int N, typename T>
-inline
-TableBase<N,T>::~TableBase () = default;
-
-
-
 template <int N, typename T>
 inline
 TableBase<N,T> &
index 84d20d893c615724e436a1fff8ecfa01503e2d78..65e170673c678370bb05945b29bcc588e21b1176 100644 (file)
@@ -241,7 +241,7 @@ public:
   /**
    * Default constructor. Some compilers require this for some reasons.
    */
-  StraightBoundary () = default;
+  StraightBoundary ();
 
   /**
    * Let the new point be the arithmetic mean of the two vertices of the line.
index f3162c43d78de904954ddda67dc7c8efec4c3807..5c6bc666e8f050b61f62f4d336dbb64c7dc3f686 100644 (file)
@@ -136,7 +136,7 @@ namespace PETScWrappers
       /**
        * Destructor. Clears memory
        */
-      ~BlockVector ();
+      ~BlockVector () = default;
 
       /**
        * Copy operator: fill all components of the vector that are locally
@@ -356,8 +356,6 @@ namespace PETScWrappers
       return *this;
     }
 
-    inline
-    BlockVector::~BlockVector () = default;
 
 
     inline
index 50cf9e7f7c943ad0b2fab06ba4db9c12aa21a438..36750f2dcff0ecd2731dffd91cc42aef1b3ac102 100644 (file)
@@ -68,7 +68,7 @@ public:
    * of any derived class is called whenever a pointer-to-base-class object is
    * destroyed.
    */
-  virtual ~PointerMatrixBase ();
+  virtual ~PointerMatrixBase () = default;
 
   /**
    * Reset the object to its original state.
@@ -581,14 +581,6 @@ new_pointer_matrix_base(const TridiagonalMatrix<numberm> &matrix, const Vector<n
 /*@}*/
 //---------------------------------------------------------------------------
 
-template <typename VectorType>
-inline
-PointerMatrixBase<VectorType>::~PointerMatrixBase () = default;
-
-
-
-//----------------------------------------------------------------------//
-
 
 template <typename MatrixType, typename VectorType>
 PointerMatrix<MatrixType, VectorType>::PointerMatrix (const MatrixType *M)
index f74ed6874bc3f9bf15767807a920a5db328ca914..aa83f867b4bb92c17024b757a8fec48a1f28cb61 100644 (file)
@@ -155,7 +155,7 @@ public:
   /**
    * Destructor.
    */
-  ~PreconditionBlock();
+  ~PreconditionBlock() = default;
 
   /**
    * Initialize matrix and block size.  We store the matrix and the block size
index d2c8f8c897e69df3c8bbbc9b4ff819d9dfc4ad57..6d0a65ec47a4d0ee0c0f82daecc6553ecfeb06ec 100644 (file)
@@ -55,9 +55,6 @@ PreconditionBlock<MatrixType,inverse_type>::PreconditionBlock (bool store)
 {}
 
 
-template <typename MatrixType, typename inverse_type>
-PreconditionBlock<MatrixType,inverse_type>::~PreconditionBlock () = default;
-
 
 template <typename MatrixType, typename inverse_type>
 void PreconditionBlock<MatrixType,inverse_type>::clear ()
index 8fca65a3049a38b356a40788fe2d78a39fefe9a8..3d784e5867a1a4176c2888b1e4b454cee4c8c480 100644 (file)
@@ -120,7 +120,7 @@ public:
   /**
    * Virtual destructor.
    */
-  virtual ~SolverCG ();
+  virtual ~SolverCG () = default;
 
   /**
    * Solve the linear system $Ax=b$ for x.
@@ -248,11 +248,6 @@ SolverCG<VectorType>::SolverCG (SolverControl        &cn,
 
 
 
-template <typename VectorType>
-SolverCG<VectorType>::~SolverCG () = default;
-
-
-
 template <typename VectorType>
 void
 SolverCG<VectorType>::print_vectors(const unsigned int,
index b2acd5ece4ee2ac2d891557d9638bbaf565db3aa..91a5967ef01c5231591bf06d197c62c5987553bb 100644 (file)
@@ -182,7 +182,7 @@ namespace MatrixFreeOperators
     /**
      * Virtual destructor.
      */
-    virtual ~Base();
+    virtual ~Base() = default;
 
     /**
      * Release all memory and return to a state just like after having called
@@ -882,11 +882,6 @@ namespace MatrixFreeOperators
   }
 
   //----------------- Base operator -----------------------------
-  template <int dim, typename VectorType>
-  Base<dim,VectorType>::~Base () = default;
-
-
-
   template <int dim, typename VectorType>
   Base<dim,VectorType>::Base ()
     :
index 6345390c3487612cb90ccccb8069b5f991b2385e..0a50b67336c23de01f61b80b7b7c0a2c8fa89902 100644 (file)
@@ -57,7 +57,7 @@ FunctionParser<dim>::FunctionParser(const unsigned int n_components,
 
 
 
-// We deliberately delay the definition of the default constructor
+// We deliberately delay the definition of the default destructor
 // so that we don't need to include the definition of mu::Parser
 // in the header file.
 template <int dim>
index aa13cfb072ace865a82b5efb21570d988c18fe00..250dd5dbae7604e6d39932a3c876f7bdd79b34b8 100644 (file)
@@ -180,6 +180,12 @@ get_line_support_points (const unsigned int n_intermediate_points) const
 
 /* -------------------------- StraightBoundary --------------------- */
 
+// At least clang < 3.9.0 complains if we move this definition to its
+// declaration when a 'const StraightBoundary' object is built.
+template <int dim, int spacedim>
+StraightBoundary<dim, spacedim>::StraightBoundary () = default;
+
+
 
 template <int dim, int spacedim>
 Point<spacedim>

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.