Reference documentation for deal.II version 9.0.0
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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_q_collection_h
17 #define dealii_q_collection_h
18 
19 #include <deal.II/base/config.h>
20 #include <deal.II/base/subscriptor.h>
21 #include <deal.II/base/quadrature.h>
22 #include <deal.II/base/memory_consumption.h>
23 #include <deal.II/fe/fe.h>
24 
25 #include <vector>
26 #include <memory>
27 
28 DEAL_II_NAMESPACE_OPEN
29 
30 namespace hp
31 {
45  template <int dim>
46  class QCollection : public Subscriptor
47  {
48  public:
53  QCollection () = default;
54 
61  explicit QCollection (const Quadrature<dim> &quadrature);
62 
87  void push_back (const Quadrature<dim> &new_quadrature);
88 
95  const Quadrature<dim> &
96  operator[] (const unsigned int index) const;
97 
101  unsigned int size () const;
102 
109  unsigned int max_n_quadrature_points () const;
110 
115  std::size_t memory_consumption () const;
116 
121 
122  private:
127  std::vector<std::shared_ptr<const Quadrature<dim> > > quadratures;
128  };
129 
130 
131 
132  /* --------------- inline functions ------------------- */
133 
134  template <int dim>
135  inline
136  unsigned int
138  {
139  return quadratures.size();
140  }
141 
142 
143 
144  template <int dim>
145  inline
146  unsigned int
148  {
149  Assert (quadratures.size() > 0,
150  ExcMessage ("You can't call this function for an empty collection"));
151 
152  unsigned int m = 0;
153  for (unsigned int i=0; i<quadratures.size(); ++i)
154  if (quadratures[i]->size() > m)
155  m = quadratures[i]->size();
156 
157  return m;
158  }
159 
160 
161 
162  template <int dim>
163  inline
164  const Quadrature<dim> &
165  QCollection<dim>::operator[] (const unsigned int index) const
166  {
167  Assert (index < quadratures.size (),
168  ExcIndexRange (index, 0, quadratures.size ()));
169  return *quadratures[index];
170  }
171 
172 
173 
174  template <int dim>
175  inline
177  {
178  quadratures.push_back (std::make_shared<const Quadrature<dim> >(quadrature));
179  }
180 
181 
182 
183  template <int dim>
184  inline
185  std::size_t
187  {
188  return (sizeof(*this) +
190  }
191 
192 
193  template <int dim>
194  inline
195  void
197  {
198  quadratures.push_back (std::make_shared<const Quadrature<dim> >(new_quadrature));
199  }
200 
201 } // namespace hp
202 
203 
204 DEAL_II_NAMESPACE_CLOSE
205 
206 #endif
std::size_t memory_consumption() const
Definition: q_collection.h:186
QCollection()=default
std::vector< std::shared_ptr< const Quadrature< dim > > > quadratures
Definition: q_collection.h:127
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
const Quadrature< dim > & operator[](const unsigned int index) const
Definition: q_collection.h:165
static ::ExceptionBase & ExcMessage(std::string arg1)
unsigned int max_n_quadrature_points() const
Definition: q_collection.h:147
#define Assert(cond, exc)
Definition: exceptions.h:1142
#define DeclException0(Exception0)
Definition: exceptions.h:323
unsigned int size() const
Definition: q_collection.h:137
Definition: hp.h:102
static ::ExceptionBase & ExcNoQuadrature()
void push_back(const Quadrature< dim > &new_quadrature)
Definition: q_collection.h:196
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)