]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Fix some TODOs. Add memory_consumption functions in several places.
authorwolf <wolf@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 13 Jul 2001 08:08:27 +0000 (08:08 +0000)
committerwolf <wolf@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 13 Jul 2001 08:08:27 +0000 (08:08 +0000)
git-svn-id: https://svn.dealii.org/trunk@4841 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/fe/fe_values.h
deal.II/deal.II/include/fe/mapping.h
deal.II/deal.II/include/fe/mapping_cartesian.h
deal.II/deal.II/include/fe/mapping_q.h
deal.II/deal.II/include/fe/mapping_q1.h
deal.II/deal.II/source/fe/fe_values.cc
deal.II/deal.II/source/fe/mapping.cc
deal.II/deal.II/source/fe/mapping_cartesian.cc
deal.II/deal.II/source/fe/mapping_q.cc
deal.II/deal.II/source/fe/mapping_q1.cc

index 32e9b974a62f59aabf92fc4ab58ef69fb944dccf..713505d48121e0e6799fe1f542976c26c94c24ce 100644 (file)
@@ -802,6 +802,13 @@ class FEFaceValuesBase : public FEValuesBase<dim>
                                      */
     const Quadrature<dim-1> & get_quadrature () const;
     
+                                    /**
+                                     * Determine an estimate for the
+                                     * memory consumption (in bytes)
+                                     * of this object.
+                                     */
+    unsigned int memory_consumption () const;
+    
   protected:
                                     /**
                                      * Store a copy of the quadrature
@@ -820,6 +827,7 @@ class FEFaceValuesBase : public FEValuesBase<dim>
 };
 
 
+
 /**
  * Finite element evaluated in quadrature points on a face.
  *
index 46eb825b24cfa7e83a09174b91630e86ae5c3cfc..9c2793ab8d2b1f288d90ff88a746f45571c9c7f2 100644 (file)
@@ -142,6 +142,14 @@ class Mapping : public Subscriptor
                                          * cell-independent fields.
                                          */
        bool first_cell;
+
+                                        /**
+                                         * Return an estimate (in
+                                         * bytes) or the memory
+                                         * consumption of this
+                                         * object.
+                                         */
+       virtual unsigned int memory_consumption () const;
     };
     
                                     /**
index a724685948db4e38dc2b1b51094dfe6bcf9694e4..43534ca29f7752c0fb4d408b642ff0ba15cbafde 100644 (file)
@@ -213,6 +213,14 @@ class MappingCartesian : public Mapping<dim>
                                          */
        InternalData (const Quadrature<dim> &quadrature);
 
+                                        /**
+                                         * Return an estimate (in
+                                         * bytes) or the memory
+                                         * consumption of this
+                                         * object.
+                                         */
+       virtual unsigned int memory_consumption () const;
+
                                         /**
                                          * Length of the cell in
                                          * different coordinate
index 5f3f5ce849726be9ed9e668a1d6bb7c2d95f4609..4433a4d04d3de89f76bf744989e40bcca8c41819 100644 (file)
@@ -168,6 +168,15 @@ class MappingQ : public MappingQ1<dim>
                                          */
        InternalData (const unsigned int n_shape_functions);
        
+
+                                        /**
+                                         * Return an estimate (in
+                                         * bytes) or the memory
+                                         * consumption of this
+                                         * object.
+                                         */
+       virtual unsigned int memory_consumption () const;
+
                                         /**
                                          * Unit normal vectors. Used
                                          * for the alternative
index 56f44931b6008e7fb5d56c7a3e7ea84e3da80823..9f267dc62c073e48ff25bf4fc1ca0b967650a39b 100644 (file)
@@ -226,6 +226,14 @@ class MappingQ1 : public Mapping<dim>
                                          */
        Tensor<1,dim> &derivative (const unsigned int qpoint,
                                   const unsigned int shape_nr);
+
+                                        /**
+                                         * Return an estimate (in
+                                         * bytes) or the memory
+                                         * consumption of this
+                                         * object.
+                                         */
+       virtual unsigned int memory_consumption () const;
        
                                         /**
                                          * Values of shape
index 824a129154eebf61d180cc150c43db8624866063..193eb7cd357b105dd98bd98ebb6faa8e38119261 100644 (file)
@@ -431,15 +431,23 @@ template <int dim>
 unsigned int
 FEValuesBase<dim>::memory_consumption () const
 {
-//TODO:[WB] check whether this is up to date  
   return (MemoryConsumption::memory_consumption (shape_values) +
          MemoryConsumption::memory_consumption (shape_gradients) +
          MemoryConsumption::memory_consumption (shape_2nd_derivatives) +
          MemoryConsumption::memory_consumption (JxW_values) +
          MemoryConsumption::memory_consumption (quadrature_points) +
+         MemoryConsumption::memory_consumption (normal_vectors) +
+         MemoryConsumption::memory_consumption (boundary_forms) +
          sizeof(update_flags) +
          MemoryConsumption::memory_consumption (present_cell) +
-         MemoryConsumption::memory_consumption (fe));
+         MemoryConsumption::memory_consumption (n_quadrature_points) +
+         MemoryConsumption::memory_consumption (dofs_per_cell) +
+         MemoryConsumption::memory_consumption (mapping) +
+         MemoryConsumption::memory_consumption (fe) +
+         MemoryConsumption::memory_consumption (mapping_data) +
+         MemoryConsumption::memory_consumption (*mapping_data) +
+         MemoryConsumption::memory_consumption (fe_data) +
+         MemoryConsumption::memory_consumption (*fe_data));
 };
 
 
@@ -567,8 +575,8 @@ template <int dim>
 unsigned int
 FEValues<dim>::memory_consumption () const
 {
-//TODO:[WB] check whether this is up to date  
-  return FEValuesBase<dim>::memory_consumption ();
+  return (FEValuesBase<dim>::memory_consumption () +
+         MemoryConsumption::memory_consumption (quadrature));
 };
 
 
@@ -616,6 +624,17 @@ FEFaceValuesBase<dim>::get_boundary_forms () const
 };
 
 
+
+template <int dim>
+unsigned int
+FEFaceValuesBase<dim>::memory_consumption () const
+{
+  return (FEValuesBase<dim>::memory_consumption () +
+         MemoryConsumption::memory_consumption (quadrature) +
+         MemoryConsumption::memory_consumption (present_face));
+};
+
+
 /*------------------------------- FEFaceValues -------------------------------*/
 
 
index 0a8fd3918d46165b68998051f570b8dc762a48b5..fde86eda639369eb597b869b9423b7a821b2f1e5 100644 (file)
@@ -69,6 +69,15 @@ Mapping<dim>::InternalDataBase::~InternalDataBase ()
 {}
 
 
+
+template <int dim>
+unsigned int
+Mapping<dim>::InternalDataBase::memory_consumption () const
+{
+  return sizeof(*this);
+}
+
+
 /*------------------------------ InternalData ------------------------------*/
 
 
index 79509f6cfc04f16a55d98f1631f610fcdf6d6da0..0827a5944442c10642e4e411d0a3e719014f4097 100644 (file)
@@ -39,6 +39,20 @@ MappingCartesian<dim>::InternalData::InternalData (const Quadrature<dim>& q)
 
 
 
+template <int dim>
+unsigned int
+MappingCartesian<dim>::InternalData::memory_consumption () const
+{
+  return (Mapping<dim>::InternalDataBase::memory_consumption() +
+         MemoryConsumption::memory_consumption (length) +
+         MemoryConsumption::memory_consumption (quadrature_points) +
+         MemoryConsumption::memory_consumption (unit_tangentials) +
+         MemoryConsumption::memory_consumption (aux));
+};
+
+
+
+
 template <int dim>
 UpdateFlags
 MappingCartesian<dim>::update_once (const UpdateFlags) const
index 582866a1162f36757f48378c982f68630a9b23ae..649292e0468b007568c58a494fe8d0484fc6ece3 100644 (file)
@@ -15,6 +15,7 @@
 #include <fe/fe_q.h>
 #include <base/quadrature.h>
 #include <base/quadrature_lib.h>
+#include <base/memory_consumption.h>
 #include <base/tensor_product_polynomials.h>
 #include <lac/full_matrix.h>
 #include <grid/tria_iterator.h>
@@ -40,6 +41,19 @@ MappingQ<dim>::InternalData::InternalData (const unsigned int n_shape_functions)
 }
 
 
+
+template<int dim>
+unsigned int
+MappingQ<dim>::InternalData::memory_consumption () const 
+{
+  return (MappingQ1<dim>::InternalData::memory_consumption () +
+         MemoryConsumption::memory_consumption (unit_normals) +
+         MemoryConsumption::memory_consumption (use_mapping_q1_on_current_cell) +
+         MemoryConsumption::memory_consumption (mapping_q1_data));
+};
+
+
+
 #if deal_II_dimension == 1
 
 // in 1d, it is irrelevant which polynomial degree to use, since all
index 3632debff4084956b58dba828f057e522a38c352..2b5837bea7064f3f8402e607f902cd284fb2d4ff 100644 (file)
@@ -13,6 +13,7 @@
 
 #include <base/tensor.h>
 #include <base/quadrature.h>
+#include <base/memory_consumption.h>
 #include <lac/full_matrix.h>
 #include <grid/tria.h>
 #include <grid/tria_iterator.h>
@@ -91,6 +92,25 @@ MappingQ1<dim>::InternalData::derivative (const unsigned int qpoint,
 
 
 
+template <int dim>
+unsigned int
+MappingQ1<dim>::InternalData::memory_consumption () const
+{
+  return (Mapping<dim>::InternalDataBase::memory_consumption() +
+         MemoryConsumption::memory_consumption (shape_values) +
+         MemoryConsumption::memory_consumption (shape_derivatives) +
+         MemoryConsumption::memory_consumption (covariant) +
+         MemoryConsumption::memory_consumption (contravariant) +
+         MemoryConsumption::memory_consumption (unit_tangentials) +
+         MemoryConsumption::memory_consumption (aux) +
+         MemoryConsumption::memory_consumption (mapping_support_points) +
+         MemoryConsumption::memory_consumption (cell_of_current_support_points) +
+         MemoryConsumption::memory_consumption (is_mapping_q1_data) +
+         MemoryConsumption::memory_consumption (n_shape_functions));
+};
+
+
+
 template <int dim>
 MappingQ1<dim>::MappingQ1 ()
 {};

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.