* dimensional polynomials for each space direction, given the index
* <i>i</i>.
*/
- // fix to avoid compiler warnings about zero length arrays
void
- compute_index(const unsigned int i,
- unsigned int (&indices)[(dim > 0 ? dim : 1)]) const;
+ compute_index(const unsigned int i,
+ std::array<unsigned int, dim> &indices) const;
/**
* TensorProductPolynomialsBubbles has a TensorProductPolynomials class
* <tt>i</tt>.
*/
void
- compute_index(const unsigned int i,
- unsigned int (&indices)[(dim > 0 ? dim : 1)]) const;
+ compute_index(const unsigned int i,
+ std::array<unsigned int, dim> &indices) const;
/**
* Given the input to the constructor, compute <tt>n_pols</tt>.
const unsigned int i,
const Point<dim> & p) const
{
- unsigned int indices[dim];
+ std::array<unsigned int, dim> indices;
compute_index(i, indices);
- double v[dim][5];
+ std::array<std::array<double, 5>, dim> v;
{
std::vector<double> tmp(5);
for (unsigned int d = 0; d < dim; ++d)
AnisotropicPolynomials<dim>::compute_derivative(const unsigned int i,
const Point<dim> & p) const
{
- unsigned int indices[dim];
+ std::array<unsigned int, dim> indices;
compute_index(i, indices);
std::vector<std::vector<double>> v(dim, std::vector<double>(order + 1));
compute_tensor_index(const unsigned int,
const unsigned int,
const unsigned int,
- unsigned int (&)[dim])
+ std::array<unsigned int, dim> &)
{
Assert(false, ExcNotImplemented());
}
+ template <>
inline void
- compute_tensor_index(const unsigned int n,
- const unsigned int,
- const unsigned int,
- unsigned int (&indices)[1])
+ compute_tensor_index<1>(const unsigned int n,
+ const unsigned int,
+ const unsigned int,
+ std::array<unsigned int, 1> &indices)
{
indices[0] = n;
}
+ template <>
inline void
- compute_tensor_index(const unsigned int n,
- const unsigned int n_pols_0,
- const unsigned int,
- unsigned int (&indices)[2])
+ compute_tensor_index<2>(const unsigned int n,
+ const unsigned int n_pols_0,
+ const unsigned int,
+ std::array<unsigned int, 2> &indices)
{
indices[0] = n % n_pols_0;
indices[1] = n / n_pols_0;
}
+ template <>
inline void
- compute_tensor_index(const unsigned int n,
- const unsigned int n_pols_0,
- const unsigned int n_pols_1,
- unsigned int (&indices)[3])
+ compute_tensor_index<3>(const unsigned int n,
+ const unsigned int n_pols_0,
+ const unsigned int n_pols_1,
+ std::array<unsigned int, 3> &indices)
{
indices[0] = n % n_pols_0;
indices[1] = (n / n_pols_0) % n_pols_1;
template <int dim, typename PolynomialType>
inline void
TensorProductPolynomials<dim, PolynomialType>::compute_index(
- const unsigned int i,
- unsigned int (&indices)[(dim > 0 ? dim : 1)]) const
+ const unsigned int i,
+ std::array<unsigned int, dim> &indices) const
{
Assert(i < Utilities::fixed_power<dim>(polynomials.size()),
ExcInternalError());
- internal::compute_tensor_index(index_map[i],
- polynomials.size(),
- polynomials.size(),
- indices);
+ internal::compute_tensor_index<dim>(index_map[i],
+ polynomials.size(),
+ polynomials.size(),
+ indices);
}
TensorProductPolynomials<dim, PolynomialType>::output_indices(
std::ostream &out) const
{
- unsigned int ix[(dim > 0) ? dim : 1];
+ std::array<unsigned int, dim> ix;
for (unsigned int i = 0; i < this->n(); ++i)
{
compute_index(i, ix);
{
Assert(dim > 0, ExcNotImplemented());
- unsigned int indices[dim];
+ std::array<unsigned int, dim> indices;
compute_index(i, indices);
double value = 1.;
const unsigned int i,
const Point<dim> & p) const
{
- unsigned int indices[dim];
+ std::array<unsigned int, dim> indices;
compute_index(i, indices);
// compute values and
// uni-directional derivatives at
// the given point in each
// co-ordinate direction
- double v[dim][2];
+ std::array<std::array<double, 2>, dim> v;
{
std::vector<double> tmp(2);
for (unsigned int d = 0; d < dim; ++d)
const unsigned int i,
const Point<dim> & p) const
{
- unsigned int indices[(dim > 0) ? dim : 1];
+ std::array<unsigned int, dim> indices;
compute_index(i, indices);
- double v[dim][3];
+ std::array<std::array<double, 3>, dim> v;
{
std::vector<double> tmp(3);
for (unsigned int d = 0; d < dim; ++d)
template <int dim>
void
AnisotropicPolynomials<dim>::compute_index(
- const unsigned int i,
- unsigned int (&indices)[(dim > 0 ? dim : 1)]) const
+ const unsigned int i,
+ std::array<unsigned int, dim> &indices) const
{
#ifdef DEBUG
unsigned int n_poly = 1;
{
}
else if (dim == 1)
- internal::compute_tensor_index(i,
- polynomials[0].size(),
- 0 /*not used*/,
- indices);
+ internal::compute_tensor_index<dim>(i,
+ polynomials[0].size(),
+ 0 /*not used*/,
+ indices);
else
- internal::compute_tensor_index(i,
- polynomials[0].size(),
- polynomials[1].size(),
- indices);
+ internal::compute_tensor_index<dim>(i,
+ polynomials[0].size(),
+ polynomials[1].size(),
+ indices);
}
AnisotropicPolynomials<dim>::compute_value(const unsigned int i,
const Point<dim> & p) const
{
- unsigned int indices[(dim > 0) ? dim : 1];
+ std::array<unsigned int, dim> indices;
compute_index(i, indices);
double value = 1.;
AnisotropicPolynomials<dim>::compute_grad(const unsigned int i,
const Point<dim> & p) const
{
- unsigned int indices[(dim > 0) ? dim : 1];
+ std::array<unsigned int, dim> indices;
compute_index(i, indices);
// compute values and
AnisotropicPolynomials<dim>::compute_grad_grad(const unsigned int i,
const Point<dim> & p) const
{
- unsigned int indices[(dim > 0) ? dim : 1];
+ std::array<unsigned int, dim> indices;
compute_index(i, indices);
std::vector<std::vector<double>> v(dim, std::vector<double>(3));
// one-dimensional indices of
// this particular tensor
// product polynomial
- unsigned int indices[(dim > 0) ? dim : 1];
+ std::array<unsigned int, dim> indices;
compute_index(i, indices);
if (update_values)