]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix warning about default constructors and default assignment operators
authorDaniel Arndt <arndtd@ornl.gov>
Sun, 3 May 2020 03:53:58 +0000 (23:53 -0400)
committerDaniel Arndt <arndtd@ornl.gov>
Sun, 3 May 2020 03:59:20 +0000 (23:59 -0400)
examples/step-32/step-32.cc
include/deal.II/meshworker/dof_info.h

index 6877affc8445a3e5cab15716f5dee341cbbacca2..19a1be8b212c4d172a72c4b2d536f83995958c70 100644 (file)
@@ -654,7 +654,6 @@ namespace Step32
       struct StokesSystem : public StokesPreconditioner<dim>
       {
         StokesSystem(const FiniteElement<dim> &stokes_fe);
-        StokesSystem(const StokesSystem<dim> &data);
 
         Vector<double> local_rhs;
       };
@@ -665,19 +664,12 @@ namespace Step32
         , local_rhs(stokes_fe.dofs_per_cell)
       {}
 
-      template <int dim>
-      StokesSystem<dim>::StokesSystem(const StokesSystem<dim> &data)
-        : StokesPreconditioner<dim>(data)
-        , local_rhs(data.local_rhs)
-      {}
-
 
 
       template <int dim>
       struct TemperatureMatrix
       {
         TemperatureMatrix(const FiniteElement<dim> &temperature_fe);
-        TemperatureMatrix(const TemperatureMatrix &data);
 
         FullMatrix<double>                   local_mass_matrix;
         FullMatrix<double>                   local_stiffness_matrix;
@@ -694,20 +686,12 @@ namespace Step32
         , local_dof_indices(temperature_fe.dofs_per_cell)
       {}
 
-      template <int dim>
-      TemperatureMatrix<dim>::TemperatureMatrix(const TemperatureMatrix &data)
-        : local_mass_matrix(data.local_mass_matrix)
-        , local_stiffness_matrix(data.local_stiffness_matrix)
-        , local_dof_indices(data.local_dof_indices)
-      {}
-
 
 
       template <int dim>
       struct TemperatureRHS
       {
         TemperatureRHS(const FiniteElement<dim> &temperature_fe);
-        TemperatureRHS(const TemperatureRHS &data);
 
         Vector<double>                       local_rhs;
         std::vector<types::global_dof_index> local_dof_indices;
@@ -722,13 +706,6 @@ namespace Step32
         , matrix_for_bc(temperature_fe.dofs_per_cell,
                         temperature_fe.dofs_per_cell)
       {}
-
-      template <int dim>
-      TemperatureRHS<dim>::TemperatureRHS(const TemperatureRHS &data)
-        : local_rhs(data.local_rhs)
-        , local_dof_indices(data.local_dof_indices)
-        , matrix_for_bc(data.matrix_for_bc)
-      {}
     } // namespace CopyData
   }   // namespace Assembly
 
index 177c7a20e0fec5de500f4ff0030de61d3a4b06d2..822801b991c353d2401fc815c3087756c1464379 100644 (file)
@@ -234,6 +234,12 @@ namespace MeshWorker
      */
     DoFInfoBox(const DoFInfoBox<dim, DOFINFO> &);
 
+    /**
+     * Copy assignment operator, taking another object as seed.
+     */
+    DoFInfoBox &
+    operator=(const DoFInfoBox<dim, DOFINFO> &);
+
     /**
      * Reset all the availability flags.
      */
@@ -453,6 +459,23 @@ namespace MeshWorker
   }
 
 
+  template <int dim, class DOFINFO>
+  inline DoFInfoBox<dim, DOFINFO> &
+  DoFInfoBox<dim, DOFINFO>::operator=(const DoFInfoBox<dim, DOFINFO> &other)
+  {
+    cell       = other.cell;
+    cell_valid = other.cell_valid;
+    for (unsigned int i : GeometryInfo<dim>::face_indices())
+      {
+        exterior[i]                = other.exterior[i];
+        interior[i]                = other.interior[i];
+        interior_face_available[i] = false;
+        exterior_face_available[i] = false;
+      }
+    return *this;
+  }
+
+
   template <int dim, class DOFINFO>
   inline void
   DoFInfoBox<dim, DOFINFO>::reset()

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.