17 #ifndef dealii_matrix_free_operators_h 18 #define dealii_matrix_free_operators_h 46 template <
typename VectorType>
51 return vector.n_blocks();
54 template <
typename VectorType>
62 template <
typename VectorType>
64 typename VectorType::BlockType &>::type
67 return vector.block(block_no);
70 template <
typename VectorType>
72 const typename VectorType::BlockType &>::type
76 return vector.block(block_no);
79 template <
typename VectorType>
87 template <
typename VectorType>
95 template <
typename VectorType>
99 vector.collect_sizes();
102 template <
typename VectorType>
189 typename VectorizedArrayType =
212 virtual ~
Base()
override =
default;
238 initialize(std::shared_ptr<
240 const std::vector<unsigned int> &selected_row_blocks =
241 std::vector<unsigned int>(),
242 const std::vector<unsigned int> &selected_column_blocks =
243 std::vector<unsigned int>());
259 initialize(std::shared_ptr<
262 const unsigned int level,
263 const std::vector<unsigned int> &selected_row_blocks =
264 std::vector<unsigned int>());
281 initialize(std::shared_ptr<
283 const std::vector<MGConstrainedDoFs> &mg_constrained_dofs,
284 const unsigned int level,
285 const std::vector<unsigned int> & selected_row_blocks =
286 std::vector<unsigned int>());
304 vmult_interface_down(VectorType &dst,
const VectorType &src)
const;
310 vmult_interface_up(VectorType &dst,
const VectorType &src)
const;
316 vmult(VectorType &dst,
const VectorType &src)
const;
322 Tvmult(VectorType &dst,
const VectorType &src)
const;
328 vmult_add(VectorType &dst,
const VectorType &src)
const;
334 Tvmult_add(VectorType &dst,
const VectorType &src)
const;
341 el(
const unsigned int row,
const unsigned int col)
const;
354 initialize_dof_vector(VectorType &vec)
const;
364 compute_diagonal() = 0;
369 std::shared_ptr<const MatrixFree<dim, value_type, VectorizedArrayType>>
370 get_matrix_free()
const;
375 const std::shared_ptr<DiagonalMatrix<VectorType>> &
376 get_matrix_diagonal_inverse()
const;
381 const std::shared_ptr<DiagonalMatrix<VectorType>> &
382 get_matrix_diagonal()
const;
391 precondition_Jacobi(VectorType & dst,
392 const VectorType &src,
401 preprocess_constraints(VectorType &dst,
const VectorType &src)
const;
408 postprocess_constraints(VectorType &dst,
const VectorType &src)
const;
415 set_constrained_entries_to_one(VectorType &dst)
const;
421 apply_add(VectorType &dst,
const VectorType &src)
const = 0;
429 Tapply_add(VectorType &dst,
const VectorType &src)
const;
434 std::shared_ptr<const MatrixFree<dim, value_type, VectorizedArrayType>>
470 mutable std::vector<std::vector<std::pair<value_type, value_type>>>
484 mult_add(VectorType & dst,
485 const VectorType &src,
496 adjust_ghost_range_if_necessary(
const VectorType &vec,
497 const bool is_row)
const;
538 template <
typename OperatorType>
567 initialize(
const OperatorType &operator_in);
572 template <
typename VectorType>
574 vmult(VectorType &dst,
const VectorType &src)
const;
579 template <
typename VectorType>
581 Tvmult(VectorType &dst,
const VectorType &src)
const;
586 template <
typename VectorType>
588 initialize_dof_vector(VectorType &vec)
const;
628 "Type of Number and of VectorizedArrayType do not match.");
640 VectorizedArrayType> &fe_eval);
652 const unsigned int n_actual_components,
653 const VectorizedArrayType * in_array,
654 VectorizedArrayType * out_array)
const;
668 apply(
const VectorizedArrayType *in_array,
669 VectorizedArrayType * out_array)
const;
705 transform_from_q_points_to_basis(
const unsigned int n_actual_components,
706 const VectorizedArrayType *in_array,
707 VectorizedArrayType *out_array)
const;
715 fill_inverse_JxW_values(
742 int n_q_points_1d = fe_degree + 1,
743 int n_components = 1,
745 typename VectorizedArrayType =
772 compute_diagonal()
override;
781 apply_add(VectorType &dst,
const VectorType &src)
const override;
790 const VectorType & src,
791 const std::pair<unsigned int, unsigned int> & cell_range)
const;
810 int n_q_points_1d = fe_degree + 1,
811 int n_components = 1,
813 typename VectorizedArrayType =
842 compute_diagonal()
override;
911 std::shared_ptr<Table<2, VectorizedArrayType>>
921 apply_add(VectorType &dst,
const VectorType &src)
const override;
930 const VectorType & src,
931 const std::pair<unsigned int, unsigned int> & cell_range)
const;
941 const std::pair<unsigned int, unsigned int> &cell_range)
const;
947 do_operation_on_cell(
950 const unsigned int cell)
const;
966 typename VectorizedArrayType>
971 VectorizedArrayType>::
972 CellwiseInverseMassMatrix(
977 VectorizedArrayType> &fe_eval)
987 typename VectorizedArrayType>
993 VectorizedArrayType>
:: 997 constexpr
unsigned int dofs_per_component_on_cell =
1000 inverse_jxw.
size() % dofs_per_component_on_cell == 0,
1002 "Expected diagonal to be a multiple of scalar dof per cells"));
1006 for (
unsigned int q = 0; q < dofs_per_component_on_cell; ++q)
1009 for (
unsigned int q = dofs_per_component_on_cell; q < inverse_jxw.
size();)
1010 for (
unsigned int i = 0; i < dofs_per_component_on_cell; ++i, ++q)
1011 inverse_jxw[q] = inverse_jxw[i];
1020 typename VectorizedArrayType>
1027 VectorizedArrayType>
::apply(
const VectorizedArrayType *in_array,
1028 VectorizedArrayType * out_array)
const 1043 typename VectorizedArrayType>
1049 VectorizedArrayType>
:: 1051 const unsigned int n_actual_components,
1052 const VectorizedArrayType * in_array,
1053 VectorizedArrayType * out_array)
const 1058 VectorizedArrayType>
:: 1060 inverse_coefficients,
1061 n_actual_components,
1072 typename VectorizedArrayType>
1078 VectorizedArrayType>
:: 1080 const VectorizedArrayType *in_array,
1081 VectorizedArrayType * out_array)
const 1086 VectorizedArrayType>
:: 1089 n_actual_components,
1097 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1100 , have_interface_matrices(false)
1105 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1113 total_size +=
data->get_vector_partitioner(selected_row)->size();
1119 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1127 total_size +=
data->get_vector_partitioner(selected_column)->size();
1133 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1143 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1146 const unsigned int col)
const 1153 return 1.0 / (*inverse_diagonal_entries)(row, row);
1158 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1161 VectorType &vec)
const 1169 .partitioners_are_compatible(
1170 *
data->get_dof_info(index).vector_partitioner))
1174 .partitioners_are_globally_compatible(
1175 *
data->get_dof_info(index).vector_partitioner),
1183 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1188 const std::vector<unsigned int> &given_row_selection,
1189 const std::vector<unsigned int> &given_column_selection)
1195 if (given_row_selection.empty())
1196 for (
unsigned int i = 0; i < data_->n_components(); ++i)
1200 for (
unsigned int i = 0; i < given_row_selection.size(); ++i)
1203 for (
unsigned int j = 0; j < given_row_selection.size(); ++j)
1205 Assert(given_row_selection[j] != given_row_selection[i],
1206 ExcMessage(
"Given row indices must be unique"));
1211 if (given_column_selection.size() == 0)
1215 for (
unsigned int i = 0; i < given_column_selection.size(); ++i)
1218 for (
unsigned int j = 0; j < given_column_selection.size(); ++j)
1220 Assert(given_column_selection[j] != given_column_selection[i],
1221 ExcMessage(
"Given column indices must be unique"));
1236 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1242 const unsigned int level,
1243 const std::vector<unsigned int> &given_row_selection)
1245 std::vector<MGConstrainedDoFs> mg_constrained_dofs_vector(
1246 1, mg_constrained_dofs);
1247 initialize(data_, mg_constrained_dofs_vector, level, given_row_selection);
1252 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1257 const std::vector<MGConstrainedDoFs> &mg_constrained_dofs,
1258 const unsigned int level,
1259 const std::vector<unsigned int> & given_row_selection)
1266 if (given_row_selection.empty())
1267 for (
unsigned int i = 0; i < data_->n_components(); ++i)
1271 for (
unsigned int i = 0; i < given_row_selection.size(); ++i)
1274 for (
unsigned int j = 0; j < given_row_selection.size(); ++j)
1276 Assert(given_row_selection[j] != given_row_selection[i],
1277 ExcMessage(
"Given row indices must be unique"));
1294 if (data_->n_macro_cells() > 0)
1300 std::vector<types::global_dof_index> interface_indices;
1301 mg_constrained_dofs[j]
1302 .get_refinement_edge_indices(level)
1303 .fill_index_vector(interface_indices);
1309 for (
const auto interface_index : interface_indices)
1310 if (locally_owned.
is_element(interface_index))
1316 data->get_vector_partitioner()->get_mpi_communicator()) > 0;
1322 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1325 VectorType &dst)
const 1329 const std::vector<unsigned int> &constrained_dofs =
1331 for (
const auto constrained_dof : constrained_dofs)
1341 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1344 const VectorType &src)
const 1354 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1358 const VectorType &src)
const 1365 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1369 const VectorType &src)
const 1376 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1379 const VectorType &src,
1380 const bool is_row)
const 1386 const unsigned int mf_component =
1390 data->get_dof_info(mf_component).vector_partitioner.get())
1397 data->get_dof_info(mf_component).vector_partitioner->local_size(),
1398 ExcMessage(
"The vector passed to the vmult() function does not have " 1399 "the correct size for compatibility with MatrixFree."));
1406 .reinit(
data->get_dof_info(mf_component).vector_partitioner);
1408 .copy_locally_owned_data_from(copy_vec);
1414 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1418 const VectorType &src)
const 1444 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1448 const VectorType &src,
1464 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1468 const VectorType &src)
const 1472 const std::vector<unsigned int> &constrained_dofs =
1474 for (
const auto constrained_dof : constrained_dofs)
1498 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1502 const VectorType &src)
const 1542 .local_element(edge_constrained_indices[j][i]) =
1552 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1556 const VectorType &src)
const 1569 VectorType src_cpy(src);
1593 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1597 const VectorType &src)
const 1607 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1618 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1622 VectorizedArrayType>>
1630 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1631 const std::shared_ptr<DiagonalMatrix<VectorType>> &
1643 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1644 const std::shared_ptr<DiagonalMatrix<VectorType>> &
1654 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1658 const VectorType &src)
const 1665 template <
int dim,
typename VectorType,
typename VectorizedArrayType>
1669 const VectorType & src,
1670 const typename Base<dim, VectorType, VectorizedArrayType>::value_type omega)
1683 template <
typename OperatorType>
1686 , mf_base_operator(nullptr)
1691 template <
typename OperatorType>
1700 template <
typename OperatorType>
1709 template <
typename OperatorType>
1710 template <
typename VectorType>
1713 const VectorType &src)
const 1715 #ifndef DEAL_II_MSVC 1718 "The vector type must be based on the same value type as this " 1729 template <
typename OperatorType>
1730 template <
typename VectorType>
1733 const VectorType &src)
const 1735 #ifndef DEAL_II_MSVC 1738 "The vector type must be based on the same value type as this " 1749 template <
typename OperatorType>
1750 template <
typename VectorType>
1753 VectorType &vec)
const 1768 typename VectorType,
1769 typename VectorizedArrayType>
1776 :
Base<dim, VectorType, VectorizedArrayType>()
1785 typename VectorType,
1786 typename VectorizedArrayType>
1796 typename Base<dim, VectorType, VectorizedArrayType>::value_type;
1801 std::make_shared<DiagonalMatrix<VectorType>>();
1803 VectorType &inverse_diagonal_vector =
1805 VectorType &diagonal_vector = this->diagonal_entries->get_vector();
1808 inverse_diagonal_vector = Number(1.);
1809 apply_add(diagonal_vector, inverse_diagonal_vector);
1812 inverse_diagonal_vector = diagonal_vector;
1814 const unsigned int local_size = inverse_diagonal_vector.local_size();
1815 for (
unsigned int i = 0; i < local_size; ++i)
1816 inverse_diagonal_vector.local_element(i) =
1817 Number(1.) / inverse_diagonal_vector.local_element(i);
1819 inverse_diagonal_vector.update_ghost_values();
1820 diagonal_vector.update_ghost_values();
1829 typename VectorType,
1830 typename VectorizedArrayType>
1838 const VectorType &src)
const 1850 typename VectorType,
1851 typename VectorizedArrayType>
1858 VectorizedArrayType>
:: 1862 typename Base<dim, VectorType, VectorizedArrayType>::value_type,
1863 VectorizedArrayType> &
data,
1865 const VectorType & src,
1866 const std::pair<unsigned int, unsigned int> &cell_range)
const 1869 typename Base<dim, VectorType, VectorizedArrayType>::value_type;
1875 VectorizedArrayType>
1877 for (
unsigned int cell = cell_range.first; cell < cell_range.second; ++cell)
1880 phi.read_dof_values(src);
1881 phi.evaluate(
true,
false,
false);
1882 for (
unsigned int q = 0; q < phi.n_q_points; ++q)
1883 phi.submit_value(phi.get_value(q), q);
1884 phi.integrate(
true,
false);
1885 phi.distribute_local_to_global(dst);
1896 typename VectorType,
1897 typename VectorizedArrayType>
1904 :
Base<dim, VectorType, VectorizedArrayType>()
1913 typename VectorType,
1914 typename VectorizedArrayType>
1933 typename VectorType,
1934 typename VectorizedArrayType>
1941 VectorizedArrayType>
:: 1954 typename VectorType,
1955 typename VectorizedArrayType>
1956 std::shared_ptr<Table<2, VectorizedArrayType>>
1974 typename VectorType,
1975 typename VectorizedArrayType>
1985 typename Base<dim, VectorType, VectorizedArrayType>::value_type;
1990 std::make_shared<DiagonalMatrix<VectorType>>();
1992 VectorType &inverse_diagonal_vector =
1994 VectorType &diagonal_vector = this->diagonal_entries->get_vector();
2004 inverse_diagonal_vector = diagonal_vector;
2006 for (
unsigned int i = 0; i < inverse_diagonal_vector.local_size(); ++i)
2007 if (std::abs(inverse_diagonal_vector.local_element(i)) >
2009 inverse_diagonal_vector.local_element(i) =
2010 1. / inverse_diagonal_vector.local_element(i);
2012 inverse_diagonal_vector.local_element(i) = 1.;
2014 inverse_diagonal_vector.update_ghost_values();
2015 diagonal_vector.update_ghost_values();
2024 typename VectorType,
2025 typename VectorizedArrayType>
2033 const VectorType &src)
const 2039 namespace Implementation
2041 template <
typename VectorizedArrayType>
2045 for (
unsigned int v = 0; v < VectorizedArrayType::size(); ++v)
2059 typename VectorType,
2060 typename VectorizedArrayType>
2067 VectorizedArrayType>
:: 2074 typename Base<dim, VectorType, VectorizedArrayType>::value_type> &phi,
2075 const unsigned int cell)
const 2077 phi.evaluate(
false,
true,
false);
2082 ExcMessage(
"The number of columns in the coefficient table must " 2083 "be either 1 or the number of quadrature points " +
2085 ", but the given value was " +
2088 for (
unsigned int q = 0; q < phi.n_q_points; ++q)
2092 ExcMessage(
"Coefficient must be non-negative"));
2094 phi.get_gradient(q),
2100 ExcMessage(
"Coefficient must be non-negative"));
2101 const VectorizedArrayType coefficient =
2102 (*scalar_coefficient)(cell, 0);
2103 for (
unsigned int q = 0; q < phi.n_q_points; ++q)
2104 phi.submit_gradient(coefficient * phi.get_gradient(q), q);
2109 for (
unsigned int q = 0; q < phi.n_q_points; ++q)
2111 phi.submit_gradient(phi.get_gradient(q), q);
2114 phi.integrate(
false,
true);
2123 typename VectorType,
2124 typename VectorizedArrayType>
2131 VectorizedArrayType>
:: 2135 typename Base<dim, VectorType, VectorizedArrayType>::value_type,
2136 VectorizedArrayType> &
data,
2138 const VectorType & src,
2139 const std::pair<unsigned int, unsigned int> &cell_range)
const 2142 typename Base<dim, VectorType, VectorizedArrayType>::value_type;
2145 for (
unsigned int cell = cell_range.first; cell < cell_range.second; ++cell)
2159 typename VectorType,
2160 typename VectorizedArrayType>
2167 VectorizedArrayType>
:: 2171 typename Base<dim, VectorType, VectorizedArrayType>::value_type,
2172 VectorizedArrayType> &
data,
2175 const std::pair<unsigned int, unsigned int> &cell_range)
const 2178 typename Base<dim, VectorType, VectorizedArrayType>::value_type;
2182 for (
unsigned int cell = cell_range.first; cell < cell_range.second; ++cell)
2197 local_diagonal_vector[i];
typename OperatorType::value_type value_type
VectorizedArrayType JxW(const unsigned int q_index) const
virtual void apply_add(VectorType &dst, const VectorType &src) const =0
void fill_inverse_JxW_values(AlignedVector< VectorizedArrayType > &inverse_jxw) const
static const unsigned int invalid_unsigned_int
const std::shared_ptr< DiagonalMatrix< VectorType > > & get_matrix_diagonal_inverse() const
#define AssertDimension(dim1, dim2)
virtual std::size_t memory_consumption() const
const FEEvaluationBase< dim, n_components, Number, false, VectorizedArrayType > & fe_eval
types::global_dof_index size_type
void local_diagonal_cell(const MatrixFree< dim, value_type, VectorizedArrayType > &data, VectorType &dst, const VectorType &, const std::pair< unsigned int, unsigned int > &cell_range) const
void precondition_Jacobi(VectorType &dst, const VectorType &src, const value_type omega) const
std::shared_ptr< DiagonalMatrix< VectorType > > diagonal_entries
void read_dof_values(const VectorType &src, const unsigned int first_index=0)
#define AssertIndexRange(index, range)
std::vector< std::vector< std::pair< value_type, value_type > > > edge_constrained_values
static ::ExceptionBase & ExcNotInitialized()
std::shared_ptr< const MatrixFree< dim, value_type, VectorizedArrayType > > data
void Tvmult(VectorType &dst, const VectorType &src) const
#define AssertThrow(cond, exc)
std::enable_if< IsBlockVector< VectorType >::value, void >::type collect_sizes(VectorType &vector)
typename VectorType::value_type value_type
typename OperatorType::size_type size_type
virtual void compute_diagonal() override
SymmetricTensor< 2, dim, Number > epsilon(const Tensor< 2, dim, Number > &Grad_u)
constexpr T pow(const T base, const int iexp)
std::shared_ptr< const MatrixFree< dim, value_type, VectorizedArrayType > > get_matrix_free() const
void initialize_dof_vector(VectorType &vec) const
void vmult_add(VectorType &dst, const VectorType &src) const
const std::shared_ptr< DiagonalMatrix< VectorType > > & get_matrix_diagonal() const
static ::ExceptionBase & ExcMessage(std::string arg1)
void vmult_interface_up(VectorType &dst, const VectorType &src) const
std::enable_if< IsBlockVector< VectorType >::value, typename VectorType::BlockType & >::type subblock(VectorType &vector, unsigned int block_no)
#define Assert(cond, exc)
void local_apply_cell(const MatrixFree< dim, value_type, VectorizedArrayType > &data, VectorType &dst, const VectorType &src, const std::pair< unsigned int, unsigned int > &cell_range) const
size_type index_within_set(const size_type global_index) const
void postprocess_constraints(VectorType &dst, const VectorType &src) const
virtual void Tapply_add(VectorType &dst, const VectorType &src) const
std::shared_ptr< Table< 2, VectorizedArrayType > > scalar_coefficient
#define DEAL_II_NAMESPACE_CLOSE
auto apply(F &&fn, Tuple &&t) -> decltype(apply_impl(std::forward< F >(fn), std::forward< Tuple >(t), std_cxx14::make_index_sequence< std::tuple_size< typename std::remove_reference< Tuple >::type >::value >()))
const unsigned int dofs_per_component
std::vector< unsigned int > selected_rows
bool have_interface_matrices
void reinit(const unsigned int cell_batch_index)
std::enable_if< IsBlockVector< VectorType >::value, unsigned int >::type n_blocks(const VectorType &vector)
void set_constrained_entries_to_one(VectorType &dst) const
virtual void clear() override
void initialize_dof_vector(VectorType &vec) const
void vmult(VectorType &dst, const VectorType &src) const
virtual void apply_add(VectorType &dst, const VectorType &src) const override
virtual void apply_add(VectorType &dst, const VectorType &src) const override
void Tvmult_add(VectorType &dst, const VectorType &src) const
static constexpr unsigned int n_components
static constexpr unsigned int static_dofs_per_cell
typename VectorType::size_type size_type
void mult_add(VectorType &dst, const VectorType &src, const bool transpose) const
std::shared_ptr< Table< 2, VectorizedArrayType > > get_coefficient()
void do_operation_on_cell(FEEvaluation< dim, fe_degree, n_q_points_1d, n_components, value_type > &phi, const unsigned int cell) const
const VectorizedArrayType * begin_dof_values() const
const internal::MatrixFreeFunctions::ShapeInfo< VectorizedArrayType > & get_shape_info() const
#define DEAL_II_NAMESPACE_OPEN
std::vector< unsigned int > selected_columns
void vmult_interface_down(VectorType &dst, const VectorType &src) const
void distribute_local_to_global(VectorType &dst, const unsigned int first_index=0, const std::bitset< VectorizedArrayType::size()> &mask=std::bitset< VectorizedArrayType::size()>().flip()) const
static ::ExceptionBase & ExcNotImplemented()
value_type el(const unsigned int row, const unsigned int col) const
std::vector< UnivariateShapeData< Number > > data
bool is_element(const size_type index) const
void initialize(const OperatorType &operator_in)
SmartPointer< const OperatorType > mf_base_operator
void Tvmult(VectorType &dst, const VectorType &src) const
void apply(const AlignedVector< VectorizedArrayType > &inverse_coefficient, const unsigned int n_actual_components, const VectorizedArrayType *in_array, VectorizedArrayType *out_array) const
std::shared_ptr< DiagonalMatrix< VectorType > > inverse_diagonal_entries
void vmult(VectorType &dst, const VectorType &src) const
bool non_negative(const VectorizedArrayType &n)
virtual void compute_diagonal() override
T max(const T &t, const MPI_Comm &mpi_communicator)
void local_apply_cell(const MatrixFree< dim, value_type, VectorizedArrayType > &data, VectorType &dst, const VectorType &src, const std::pair< unsigned int, unsigned int > &cell_range) const
std::vector< std::vector< unsigned int > > edge_constrained_indices
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
void initialize(std::shared_ptr< const MatrixFree< dim, value_type, VectorizedArrayType >> data, const std::vector< unsigned int > &selected_row_blocks=std::vector< unsigned int >(), const std::vector< unsigned int > &selected_column_blocks=std::vector< unsigned int >())
void preprocess_constraints(VectorType &dst, const VectorType &src) const
void transform_from_q_points_to_basis(const unsigned int n_actual_components, const VectorizedArrayType *in_array, VectorizedArrayType *out_array) const
void set_coefficient(const std::shared_ptr< Table< 2, VectorizedArrayType >> &scalar_coefficient)
static ::ExceptionBase & ExcInternalError()
void adjust_ghost_range_if_necessary(const VectorType &vec, const bool is_row) const