34 for (
unsigned int i = 0; i < N; ++i)
42 for (
unsigned int i = 0; i < N; ++i)
43 for (
unsigned int j = 0; j < N; ++j)
54 for (
unsigned int i = 0; i < N; ++i)
55 for (
unsigned int j = 0; j < N; ++j)
56 for (
unsigned int k = 0; k < N; ++k)
66 template <
int dim,
int spacedim>
72 const unsigned int component)
74 std::complex<double>
sum = 0;
75 for (
unsigned int q = 0; q < quadrature.
size(); ++q)
78 sum +=
std::exp(std::complex<double>(0, 1) * (k_vector * x_q)) *
91 template <
int spacedim>
94 const std::vector<unsigned int> & n_coefficients_per_direction,
98 const unsigned int fe,
99 const unsigned int component,
100 std::vector<
FullMatrix<std::complex<double>>> &fourier_transform_matrices)
104 if (fourier_transform_matrices[fe].m() == 0)
106 fourier_transform_matrices[fe].reinit(
107 n_coefficients_per_direction[fe],
108 fe_collection[fe].n_dofs_per_cell());
110 for (
unsigned int k = 0; k < n_coefficients_per_direction[fe]; ++k)
111 for (
unsigned int j = 0; j < fe_collection[fe].n_dofs_per_cell(); ++j)
112 fourier_transform_matrices[fe](k, j) = integrate(
113 fe_collection[fe], q_collection[fe], k_vectors(k), j, component);
117 template <
int spacedim>
120 const std::vector<unsigned int> & n_coefficients_per_direction,
124 const unsigned int fe,
125 const unsigned int component,
126 std::vector<
FullMatrix<std::complex<double>>> &fourier_transform_matrices)
130 if (fourier_transform_matrices[fe].m() == 0)
132 fourier_transform_matrices[fe].reinit(
133 Utilities::fixed_power<2>(n_coefficients_per_direction[fe]),
134 fe_collection[fe].n_dofs_per_cell());
137 for (
unsigned int k1 = 0; k1 < n_coefficients_per_direction[fe]; ++k1)
138 for (
unsigned int k2 = 0; k2 < n_coefficients_per_direction[fe];
140 for (
unsigned int j = 0; j < fe_collection[fe].n_dofs_per_cell();
142 fourier_transform_matrices[fe](k, j) =
143 integrate(fe_collection[fe],
151 template <
int spacedim>
154 const std::vector<unsigned int> & n_coefficients_per_direction,
158 const unsigned int fe,
159 const unsigned int component,
160 std::vector<
FullMatrix<std::complex<double>>> &fourier_transform_matrices)
164 if (fourier_transform_matrices[fe].m() == 0)
166 fourier_transform_matrices[fe].reinit(
167 Utilities::fixed_power<3>(n_coefficients_per_direction[fe]),
168 fe_collection[fe].n_dofs_per_cell());
171 for (
unsigned int k1 = 0; k1 < n_coefficients_per_direction[fe]; ++k1)
172 for (
unsigned int k2 = 0; k2 < n_coefficients_per_direction[fe]; ++k2)
173 for (
unsigned int k3 = 0; k3 < n_coefficients_per_direction[fe];
175 for (
unsigned int j = 0; j < fe_collection[fe].n_dofs_per_cell();
177 fourier_transform_matrices[fe](k, j) =
178 integrate(fe_collection[fe],
180 k_vectors(k1, k2, k3),
191 template <
int dim,
int spacedim>
193 const std::vector<unsigned int> & n_coefficients_per_direction,
196 const unsigned int component_)
197 : n_coefficients_per_direction(n_coefficients_per_direction)
198 , fe_collection(&fe_collection)
199 , q_collection(q_collection)
200 , fourier_transform_matrices(fe_collection.size())
201 , component(component_ !=
numbers::invalid_unsigned_int ? component_ : 0)
205 ExcMessage(
"All parameters are supposed to have the same size."));
211 "For vector-valued problems, you need to explicitly specify for "
212 "which vector component you will want to do a Fourier decomposition "
213 "by setting the 'component' argument of this constructor."));
217 const unsigned int max_n_coefficients_per_direction =
220 set_k_vectors(
k_vectors, max_n_coefficients_per_direction);
228 template <
int dim,
int spacedim>
244 template <
int dim,
int spacedim>
249 for (
unsigned int fe = 0; fe < fe_collection->
size(); ++fe)
251 ensure_existence(n_coefficients_per_direction,
257 fourier_transform_matrices);
265 template <
int dim,
int spacedim>
268 const unsigned int index)
const
270 return n_coefficients_per_direction[index];
275 template <
int dim,
int spacedim>
276 template <
typename Number>
280 const unsigned int cell_active_fe_index,
283 for (
unsigned int d = 0; d < dim; ++d)
285 n_coefficients_per_direction[cell_active_fe_index]);
287 ensure_existence(n_coefficients_per_direction,
291 cell_active_fe_index,
293 fourier_transform_matrices);
296 fourier_transform_matrices[cell_active_fe_index];
298 unrolled_coefficients.resize(Utilities::fixed_power<dim>(
299 n_coefficients_per_direction[cell_active_fe_index]));
300 std::fill(unrolled_coefficients.begin(),
301 unrolled_coefficients.end(),
309 for (
unsigned int i = 0; i < unrolled_coefficients.size(); ++i)
310 for (
unsigned int j = 0; j < local_dof_values.
size(); ++j)
311 unrolled_coefficients[i] += matrix[i][j] * local_dof_values[j];
313 fourier_coefficients.fill(unrolled_coefficients.begin());
319#include "fe_series_fourier.inst"
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()
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
const Point< dim > & point(const unsigned int i) const
double weight(const unsigned int i) const
unsigned int size() const
unsigned int size() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
Task< RT > new_task(const std::function< RT()> &function)
T sum(const T &t, const MPI_Comm mpi_communicator)
static constexpr double PI
static const unsigned int invalid_unsigned_int
::VectorizedArray< Number, width > exp(const ::VectorizedArray< Number, width > &)