16 #ifndef dealii_fe_series_h
17 #define dealii_fe_series_h
90 template <
int dim,
int spacedim = dim>
127 template <
typename Number>
129 calculate(const ::Vector<Number> &local_dof_values,
130 const unsigned int cell_active_fe_index,
162 template <
class Archive>
170 template <
class Archive>
258 template <
int dim,
int spacedim = dim>
294 template <
typename Number>
296 calculate(const ::Vector<Number> &local_dof_values,
297 const unsigned int cell_active_fe_index,
329 template <
class Archive>
337 template <
class Archive>
394 template <
int dim,
typename CoefficientType>
395 std::pair<std::vector<unsigned int>, std::vector<double>>
397 const std::function<std::pair<bool, unsigned int>(
400 const double smallest_abs_coefficient = 1
e-10);
407 std::pair<double, double>
422 namespace FESeriesImplementation
424 template <
int dim,
typename CoefficientType>
431 std::map<
unsigned int, std::vector<CoefficientType>> &pred_to_values)
433 const std::pair<bool, unsigned int> pred_pair = predicate(ind);
435 if (pred_pair.first ==
false)
438 const unsigned int pred_value = pred_pair.second;
439 const CoefficientType &coeff_value = coefficients(ind);
442 pred_to_values[pred_value].push_back(coeff_value);
447 template <
typename CoefficientType>
453 std::map<
unsigned int, std::vector<CoefficientType>> & pred_to_values)
455 for (
unsigned int i = 0; i < coefficients.
size(0); i++)
458 fill_map_index(coefficients, ind, predicate, pred_to_values);
464 template <
typename CoefficientType>
470 std::map<
unsigned int, std::vector<CoefficientType>> & pred_to_values)
472 for (
unsigned int i = 0; i < coefficients.
size(0); i++)
473 for (
unsigned int j = 0; j < coefficients.
size(1); j++)
476 fill_map_index(coefficients, ind, predicate, pred_to_values);
482 template <
typename CoefficientType>
488 std::map<
unsigned int, std::vector<CoefficientType>> & pred_to_values)
490 for (
unsigned int i = 0; i < coefficients.
size(0); i++)
491 for (
unsigned int j = 0; j < coefficients.
size(1); j++)
492 for (
unsigned int k = 0; k < coefficients.
size(2); k++)
495 fill_map_index(coefficients, ind, predicate, pred_to_values);
501 template <
typename Number>
503 complex_mean_value(
const Number &
value)
510 template <
typename Number>
512 complex_mean_value(
const std::complex<Number> &
value)
516 "FESeries::process_coefficients() can not be used with "
517 "complex-valued coefficients and VectorTools::mean norm."));
525 template <
int dim,
typename CoefficientType>
526 std::pair<std::vector<unsigned int>, std::vector<double>>
532 const double smallest_abs_coefficient)
534 Assert(smallest_abs_coefficient >= 0.,
535 ExcMessage(
"smallest_abs_coefficient should be non-negative."));
537 std::vector<unsigned int> predicate_values;
538 std::vector<double> norm_values;
543 std::map<unsigned int, std::vector<CoefficientType>> pred_to_values;
544 internal::FESeriesImplementation::fill_map(coefficients,
549 for (
const auto &pred_to_value : pred_to_values)
551 Vector<CoefficientType> values(pred_to_value.second.cbegin(),
552 pred_to_value.second.cend());
554 double norm_value = 0;
559 norm_value = values.l2_norm();
564 norm_value = values.l1_norm();
569 norm_value = values.linfty_norm();
574 norm_value = internal::FESeriesImplementation::complex_mean_value(
575 values.mean_value());
584 if (std::abs(norm_value) > smallest_abs_coefficient)
586 predicate_values.push_back(pred_to_value.first);
587 norm_values.push_back(norm_value);
591 return std::make_pair(predicate_values, norm_values);
596 template <
int dim,
int spacedim>
597 template <
class Archive>
606 ar &n_coefficients_per_direction;
609 unsigned int size = fe_collection->size();
611 for (
unsigned int i = 0; i < size; ++i)
612 ar &(*fe_collection)[i].get_name();
615 size = q_collection.size();
617 for (
unsigned int i = 0; i < size; ++i)
621 ar &fourier_transform_matrices;
626 template <
int dim,
int spacedim>
627 template <
class Archive>
636 std::vector<unsigned int> compare_coefficients;
637 ar & compare_coefficients;
638 Assert(compare_coefficients == n_coefficients_per_direction,
639 ExcMessage(
"A different number of coefficients vector has been used "
640 "to generate the transformation matrices you are about "
648 for (
unsigned int i = 0; i < size; ++i)
651 Assert(name.compare((*fe_collection)[i].get_name()) == 0,
652 ExcMessage(
"A different FECollection has been used to generate "
653 "the transformation matrices you are about to load!"));
660 for (
unsigned int i = 0; i < size; ++i)
663 Assert(quadrature == q_collection[i],
664 ExcMessage(
"A different QCollection has been used to generate "
665 "the transformation matrices you are about to load!"));
669 ar &fourier_transform_matrices;
674 template <
int dim,
int spacedim>
675 template <
class Archive>
684 ar &n_coefficients_per_direction;
687 unsigned int size = fe_collection->size();
689 for (
unsigned int i = 0; i < size; ++i)
690 ar &(*fe_collection)[i].get_name();
693 size = q_collection.size();
695 for (
unsigned int i = 0; i < size; ++i)
699 ar &legendre_transform_matrices;
704 template <
int dim,
int spacedim>
705 template <
class Archive>
714 std::vector<unsigned int> compare_coefficients;
715 ar & compare_coefficients;
716 Assert(compare_coefficients == n_coefficients_per_direction,
717 ExcMessage(
"A different number of coefficients vector has been used "
718 "to generate the transformation matrices you are about "
726 for (
unsigned int i = 0; i < size; ++i)
729 Assert(name.compare((*fe_collection)[i].get_name()) == 0,
730 ExcMessage(
"A different FECollection has been used to generate "
731 "the transformation matrices you are about to load!"));
738 for (
unsigned int i = 0; i < size; ++i)
741 Assert(quadrature == q_collection[i],
742 ExcMessage(
"A different QCollection has been used to generate "
743 "the transformation matrices you are about to load!"));
747 ar &legendre_transform_matrices;
755 #endif // dealii_fe_series_h