Skip to content

Commit

Permalink
Implement FE_Poly::memory_consumption so that memory_consumption work…
Browse files Browse the repository at this point in the history
…s for FE_Q and FE_DGQ
  • Loading branch information
peterrum committed Dec 21, 2019
1 parent 150f663 commit 773fa71
Show file tree
Hide file tree
Showing 22 changed files with 187 additions and 28 deletions.
4 changes: 4 additions & 0 deletions doc/news/changes/minor/20191222PeterMunch
@@ -0,0 +1,4 @@
New: The method memory_consumption has been implemented in and added
to a variety of classes, e.g. FE_DGQ.
<br>
(Peter Munch, 2019/12/22)
6 changes: 6 additions & 0 deletions include/deal.II/base/polynomial.h
Expand Up @@ -225,6 +225,12 @@ namespace Polynomials
void
serialize(Archive &ar, const unsigned int version);

/**
* Return an estimate (in bytes) for the memory consumption of this object.
*/
virtual std::size_t
memory_consumption() const;

protected:
/**
* This function performs the actual scaling.
Expand Down
6 changes: 6 additions & 0 deletions include/deal.II/base/polynomials_piecewise.h
Expand Up @@ -132,6 +132,12 @@ namespace Polynomials
void
serialize(Archive &ar, const unsigned int version);

/**
* Return an estimate (in bytes) for the memory consumption of this object.
*/
virtual std::size_t
memory_consumption() const;

protected:
/**
* Underlying polynomial object that is scaled to a subinterval and
Expand Down
6 changes: 6 additions & 0 deletions include/deal.II/base/scalar_polynomials_base.h
Expand Up @@ -138,6 +138,12 @@ class ScalarPolynomialsBase
virtual std::string
name() const = 0;

/**
* Return an estimate (in bytes) for the memory consumption of this object.
*/
virtual std::size_t
memory_consumption() const;

private:
/**
* The highest polynomial degree of this functions represented by this object.
Expand Down
6 changes: 6 additions & 0 deletions include/deal.II/base/tensor_product_polynomials.h
Expand Up @@ -209,6 +209,12 @@ class TensorProductPolynomials : public ScalarPolynomialsBase<dim>
virtual std::unique_ptr<ScalarPolynomialsBase<dim>>
clone() const override;

/**
* Return an estimate (in bytes) for the memory consumption of this object.
*/
virtual std::size_t
memory_consumption() const override;

protected:
/**
* Copy of the vector <tt>pols</tt> of polynomials given to the constructor.
Expand Down
2 changes: 1 addition & 1 deletion include/deal.II/fe/fe.h
Expand Up @@ -715,7 +715,7 @@ class FiniteElement : public Subscriptor, public FiniteElementData<dim>
UpdateFlags update_each;

/**
* Return an estimate (in bytes) or the memory consumption of this object.
* Return an estimate (in bytes) for the memory consumption of this object.
*/
virtual std::size_t
memory_consumption() const;
Expand Down
11 changes: 0 additions & 11 deletions include/deal.II/fe/fe_dgq.h
Expand Up @@ -320,17 +320,6 @@ class FE_DGQ : public FE_Poly<TensorProductPolynomials<dim>, dim, spacedim>
const std::vector<Vector<double>> &support_point_values,
std::vector<double> & nodal_values) const override;

/**
* Determine an estimate for the memory consumption (in bytes) of this
* object.
*
* This function is made virtual, since finite element objects are usually
* accessed through pointers to their base class, rather than the class
* itself.
*/
virtual std::size_t
memory_consumption() const override;

virtual std::unique_ptr<FiniteElement<dim, spacedim>>
clone() const override;

Expand Down
6 changes: 6 additions & 0 deletions include/deal.II/fe/fe_poly.h
Expand Up @@ -227,6 +227,12 @@ class FE_Poly : public FiniteElement<dim, spacedim>
const Point<dim> & p,
const unsigned int component) const override;

/**
* Return an estimate (in bytes) for the memory consumption of this object.
*/
virtual std::size_t
memory_consumption() const override;

protected:
/*
* NOTE: The following function has its definition inlined into the class
Expand Down
10 changes: 10 additions & 0 deletions include/deal.II/fe/fe_poly.templates.h
Expand Up @@ -581,6 +581,16 @@ FE_Poly<PolynomialType, dim, spacedim>::get_poly_space_numbering_inverse() const



template <class PolynomialType, int dim, int spacedim>
std::size_t
FE_Poly<PolynomialType, dim, spacedim>::memory_consumption() const
{
return FiniteElement<dim, spacedim>::memory_consumption() +
poly_space.memory_consumption();
}



DEAL_II_NAMESPACE_CLOSE

#endif
2 changes: 1 addition & 1 deletion include/deal.II/fe/mapping.h
Expand Up @@ -630,7 +630,7 @@ class Mapping : public Subscriptor
UpdateFlags update_each;

/**
* Return an estimate (in bytes) or the memory consumption of this object.
* Return an estimate (in bytes) for the memory consumption of this object.
*/
virtual std::size_t
memory_consumption() const;
Expand Down
2 changes: 1 addition & 1 deletion include/deal.II/fe/mapping_cartesian.h
Expand Up @@ -182,7 +182,7 @@ class MappingCartesian : public Mapping<dim, spacedim>
InternalData(const Quadrature<dim> &quadrature);

/**
* Return an estimate (in bytes) or the memory consumption of this object.
* Return an estimate (in bytes) for the memory consumption of this object.
*/
virtual std::size_t
memory_consumption() const override;
Expand Down
2 changes: 1 addition & 1 deletion include/deal.II/fe/mapping_fe_field.h
Expand Up @@ -355,7 +355,7 @@ class MappingFEField : public Mapping<dim, spacedim>
fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr);

/**
* Return an estimate (in bytes) or the memory consumption of this object.
* Return an estimate (in bytes) for the memory consumption of this object.
*/
virtual std::size_t
memory_consumption() const override;
Expand Down
2 changes: 1 addition & 1 deletion include/deal.II/fe/mapping_manifold.h
Expand Up @@ -211,7 +211,7 @@ class MappingManifold : public Mapping<dim, spacedim>
const typename Triangulation<dim, spacedim>::cell_iterator &cell) const;

/**
* Return an estimate (in bytes) or the memory consumption of this object.
* Return an estimate (in bytes) for the memory consumption of this object.
*/
virtual std::size_t
memory_consumption() const override;
Expand Down
2 changes: 1 addition & 1 deletion include/deal.II/fe/mapping_q.h
Expand Up @@ -237,7 +237,7 @@ class MappingQ : public Mapping<dim, spacedim>


/**
* Return an estimate (in bytes) or the memory consumption of this object.
* Return an estimate (in bytes) for the memory consumption of this object.
*/
virtual std::size_t
memory_consumption() const override;
Expand Down
2 changes: 1 addition & 1 deletion include/deal.II/fe/mapping_q_generic.h
Expand Up @@ -364,7 +364,7 @@ class MappingQGeneric : public Mapping<dim, spacedim>
fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr);

/**
* Return an estimate (in bytes) or the memory consumption of this object.
* Return an estimate (in bytes) for the memory consumption of this object.
*/
virtual std::size_t
memory_consumption() const override;
Expand Down
12 changes: 12 additions & 0 deletions source/base/polynomial.cc
Expand Up @@ -14,6 +14,7 @@
// ---------------------------------------------------------------------

#include <deal.II/base/exceptions.h>
#include <deal.II/base/memory_consumption.h>
#include <deal.II/base/point.h>
#include <deal.II/base/polynomial.h>
#include <deal.II/base/quadrature_lib.h>
Expand Down Expand Up @@ -660,6 +661,17 @@ namespace Polynomials
}


template <typename number>
std::size_t
Polynomial<number>::memory_consumption() const
{
return (MemoryConsumption::memory_consumption(coefficients) +
MemoryConsumption::memory_consumption(in_lagrange_product_form) +
MemoryConsumption::memory_consumption(lagrange_support_points) +
MemoryConsumption::memory_consumption(lagrange_weight));
}



// ------------------ class Monomial -------------------------- //

Expand Down
13 changes: 13 additions & 0 deletions source/base/polynomials_piecewise.cc
Expand Up @@ -13,6 +13,7 @@
//
// ---------------------------------------------------------------------

#include <deal.II/base/memory_consumption.h>
#include <deal.II/base/polynomials_piecewise.h>


Expand Down Expand Up @@ -117,6 +118,18 @@ namespace Polynomials



template <typename number>
std::size_t
PiecewisePolynomial<number>::memory_consumption() const
{
return (polynomial.memory_consumption() +
MemoryConsumption::memory_consumption(n_intervals) +
MemoryConsumption::memory_consumption(interval) +
MemoryConsumption::memory_consumption(spans_two_intervals));
}



std::vector<PiecewisePolynomial<double>>
generate_complete_Lagrange_basis_on_subdivisions(
const unsigned int n_subdivisions,
Expand Down
10 changes: 10 additions & 0 deletions source/base/scalar_polynomials_base.cc
Expand Up @@ -34,6 +34,16 @@ ScalarPolynomialsBase<dim>::ScalarPolynomialsBase(



template <int dim>
std::size_t
ScalarPolynomialsBase<dim>::memory_consumption() const
{
Assert(false, ExcNotImplemented());
return 0;
}



template class ScalarPolynomialsBase<0>;
template class ScalarPolynomialsBase<1>;
template class ScalarPolynomialsBase<2>;
Expand Down
12 changes: 12 additions & 0 deletions source/base/tensor_product_polynomials.cc
Expand Up @@ -14,6 +14,7 @@
// ---------------------------------------------------------------------

#include <deal.II/base/exceptions.h>
#include <deal.II/base/memory_consumption.h>
#include <deal.II/base/polynomials_piecewise.h>
#include <deal.II/base/std_cxx14/memory.h>
#include <deal.II/base/tensor_product_polynomials.h>
Expand Down Expand Up @@ -435,6 +436,17 @@ TensorProductPolynomials<dim, PolynomialType>::clone() const



template <int dim, typename PolynomialType>
std::size_t
TensorProductPolynomials<dim, PolynomialType>::memory_consumption() const
{
return (MemoryConsumption::memory_consumption(polynomials) +
MemoryConsumption::memory_consumption(index_map) +
MemoryConsumption::memory_consumption(index_map_inverse));
}



/* ------------------- AnisotropicPolynomials -------------- */


Expand Down
10 changes: 0 additions & 10 deletions source/fe/fe_dgq.cc
Expand Up @@ -819,16 +819,6 @@ FE_DGQ<dim, spacedim>::get_constant_modes() const



template <int dim, int spacedim>
std::size_t
FE_DGQ<dim, spacedim>::memory_consumption() const
{
Assert(false, ExcNotImplemented());
return 0;
}



// ------------------------------ FE_DGQArbitraryNodes -----------------------

template <int dim, int spacedim>
Expand Down
67 changes: 67 additions & 0 deletions tests/fe/memory_consumption_01.cc
@@ -0,0 +1,67 @@
// ---------------------------------------------------------------------
//
// Copyright (C) 2019 by the deal.II authors
//
// This file is part of the deal.II library.
//
// The deal.II library is free software; you can use it, redistribute
// it, and/or modify it under the terms of the GNU Lesser General
// Public License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// The full text of the license can be found in the file LICENSE.md at
// the top level directory of deal.II.
//
// ---------------------------------------------------------------------


// Test memory_consumption of FE_Q and FE_DGQ.


#include <deal.II/fe/fe_dgq.h>
#include <deal.II/fe/fe_q.h>

#include "../tests.h"

using namespace dealii;

template <int dim>
void
test(const FiniteElement<dim> &fe)
{
deallog << fe.get_name() << " " << fe.memory_consumption() << std::endl;
}

template <int dim>
void
tests()
{
// FE_Q
for (unsigned int i = 1; i < 4; i++)
test(FE_Q<dim>(i));

// FE_DGQ
for (unsigned int i = 0; i < 4; i++)
test(FE_DGQ<dim>(i));
}

int
main()
{
initlog();

{
deallog.push("1d");
tests<1>();
deallog.pop();
}
{
deallog.push("2d");
tests<2>();
deallog.pop();
}
{
deallog.push("3d");
tests<3>();
deallog.pop();
}
}
22 changes: 22 additions & 0 deletions tests/fe/memory_consumption_01.output
@@ -0,0 +1,22 @@

DEAL:1d::FE_Q<1>(1) 1770
DEAL:1d::FE_Q<1>(2) 1979
DEAL:1d::FE_Q<1>(3) 2204
DEAL:1d::FE_DGQ<1>(0) 1557
DEAL:1d::FE_DGQ<1>(1) 1750
DEAL:1d::FE_DGQ<1>(2) 1959
DEAL:1d::FE_DGQ<1>(3) 2184
DEAL:2d::FE_Q<2>(1) 4062
DEAL:2d::FE_Q<2>(2) 4827
DEAL:2d::FE_Q<2>(3) 5880
DEAL:2d::FE_DGQ<2>(0) 3573
DEAL:2d::FE_DGQ<2>(1) 4006
DEAL:2d::FE_DGQ<2>(2) 4695
DEAL:2d::FE_DGQ<2>(3) 5640
DEAL:3d::FE_Q<3>(1) 10678
DEAL:3d::FE_Q<3>(2) 14499
DEAL:3d::FE_Q<3>(3) 23432
DEAL:3d::FE_DGQ<3>(0) 9525
DEAL:3d::FE_DGQ<3>(1) 10438
DEAL:3d::FE_DGQ<3>(2) 12807
DEAL:3d::FE_DGQ<3>(3) 17352

0 comments on commit 773fa71

Please sign in to comment.