*
* @tparam Number Arithmetic type of the underlying array elements.
*
- * @tparam size Compile-time array lengths. By default at -1, which means that
- * the run-time information stored in the matrices passed to the reinit()
- * function is used.
+ * @tparam n_rows_1d Compile-time number of rows of 1D matrices (only
+ * valid if the number of rows and columns coincide for each
+ * dimension). By default at -1, which means that the number of rows
+ * is determined at run-time by means of the matrices passed to the
+ * reinit() function.
*
* @author Martin Kronbichler and Julius Witte, 2017
*/
-template <int dim, typename Number, int size = -1>
+template <int dim, typename Number, int n_rows_1d = -1>
class TensorProductMatrixSymmetricSumBase
{
public:
using value_type = Number;
/**
- * Compile-time array lengths given by the template parameter size.
+ * The static number of rows of the 1D matrices. For more details,
+ * see the description of the template parameter <tt>n_rows_1d</tt>.
*/
- static constexpr int static_size = size;
+ static constexpr int n_rows_1d_static = n_rows_1d;
/**
* Return the number of rows of the tensor product matrix
*
* Note that this class allows for two modes of usage. The first is a use case
* with run time constants for the matrix dimensions that is achieved by
- * setting the optional template parameter for the size to -1. The second mode
- * of usage that is faster allows to set the template parameter as a compile
- * time constant, giving significantly faster code in particular for small
- * sizes of the matrix.
+ * setting the optional template parameter <tt>n_rows_1d</tt> to -1. The second
+ * mode of usage that is faster allows to set the template parameter as a
+ * compile time constant, giving significantly faster code in particular for
+ * small sizes of the matrix.
*
* @tparam dim Dimension of the problem. Currently, 1D, 2D, and 3D codes are
* implemented.
* to perform LAPACK calculations for each vectorization lane, i.e. for the
* supported float and double numbers.
*
- * @tparam size Compile-time array lengths. By default at -1, which means that
- * the run-time information stored in the matrices passed to the reinit()
- * function is used.
+ * @tparam n_rows_1d Compile-time number of rows of 1D matrices (only
+ * valid if the number of rows and columns coincide for each
+ * dimension). By default at -1, which means that the number of rows
+ * is determined at run-time by means of the matrices passed to the
+ * reinit() function.
*
* @author Martin Kronbichler and Julius Witte, 2017
*/
-template <int dim, typename Number, int size = -1>
+template <int dim, typename Number, int n_rows_1d = -1>
class TensorProductMatrixSymmetricSum
- : public TensorProductMatrixSymmetricSumBase<dim, Number, size>
+ : public TensorProductMatrixSymmetricSumBase<dim, Number, n_rows_1d>
{
public:
/**
*
* @author Martin Kronbichler and Julius Witte, 2017
*/
-template <int dim, typename Number, int size>
-class TensorProductMatrixSymmetricSum<dim, VectorizedArray<Number>, size>
+template <int dim, typename Number, int n_rows_1d>
+class TensorProductMatrixSymmetricSum<dim, VectorizedArray<Number>, n_rows_1d>
: public TensorProductMatrixSymmetricSumBase<dim,
VectorizedArray<Number>,
- size>
+ n_rows_1d>
{
public:
/**
} // namespace internal
-template <int dim, typename Number, int size>
+template <int dim, typename Number, int n_rows_1d>
inline unsigned int
-TensorProductMatrixSymmetricSumBase<dim, Number, size>::m() const
+TensorProductMatrixSymmetricSumBase<dim, Number, n_rows_1d>::m() const
{
unsigned int m = mass_matrix[0].n_rows();
for (unsigned int d = 1; d < dim; ++d)
-template <int dim, typename Number, int size>
+template <int dim, typename Number, int n_rows_1d>
inline unsigned int
-TensorProductMatrixSymmetricSumBase<dim, Number, size>::n() const
+TensorProductMatrixSymmetricSumBase<dim, Number, n_rows_1d>::n() const
{
unsigned int n = mass_matrix[0].n_cols();
for (unsigned int d = 1; d < dim; ++d)
-template <int dim, typename Number, int size>
+template <int dim, typename Number, int n_rows_1d>
inline void
-TensorProductMatrixSymmetricSumBase<dim, Number, size>::vmult(
+TensorProductMatrixSymmetricSumBase<dim, Number, n_rows_1d>::vmult(
const ArrayView<Number> & dst_view,
const ArrayView<const Number> &src_view) const
{
AssertDimension(dst_view.size(), this->m());
AssertDimension(src_view.size(), this->n());
std::lock_guard<std::mutex> lock(this->mutex);
- const unsigned int n =
- Utilities::fixed_power<dim>(size > 0 ? size : eigenvalues[0].size());
+ const unsigned int n = Utilities::fixed_power<dim>(
+ n_rows_1d > 0 ? n_rows_1d : eigenvalues[0].size());
tmp_array.resize_fast(n * 2);
- constexpr int kernel_size = size > 0 ? size : 0;
+ constexpr int kernel_size = n_rows_1d > 0 ? n_rows_1d : 0;
internal::EvaluatorTensorProduct<internal::evaluate_general,
dim,
kernel_size,
-template <int dim, typename Number, int size>
+template <int dim, typename Number, int n_rows_1d>
inline void
-TensorProductMatrixSymmetricSumBase<dim, Number, size>::apply_inverse(
+TensorProductMatrixSymmetricSumBase<dim, Number, n_rows_1d>::apply_inverse(
const ArrayView<Number> & dst_view,
const ArrayView<const Number> &src_view) const
{
AssertDimension(dst_view.size(), this->n());
AssertDimension(src_view.size(), this->m());
std::lock_guard<std::mutex> lock(this->mutex);
- const unsigned int n = size > 0 ? size : eigenvalues[0].size();
+ const unsigned int n = n_rows_1d > 0 ? n_rows_1d : eigenvalues[0].size();
tmp_array.resize_fast(Utilities::fixed_power<dim>(n));
- constexpr int kernel_size = size > 0 ? size : 0;
+ constexpr int kernel_size = n_rows_1d > 0 ? n_rows_1d : 0;
internal::EvaluatorTensorProduct<internal::evaluate_general,
dim,
kernel_size,
//---------------------- TensorProductMatrixSymmetricSum ----------------------
-template <int dim, typename Number, int size>
-inline TensorProductMatrixSymmetricSum<dim, Number, size>::
+template <int dim, typename Number, int n_rows_1d>
+inline TensorProductMatrixSymmetricSum<dim, Number, n_rows_1d>::
TensorProductMatrixSymmetricSum(
const std::array<Table<2, Number>, dim> &mass_matrix,
const std::array<Table<2, Number>, dim> &derivative_matrix)
-template <int dim, typename Number, int size>
-inline TensorProductMatrixSymmetricSum<dim, Number, size>::
+template <int dim, typename Number, int n_rows_1d>
+inline TensorProductMatrixSymmetricSum<dim, Number, n_rows_1d>::
TensorProductMatrixSymmetricSum(
const std::array<FullMatrix<Number>, dim> &mass_matrix,
const std::array<FullMatrix<Number>, dim> &derivative_matrix)
-template <int dim, typename Number, int size>
-inline TensorProductMatrixSymmetricSum<dim, Number, size>::
+template <int dim, typename Number, int n_rows_1d>
+inline TensorProductMatrixSymmetricSum<dim, Number, n_rows_1d>::
TensorProductMatrixSymmetricSum(const Table<2, Number> &mass_matrix,
const Table<2, Number> &derivative_matrix)
{
-template <int dim, typename Number, int size>
+template <int dim, typename Number, int n_rows_1d>
template <typename MatrixArray>
inline void
-TensorProductMatrixSymmetricSum<dim, Number, size>::reinit_impl(
+TensorProductMatrixSymmetricSum<dim, Number, n_rows_1d>::reinit_impl(
MatrixArray &&mass_matrices_,
MatrixArray &&derivative_matrices_)
{
for (int dir = 0; dir < dim; ++dir)
{
- Assert(size == -1 || (size > 0 && static_cast<unsigned int>(size) ==
- mass_matrices[dir].n_rows()),
- ExcDimensionMismatch(size, mass_matrices[dir].n_rows()));
+ Assert(n_rows_1d == -1 ||
+ (n_rows_1d > 0 && static_cast<unsigned int>(n_rows_1d) ==
+ mass_matrices[dir].n_rows()),
+ ExcDimensionMismatch(n_rows_1d, mass_matrices[dir].n_rows()));
AssertDimension(mass_matrices[dir].n_rows(), mass_matrices[dir].n_cols());
AssertDimension(mass_matrices[dir].n_rows(),
derivative_matrices[dir].n_rows());
-template <int dim, typename Number, int size>
+template <int dim, typename Number, int n_rows_1d>
inline void
-TensorProductMatrixSymmetricSum<dim, Number, size>::reinit(
+TensorProductMatrixSymmetricSum<dim, Number, n_rows_1d>::reinit(
const std::array<Table<2, Number>, dim> &mass_matrix,
const std::array<Table<2, Number>, dim> &derivative_matrix)
{
-template <int dim, typename Number, int size>
+template <int dim, typename Number, int n_rows_1d>
inline void
-TensorProductMatrixSymmetricSum<dim, Number, size>::reinit(
+TensorProductMatrixSymmetricSum<dim, Number, n_rows_1d>::reinit(
const std::array<FullMatrix<Number>, dim> &mass_matrix,
const std::array<FullMatrix<Number>, dim> &derivative_matrix)
{
-template <int dim, typename Number, int size>
+template <int dim, typename Number, int n_rows_1d>
inline void
-TensorProductMatrixSymmetricSum<dim, Number, size>::reinit(
+TensorProductMatrixSymmetricSum<dim, Number, n_rows_1d>::reinit(
const Table<2, Number> &mass_matrix,
const Table<2, Number> &derivative_matrix)
{
//------------- vectorized spec.: TensorProductMatrixSymmetricSum -------------
-template <int dim, typename Number, int size>
-inline TensorProductMatrixSymmetricSum<dim, VectorizedArray<Number>, size>::
+template <int dim, typename Number, int n_rows_1d>
+inline TensorProductMatrixSymmetricSum<dim,
+ VectorizedArray<Number>,
+ n_rows_1d>::
TensorProductMatrixSymmetricSum(
const std::array<Table<2, VectorizedArray<Number>>, dim> &mass_matrix,
const std::array<Table<2, VectorizedArray<Number>>, dim> &derivative_matrix)
-template <int dim, typename Number, int size>
-inline TensorProductMatrixSymmetricSum<dim, VectorizedArray<Number>, size>::
+template <int dim, typename Number, int n_rows_1d>
+inline TensorProductMatrixSymmetricSum<dim,
+ VectorizedArray<Number>,
+ n_rows_1d>::
TensorProductMatrixSymmetricSum(
const Table<2, VectorizedArray<Number>> &mass_matrix,
const Table<2, VectorizedArray<Number>> &derivative_matrix)
-template <int dim, typename Number, int size>
+template <int dim, typename Number, int n_rows_1d>
template <typename MatrixArray>
inline void
-TensorProductMatrixSymmetricSum<dim, VectorizedArray<Number>, size>::
+TensorProductMatrixSymmetricSum<dim, VectorizedArray<Number>, n_rows_1d>::
reinit_impl(MatrixArray &&mass_matrices_, MatrixArray &&derivative_matrices_)
{
auto &&mass_matrix = std::forward<MatrixArray>(mass_matrices_);
this->derivative_matrix = derivative_matrix;
constexpr unsigned int macro_size = VectorizedArray<Number>::n_array_elements;
- std::size_t n_rows_max = (size > 0) ? size : 0;
- if (size == -1)
+ std::size_t n_rows_max = (n_rows_1d > 0) ? n_rows_1d : 0;
+ if (n_rows_1d == -1)
for (unsigned int d = 0; d < dim; ++d)
n_rows_max = std::max(n_rows_max, mass_matrix[d].n_rows());
const std::size_t nm_flat_size_max = n_rows_max * n_rows_max * macro_size;
std::array<unsigned int, macro_size> offsets_n;
for (int dir = 0; dir < dim; ++dir)
{
- Assert(size == -1 || (size > 0 && static_cast<unsigned int>(size) ==
- mass_matrix[dir].n_rows()),
- ExcDimensionMismatch(size, mass_matrix[dir].n_rows()));
+ Assert(n_rows_1d == -1 ||
+ (n_rows_1d > 0 && static_cast<unsigned int>(n_rows_1d) ==
+ mass_matrix[dir].n_rows()),
+ ExcDimensionMismatch(n_rows_1d, mass_matrix[dir].n_rows()));
AssertDimension(mass_matrix[dir].n_rows(), mass_matrix[dir].n_cols());
AssertDimension(mass_matrix[dir].n_rows(),
derivative_matrix[dir].n_rows());
-template <int dim, typename Number, int size>
+template <int dim, typename Number, int n_rows_1d>
inline void
-TensorProductMatrixSymmetricSum<dim, VectorizedArray<Number>, size>::reinit(
- const std::array<Table<2, VectorizedArray<Number>>, dim> &mass_matrix,
- const std::array<Table<2, VectorizedArray<Number>>, dim> &derivative_matrix)
+TensorProductMatrixSymmetricSum<dim, VectorizedArray<Number>, n_rows_1d>::
+ reinit(
+ const std::array<Table<2, VectorizedArray<Number>>, dim> &mass_matrix,
+ const std::array<Table<2, VectorizedArray<Number>>, dim> &derivative_matrix)
{
reinit_impl(mass_matrix, derivative_matrix);
}
-template <int dim, typename Number, int size>
+template <int dim, typename Number, int n_rows_1d>
inline void
-TensorProductMatrixSymmetricSum<dim, VectorizedArray<Number>, size>::reinit(
- const Table<2, VectorizedArray<Number>> &mass_matrix,
- const Table<2, VectorizedArray<Number>> &derivative_matrix)
+TensorProductMatrixSymmetricSum<dim, VectorizedArray<Number>, n_rows_1d>::
+ reinit(const Table<2, VectorizedArray<Number>> &mass_matrix,
+ const Table<2, VectorizedArray<Number>> &derivative_matrix)
{
std::array<Table<2, VectorizedArray<Number>>, dim> mass_matrices;
std::array<Table<2, VectorizedArray<Number>>, dim> derivative_matrices;