16 #ifndef dealii_fe_series_H 17 #define dealii_fe_series_H 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> 38 DEAL_II_NAMESPACE_OPEN
97 Fourier(
const unsigned int size_in_each_direction,
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);
198 Legendre(
const unsigned int size_in_each_direction,
207 void calculate(const ::Vector<double> &local_dof_values,
208 const unsigned int cell_active_fe_index,
215 const unsigned int N;
260 template <
int dim,
typename T>
261 std::pair<std::vector<unsigned int>,std::vector<double> >
263 const std::function<std::pair<bool,unsigned int>(
const TableIndices<dim> &)> &predicate,
274 const std::vector<double> &y);
286 template <
int dim,
typename T>
289 const std::function<std::pair<bool,unsigned int>(
const TableIndices<dim> &)> &predicate,
290 std::map<
unsigned int, std::vector<T> > &pred_to_values)
292 const std::pair<bool,unsigned int> pred_pair = predicate(ind);
294 if (pred_pair.first ==
false)
297 const unsigned int &pred_value = pred_pair.second;
298 const T &coeff_value = coefficients(ind);
301 pred_to_values[pred_value].push_back(coeff_value);
304 template <
typename T>
306 const std::function<std::pair<bool,unsigned int>(
const TableIndices<1> &)> &predicate,
307 std::map<
unsigned int, std::vector<T> > &pred_to_values)
309 for (
unsigned int i = 0; i < coefficients.
size(0); i++)
312 fill_map_index(coefficients,ind,predicate,pred_to_values);
317 template <
typename T>
319 const std::function<std::pair<bool,unsigned int>(
const TableIndices<2> &)> &predicate,
320 std::map<
unsigned int, std::vector<T> > &pred_to_values)
322 for (
unsigned int i = 0; i < coefficients.
size(0); i++)
323 for (
unsigned int j = 0; j < coefficients.
size(1); j++)
326 fill_map_index(coefficients,ind,predicate,pred_to_values);
331 template <
typename T>
333 const std::function<std::pair<bool,unsigned int>(
const TableIndices<3> &)> &predicate,
334 std::map<
unsigned int, std::vector<T> > &pred_to_values)
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++)
341 fill_map_index(coefficients,ind,predicate,pred_to_values);
346 template <
typename T>
347 double complex_mean_value(
const T &value)
352 template <
typename T>
353 double complex_mean_value(
const std::complex<T> &value)
356 "complex-valued coefficients and VectorTools::mean norm."));
357 return std::abs(value);
363 template <
int dim,
typename T>
364 std::pair<std::vector<unsigned int>,std::vector<double> >
366 const std::function<std::pair<bool,unsigned int>(
const TableIndices<dim> &)> &predicate,
369 std::vector<unsigned int> predicate_values;
370 std::vector<double> norm_values;
375 std::map<unsigned int, std::vector<T> > pred_to_values;
376 fill_map(coefficients,predicate,pred_to_values);
379 for (
typename std::map<
unsigned int, std::vector<T> >::const_iterator it = pred_to_values.begin();
380 it != pred_to_values.end(); ++it)
382 predicate_values.push_back(it->first);
383 Vector<T> values(it->second.begin(),
390 norm_values.push_back(values.l2_norm());
395 norm_values.push_back(values.l1_norm());
400 norm_values.push_back(values.linfty_norm());
405 norm_values.push_back(complex_mean_value(values.mean_value()));
414 return std::make_pair(predicate_values,norm_values);
420 DEAL_II_NAMESPACE_CLOSE
422 #endif // dealii_fe_series_H
void ensure_existence(const unsigned int fe_index)
SmartPointer< const hp::FECollection< dim > > fe_collection
Legendre(const unsigned int size_in_each_direction, const hp::FECollection< dim > &fe_collection, const hp::QCollection< dim > &q_collection)
std::pair< double, double > linear_regression(const std::vector< double > &x, const std::vector< double > &y)
SmartPointer< const hp::QCollection< dim > > q_collection
#define AssertThrow(cond, exc)
void calculate(const ::Vector< double > &local_dof_values, const unsigned int cell_active_fe_index, Table< dim, double > &legendre_coefficients)
Fourier(const unsigned int size_in_each_direction, const hp::FECollection< dim > &fe_collection, const hp::QCollection< dim > &q_collection)
std::vector< std::complex< double > > unrolled_coefficients
SmartPointer< const hp::QCollection< dim > > q_collection
static ::ExceptionBase & ExcMessage(std::string arg1)
SmartPointer< const hp::FECollection< dim > > fe_collection
void calculate(const ::Vector< double > &local_dof_values, const unsigned int cell_active_fe_index, Table< dim, std::complex< double > > &fourier_coefficients)
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
static ::ExceptionBase & ExcNotImplemented()
Table< dim, Tensor< 1, dim > > k_vectors
std::vector< FullMatrix< std::complex< double > > > fourier_transform_matrices
std::vector< FullMatrix< double > > legendre_transform_matrices