33 <<
"x[" << arg1 <<
"] = " << arg2 <<
" is not in [0,1]");
44 for (
unsigned int d = 0;
d < dim; ++
d)
46 const double x = 2.0 * (x_q[
d] - 0.5);
47 Assert((x_q[d] <= 1.0) && (x_q[d] >= 0.), ExcLegendre(d, x_q[d]));
48 const unsigned int ind = indices[
d];
64 for (
unsigned int d = 0;
d < dim; ++
d)
65 res *= (0.5 + indices[d]);
72 template <
int dim,
int spacedim>
77 const unsigned int dof,
78 const unsigned int component)
81 for (
unsigned int q = 0; q < quadrature.
size(); ++q)
84 sum += Lh(x_q, indices) *
88 return sum * multiplier(indices);
97 template <
int spacedim>
100 const std::vector<unsigned int> & n_coefficients_per_direction,
103 const unsigned int fe,
104 const unsigned int component,
109 if (legendre_transform_matrices[fe].m() == 0)
111 legendre_transform_matrices[fe].reinit(
112 n_coefficients_per_direction[fe],
113 fe_collection[fe].n_dofs_per_cell());
115 for (
unsigned int k = 0; k < n_coefficients_per_direction[fe]; ++k)
116 for (
unsigned int j = 0; j < fe_collection[fe].n_dofs_per_cell(); ++j)
117 legendre_transform_matrices[fe](k, j) =
118 integrate(fe_collection[fe],
126 template <
int spacedim>
129 const std::vector<unsigned int> & n_coefficients_per_direction,
132 const unsigned int fe,
133 const unsigned int component,
138 if (legendre_transform_matrices[fe].m() == 0)
140 legendre_transform_matrices[fe].reinit(
141 Utilities::fixed_power<2>(n_coefficients_per_direction[fe]),
142 fe_collection[fe].n_dofs_per_cell());
145 for (
unsigned int k1 = 0; k1 < n_coefficients_per_direction[fe]; ++k1)
146 for (
unsigned int k2 = 0; k2 < n_coefficients_per_direction[fe];
148 for (
unsigned int j = 0; j < fe_collection[fe].n_dofs_per_cell();
150 legendre_transform_matrices[fe](k, j) =
151 integrate(fe_collection[fe],
159 template <
int spacedim>
162 const std::vector<unsigned int> & n_coefficients_per_direction,
165 const unsigned int fe,
166 const unsigned int component,
171 if (legendre_transform_matrices[fe].m() == 0)
173 legendre_transform_matrices[fe].reinit(
174 Utilities::fixed_power<3>(n_coefficients_per_direction[fe]),
175 fe_collection[fe].n_dofs_per_cell());
178 for (
unsigned int k1 = 0; k1 < n_coefficients_per_direction[fe]; ++k1)
179 for (
unsigned int k2 = 0; k2 < n_coefficients_per_direction[fe]; ++k2)
180 for (
unsigned int k3 = 0; k3 < n_coefficients_per_direction[fe];
182 for (
unsigned int j = 0; j < fe_collection[fe].n_dofs_per_cell();
184 legendre_transform_matrices[fe](k, j) =
185 integrate(fe_collection[fe],
198 template <
int dim,
int spacedim>
200 const std::vector<unsigned int> & n_coefficients_per_direction,
203 const unsigned int component_)
204 : n_coefficients_per_direction(n_coefficients_per_direction)
205 , fe_collection(&fe_collection)
206 , q_collection(q_collection)
207 , legendre_transform_matrices(fe_collection.size())
208 , component(component_ !=
numbers::invalid_unsigned_int ? component_ : 0)
212 ExcMessage(
"All parameters are supposed to have the same size."));
218 "For vector-valued problems, you need to explicitly specify for "
219 "which vector component you will want to do a Legendre decomposition "
220 "by setting the 'component' argument of this constructor."));
225 const unsigned int max_n_coefficients_per_direction =
229 Utilities::fixed_power<dim>(max_n_coefficients_per_direction));
234 template <
int dim,
int spacedim>
240 (n_coefficients_per_direction == legendre.n_coefficients_per_direction) &&
241 (*fe_collection == *(legendre.fe_collection)) &&
242 (q_collection == legendre.q_collection) &&
243 (legendre_transform_matrices == legendre.legendre_transform_matrices) &&
244 (component == legendre.component));
249 template <
int dim,
int spacedim>
254 for (
unsigned int fe = 0; fe < fe_collection->
size(); ++fe)
256 ensure_existence(n_coefficients_per_direction,
261 legendre_transform_matrices);
269 template <
int dim,
int spacedim>
272 const unsigned int index)
const
274 return n_coefficients_per_direction[index];
279 template <
int dim,
int spacedim>
280 template <
typename Number>
283 const ::Vector<Number> &local_dof_values,
284 const unsigned int cell_active_fe_index,
287 for (
unsigned int d = 0; d < dim; ++d)
289 n_coefficients_per_direction[cell_active_fe_index]);
291 ensure_existence(n_coefficients_per_direction,
294 cell_active_fe_index,
296 legendre_transform_matrices);
299 legendre_transform_matrices[cell_active_fe_index];
301 unrolled_coefficients.resize(Utilities::fixed_power<dim>(
302 n_coefficients_per_direction[cell_active_fe_index]));
303 std::fill(unrolled_coefficients.begin(),
304 unrolled_coefficients.end(),
309 Assert(local_dof_values.size() == matrix.n(),
312 for (
unsigned int i = 0; i < unrolled_coefficients.size(); ++i)
313 for (
unsigned int j = 0; j < local_dof_values.size(); ++j)
314 unrolled_coefficients[i] += matrix[i][j] * local_dof_values[j];
316 legendre_coefficients.fill(unrolled_coefficients.begin());
322#include "fe_series_legendre.inst"
const std::vector< unsigned int > n_coefficients_per_direction
bool operator==(const Legendre< dim, spacedim > &legendre) const
const unsigned int component
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 calculate(const ::Vector< Number > &local_dof_values, const unsigned int cell_active_fe_index, Table< dim, CoefficientType > &legendre_coefficients)
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 DeclException2(Exception2, type1, type2, outsequence)
#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)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
T sum(const T &t, const MPI_Comm mpi_communicator)
static const unsigned int invalid_unsigned_int
double legendre(unsigned int l, double x)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)