]> https://gitweb.dealii.org/ - dealii.git/commitdiff
RCM: Update dofs, part 1: dof_accessor
authorMatthias Maier <tamiko@43-1.org>
Fri, 25 May 2018 18:09:26 +0000 (13:09 -0500)
committerMatthias Maier <tamiko@43-1.org>
Wed, 6 Jun 2018 15:19:38 +0000 (10:19 -0500)
include/deal.II/dofs/dof_accessor.h
include/deal.II/dofs/dof_accessor.templates.h

index 4e79a35b2c9dc112f3dae3e00366c9fc6e56d171..cf4a50abc50e69eeee4afb676c814026968198e8 100644 (file)
@@ -33,7 +33,8 @@ template <typename number>
 class FullMatrix;
 template <typename number>
 class Vector;
-class ConstraintMatrix;
+template <typename number>
+class AffineConstraints;
 
 template <typename Accessor>
 class TriaRawIterator;
@@ -1555,16 +1556,18 @@ public:
    * or a BlockVector<double>, or a PETSc or Trilinos vector if deal.II is
    * compiled to support these libraries. It is in the responsibility of the
    * caller to assure that the types of the numbers stored in input and output
-   * vectors are compatible and with similar accuracy. The ConstraintMatrix
-   * passed as an argument to this function makes sure that constraints are
-   * correctly distributed when the dof values are calculated.
+   * vectors are compatible and with similar accuracy. The
+   * AffineConstraints object passed as an argument to this function makes
+   * sure that constraints are correctly distributed when the dof values
+   * are calculated.
    */
   template <class InputVector, typename ForwardIterator>
   void
-  get_dof_values(const ConstraintMatrix &constraints,
-                 const InputVector &     values,
-                 ForwardIterator         local_values_begin,
-                 ForwardIterator         local_values_end) const;
+  get_dof_values(
+    const AffineConstraints<typename InputVector::value_type> &constraints,
+    const InputVector &                                        values,
+    ForwardIterator local_values_begin,
+    ForwardIterator local_values_end) const;
 
   /**
    * This function is the counterpart to get_dof_values(): it takes a vector
@@ -1723,15 +1726,16 @@ public:
    *
    * The elements are <em>added</em> up to the elements in the global vector,
    * rather than just set, since this is usually what one wants. Moreover, the
-   * ConstraintMatrix passed to this function makes sure that also constraints
-   * are eliminated in this process.
+   * AffineConstraints object passed to this function makes sure that also
+   * constraints are eliminated in this process.
    */
   template <typename ForwardIterator, typename OutputVector>
   void
-  distribute_local_to_global(const ConstraintMatrix &constraints,
-                             ForwardIterator         local_source_begin,
-                             ForwardIterator         local_source_end,
-                             OutputVector &          global_destination) const;
+  distribute_local_to_global(
+    const AffineConstraints<typename OutputVector::value_type> &constraints,
+    ForwardIterator local_source_begin,
+    ForwardIterator local_source_end,
+    OutputVector &  global_destination) const;
 
   /**
    * This function does much the same as the
index 5267254b031343824a6785abc1e5366187ce592d..633b541baddcd788c4e52d59445dbaf13060e519 100644 (file)
@@ -31,7 +31,7 @@
 #include <deal.II/hp/dof_faces.h>
 #include <deal.II/hp/dof_level.h>
 
-#include <deal.II/lac/constraint_matrix.h>
+#include <deal.II/lac/affine_constraints.h>
 #include <deal.II/lac/read_write_vector.h>
 
 #include <limits>
@@ -3096,11 +3096,11 @@ namespace internal
       static void
       distribute_local_to_global(
         const DoFCellAccessor<dealii::DoFHandler<dim, spacedim>,
-                              level_dof_access> &accessor,
-        const ConstraintMatrix &                 constraints,
-        ForwardIterator                          local_source_begin,
-        ForwardIterator                          local_source_end,
-        OutputVector &                           global_destination)
+                              level_dof_access> &                   accessor,
+        const AffineConstraints<typename OutputVector::value_type> &constraints,
+        ForwardIterator local_source_begin,
+        ForwardIterator local_source_end,
+        OutputVector &  global_destination)
       {
         Assert(
           accessor.dof_handler != nullptr,
@@ -3136,11 +3136,11 @@ namespace internal
       static void
       distribute_local_to_global(
         const DoFCellAccessor<dealii::hp::DoFHandler<dim, spacedim>,
-                              level_dof_access> &accessor,
-        const ConstraintMatrix &                 constraints,
-        ForwardIterator                          local_source_begin,
-        ForwardIterator                          local_source_end,
-        OutputVector &                           global_destination)
+                              level_dof_access> &                   accessor,
+        const AffineConstraints<typename OutputVector::value_type> &constraints,
+        ForwardIterator local_source_begin,
+        ForwardIterator local_source_end,
+        OutputVector &  global_destination)
       {
         Assert(
           accessor.dof_handler != nullptr,
@@ -3352,10 +3352,10 @@ namespace internal
 
         const unsigned int n_dofs = local_matrix.size();
 
-        // TODO[WB/MK]: This function could me made more efficient because it
-        // allocates memory, which could be avoided by passing in another
-        // argument as a scratch array. Comment(GK) Do not bother and leave this
-        // to ConstraintMatrix or MeshWorker::Assembler
+        // TODO[WB/MK]: This function could me made more efficient because
+        // it allocates memory, which could be avoided by passing in
+        // another argument as a scratch array. Comment(GK) Do not bother
+        // and leave this to AffineConstraints or MeshWorker::Assembler
 
         // get indices of dofs
         std::vector<types::global_dof_index> dofs(n_dofs);
@@ -3627,10 +3627,10 @@ template <typename DoFHandlerType, bool level_dof_access>
 template <class InputVector, typename ForwardIterator>
 inline void
 DoFCellAccessor<DoFHandlerType, level_dof_access>::get_dof_values(
-  const ConstraintMatrix &constraints,
-  const InputVector &     values,
-  ForwardIterator         local_values_begin,
-  ForwardIterator         local_values_end) const
+  const AffineConstraints<typename InputVector::value_type> &constraints,
+  const InputVector &                                        values,
+  ForwardIterator                                            local_values_begin,
+  ForwardIterator local_values_end) const
 {
   Assert(this->is_artificial() == false,
          ExcMessage("Can't ask for DoF indices on artificial cells."));
@@ -3793,10 +3793,10 @@ template <typename DoFHandlerType, bool level_dof_access>
 template <typename ForwardIterator, typename OutputVector>
 inline void
 DoFCellAccessor<DoFHandlerType, level_dof_access>::distribute_local_to_global(
-  const ConstraintMatrix &constraints,
-  ForwardIterator         local_source_begin,
-  ForwardIterator         local_source_end,
-  OutputVector &          global_destination) const
+  const AffineConstraints<typename OutputVector::value_type> &constraints,
+  ForwardIterator local_source_begin,
+  ForwardIterator local_source_end,
+  OutputVector &  global_destination) const
 {
   dealii::internal::DoFCellAccessorImplementation::Implementation::
     distribute_local_to_global(*this,

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.