16#ifndef dealii_fe_series_h
17#define dealii_fe_series_h
87 template <
int dim,
int spacedim = dim>
135 template <
typename Number>
138 const unsigned int cell_active_fe_index,
170 template <
class Archive>
178 template <
class Archive>
271 template <
int dim,
int spacedim = dim>
318 template <
typename Number>
320 calculate(const ::Vector<Number> &local_dof_values,
321 const unsigned int cell_active_fe_index,
353 template <
class Archive>
361 template <
class Archive>
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>(
430 const double smallest_abs_coefficient = 1
e-10);
437 std::pair<double, double>
452 namespace FESeriesImplementation
454 template <
int dim,
typename CoefficientType>
461 std::map<
unsigned int, std::vector<CoefficientType>> &pred_to_values)
463 const std::pair<bool, unsigned int> pred_pair = predicate(ind);
465 if (pred_pair.first ==
false)
468 const unsigned int pred_value = pred_pair.second;
469 const CoefficientType &coeff_value = coefficients(ind);
472 pred_to_values[pred_value].push_back(coeff_value);
477 template <
typename CoefficientType>
483 std::map<
unsigned int, std::vector<CoefficientType>> & pred_to_values)
485 for (
unsigned int i = 0; i < coefficients.size(0); i++)
488 fill_map_index(coefficients, ind, predicate, pred_to_values);
494 template <
typename CoefficientType>
500 std::map<
unsigned int, std::vector<CoefficientType>> & pred_to_values)
502 for (
unsigned int i = 0; i < coefficients.size(0); i++)
503 for (
unsigned int j = 0; j < coefficients.size(1); j++)
506 fill_map_index(coefficients, ind, predicate, pred_to_values);
512 template <
typename CoefficientType>
518 std::map<
unsigned int, std::vector<CoefficientType>> & pred_to_values)
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++)
525 fill_map_index(coefficients, ind, predicate, pred_to_values);
531 template <
typename Number>
533 complex_mean_value(
const Number &value)
540 template <
typename Number>
542 complex_mean_value(
const std::complex<Number> &value)
546 "FESeries::process_coefficients() can not be used with "
547 "complex-valued coefficients and VectorTools::mean norm."));
555template <
int dim,
typename CoefficientType>
556std::pair<std::vector<unsigned int>, std::vector<double>>
562 const double smallest_abs_coefficient)
564 Assert(smallest_abs_coefficient >= 0.,
565 ExcMessage(
"smallest_abs_coefficient should be non-negative."));
567 std::vector<unsigned int> predicate_values;
568 std::vector<double> norm_values;
573 std::map<unsigned int, std::vector<CoefficientType>> pred_to_values;
574 internal::FESeriesImplementation::fill_map(coefficients,
579 for (
const auto &pred_to_value : pred_to_values)
582 pred_to_value.second.cend());
584 double norm_value = 0;
589 norm_value =
values.l2_norm();
594 norm_value =
values.l1_norm();
599 norm_value =
values.linfty_norm();
604 norm_value = internal::FESeriesImplementation::complex_mean_value(
614 if (
std::abs(norm_value) > smallest_abs_coefficient)
616 predicate_values.push_back(pred_to_value.first);
617 norm_values.push_back(norm_value);
621 return std::make_pair(predicate_values, norm_values);
626template <
int dim,
int spacedim>
627template <
class Archive>
636 ar &n_coefficients_per_direction;
639 unsigned int size = fe_collection->size();
641 for (
unsigned int i = 0; i < size; ++i)
642 ar &(*fe_collection)[i].get_name();
645 size = q_collection.size();
647 for (
unsigned int i = 0; i < size; ++i)
651 ar &fourier_transform_matrices;
656template <
int dim,
int spacedim>
657template <
class Archive>
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 "
678 for (
unsigned int i = 0; i < size; ++i)
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!"));
690 for (
unsigned int i = 0; i < size; ++i)
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!"));
699 ar &fourier_transform_matrices;
704template <
int dim,
int spacedim>
705template <
class Archive>
714 ar &n_coefficients_per_direction;
717 unsigned int size = fe_collection->size();
719 for (
unsigned int i = 0; i < size; ++i)
720 ar &(*fe_collection)[i].get_name();
723 size = q_collection.size();
725 for (
unsigned int i = 0; i < size; ++i)
729 ar &legendre_transform_matrices;
734template <
int dim,
int spacedim>
735template <
class Archive>
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 "
756 for (
unsigned int i = 0; i < size; ++i)
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!"));
768 for (
unsigned int i = 0; i < size; ++i)
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!"));
777 ar &legendre_transform_matrices;
SmartPointer< const hp::FECollection< dim, spacedim > > fe_collection
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
std::vector< CoefficientType > unrolled_coefficients
const unsigned int component
typename std::complex< double > CoefficientType
unsigned int get_n_coefficients_per_direction(const unsigned int index) const
const std::vector< unsigned int > n_coefficients_per_direction
bool operator==(const Fourier< dim, spacedim > &fourier) const
std::vector< FullMatrix< CoefficientType > > fourier_transform_matrices
const hp::QCollection< dim > q_collection
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
void save_transformation_matrices(Archive &ar, const unsigned int version)
bool operator==(const Legendre< dim, spacedim > &legendre) const
const unsigned int component
std::vector< FullMatrix< CoefficientType > > legendre_transform_matrices
const hp::QCollection< dim > q_collection
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
SmartPointer< const hp::FECollection< dim, spacedim > > fe_collection
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)
#define DEAL_II_DEPRECATED
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
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)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
static const unsigned int invalid_unsigned_int
double legendre(unsigned int l, double x)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)