Reference documentation for deal.II version 9.1.0-pre
q_collection.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2005 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_q_collection_h
17 #define dealii_q_collection_h
18 
19 #include <deal.II/base/config.h>
20 
21 #include <deal.II/base/memory_consumption.h>
22 #include <deal.II/base/quadrature.h>
23 #include <deal.II/base/subscriptor.h>
24 
25 #include <deal.II/fe/fe.h>
26 
27 #include <memory>
28 #include <vector>
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 namespace hp
33 {
47  template <int dim>
48  class QCollection : public Subscriptor
49  {
50  public:
55  QCollection() = default;
56 
63  explicit QCollection(const Quadrature<dim> &quadrature);
64 
89  void
90  push_back(const Quadrature<dim> &new_quadrature);
91 
98  const Quadrature<dim> &operator[](const unsigned int index) const;
99 
103  unsigned int
104  size() const;
105 
112  unsigned int
113  max_n_quadrature_points() const;
114 
119  std::size_t
120  memory_consumption() const;
121 
126 
127  private:
132  std::vector<std::shared_ptr<const Quadrature<dim>>> quadratures;
133  };
134 
135 
136 
137  /* --------------- inline functions ------------------- */
138 
139  template <int dim>
140  inline unsigned int
142  {
143  return quadratures.size();
144  }
145 
146 
147 
148  template <int dim>
149  inline unsigned int
151  {
152  Assert(quadratures.size() > 0,
153  ExcMessage("You can't call this function for an empty collection"));
154 
155  unsigned int m = 0;
156  for (unsigned int i = 0; i < quadratures.size(); ++i)
157  if (quadratures[i]->size() > m)
158  m = quadratures[i]->size();
159 
160  return m;
161  }
162 
163 
164 
165  template <int dim>
166  inline const Quadrature<dim> &QCollection<dim>::
167  operator[](const unsigned int index) const
168  {
169  Assert(index < quadratures.size(),
170  ExcIndexRange(index, 0, quadratures.size()));
171  return *quadratures[index];
172  }
173 
174 
175 
176  template <int dim>
178  {
179  quadratures.push_back(std::make_shared<const Quadrature<dim>>(quadrature));
180  }
181 
182 
183 
184  template <int dim>
185  inline std::size_t
187  {
188  return (sizeof(*this) + MemoryConsumption::memory_consumption(quadratures));
189  }
190 
191 
192  template <int dim>
193  inline void
195  {
196  quadratures.push_back(
197  std::make_shared<const Quadrature<dim>>(new_quadrature));
198  }
199 
200 } // namespace hp
201 
202 
203 DEAL_II_NAMESPACE_CLOSE
204 
205 #endif
QCollection()=default
std::vector< std::shared_ptr< const Quadrature< dim > > > quadratures
Definition: q_collection.h:132
static::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static::ExceptionBase & ExcMessage(std::string arg1)
unsigned int size() const
Definition: q_collection.h:141
#define Assert(cond, exc)
Definition: exceptions.h:1227
#define DeclException0(Exception0)
Definition: exceptions.h:385
const Quadrature< dim > & operator[](const unsigned int index) const
Definition: q_collection.h:167
unsigned int max_n_quadrature_points() const
Definition: q_collection.h:150
Definition: hp.h:102
static::ExceptionBase & ExcNoQuadrature()
void push_back(const Quadrature< dim > &new_quadrature)
Definition: q_collection.h:194
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
std::size_t memory_consumption() const
Definition: q_collection.h:186