Reference documentation for deal.II version 9.0.0
fe_series.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2016 - 2017 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_fe_series_H
17 #define dealii_fe_series_H
18 
19 
20 
21 #include <deal.II/base/config.h>
22 #include <deal.II/base/subscriptor.h>
23 #include <deal.II/base/exceptions.h>
24 #include <deal.II/base/table.h>
25 #include <deal.II/base/table_indices.h>
26 #include <deal.II/base/tensor.h>
27 #include <deal.II/hp/fe_collection.h>
28 #include <deal.II/hp/q_collection.h>
29 #include <deal.II/lac/full_matrix.h>
30 #include <deal.II/lac/vector.h>
31 #include <deal.II/numerics/vector_tools.h>
32 
33 #include <vector>
34 #include <string>
35 #include <memory>
36 
37 
38 DEAL_II_NAMESPACE_OPEN
39 
40 
43 
44 
53 namespace FESeries
54 {
86  template <int dim>
87  class Fourier : public Subscriptor
88  {
89  public:
97  Fourier(const unsigned int size_in_each_direction,
100 
106  void calculate(const ::Vector<double> &local_dof_values,
107  const unsigned int cell_active_fe_index,
108  Table<dim,std::complex<double> > &fourier_coefficients);
109 
110  private:
115 
120 
125  void ensure_existence(const unsigned int fe_index);
126 
131 
135  std::vector<FullMatrix<std::complex<double> > > fourier_transform_matrices;
136 
140  std::vector<std::complex<double> > unrolled_coefficients;
141 
142  };
143 
187  template <int dim>
188  class Legendre : public Subscriptor
189  {
190  public:
198  Legendre(const unsigned int size_in_each_direction,
201 
207  void calculate(const ::Vector<double> &local_dof_values,
208  const unsigned int cell_active_fe_index,
209  Table<dim,double> &legendre_coefficients);
210 
211  private:
215  const unsigned int N;
216 
221 
226 
231  void ensure_existence(const unsigned int fe_index);
232 
236  std::vector<FullMatrix<double> > legendre_transform_matrices;
237 
241  std::vector<double> unrolled_coefficients;
242 
243  };
244 
245 
260  template <int dim, typename T>
261  std::pair<std::vector<unsigned int>,std::vector<double> >
262  process_coefficients(const Table<dim,T> &coefficients,
263  const std::function<std::pair<bool,unsigned int>(const TableIndices<dim> &)> &predicate,
264  const VectorTools::NormType norm);
265 
266 
267 
273  std::pair<double,double> linear_regression(const std::vector<double> &x,
274  const std::vector<double> &y);
275 
276 }
277 
280 #ifndef DOXYGEN
281 
282 // ------------------- inline and template functions ----------------
283 
284 namespace
285 {
286  template <int dim,typename T>
287  void fill_map_index(const Table<dim,T> &coefficients,
288  const TableIndices<dim> &ind,
289  const std::function<std::pair<bool,unsigned int>(const TableIndices<dim> &)> &predicate,
290  std::map<unsigned int, std::vector<T> > &pred_to_values)
291  {
292  const std::pair<bool,unsigned int> pred_pair = predicate(ind);
293  // don't add a value if predicate is false
294  if (pred_pair.first == false)
295  return;
296 
297  const unsigned int &pred_value = pred_pair.second;
298  const T &coeff_value = coefficients(ind);
299  // If pred_value is not in the pred_to_values map, the element will be created.
300  // Otherwise a reference to the existing element is returned.
301  pred_to_values[pred_value].push_back(coeff_value);
302  }
303 
304  template <typename T>
305  void fill_map(const Table<1,T> &coefficients,
306  const std::function<std::pair<bool,unsigned int>(const TableIndices<1> &)> &predicate,
307  std::map<unsigned int, std::vector<T> > &pred_to_values)
308  {
309  for (unsigned int i = 0; i < coefficients.size(0); i++)
310  {
311  const TableIndices<1> ind(i);
312  fill_map_index(coefficients,ind,predicate,pred_to_values);
313  }
314 
315  }
316 
317  template <typename T>
318  void fill_map(const Table<2,T> &coefficients,
319  const std::function<std::pair<bool,unsigned int>(const TableIndices<2> &)> &predicate,
320  std::map<unsigned int, std::vector<T> > &pred_to_values)
321  {
322  for (unsigned int i = 0; i < coefficients.size(0); i++)
323  for (unsigned int j = 0; j < coefficients.size(1); j++)
324  {
325  const TableIndices<2> ind(i,j);
326  fill_map_index(coefficients,ind,predicate,pred_to_values);
327  }
328 
329  }
330 
331  template <typename T>
332  void fill_map(const Table<3,T> &coefficients,
333  const std::function<std::pair<bool,unsigned int>(const TableIndices<3> &)> &predicate,
334  std::map<unsigned int, std::vector<T> > &pred_to_values)
335  {
336  for (unsigned int i = 0; i < coefficients.size(0); i++)
337  for (unsigned int j = 0; j < coefficients.size(1); j++)
338  for (unsigned int k = 0; k < coefficients.size(2); k++)
339  {
340  const TableIndices<3> ind(i,j,k);
341  fill_map_index(coefficients,ind,predicate,pred_to_values);
342  }
343  }
344 
345 
346  template <typename T>
347  double complex_mean_value(const T &value)
348  {
349  return value;
350  }
351 
352  template <typename T>
353  double complex_mean_value(const std::complex<T> &value)
354  {
355  AssertThrow(false, ExcMessage("FESeries::process_coefficients() can not be used with"
356  "complex-valued coefficients and VectorTools::mean norm."));
357  return std::abs(value);
358  }
359 
360 }
361 
362 
363 template <int dim, typename T>
364 std::pair<std::vector<unsigned int>,std::vector<double> >
365 FESeries::process_coefficients(const Table<dim,T> &coefficients,
366  const std::function<std::pair<bool,unsigned int>(const TableIndices<dim> &)> &predicate,
367  const VectorTools::NormType norm)
368 {
369  std::vector<unsigned int> predicate_values;
370  std::vector<double> norm_values;
371 
372  // first, parse all table elements into a map of predicate values and coefficients.
373  // We could have stored (predicate values ->TableIndicies) map, but its
374  // processing would have been much harder later on.
375  std::map<unsigned int, std::vector<T> > pred_to_values;
376  fill_map(coefficients,predicate,pred_to_values);
377 
378  // now go through the map and populate the @p norm_values based on @p norm:
379  for (typename std::map<unsigned int, std::vector<T> >::const_iterator it = pred_to_values.begin();
380  it != pred_to_values.end(); ++it)
381  {
382  predicate_values.push_back(it->first);
383  Vector<T> values(it->second.begin(),
384  it->second.end());
385 
386  switch (norm)
387  {
389  {
390  norm_values.push_back(values.l2_norm());
391  break;
392  }
394  {
395  norm_values.push_back(values.l1_norm());
396  break;
397  }
399  {
400  norm_values.push_back(values.linfty_norm());
401  break;
402  }
403  case VectorTools::mean:
404  {
405  norm_values.push_back(complex_mean_value(values.mean_value()));
406  break;
407  }
408  default:
409  AssertThrow(false, ExcNotImplemented());
410  break;
411  }
412  }
413 
414  return std::make_pair(predicate_values,norm_values);
415 }
416 
417 
418 #endif // DOXYGEN
419 
420 DEAL_II_NAMESPACE_CLOSE
421 
422 #endif // dealii_fe_series_H
void ensure_existence(const unsigned int fe_index)
SmartPointer< const hp::FECollection< dim > > fe_collection
Definition: fe_series.h:114
Legendre(const unsigned int size_in_each_direction, const hp::FECollection< dim > &fe_collection, const hp::QCollection< dim > &q_collection)
Definition: fe_series.cc:245
std::pair< double, double > linear_regression(const std::vector< double > &x, const std::vector< double > &y)
Definition: fe_series.cc:372
SmartPointer< const hp::QCollection< dim > > q_collection
Definition: fe_series.h:119
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
void calculate(const ::Vector< double > &local_dof_values, const unsigned int cell_active_fe_index, Table< dim, double > &legendre_coefficients)
Definition: fe_series.cc:259
Fourier(const unsigned int size_in_each_direction, const hp::FECollection< dim > &fe_collection, const hp::QCollection< dim > &q_collection)
Definition: fe_series.cc:74
std::vector< std::complex< double > > unrolled_coefficients
Definition: fe_series.h:140
SmartPointer< const hp::QCollection< dim > > q_collection
Definition: fe_series.h:225
static ::ExceptionBase & ExcMessage(std::string arg1)
SmartPointer< const hp::FECollection< dim > > fe_collection
Definition: fe_series.h:220
void calculate(const ::Vector< double > &local_dof_values, const unsigned int cell_active_fe_index, Table< dim, std::complex< double > > &fourier_coefficients)
Definition: fe_series.cc:87
void ensure_existence(const unsigned int fe_index)
size_type size(const unsigned int i) const
std::pair< std::vector< unsigned int >, std::vector< double > > process_coefficients(const Table< dim, T > &coefficients, const std::function< std::pair< bool, unsigned int >(const TableIndices< dim > &)> &predicate, const VectorTools::NormType norm)
std::vector< double > unrolled_coefficients
Definition: fe_series.h:241
const unsigned int N
Definition: fe_series.h:215
static ::ExceptionBase & ExcNotImplemented()
Definition: table.h:34
Table< dim, Tensor< 1, dim > > k_vectors
Definition: fe_series.h:130
std::vector< FullMatrix< std::complex< double > > > fourier_transform_matrices
Definition: fe_series.h:135
std::vector< FullMatrix< double > > legendre_transform_matrices
Definition: fe_series.h:236