Reference documentation for deal.II version 9.3.2
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
fe_series.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2016 - 2021 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_fe_series_h
17 #define dealii_fe_series_h
18 
19 
20 
21 #include <deal.II/base/config.h>
22 
25 #include <deal.II/base/table.h>
27 #include <deal.II/base/tensor.h>
28 
31 
33 #include <deal.II/lac/vector.h>
34 
36 
37 #include <memory>
38 #include <string>
39 #include <vector>
40 
41 
43 
44 
47 
48 
55 namespace FESeries
56 {
87  template <int dim, int spacedim = dim>
88  class Fourier : public Subscriptor
89  {
90  public:
91  using CoefficientType = typename std::complex<double>;
92 
111  Fourier(const std::vector<unsigned int> & n_coefficients_per_direction,
114  const unsigned int component = numbers::invalid_unsigned_int);
115 
126  Fourier(const unsigned int n_coefficients_per_direction,
129 
135  template <typename Number>
136  void
137  calculate(const ::Vector<Number> &local_dof_values,
138  const unsigned int cell_active_fe_index,
139  Table<dim, CoefficientType> & fourier_coefficients);
140 
145  unsigned int
146  get_n_coefficients_per_direction(const unsigned int index) const;
147 
158  void
160 
170  template <class Archive>
171  void
172  save_transformation_matrices(Archive &ar, const unsigned int version);
173 
178  template <class Archive>
179  void
180  load_transformation_matrices(Archive &ar, const unsigned int version);
181 
185  bool
186  operator==(const Fourier<dim, spacedim> &fourier) const;
187 
188  private:
193  const std::vector<unsigned int> n_coefficients_per_direction;
194 
199 
204 
209 
213  std::vector<FullMatrix<CoefficientType>> fourier_transform_matrices;
214 
218  std::vector<CoefficientType> unrolled_coefficients;
219 
224  const unsigned int component;
225  };
226 
227 
228 
271  template <int dim, int spacedim = dim>
272  class Legendre : public Subscriptor
273  {
274  public:
276 
295  Legendre(const std::vector<unsigned int> &n_coefficients_per_direction,
298  const unsigned int component = numbers::invalid_unsigned_int);
299 
309  Legendre(const unsigned int n_coefficients_per_direction,
312 
318  template <typename Number>
319  void
320  calculate(const ::Vector<Number> &local_dof_values,
321  const unsigned int cell_active_fe_index,
322  Table<dim, CoefficientType> & legendre_coefficients);
323 
328  unsigned int
329  get_n_coefficients_per_direction(const unsigned int index) const;
330 
341  void
343 
353  template <class Archive>
354  void
355  save_transformation_matrices(Archive &ar, const unsigned int version);
356 
361  template <class Archive>
362  void
363  load_transformation_matrices(Archive &ar, const unsigned int version);
364 
368  bool
370 
371  private:
376  const std::vector<unsigned int> n_coefficients_per_direction;
377 
382 
387 
391  std::vector<FullMatrix<CoefficientType>> legendre_transform_matrices;
392 
396  std::vector<CoefficientType> unrolled_coefficients;
397 
402  const unsigned int component;
403  };
404 
405 
406 
424  template <int dim, typename CoefficientType>
425  std::pair<std::vector<unsigned int>, std::vector<double>>
427  const std::function<std::pair<bool, unsigned int>(
428  const TableIndices<dim> &)> & predicate,
429  const VectorTools::NormType norm_type,
430  const double smallest_abs_coefficient = 1e-10);
431 
437  std::pair<double, double>
438  linear_regression(const std::vector<double> &x, const std::vector<double> &y);
439 
440 } // namespace FESeries
441 
446 #ifndef DOXYGEN
447 
448 // ------------------- inline and template functions ----------------
449 
450 namespace internal
451 {
452  namespace FESeriesImplementation
453  {
454  template <int dim, typename CoefficientType>
455  void
456  fill_map_index(
457  const Table<dim, CoefficientType> &coefficients,
458  const TableIndices<dim> & ind,
459  const std::function<
460  std::pair<bool, unsigned int>(const TableIndices<dim> &)> &predicate,
461  std::map<unsigned int, std::vector<CoefficientType>> &pred_to_values)
462  {
463  const std::pair<bool, unsigned int> pred_pair = predicate(ind);
464  // don't add a value if predicate is false
465  if (pred_pair.first == false)
466  return;
467 
468  const unsigned int pred_value = pred_pair.second;
469  const CoefficientType &coeff_value = coefficients(ind);
470  // If pred_value is not in the pred_to_values map, the element will be
471  // created. Otherwise a reference to the existing element is returned.
472  pred_to_values[pred_value].push_back(coeff_value);
473  }
474 
475 
476 
477  template <typename CoefficientType>
478  void
479  fill_map(
480  const Table<1, CoefficientType> &coefficients,
481  const std::function<
482  std::pair<bool, unsigned int>(const TableIndices<1> &)> &predicate,
483  std::map<unsigned int, std::vector<CoefficientType>> & pred_to_values)
484  {
485  for (unsigned int i = 0; i < coefficients.size(0); i++)
486  {
487  const TableIndices<1> ind(i);
488  fill_map_index(coefficients, ind, predicate, pred_to_values);
489  }
490  }
491 
492 
493 
494  template <typename CoefficientType>
495  void
496  fill_map(
497  const Table<2, CoefficientType> &coefficients,
498  const std::function<
499  std::pair<bool, unsigned int>(const TableIndices<2> &)> &predicate,
500  std::map<unsigned int, std::vector<CoefficientType>> & pred_to_values)
501  {
502  for (unsigned int i = 0; i < coefficients.size(0); i++)
503  for (unsigned int j = 0; j < coefficients.size(1); j++)
504  {
505  const TableIndices<2> ind(i, j);
506  fill_map_index(coefficients, ind, predicate, pred_to_values);
507  }
508  }
509 
510 
511 
512  template <typename CoefficientType>
513  void
514  fill_map(
515  const Table<3, CoefficientType> &coefficients,
516  const std::function<
517  std::pair<bool, unsigned int>(const TableIndices<3> &)> &predicate,
518  std::map<unsigned int, std::vector<CoefficientType>> & pred_to_values)
519  {
520  for (unsigned int i = 0; i < coefficients.size(0); i++)
521  for (unsigned int j = 0; j < coefficients.size(1); j++)
522  for (unsigned int k = 0; k < coefficients.size(2); k++)
523  {
524  const TableIndices<3> ind(i, j, k);
525  fill_map_index(coefficients, ind, predicate, pred_to_values);
526  }
527  }
528 
529 
530 
531  template <typename Number>
532  double
533  complex_mean_value(const Number &value)
534  {
535  return value;
536  }
537 
538 
539 
540  template <typename Number>
541  double
542  complex_mean_value(const std::complex<Number> &value)
543  {
544  AssertThrow(false,
545  ExcMessage(
546  "FESeries::process_coefficients() can not be used with "
547  "complex-valued coefficients and VectorTools::mean norm."));
548  return std::abs(value);
549  }
550  } // namespace FESeriesImplementation
551 } // namespace internal
552 
553 
554 
555 template <int dim, typename CoefficientType>
556 std::pair<std::vector<unsigned int>, std::vector<double>>
558  const Table<dim, CoefficientType> &coefficients,
559  const std::function<std::pair<bool, unsigned int>(const TableIndices<dim> &)>
560  & predicate,
561  const VectorTools::NormType norm_type,
562  const double smallest_abs_coefficient)
563 {
564  Assert(smallest_abs_coefficient >= 0.,
565  ExcMessage("smallest_abs_coefficient should be non-negative."));
566 
567  std::vector<unsigned int> predicate_values;
568  std::vector<double> norm_values;
569 
570  // first, parse all table elements into a map of predicate values and
571  // coefficients. We could have stored (predicate values ->TableIndicies) map,
572  // but its processing would have been much harder later on.
573  std::map<unsigned int, std::vector<CoefficientType>> pred_to_values;
574  internal::FESeriesImplementation::fill_map(coefficients,
575  predicate,
576  pred_to_values);
577 
578  // now go through the map and populate the @p norm_values based on @p norm:
579  for (const auto &pred_to_value : pred_to_values)
580  {
581  Vector<CoefficientType> values(pred_to_value.second.cbegin(),
582  pred_to_value.second.cend());
583 
584  double norm_value = 0;
585  switch (norm_type)
586  {
588  {
589  norm_value = values.l2_norm();
590  break;
591  }
593  {
594  norm_value = values.l1_norm();
595  break;
596  }
598  {
599  norm_value = values.linfty_norm();
600  break;
601  }
602  case VectorTools::mean:
603  {
604  norm_value = internal::FESeriesImplementation::complex_mean_value(
605  values.mean_value());
606  break;
607  }
608  default:
609  AssertThrow(false, ExcNotImplemented());
610  break;
611  }
612 
613  // will use all non-zero coefficients
614  if (std::abs(norm_value) > smallest_abs_coefficient)
615  {
616  predicate_values.push_back(pred_to_value.first);
617  norm_values.push_back(norm_value);
618  }
619  }
620 
621  return std::make_pair(predicate_values, norm_values);
622 }
623 
624 
625 
626 template <int dim, int spacedim>
627 template <class Archive>
628 inline void
630  Archive &ar,
631  const unsigned int /*version*/)
632 {
633  // Store information about those resources which have been used to generate
634  // the transformation matrices.
635  // mode vector
636  ar &n_coefficients_per_direction;
637 
638  // finite element collection
639  unsigned int size = fe_collection->size();
640  ar & size;
641  for (unsigned int i = 0; i < size; ++i)
642  ar &(*fe_collection)[i].get_name();
643 
644  // quadrature collection
645  size = q_collection.size();
646  ar &size;
647  for (unsigned int i = 0; i < size; ++i)
648  ar &q_collection[i];
649 
650  // Store the actual transform matrices.
651  ar &fourier_transform_matrices;
652 }
653 
654 
655 
656 template <int dim, int spacedim>
657 template <class Archive>
658 inline void
660  Archive &ar,
661  const unsigned int /*version*/)
662 {
663  // Check whether the currently registered resources are compatible with
664  // the transformation matrices to load.
665  // mode vector
666  std::vector<unsigned int> compare_coefficients;
667  ar & compare_coefficients;
668  Assert(compare_coefficients == n_coefficients_per_direction,
669  ExcMessage("A different number of coefficients vector has been used "
670  "to generate the transformation matrices you are about "
671  "to load!"));
672 
673  // finite element collection
674  unsigned int size;
675  ar & size;
676  AssertDimension(size, fe_collection->size());
677  std::string name;
678  for (unsigned int i = 0; i < size; ++i)
679  {
680  ar &name;
681  Assert(name.compare((*fe_collection)[i].get_name()) == 0,
682  ExcMessage("A different FECollection has been used to generate "
683  "the transformation matrices you are about to load!"));
684  }
685 
686  // quadrature collection
687  ar &size;
688  AssertDimension(size, q_collection.size());
689  Quadrature<dim> quadrature;
690  for (unsigned int i = 0; i < size; ++i)
691  {
692  ar &quadrature;
693  Assert(quadrature == q_collection[i],
694  ExcMessage("A different QCollection has been used to generate "
695  "the transformation matrices you are about to load!"));
696  }
697 
698  // Restore the transform matrices since all prerequisites are fulfilled.
699  ar &fourier_transform_matrices;
700 }
701 
702 
703 
704 template <int dim, int spacedim>
705 template <class Archive>
706 inline void
708  Archive &ar,
709  const unsigned int /*version*/)
710 {
711  // Store information about those resources which have been used to generate
712  // the transformation matrices.
713  // mode vector
714  ar &n_coefficients_per_direction;
715 
716  // finite element collection
717  unsigned int size = fe_collection->size();
718  ar & size;
719  for (unsigned int i = 0; i < size; ++i)
720  ar &(*fe_collection)[i].get_name();
721 
722  // quadrature collection
723  size = q_collection.size();
724  ar &size;
725  for (unsigned int i = 0; i < size; ++i)
726  ar &q_collection[i];
727 
728  // Store the actual transform matrices.
729  ar &legendre_transform_matrices;
730 }
731 
732 
733 
734 template <int dim, int spacedim>
735 template <class Archive>
736 inline void
738  Archive &ar,
739  const unsigned int /*version*/)
740 {
741  // Check whether the currently registered resources are compatible with
742  // the transformation matrices to load.
743  // mode vector
744  std::vector<unsigned int> compare_coefficients;
745  ar & compare_coefficients;
746  Assert(compare_coefficients == n_coefficients_per_direction,
747  ExcMessage("A different number of coefficients vector has been used "
748  "to generate the transformation matrices you are about "
749  "to load!"));
750 
751  // finite element collection
752  unsigned int size;
753  ar & size;
754  AssertDimension(size, fe_collection->size());
755  std::string name;
756  for (unsigned int i = 0; i < size; ++i)
757  {
758  ar &name;
759  Assert(name.compare((*fe_collection)[i].get_name()) == 0,
760  ExcMessage("A different FECollection has been used to generate "
761  "the transformation matrices you are about to load!"));
762  }
763 
764  // quadrature collection
765  ar &size;
766  AssertDimension(size, q_collection.size());
767  Quadrature<dim> quadrature;
768  for (unsigned int i = 0; i < size; ++i)
769  {
770  ar &quadrature;
771  Assert(quadrature == q_collection[i],
772  ExcMessage("A different QCollection has been used to generate "
773  "the transformation matrices you are about to load!"));
774  }
775 
776  // Restore the transform matrices since all prerequisites are fulfilled.
777  ar &legendre_transform_matrices;
778 }
779 
780 
781 #endif // DOXYGEN
782 
784 
785 #endif // dealii_fe_series_h
SmartPointer< const hp::FECollection< dim, spacedim > > fe_collection
Definition: fe_series.h:198
void calculate(const ::Vector< Number > &local_dof_values, const unsigned int cell_active_fe_index, Table< dim, CoefficientType > &fourier_coefficients)
Table< dim, Tensor< 1, dim > > k_vectors
Definition: fe_series.h:208
std::vector< CoefficientType > unrolled_coefficients
Definition: fe_series.h:218
const unsigned int component
Definition: fe_series.h:224
typename std::complex< double > CoefficientType
Definition: fe_series.h:91
unsigned int get_n_coefficients_per_direction(const unsigned int index) const
const std::vector< unsigned int > n_coefficients_per_direction
Definition: fe_series.h:193
bool operator==(const Fourier< dim, spacedim > &fourier) const
std::vector< FullMatrix< CoefficientType > > fourier_transform_matrices
Definition: fe_series.h:213
const hp::QCollection< dim > q_collection
Definition: fe_series.h:203
Fourier(const std::vector< unsigned int > &n_coefficients_per_direction, const hp::FECollection< dim, spacedim > &fe_collection, const hp::QCollection< dim > &q_collection, const unsigned int component=numbers::invalid_unsigned_int)
void precalculate_all_transformation_matrices()
void load_transformation_matrices(Archive &ar, const unsigned int version)
void save_transformation_matrices(Archive &ar, const unsigned int version)
const std::vector< unsigned int > n_coefficients_per_direction
Definition: fe_series.h:376
void save_transformation_matrices(Archive &ar, const unsigned int version)
bool operator==(const Legendre< dim, spacedim > &legendre) const
const unsigned int component
Definition: fe_series.h:402
std::vector< FullMatrix< CoefficientType > > legendre_transform_matrices
Definition: fe_series.h:391
const hp::QCollection< dim > q_collection
Definition: fe_series.h:386
Legendre(const std::vector< unsigned int > &n_coefficients_per_direction, const hp::FECollection< dim, spacedim > &fe_collection, const hp::QCollection< dim > &q_collection, const unsigned int component=numbers::invalid_unsigned_int)
void precalculate_all_transformation_matrices()
unsigned int get_n_coefficients_per_direction(const unsigned int index) const
std::vector< CoefficientType > unrolled_coefficients
Definition: fe_series.h:396
SmartPointer< const hp::FECollection< dim, spacedim > > fe_collection
Definition: fe_series.h:381
void load_transformation_matrices(Archive &ar, const unsigned int version)
void calculate(const ::Vector< Number > &local_dof_values, const unsigned int cell_active_fe_index, Table< dim, CoefficientType > &legendre_coefficients)
Definition: vector.h:110
#define DEAL_II_DEPRECATED
Definition: config.h:159
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:396
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:397
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcNotImplemented()
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1622
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1575
std::pair< std::vector< unsigned int >, std::vector< double > > process_coefficients(const Table< dim, CoefficientType > &coefficients, const std::function< std::pair< bool, unsigned int >(const TableIndices< dim > &)> &predicate, const VectorTools::NormType norm_type, const double smallest_abs_coefficient=1e-10)
std::pair< double, double > linear_regression(const std::vector< double > &x, const std::vector< double > &y)
Definition: fe_series.cc:30
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
static const unsigned int invalid_unsigned_int
Definition: types.h:196
double legendre(unsigned int l, double x)
Definition: cmath.h:75
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)