ArrayView offers the same functionality but in a more general way.
This PR reimplements VectorSlice as a class inheriting from ArrayView,
which automatically gives us all of the right conversions.
--- /dev/null
+Changed: VectorSlice has been deprecated in favor of the more general ArrayView
+class. All public interfaces now take ArrayView arguments instead of VectorSlice
+arguments. Since VectorSlice now inherits from ArrayView this should be
+compatible with the overwhelming majority of use cases.
+<br>
+(David Wells, 2019/01/21)
DEAL_II_NAMESPACE_OPEN
-
/**
* Filter a range out of any object having a random access <tt>operator[]
* (unsigned int)</tt> and a function <tt>size() const</tt>.
*
* @ingroup data
* @author Guido Kanschat, 2004
+ *
+ * @deprecated Use the more general ArrayView class instead.
*/
template <typename VectorType>
-class VectorSlice
+class DEAL_II_DEPRECATED VectorSlice
+ : public ArrayView<
+ typename std::conditional<std::is_const<VectorType>::value,
+ const typename VectorType::value_type,
+ typename VectorType::value_type>::type>
{
public:
/**
*/
VectorSlice(VectorType &v, unsigned int start, unsigned int length);
+protected:
/**
- * Conversion operator to an ArrayView object that represents an array of
- * non-const elements pointing to the same location as the current object.
- */
- operator ArrayView<typename VectorType::value_type *>();
-
- /**
- * Conversion operator to an ArrayView object that represents an array of
- * const elements pointing to the same location as the current object.
- */
- operator ArrayView<const typename VectorType::value_type *>() const;
-
- /**
- * Return the length of the slice using the same interface as
- * <tt>std::vector</tt>.
- */
- unsigned int
- size() const;
-
- /**
- * Return a reference to the $i$th element of the range represented by the
- * current object.
+ * Alias for the base class name.
*/
- typename VectorType::reference operator[](unsigned int i);
-
- /**
- * Return a @p const reference to the $i$th element of the range represented
- * by the current object.
- */
- typename VectorType::const_reference operator[](unsigned int i) const;
-
- /**
- * Standard-conforming iterator function.
- */
- typename VectorType::iterator
- begin();
-
- /**
- * Standard-conforming iterator function.
- */
- typename VectorType::const_iterator
- begin() const;
-
- /**
- * Standard-conforming iterator function.
- */
- typename VectorType::iterator
- end();
-
- /**
- * Standard-conforming iterator function.
- */
- typename VectorType::const_iterator
- end() const;
-
-private:
- /**
- * The vector we extract from.
- */
- VectorType &v;
- /**
- * The start index of the slice.
- */
- const unsigned int start;
- /**
- * The length of the slice.
- */
- const unsigned int length;
+ using ArrayViewType =
+ ArrayView<typename std::conditional<std::is_const<VectorType>::value,
+ const typename VectorType::value_type,
+ typename VectorType::value_type>::type>;
};
template <typename VectorType>
inline VectorSlice<VectorType>::VectorSlice(VectorType &v)
- : v(v)
- , start(0)
- , length(v.size())
+ : ArrayViewType(make_array_view(std::begin(v), std::end(v)))
{}
+
template <typename VectorType>
inline VectorSlice<VectorType>::VectorSlice(VectorType & v,
unsigned int start,
unsigned int length)
- : v(v)
- , start(start)
- , length(length)
+ : ArrayViewType(
+ make_array_view(std::begin(v) + start, std::begin(v) + start + length))
{
Assert((start + length <= v.size()),
ExcIndexRange(length, 0, v.size() - start + 1));
}
-
-template <typename VectorType>
-inline unsigned int
-VectorSlice<VectorType>::size() const
-{
- return length;
-}
-
-
-template <typename VectorType>
-VectorSlice<VectorType>::operator ArrayView<typename VectorType::value_type *>()
-{
- return ArrayView<typename VectorType::value_type *>(&v[start], length);
-}
-
-
-template <typename VectorType>
-VectorSlice<VectorType>::
-operator ArrayView<const typename VectorType::value_type *>() const
-{
- return ArrayView<const typename VectorType::value_type *>(&v[start], length);
-}
-
-
-template <typename VectorType>
-inline typename VectorType::reference VectorSlice<VectorType>::
- operator[](unsigned int i)
-{
- Assert((i < length), ExcIndexRange(i, 0, length));
-
- return v[start + i];
-}
-
-
-template <typename VectorType>
-inline typename VectorType::const_reference VectorSlice<VectorType>::
- operator[](unsigned int i) const
-{
- Assert((i < length), ExcIndexRange(i, 0, length));
-
- return v[start + i];
-}
-
-
-template <typename VectorType>
-inline typename VectorType::const_iterator
-VectorSlice<VectorType>::begin() const
-{
- return v.begin() + start;
-}
-
-
-template <typename VectorType>
-inline typename VectorType::iterator
-VectorSlice<VectorType>::begin()
-{
- return v.begin() + start;
-}
-
-
-template <typename VectorType>
-inline typename VectorType::const_iterator
-VectorSlice<VectorType>::end() const
-{
- return v.begin() + start + length;
-}
-
-
-template <typename VectorType>
-inline typename VectorType::iterator
-VectorSlice<VectorType>::end()
-{
- return v.begin() + start + length;
-}
-
DEAL_II_NAMESPACE_CLOSE
#endif
template <class InputVector>
void
get_function_values(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
- std::vector<typename InputVector::value_type> &values) const;
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> &indices,
+ std::vector<typename InputVector::value_type> & values) const;
/**
* Generate vector function values from an arbitrary vector.
template <class InputVector>
void
get_function_values(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> & indices,
std::vector<Vector<typename InputVector::value_type>> &values) const;
template <class InputVector>
void
get_function_values(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
- VectorSlice<std::vector<std::vector<typename InputVector::value_type>>>
- values,
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> & indices,
+ ArrayView<std::vector<typename InputVector::value_type>> values,
const bool quadrature_points_fastest) const;
//@}
template <class InputVector>
void
get_function_gradients(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> &indices,
std::vector<Tensor<1, spacedim, typename InputVector::value_type>>
&gradients) const;
template <class InputVector>
void
get_function_gradients(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
- VectorSlice<std::vector<
- std::vector<Tensor<1, spacedim, typename InputVector::value_type>>>>
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> &indices,
+ ArrayView<
+ std::vector<Tensor<1, spacedim, typename InputVector::value_type>>>
gradients,
bool quadrature_points_fastest = false) const;
template <class InputVector>
void
get_function_hessians(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> &indices,
std::vector<Tensor<2, spacedim, typename InputVector::value_type>>
&hessians) const;
template <class InputVector>
void
get_function_hessians(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
- VectorSlice<std::vector<
- std::vector<Tensor<2, spacedim, typename InputVector::value_type>>>>
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> &indices,
+ ArrayView<
+ std::vector<Tensor<2, spacedim, typename InputVector::value_type>>>
hessians,
bool quadrature_points_fastest = false) const;
template <class InputVector>
void
get_function_laplacians(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
- std::vector<typename InputVector::value_type> &laplacians) const;
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> &indices,
+ std::vector<typename InputVector::value_type> & laplacians) const;
/**
* Access to the second derivatives of a function with more flexibility. See
template <class InputVector>
void
get_function_laplacians(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> & indices,
std::vector<Vector<typename InputVector::value_type>> &laplacians) const;
/**
template <class InputVector>
void
get_function_laplacians(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
- std::vector<std::vector<typename InputVector::value_type>> & laplacians,
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> & indices,
+ std::vector<std::vector<typename InputVector::value_type>> &laplacians,
bool quadrature_points_fastest = false) const;
//@}
template <class InputVector>
void
get_function_third_derivatives(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> &indices,
std::vector<Tensor<3, spacedim, typename InputVector::value_type>>
&third_derivatives) const;
template <class InputVector>
void
get_function_third_derivatives(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
- VectorSlice<std::vector<
- std::vector<Tensor<3, spacedim, typename InputVector::value_type>>>>
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> &indices,
+ ArrayView<
+ std::vector<Tensor<3, spacedim, typename InputVector::value_type>>>
third_derivatives,
bool quadrature_points_fastest = false) const;
//@}
*/
template <int dim>
void
- cell_matrix(
- FullMatrix<double> & M,
- const FEValuesBase<dim> & fe,
- const FEValuesBase<dim> & fetest,
- const VectorSlice<const std::vector<std::vector<double>>> &velocity,
- const double factor = 1.)
+ cell_matrix(FullMatrix<double> & M,
+ const FEValuesBase<dim> & fe,
+ const FEValuesBase<dim> & fetest,
+ const ArrayView<const std::vector<double>> &velocity,
+ const double factor = 1.)
{
const unsigned int n_dofs = fe.dofs_per_cell;
const unsigned int t_dofs = fetest.dofs_per_cell;
*/
template <int dim>
inline void
- cell_residual(
- Vector<double> & result,
- const FEValuesBase<dim> & fe,
- const std::vector<Tensor<1, dim>> & input,
- const VectorSlice<const std::vector<std::vector<double>>> &velocity,
- double factor = 1.)
+ cell_residual(Vector<double> & result,
+ const FEValuesBase<dim> & fe,
+ const std::vector<Tensor<1, dim>> & input,
+ const ArrayView<const std::vector<double>> &velocity,
+ double factor = 1.)
{
const unsigned int nq = fe.n_quadrature_points;
const unsigned int n_dofs = fe.dofs_per_cell;
*/
template <int dim>
inline void
- cell_residual(
- Vector<double> & result,
- const FEValuesBase<dim> & fe,
- const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &input,
- const VectorSlice<const std::vector<std::vector<double>>> &velocity,
- double factor = 1.)
+ cell_residual(Vector<double> & result,
+ const FEValuesBase<dim> & fe,
+ const ArrayView<const std::vector<Tensor<1, dim>>> &input,
+ const ArrayView<const std::vector<double>> & velocity,
+ double factor = 1.)
{
const unsigned int nq = fe.n_quadrature_points;
const unsigned int n_dofs = fe.dofs_per_cell;
*/
template <int dim>
inline void
- cell_residual(
- Vector<double> & result,
- const FEValuesBase<dim> & fe,
- const std::vector<double> & input,
- const VectorSlice<const std::vector<std::vector<double>>> &velocity,
- double factor = 1.)
+ cell_residual(Vector<double> & result,
+ const FEValuesBase<dim> & fe,
+ const std::vector<double> & input,
+ const ArrayView<const std::vector<double>> &velocity,
+ double factor = 1.)
{
const unsigned int nq = fe.n_quadrature_points;
const unsigned int n_dofs = fe.dofs_per_cell;
*/
template <int dim>
inline void
- cell_residual(
- Vector<double> & result,
- const FEValuesBase<dim> & fe,
- const VectorSlice<const std::vector<std::vector<double>>> &input,
- const VectorSlice<const std::vector<std::vector<double>>> &velocity,
- double factor = 1.)
+ cell_residual(Vector<double> & result,
+ const FEValuesBase<dim> & fe,
+ const ArrayView<const std::vector<double>> &input,
+ const ArrayView<const std::vector<double>> &velocity,
+ double factor = 1.)
{
const unsigned int nq = fe.n_quadrature_points;
const unsigned int n_dofs = fe.dofs_per_cell;
* u_i v_j \, ds
* @f]
*
- * The <tt>velocity</tt> is provided as a VectorSlice, having <tt>dim</tt>
+ * The <tt>velocity</tt> is provided as an ArrayView, having <tt>dim</tt>
* vectors, one for each velocity component. Each of the vectors must
* either have only a single entry, if the advection velocity is constant,
* or have an entry for each quadrature point.
*/
template <int dim>
void
- upwind_value_matrix(
- FullMatrix<double> & M,
- const FEValuesBase<dim> & fe,
- const FEValuesBase<dim> & fetest,
- const VectorSlice<const std::vector<std::vector<double>>> &velocity,
- double factor = 1.)
+ upwind_value_matrix(FullMatrix<double> & M,
+ const FEValuesBase<dim> & fe,
+ const FEValuesBase<dim> & fetest,
+ const ArrayView<const std::vector<double>> &velocity,
+ double factor = 1.)
{
const unsigned int n_dofs = fe.dofs_per_cell;
const unsigned int t_dofs = fetest.dofs_per_cell;
* argument `input` on the outflow boundary. On the inflow boundary, it is
* the inhomogenous boundary value in the argument `data`.
*
- * The <tt>velocity</tt> is provided as a VectorSlice, having <tt>dim</tt>
+ * The <tt>velocity</tt> is provided as an ArrayView, having <tt>dim</tt>
* vectors, one for each velocity component. Each of the vectors must
* either have only a single entry, if the advection velocity is constant,
* or have an entry for each quadrature point.
*/
template <int dim>
inline void
- upwind_value_residual(
- Vector<double> & result,
- const FEValuesBase<dim> & fe,
- const std::vector<double> & input,
- const std::vector<double> & data,
- const VectorSlice<const std::vector<std::vector<double>>> &velocity,
- double factor = 1.)
+ upwind_value_residual(Vector<double> & result,
+ const FEValuesBase<dim> & fe,
+ const std::vector<double> & input,
+ const std::vector<double> & data,
+ const ArrayView<const std::vector<double>> &velocity,
+ double factor = 1.)
{
const unsigned int n_dofs = fe.dofs_per_cell;
* argument `input` on the outflow boundary. On the inflow boundary, it is
* the inhomogenous boundary value in the argument `data`.
*
- * The <tt>velocity</tt> is provided as a VectorSlice, having <tt>dim</tt>
+ * The <tt>velocity</tt> is provided as an ArrayView, having <tt>dim</tt>
* vectors, one for each velocity component. Each of the vectors must
* either have only a single entry, if the advection velocity is constant,
* or have an entry for each quadrature point.
*/
template <int dim>
inline void
- upwind_value_residual(
- Vector<double> & result,
- const FEValuesBase<dim> & fe,
- const VectorSlice<const std::vector<std::vector<double>>> &input,
- const VectorSlice<const std::vector<std::vector<double>>> &data,
- const VectorSlice<const std::vector<std::vector<double>>> &velocity,
- double factor = 1.)
+ upwind_value_residual(Vector<double> & result,
+ const FEValuesBase<dim> & fe,
+ const ArrayView<const std::vector<double>> &input,
+ const ArrayView<const std::vector<double>> &data,
+ const ArrayView<const std::vector<double>> &velocity,
+ double factor = 1.)
{
const unsigned int n_dofs = fe.dofs_per_cell;
const unsigned int n_comp = fe.get_fe().n_components();
* \,ds
* @f]
*
- * The <tt>velocity</tt> is provided as a VectorSlice, having <tt>dim</tt>
+ * The <tt>velocity</tt> is provided as an ArrayView, having <tt>dim</tt>
* vectors, one for each velocity component. Each of the vectors must
* either have only a single entry, if the advection velocity is constant,
* or have an entry for each quadrature point.
*/
template <int dim>
void
- upwind_value_matrix(
- FullMatrix<double> & M11,
- FullMatrix<double> & M12,
- FullMatrix<double> & M21,
- FullMatrix<double> & M22,
- const FEValuesBase<dim> & fe1,
- const FEValuesBase<dim> & fe2,
- const FEValuesBase<dim> & fetest1,
- const FEValuesBase<dim> & fetest2,
- const VectorSlice<const std::vector<std::vector<double>>> &velocity,
- const double factor = 1.)
+ upwind_value_matrix(FullMatrix<double> & M11,
+ FullMatrix<double> & M12,
+ FullMatrix<double> & M21,
+ FullMatrix<double> & M22,
+ const FEValuesBase<dim> & fe1,
+ const FEValuesBase<dim> & fe2,
+ const FEValuesBase<dim> & fetest1,
+ const FEValuesBase<dim> & fetest2,
+ const ArrayView<const std::vector<double>> &velocity,
+ const double factor = 1.)
{
const unsigned int n1 = fe1.dofs_per_cell;
// Multiply the quadrature point
* \,ds
* @f]
*
- * The <tt>velocity</tt> is provided as a VectorSlice, having <tt>dim</tt>
+ * The <tt>velocity</tt> is provided as an ArrayView, having <tt>dim</tt>
* vectors, one for each velocity component. Each of the vectors must
* either have only a single entry, if the advection velocity is constant,
* or have an entry for each quadrature point.
*/
template <int dim>
void
- upwind_face_residual(
- Vector<double> & result1,
- Vector<double> & result2,
- const FEValuesBase<dim> & fe1,
- const FEValuesBase<dim> & fe2,
- const std::vector<double> & input1,
- const std::vector<double> & input2,
- const VectorSlice<const std::vector<std::vector<double>>> &velocity,
- const double factor = 1.)
+ upwind_face_residual(Vector<double> & result1,
+ Vector<double> & result2,
+ const FEValuesBase<dim> & fe1,
+ const FEValuesBase<dim> & fe2,
+ const std::vector<double> & input1,
+ const std::vector<double> & input2,
+ const ArrayView<const std::vector<double>> &velocity,
+ const double factor = 1.)
{
Assert(fe1.get_fe().n_components() == 1,
ExcDimensionMismatch(fe1.get_fe().n_components(), 1));
* \,ds
* @f]
*
- * The <tt>velocity</tt> is provided as a VectorSlice, having <tt>dim</tt>
+ * The <tt>velocity</tt> is provided as an ArrayView, having <tt>dim</tt>
* vectors, one for each velocity component. Each of the vectors must
* either have only a single entry, if the advection velocity is constant,
* or have an entry for each quadrature point.
*/
template <int dim>
void
- upwind_face_residual(
- Vector<double> & result1,
- Vector<double> & result2,
- const FEValuesBase<dim> & fe1,
- const FEValuesBase<dim> & fe2,
- const VectorSlice<const std::vector<std::vector<double>>> &input1,
- const VectorSlice<const std::vector<std::vector<double>>> &input2,
- const VectorSlice<const std::vector<std::vector<double>>> &velocity,
- const double factor = 1.)
+ upwind_face_residual(Vector<double> & result1,
+ Vector<double> & result2,
+ const FEValuesBase<dim> & fe1,
+ const FEValuesBase<dim> & fe2,
+ const ArrayView<const std::vector<double>> &input1,
+ const ArrayView<const std::vector<double>> &input2,
+ const ArrayView<const std::vector<double>> &velocity,
+ const double factor = 1.)
{
const unsigned int n_comp = fe1.get_fe().n_components();
const unsigned int n1 = fe1.dofs_per_cell;
*/
template <int dim, typename number>
void
- cell_residual(
- Vector<number> & result,
- const FEValuesBase<dim> & fetest,
- const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &input,
- const double factor = 1.)
+ cell_residual(Vector<number> & result,
+ const FEValuesBase<dim> & fetest,
+ const ArrayView<const std::vector<Tensor<1, dim>>> &input,
+ const double factor = 1.)
{
AssertDimension(fetest.get_fe().n_components(), 1);
AssertVectorVectorDimension(input, dim, fetest.n_quadrature_points);
*/
template <int dim, typename number>
void
- cell_residual(
- Vector<number> & result,
- const FEValuesBase<dim> & fetest,
- const VectorSlice<const std::vector<std::vector<double>>> &input,
- const double factor = 1.)
+ cell_residual(Vector<number> & result,
+ const FEValuesBase<dim> & fetest,
+ const ArrayView<const std::vector<double>> &input,
+ const double factor = 1.)
{
AssertDimension(fetest.get_fe().n_components(), 1);
AssertVectorVectorDimension(input, dim, fetest.n_quadrature_points);
*/
template <int dim, typename number>
void
- u_dot_n_residual(
- Vector<number> & result,
- const FEValuesBase<dim> & fe,
- const FEValuesBase<dim> & fetest,
- const VectorSlice<const std::vector<std::vector<double>>> &data,
- double factor = 1.)
+ u_dot_n_residual(Vector<number> & result,
+ const FEValuesBase<dim> & fe,
+ const FEValuesBase<dim> & fetest,
+ const ArrayView<const std::vector<double>> &data,
+ double factor = 1.)
{
const unsigned int t_dofs = fetest.dofs_per_cell;
*/
template <int dim, typename number>
DEAL_II_DEPRECATED void
- grad_div_residual(
- Vector<number> & result,
- const FEValuesBase<dim> & fetest,
- const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &input,
- const double factor = 1.);
+ grad_div_residual(Vector<number> & result,
+ const FEValuesBase<dim> &fetest,
+ const ArrayView<const std::vector<Tensor<1, dim>>> &input,
+ const double factor = 1.);
template <int dim, typename number>
void
- grad_div_residual(
- Vector<number> & result,
- const FEValuesBase<dim> & fetest,
- const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &input,
- const double factor)
+ grad_div_residual(Vector<number> & result,
+ const FEValuesBase<dim> &fetest,
+ const ArrayView<const std::vector<Tensor<1, dim>>> &input,
+ const double factor)
{
GradDiv::cell_residual(result, fetest, input, factor);
}
*/
template <int dim>
double
- norm(const FEValuesBase<dim> & fe,
- const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &Du)
+ norm(const FEValuesBase<dim> & fe,
+ const ArrayView<const std::vector<Tensor<1, dim>>> &Du)
{
AssertDimension(fe.get_fe().n_components(), dim);
AssertVectorVectorDimension(Du, dim, fe.n_quadrature_points);
*/
template <int dim, typename number>
inline void
- cell_residual(
- Vector<number> & result,
- const FEValuesBase<dim> & fe,
- const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &input,
- double factor = 1.)
+ cell_residual(Vector<number> & result,
+ const FEValuesBase<dim> & fe,
+ const ArrayView<const std::vector<Tensor<1, dim>>> &input,
+ double factor = 1.)
{
const unsigned int nq = fe.n_quadrature_points;
const unsigned int n_dofs = fe.dofs_per_cell;
*/
template <int dim, typename number>
void
- nitsche_residual(
- Vector<number> & result,
- const FEValuesBase<dim> & fe,
- const VectorSlice<const std::vector<std::vector<double>>> & input,
- const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &Dinput,
- const VectorSlice<const std::vector<std::vector<double>>> & data,
- double penalty,
- double factor = 1.)
+ nitsche_residual(Vector<number> & result,
+ const FEValuesBase<dim> & fe,
+ const ArrayView<const std::vector<double>> & input,
+ const ArrayView<const std::vector<Tensor<1, dim>>> &Dinput,
+ const ArrayView<const std::vector<double>> & data,
+ double penalty,
+ double factor = 1.)
{
const unsigned int n_dofs = fe.dofs_per_cell;
AssertVectorVectorDimension(input, dim, fe.n_quadrature_points);
template <int dim, typename number>
inline void
nitsche_tangential_residual(
- Vector<number> & result,
- const FEValuesBase<dim> & fe,
- const VectorSlice<const std::vector<std::vector<double>>> & input,
- const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &Dinput,
- const VectorSlice<const std::vector<std::vector<double>>> & data,
- double penalty,
- double factor = 1.)
+ Vector<number> & result,
+ const FEValuesBase<dim> & fe,
+ const ArrayView<const std::vector<double>> & input,
+ const ArrayView<const std::vector<Tensor<1, dim>>> &Dinput,
+ const ArrayView<const std::vector<double>> & data,
+ double penalty,
+ double factor = 1.)
{
const unsigned int n_dofs = fe.dofs_per_cell;
AssertVectorVectorDimension(input, dim, fe.n_quadrature_points);
template <int dim, typename number>
void
nitsche_residual_homogeneous(
- Vector<number> & result,
- const FEValuesBase<dim> & fe,
- const VectorSlice<const std::vector<std::vector<double>>> & input,
- const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &Dinput,
- double penalty,
- double factor = 1.)
+ Vector<number> & result,
+ const FEValuesBase<dim> & fe,
+ const ArrayView<const std::vector<double>> & input,
+ const ArrayView<const std::vector<Tensor<1, dim>>> &Dinput,
+ double penalty,
+ double factor = 1.)
{
const unsigned int n_dofs = fe.dofs_per_cell;
AssertVectorVectorDimension(input, dim, fe.n_quadrature_points);
*/
template <int dim, typename number>
void
- ip_residual(
- Vector<number> & result1,
- Vector<number> & result2,
- const FEValuesBase<dim> & fe1,
- const FEValuesBase<dim> & fe2,
- const VectorSlice<const std::vector<std::vector<double>>> &input1,
- const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>>
- & Dinput1,
- const VectorSlice<const std::vector<std::vector<double>>> &input2,
- const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>>
- & Dinput2,
- double pen,
- double int_factor = 1.,
- double ext_factor = -1.)
+ ip_residual(Vector<number> & result1,
+ Vector<number> & result2,
+ const FEValuesBase<dim> & fe1,
+ const FEValuesBase<dim> & fe2,
+ const ArrayView<const std::vector<double>> & input1,
+ const ArrayView<const std::vector<Tensor<1, dim>>> &Dinput1,
+ const ArrayView<const std::vector<double>> & input2,
+ const ArrayView<const std::vector<Tensor<1, dim>>> &Dinput2,
+ double pen,
+ double int_factor = 1.,
+ double ext_factor = -1.)
{
const unsigned int n1 = fe1.dofs_per_cell;
*/
template <int dim, typename number>
void
- cell_residual(
- Vector<number> & result,
- const FEValuesBase<dim> & fetest,
- const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &input,
- const double factor = 1.)
+ cell_residual(Vector<number> & result,
+ const FEValuesBase<dim> & fetest,
+ const ArrayView<const std::vector<Tensor<1, dim>>> &input,
+ const double factor = 1.)
{
const unsigned int n_dofs = fetest.dofs_per_cell;
*/
template <int dim>
void
- nitsche_residual(
- Vector<double> & result,
- const FEValuesBase<dim> & fe,
- const VectorSlice<const std::vector<std::vector<double>>> & input,
- const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &Dinput,
- const VectorSlice<const std::vector<std::vector<double>>> & data,
- double penalty,
- double factor = 1.)
+ nitsche_residual(Vector<double> & result,
+ const FEValuesBase<dim> & fe,
+ const ArrayView<const std::vector<double>> & input,
+ const ArrayView<const std::vector<Tensor<1, dim>>> &Dinput,
+ const ArrayView<const std::vector<double>> & data,
+ double penalty,
+ double factor = 1.)
{
const unsigned int n_dofs = fe.dofs_per_cell;
AssertDimension(fe.get_fe().n_components(), dim)
*/
template <int dim>
void
- ip_residual(
- Vector<double> & result1,
- Vector<double> & result2,
- const FEValuesBase<dim> & fe1,
- const FEValuesBase<dim> & fe2,
- const VectorSlice<const std::vector<std::vector<double>>> &input1,
- const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>>
- & Dinput1,
- const VectorSlice<const std::vector<std::vector<double>>> &input2,
- const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>>
- & Dinput2,
- double pen,
- double int_factor = 1.,
- double ext_factor = -1.)
+ ip_residual(Vector<double> & result1,
+ Vector<double> & result2,
+ const FEValuesBase<dim> & fe1,
+ const FEValuesBase<dim> & fe2,
+ const ArrayView<const std::vector<double>> & input1,
+ const ArrayView<const std::vector<Tensor<1, dim>>> &Dinput1,
+ const ArrayView<const std::vector<double>> & input2,
+ const ArrayView<const std::vector<Tensor<1, dim>>> &Dinput2,
+ double pen,
+ double int_factor = 1.,
+ double ext_factor = -1.)
{
const unsigned int n1 = fe1.dofs_per_cell;
*/
template <int dim, typename number>
void
- L2(Vector<number> & result,
- const FEValuesBase<dim> & fe,
- const VectorSlice<const std::vector<std::vector<double>>> &input,
- const double factor = 1.)
+ L2(Vector<number> & result,
+ const FEValuesBase<dim> & fe,
+ const ArrayView<const std::vector<double>> &input,
+ const double factor = 1.)
{
const unsigned int n_dofs = fe.dofs_per_cell;
const unsigned int n_components = input.size();
*/
template <int dim>
inline void
- cell_residual(
- Vector<double> & result,
- const FEValuesBase<dim> & fe,
- const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &input,
- double factor = 1.)
+ cell_residual(Vector<double> & result,
+ const FEValuesBase<dim> & fe,
+ const ArrayView<const std::vector<Tensor<1, dim>>> &input,
+ double factor = 1.)
{
const unsigned int nq = fe.n_quadrature_points;
const unsigned int n_dofs = fe.dofs_per_cell;
*/
template <int dim>
void
- nitsche_residual(
- Vector<double> & result,
- const FEValuesBase<dim> & fe,
- const VectorSlice<const std::vector<std::vector<double>>> & input,
- const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>> &Dinput,
- const VectorSlice<const std::vector<std::vector<double>>> & data,
- double penalty,
- double factor = 1.)
+ nitsche_residual(Vector<double> & result,
+ const FEValuesBase<dim> & fe,
+ const ArrayView<const std::vector<double>> & input,
+ const ArrayView<const std::vector<Tensor<1, dim>>> &Dinput,
+ const ArrayView<const std::vector<double>> & data,
+ double penalty,
+ double factor = 1.)
{
const unsigned int n_dofs = fe.dofs_per_cell;
const unsigned int n_comp = fe.get_fe().n_components();
*/
template <int dim>
void
- ip_residual(
- Vector<double> & result1,
- Vector<double> & result2,
- const FEValuesBase<dim> & fe1,
- const FEValuesBase<dim> & fe2,
- const VectorSlice<const std::vector<std::vector<double>>> &input1,
- const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>>
- & Dinput1,
- const VectorSlice<const std::vector<std::vector<double>>> &input2,
- const VectorSlice<const std::vector<std::vector<Tensor<1, dim>>>>
- & Dinput2,
- double pen,
- double int_factor = 1.,
- double ext_factor = -1.)
+ ip_residual(Vector<double> & result1,
+ Vector<double> & result2,
+ const FEValuesBase<dim> & fe1,
+ const FEValuesBase<dim> & fe2,
+ const ArrayView<const std::vector<double>> & input1,
+ const ArrayView<const std::vector<Tensor<1, dim>>> &Dinput1,
+ const ArrayView<const std::vector<double>> & input2,
+ const ArrayView<const std::vector<Tensor<1, dim>>> &Dinput2,
+ double pen,
+ double int_factor = 1.,
+ double ext_factor = -1.)
{
const unsigned int n_comp = fe1.get_fe().n_components();
const unsigned int n1 = fe1.dofs_per_cell;
namespace Patches
{
template <int dim>
- inline void points_and_values(
- Table<2, double> & result,
- const FEValuesBase<dim> & fe,
- const VectorSlice<const std::vector<std::vector<double>>> &input)
+ inline void
+ points_and_values(Table<2, double> & result,
+ const FEValuesBase<dim> & fe,
+ const ArrayView<const std::vector<double>> &input)
{
const unsigned int n_comp = fe.get_fe().n_components();
AssertVectorVectorDimension(input, n_comp, fe.n_quadrature_points);
const size_type chunk_size);
/**
- * Same as above, but with a VectorSlice argument instead.
+ * Same as above, but with an ArrayView argument instead.
*/
void
- reinit(const size_type m,
- const size_type n,
- const VectorSlice<const std::vector<size_type>> &row_lengths,
- const size_type chunk_size);
+ reinit(const size_type m,
+ const size_type n,
+ const ArrayView<const size_type> &row_lengths,
+ const size_type chunk_size);
/**
* This function compresses the sparsity structure that this object
#include <deal.II/base/config.h>
+#include <deal.II/base/array_view.h>
#include <deal.II/base/exceptions.h>
#include <deal.II/base/linear_index_iterator.h>
#include <deal.II/base/std_cxx14/memory.h>
class SparseLUDecomposition;
template <typename number>
class SparseILU;
-template <typename VectorType>
-class VectorSlice;
namespace ChunkSparsityPatternIterators
{
/**
- * Same as above, but with a VectorSlice argument instead.
+ * Same as above, but with an ArrayView argument instead.
*/
void
- reinit(const size_type m,
- const size_type n,
- const VectorSlice<const std::vector<unsigned int>> &row_lengths);
+ reinit(const size_type m,
+ const size_type n,
+ const ArrayView<const unsigned int> &row_lengths);
/**
* This function compresses the sparsity structure that this object
template <class InputVector>
void
FEValuesBase<dim, spacedim>::get_function_values(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
- std::vector<typename InputVector::value_type> & values) const
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> &indices,
+ std::vector<typename InputVector::value_type> & values) const
{
using Number = typename InputVector::value_type;
Assert(this->update_flags & update_values,
template <class InputVector>
void
FEValuesBase<dim, spacedim>::get_function_values(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
- std::vector<Vector<typename InputVector::value_type>> & values) const
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> & indices,
+ std::vector<Vector<typename InputVector::value_type>> &values) const
{
using Number = typename InputVector::value_type;
// Size of indices must be a multiple of dofs_per_cell such that an integer
template <class InputVector>
void
FEValuesBase<dim, spacedim>::get_function_values(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
- VectorSlice<std::vector<std::vector<typename InputVector::value_type>>>
- values,
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> & indices,
+ ArrayView<std::vector<typename InputVector::value_type>> values,
bool quadrature_points_fastest) const
{
using Number = typename InputVector::value_type;
template <class InputVector>
void
FEValuesBase<dim, spacedim>::get_function_gradients(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> &indices,
std::vector<Tensor<1, spacedim, typename InputVector::value_type>> &gradients)
const
{
template <class InputVector>
void
FEValuesBase<dim, spacedim>::get_function_gradients(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
- VectorSlice<std::vector<
- std::vector<Tensor<1, spacedim, typename InputVector::value_type>>>>
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> &indices,
+ ArrayView<std::vector<Tensor<1, spacedim, typename InputVector::value_type>>>
gradients,
bool quadrature_points_fastest) const
{
template <class InputVector>
void
FEValuesBase<dim, spacedim>::get_function_hessians(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> &indices,
std::vector<Tensor<2, spacedim, typename InputVector::value_type>> &hessians)
const
{
template <class InputVector>
void
FEValuesBase<dim, spacedim>::get_function_hessians(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
- VectorSlice<std::vector<
- std::vector<Tensor<2, spacedim, typename InputVector::value_type>>>>
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> &indices,
+ ArrayView<std::vector<Tensor<2, spacedim, typename InputVector::value_type>>>
hessians,
bool quadrature_points_fastest) const
{
template <class InputVector>
void
FEValuesBase<dim, spacedim>::get_function_laplacians(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
- std::vector<typename InputVector::value_type> &laplacians) const
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> &indices,
+ std::vector<typename InputVector::value_type> & laplacians) const
{
using Number = typename InputVector::value_type;
Assert(this->update_flags & update_hessians,
template <class InputVector>
void
FEValuesBase<dim, spacedim>::get_function_laplacians(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> & indices,
std::vector<Vector<typename InputVector::value_type>> &laplacians) const
{
using Number = typename InputVector::value_type;
template <class InputVector>
void
FEValuesBase<dim, spacedim>::get_function_laplacians(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
- std::vector<std::vector<typename InputVector::value_type>> & laplacians,
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> & indices,
+ std::vector<std::vector<typename InputVector::value_type>> &laplacians,
bool quadrature_points_fastest) const
{
using Number = typename InputVector::value_type;
template <class InputVector>
void
FEValuesBase<dim, spacedim>::get_function_third_derivatives(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> &indices,
std::vector<Tensor<3, spacedim, typename InputVector::value_type>>
&third_derivatives) const
{
template <class InputVector>
void
FEValuesBase<dim, spacedim>::get_function_third_derivatives(
- const InputVector & fe_function,
- const VectorSlice<const std::vector<types::global_dof_index>> &indices,
- VectorSlice<std::vector<
- std::vector<Tensor<3, spacedim, typename InputVector::value_type>>>>
+ const InputVector & fe_function,
+ const ArrayView<const types::global_dof_index> &indices,
+ ArrayView<std::vector<Tensor<3, spacedim, typename InputVector::value_type>>>
third_derivatives,
bool quadrature_points_fastest) const
{
get_function_values<VEC>(const VEC &, std::vector<VEC::value_type> &)
const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
- get_function_values<VEC>(
- const VEC &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
- std::vector<VEC::value_type> &) const;
+ get_function_values<VEC>(const VEC &,
+ const ArrayView<const types::global_dof_index> &,
+ std::vector<VEC::value_type> &) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_values<VEC>(const VEC &,
std::vector<Vector<VEC::value_type>> &) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
- get_function_values<VEC>(
- const VEC &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
- std::vector<Vector<VEC::value_type>> &) const;
+ get_function_values<VEC>(const VEC &,
+ const ArrayView<const types::global_dof_index> &,
+ std::vector<Vector<VEC::value_type>> &) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
- get_function_values<VEC>(
- const VEC &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
- VectorSlice<std::vector<std::vector<VEC::value_type>>>,
- bool) const;
+ get_function_values<VEC>(const VEC &,
+ const ArrayView<const types::global_dof_index> &,
+ ArrayView<std::vector<VEC::value_type>>,
+ bool) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_gradients<VEC>(
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_gradients<VEC>(
const VEC &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
+ const ArrayView<const types::global_dof_index> &,
std::vector<dealii::Tensor<1, deal_II_space_dimension, VEC::value_type>>
&) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_gradients<VEC>(
const VEC &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
- VectorSlice<std::vector<std::vector<
- dealii::Tensor<1, deal_II_space_dimension, VEC::value_type>>>>,
+ const ArrayView<const types::global_dof_index> &,
+ ArrayView<std::vector<
+ dealii::Tensor<1, deal_II_space_dimension, VEC::value_type>>>,
bool) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_hessians<VEC>(
const VEC &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
+ const ArrayView<const types::global_dof_index> &,
std::vector<dealii::Tensor<2, deal_II_space_dimension, VEC::value_type>>
&) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_hessians<VEC>(
const VEC &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
- VectorSlice<std::vector<std::vector<
- dealii::Tensor<2, deal_II_space_dimension, VEC::value_type>>>>,
+ const ArrayView<const types::global_dof_index> &,
+ ArrayView<std::vector<
+ dealii::Tensor<2, deal_II_space_dimension, VEC::value_type>>>,
bool) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_laplacians<VEC>(
const VEC &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
+ const ArrayView<const types::global_dof_index> &,
std::vector<VEC::value_type> &) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_laplacians<VEC>(
const VEC &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
+ const ArrayView<const types::global_dof_index> &,
std::vector<Vector<VEC::value_type>> &) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_laplacians<VEC>(
const VEC &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
+ const ArrayView<const types::global_dof_index> &,
std::vector<std::vector<VEC::value_type>> &,
bool) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_third_derivatives<VEC>(
const VEC &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
+ const ArrayView<const types::global_dof_index> &,
std::vector<dealii::Tensor<3, deal_II_space_dimension, VEC::value_type>>
&) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_third_derivatives<VEC>(
const VEC &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
- VectorSlice<std::vector<std::vector<
- dealii::Tensor<3, deal_II_space_dimension, VEC::value_type>>>>,
+ const ArrayView<const types::global_dof_index> &,
+ ArrayView<std::vector<
+ dealii::Tensor<3, deal_II_space_dimension, VEC::value_type>>>,
bool) const;
#endif
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_values<IndexSet>(
const IndexSet &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
+ const ArrayView<const types::global_dof_index> &,
std::vector<IndexSet::value_type> &) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_values<IndexSet>(
const IndexSet &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
+ const ArrayView<const types::global_dof_index> &,
std::vector<Vector<IndexSet::value_type>> &) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_values<IndexSet>(
const IndexSet &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
- VectorSlice<std::vector<std::vector<IndexSet::value_type>>>,
+ const ArrayView<const types::global_dof_index> &,
+ ArrayView<std::vector<IndexSet::value_type>>,
bool) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_gradients<IndexSet>(
const IndexSet &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
+ const ArrayView<const types::global_dof_index> &,
std::vector<
dealii::Tensor<1, deal_II_space_dimension, IndexSet::value_type>> &)
const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_gradients<IndexSet>(
const IndexSet &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
- VectorSlice<std::vector<std::vector<
- dealii::Tensor<1, deal_II_space_dimension, IndexSet::value_type>>>>,
+ const ArrayView<const types::global_dof_index> &,
+ ArrayView<std::vector<
+ dealii::Tensor<1, deal_II_space_dimension, IndexSet::value_type>>>,
bool) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_hessians<IndexSet>(
const IndexSet &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
+ const ArrayView<const types::global_dof_index> &,
std::vector<
dealii::Tensor<2, deal_II_space_dimension, IndexSet::value_type>> &)
const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_hessians<IndexSet>(
const IndexSet &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
- VectorSlice<std::vector<std::vector<
- dealii::Tensor<2, deal_II_space_dimension, IndexSet::value_type>>>>,
+ const ArrayView<const types::global_dof_index> &,
+ ArrayView<std::vector<
+ dealii::Tensor<2, deal_II_space_dimension, IndexSet::value_type>>>,
bool) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_laplacians<IndexSet>(
const IndexSet &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
+ const ArrayView<const types::global_dof_index> &,
std::vector<IndexSet::value_type> &) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_laplacians<IndexSet>(
const IndexSet &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
+ const ArrayView<const types::global_dof_index> &,
std::vector<Vector<IndexSet::value_type>> &) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_laplacians<IndexSet>(
const IndexSet &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
+ const ArrayView<const types::global_dof_index> &,
std::vector<std::vector<IndexSet::value_type>> &,
bool) const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_third_derivatives<IndexSet>(
const IndexSet &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
+ const ArrayView<const types::global_dof_index> &,
std::vector<
dealii::Tensor<3, deal_II_space_dimension, IndexSet::value_type>> &)
const;
template void FEValuesBase<deal_II_dimension, deal_II_space_dimension>::
get_function_third_derivatives<IndexSet>(
const IndexSet &,
- const VectorSlice<const std::vector<types::global_dof_index>> &,
- VectorSlice<std::vector<std::vector<
- dealii::Tensor<3, deal_II_space_dimension, IndexSet::value_type>>>>,
+ const ArrayView<const types::global_dof_index> &,
+ ArrayView<std::vector<
+ dealii::Tensor<3, deal_II_space_dimension, IndexSet::value_type>>>,
bool) const;
#endif
}
row_lengths[j][0]);
else
{
- VectorSlice<const std::vector<unsigned int>> block_rows(
- row_lengths[j], start, length);
+ Assert(row_lengths[j].begin() + start + length <=
+ row_lengths[j].end(),
+ ExcInternalError());
+ ArrayView<const unsigned int> block_rows(row_lengths[j].data() +
+ start,
+ length);
block(i, j).reinit(rows.block_size(i),
cols.block_size(j),
block_rows);
void
-ChunkSparsityPattern::reinit(
- const size_type m,
- const size_type n,
- const VectorSlice<const std::vector<size_type>> &row_lengths,
- const size_type chunk_size)
+ChunkSparsityPattern::reinit(const size_type m,
+ const size_type n,
+ const ArrayView<const size_type> &row_lengths,
+ const size_type chunk_size)
{
Assert(row_lengths.size() == m, ExcInvalidNumber(m));
Assert(chunk_size > 0, ExcInvalidNumber(chunk_size));
void
-SparsityPattern::reinit(
- const size_type m,
- const size_type n,
- const VectorSlice<const std::vector<unsigned int>> &row_lengths)
+SparsityPattern::reinit(const size_type m,
+ const size_type n,
+ const ArrayView<const unsigned int> &row_lengths)
{
AssertDimension(row_lengths.size(), m);