]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Make the AffineConstraints scratch storage object-local.
authorWolfgang Bangerth <bangerth@colostate.edu>
Fri, 5 Jun 2020 13:17:46 +0000 (07:17 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Wed, 10 Jun 2020 00:03:58 +0000 (18:03 -0600)
include/deal.II/lac/affine_constraints.h
include/deal.II/lac/affine_constraints.templates.h
source/lac/affine_constraints.cc

index 7f22f9ccbe2318f200c122e10241449c78db8bb2..f0b8167cff2fbd0e7f1b98aa8e92e874cdd3452f 100644 (file)
@@ -1830,6 +1830,10 @@ private:
    */
   bool sorted;
 
+  mutable Threads::ThreadLocalStorage<
+    internal::AffineConstraints::ScratchData<number>>
+    scratch_data;
+
   /**
    * Internal function to calculate the index of line @p line_n in the vector
    * lines_cache using local_lines.
index 5c439dadaa5fd3ec20178ee0c9cba159a23660bb..c2f4eed49844659decf235389641992166ec46e3 100644 (file)
@@ -2796,11 +2796,6 @@ namespace internal
       private:
         ScratchData<number> *my_scratch_data;
       };
-
-      /**
-       * The actual data object that contains a scratch data for each thread.
-       */
-      static Threads::ThreadLocalStorage<ScratchData<number>> scratch_data;
     };
 
 
@@ -3560,9 +3555,8 @@ AffineConstraints<number>::distribute_local_to_global(
 
   const size_type n_local_dofs = local_dof_indices.size();
 
-  typename internal::AffineConstraints::AffineConstraintsData<number>::
-    ScratchDataAccessor scratch_data(
-      internal::AffineConstraints::AffineConstraintsData<number>::scratch_data);
+  typename internal::AffineConstraints::AffineConstraintsData<
+    number>::ScratchDataAccessor scratch_data(this->scratch_data);
 
   internal::AffineConstraints::GlobalRowsFromLocal<number> &global_rows =
     scratch_data->global_rows;
@@ -3712,9 +3706,8 @@ AffineConstraints<number>::distribute_local_to_global(
     }
   Assert(sorted == true, ExcMatrixNotClosed());
 
-  typename internal::AffineConstraints::AffineConstraintsData<number>::
-    ScratchDataAccessor scratch_data(
-      internal::AffineConstraints::AffineConstraintsData<number>::scratch_data);
+  typename internal::AffineConstraints::AffineConstraintsData<
+    number>::ScratchDataAccessor scratch_data(this->scratch_data);
 
   const size_type n_local_dofs = local_dof_indices.size();
   internal::AffineConstraints::GlobalRowsFromLocal<number> &global_rows =
@@ -3851,8 +3844,7 @@ AffineConstraints<number>::distribute_local_to_global(
 
   typename internal::AffineConstraints::AffineConstraintsData<
     typename MatrixType::value_type>::ScratchDataAccessor
-    scratch_data(
-      internal::AffineConstraints::AffineConstraintsData<number>::scratch_data);
+    scratch_data(this->scratch_data);
 
   internal::AffineConstraints::GlobalRowsFromLocal<number> &global_rows =
     scratch_data->global_rows;
@@ -3911,9 +3903,8 @@ AffineConstraints<number>::add_entries_local_to_global(
          ExcNotQuadratic());
 
   const size_type n_local_dofs = local_dof_indices.size();
-  typename internal::AffineConstraints::AffineConstraintsData<number>::
-    ScratchDataAccessor scratch_data(
-      internal::AffineConstraints::AffineConstraintsData<number>::scratch_data);
+  typename internal::AffineConstraints::AffineConstraintsData<
+    number>::ScratchDataAccessor scratch_data(this->scratch_data);
 
   const bool dof_mask_is_active = (dof_mask.n_rows() == n_local_dofs);
   if (dof_mask_is_active == true)
@@ -4068,9 +4059,8 @@ AffineConstraints<number>::add_entries_local_to_global(
   const size_type n_local_dofs = local_dof_indices.size();
   const size_type num_blocks   = sparsity_pattern.n_block_rows();
 
-  typename internal::AffineConstraints::AffineConstraintsData<number>::
-    ScratchDataAccessor scratch_data(
-      internal::AffineConstraints::AffineConstraintsData<number>::scratch_data);
+  typename internal::AffineConstraints::AffineConstraintsData<
+    number>::ScratchDataAccessor scratch_data(this->scratch_data);
 
   const bool dof_mask_is_active = (dof_mask.n_rows() == n_local_dofs);
   if (dof_mask_is_active == true)
index 2ddddb53250c2db80025f1775e05a219b4b8012d..42ad4994f35ccb1213e77b61c9996552f74bdd17 100644 (file)
@@ -137,35 +137,4 @@ dealii::AffineConstraints<double>::distribute<
 #  endif
 #endif
 
-/*
- * Allocate scratch data.
- *
- * We cannot use the generic template instantiation because we need to
- * provide an initializer object of type
- * internals::AffineConstraintsData<Number> that can be passed to the
- * constructor of scratch_data (it won't allow one to be constructed in
- * place).
- */
-
-namespace internal
-{
-  namespace AffineConstraints
-  {
-#define SCRATCH_INITIALIZER(number, Name)              \
-  ScratchData<number> scratch_data_initializer_##Name; \
-  template <>                                          \
-  Threads::ThreadLocalStorage<ScratchData<number>>     \
-    AffineConstraintsData<number>::scratch_data(       \
-      scratch_data_initializer_##Name)
-
-    SCRATCH_INITIALIZER(double, d);
-    SCRATCH_INITIALIZER(float, f);
-#ifdef DEAL_II_WITH_COMPLEX_VALUES
-    SCRATCH_INITIALIZER(std::complex<double>, cd);
-    SCRATCH_INITIALIZER(std::complex<float>, cf);
-#endif
-#undef SCRATCH_INITIALIZER
-  } // namespace AffineConstraints
-} // namespace internal
-
 DEAL_II_NAMESPACE_CLOSE

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.