]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Implement MGTransferGlobalCoarsening::memory_consumption() 12731/head
authorPeter Munch <peterrmuench@gmail.com>
Wed, 1 Sep 2021 20:44:37 +0000 (22:44 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Sun, 5 Sep 2021 18:15:06 +0000 (20:15 +0200)
include/deal.II/multigrid/mg_transfer_global_coarsening.h
include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h

index 64bbd4aa1ed74ab6989749ac36a6aee7dbb88696..62b285258564412e800e86481e7b6c9ae9b853be 100644 (file)
@@ -175,6 +175,12 @@ public:
    */
   void
   interpolate(VectorType &dst, const VectorType &src) const;
+
+  /**
+   * Return the memory consumption of the allocated memory in this class.
+   */
+  std::size_t
+  memory_consumption() const;
 };
 
 
@@ -277,6 +283,12 @@ public:
   interpolate(LinearAlgebra::distributed::Vector<Number> &      dst,
               const LinearAlgebra::distributed::Vector<Number> &src) const;
 
+  /**
+   * Return the memory consumption of the allocated memory in this class.
+   */
+  std::size_t
+  memory_consumption() const;
+
 private:
   /**
    * A multigrid transfer scheme. A multrigrid transfer class can have different
@@ -521,6 +533,15 @@ public:
                     MGLevelObject<VectorType> &      dst,
                     const InVector &                 src) const;
 
+  /**
+   * Return the memory consumption of the allocated memory in this class.
+   *
+   * @note Counts also the memory consumption of the underlying two-level
+   *   transfer operators.
+   */
+  std::size_t
+  memory_consumption() const;
+
 private:
   /**
    * Collection of the two-level transfer operators.
@@ -662,6 +683,23 @@ MGTransferGlobalCoarsening<dim, VectorType>::interpolate_to_mg(
     this->transfer[l].interpolate(dst[l - 1], dst[l]);
 }
 
+
+
+template <int dim, typename VectorType>
+std::size_t
+MGTransferGlobalCoarsening<dim, VectorType>::memory_consumption() const
+{
+  std::size_t size = 0;
+
+  const unsigned int min_level = transfer.min_level();
+  const unsigned int max_level = transfer.max_level();
+
+  for (unsigned int l = min_level + 1; l <= max_level; ++l)
+    size += this->transfer[l].memory_consumption();
+
+  return size;
+}
+
 #endif
 
 DEAL_II_NAMESPACE_CLOSE
index 9a84dcabee51cb8cfcd46c10a2900bf21bc31986..71c38856ef954193670391ff4fa684a13debaf0f 100644 (file)
@@ -2822,6 +2822,39 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
   return cell_transfer.run(cell_transfer_test);
 }
 
+
+
+template <int dim, typename Number>
+std::size_t
+MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
+  memory_consumption() const
+{
+  std::size_t size = 0;
+
+  for (const auto &scheme : schemes)
+    {
+      size += scheme.prolongation_matrix.memory_consumption();
+      size += scheme.prolongation_matrix_1d.memory_consumption();
+      size += scheme.restriction_matrix.memory_consumption();
+      size += scheme.restriction_matrix_1d.memory_consumption();
+    }
+
+  size += partitioner_fine->memory_consumption();
+  size += partitioner_coarse->memory_consumption();
+  size += vec_fine.memory_consumption();
+  size += vec_coarse.memory_consumption();
+  size +=
+    MemoryConsumption::memory_consumption(distribute_local_to_global_indices);
+  size +=
+    MemoryConsumption::memory_consumption(distribute_local_to_global_values);
+  size += MemoryConsumption::memory_consumption(distribute_local_to_global_ptr);
+  size += MemoryConsumption::memory_consumption(weights);
+  size += MemoryConsumption::memory_consumption(level_dof_indices_coarse);
+  size += MemoryConsumption::memory_consumption(level_dof_indices_fine);
+
+  return size;
+}
+
 DEAL_II_NAMESPACE_CLOSE
 
 #endif

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.