17 #ifndef dealii__matrix_free_operators_h 18 #define dealii__matrix_free_operators_h 21 #include <deal.II/base/exceptions.h> 22 #include <deal.II/base/vectorization.h> 23 #include <deal.II/base/subscriptor.h> 25 #include <deal.II/lac/diagonal_matrix.h> 26 #include <deal.II/lac/la_parallel_vector.h> 27 #include <deal.II/lac/vector_view.h> 28 #include <deal.II/multigrid/mg_constrained_dofs.h> 29 #include <deal.II/matrix_free/matrix_free.h> 30 #include <deal.II/matrix_free/fe_evaluation.h> 33 DEAL_II_NAMESPACE_OPEN
43 template <
typename VectorType>
44 typename std_cxx11::enable_if<IsBlockVector<VectorType>::value,
46 n_blocks(
const VectorType &vector)
48 return vector.n_blocks();
51 template <
typename VectorType>
52 typename std_cxx11::enable_if<!IsBlockVector<VectorType>::value,
54 n_blocks(
const VectorType &)
59 template <
typename VectorType>
60 typename std_cxx11::enable_if<IsBlockVector<VectorType>::value,
61 typename VectorType::BlockType &>::type
62 subblock(VectorType &vector,
unsigned int block_no)
64 return vector.block(block_no);
67 template <
typename VectorType>
68 typename std_cxx11::enable_if<IsBlockVector<VectorType>::value,
69 const typename VectorType::BlockType &>::type
70 subblock(
const VectorType &vector,
unsigned int block_no)
73 return vector.block(block_no);
76 template <
typename VectorType>
77 typename std_cxx11::enable_if<!IsBlockVector<VectorType>::value,
79 subblock(VectorType &vector,
unsigned int)
84 template <
typename VectorType>
85 typename std_cxx11::enable_if<!IsBlockVector<VectorType>::value,
86 const VectorType &>::type
87 subblock(
const VectorType &vector,
unsigned int)
92 template <
typename VectorType>
93 typename std_cxx11::enable_if<IsBlockVector<VectorType>::value,
95 collect_sizes(VectorType &vector)
97 vector.collect_sizes();
100 template <
typename VectorType>
101 typename std_cxx11::enable_if<!IsBlockVector<VectorType>::value,
103 collect_sizes(
const VectorType &)
163 template <
int dim,
typename VectorType = LinearAlgebra::distributed::Vector<
double> >
191 virtual void clear();
210 const std::vector<unsigned int> &selected_row_blocks = std::vector<unsigned int>(),
211 const std::vector<unsigned int> &selected_column_blocks = std::vector<unsigned int>());
228 const unsigned int level,
229 const std::vector<unsigned int> &selected_row_blocks = std::vector<unsigned int>());
247 const std::vector<MGConstrainedDoFs> &mg_constrained_dofs,
248 const unsigned int level,
249 const std::vector<unsigned int> &selected_row_blocks = std::vector<unsigned int>());
265 const VectorType &src)
const;
271 const VectorType &src)
const;
276 void vmult (VectorType &dst,
277 const VectorType &src)
const;
282 void Tvmult (VectorType &dst,
283 const VectorType &src)
const;
289 const VectorType &src)
const;
295 const VectorType &src)
const;
302 const unsigned int col)
const;
325 std_cxx11::shared_ptr<const MatrixFree<dim,value_type> >
331 const std_cxx11::shared_ptr<DiagonalMatrix<VectorType> > &
340 const VectorType &src,
355 const VectorType &src)
const = 0;
363 const VectorType &src)
const;
368 std_cxx11::shared_ptr<const MatrixFree<dim,value_type> >
data;
398 mutable std::vector<std::vector<std::pair<value_type,value_type> > >
412 const VectorType &src,
413 const bool transpose)
const;
423 const bool is_row)
const;
464 template <
typename OperatorType>
486 void initialize (
const OperatorType &operator_in);
491 template <
typename VectorType>
492 void vmult (VectorType &dst,
493 const VectorType &src)
const;
498 template <
typename VectorType>
499 void Tvmult (VectorType &dst,
500 const VectorType &src)
const;
530 template <
int dim,
int fe_degree,
int n_components = 1,
typename Number =
double>
549 const unsigned int n_actual_components,
583 template <
int dim,
int fe_degree,
int n_q_po
ints_1d = fe_degree+1,
int n_components = 1,
typename VectorType = LinearAlgebra::distributed::Vector<
double> >
614 const VectorType &src)
const;
621 const VectorType &src,
622 const std::pair<unsigned int,unsigned int> &cell_range)
const;
638 template <
int dim,
int fe_degree,
int n_q_po
ints_1d = fe_degree+1,
int n_components = 1,
typename VectorType = LinearAlgebra::distributed::Vector<
double> >
710 virtual void clear();
718 std_cxx11::shared_ptr< Table<2, VectorizedArray<value_type> > >
get_coefficient();
727 const VectorType &src)
const;
734 const VectorType &src,
735 const std::pair<unsigned int,unsigned int> &cell_range)
const;
742 const unsigned int &,
743 const std::pair<unsigned int,unsigned int> &cell_range)
const;
749 const unsigned int cell)
const;
761 template <
int dim,
int fe_degree,
int n_components,
typename Number>
769 for (
unsigned int i=0, c=0; i<shapes_1d.
m(); ++i)
770 for (
unsigned int j=0; j<shapes_1d.
n(); ++j, ++c)
773 const unsigned int stride = (fe_degree+2)/2;
774 inverse_shape.resize(stride*(fe_degree+1));
775 for (
unsigned int i=0; i<stride; ++i)
776 for (
unsigned int q=0; q<(fe_degree+2)/2; ++q)
778 inverse_shape[i*stride+q] =
779 0.5 * (shapes_1d(i,q) + shapes_1d(i,fe_degree-q));
780 inverse_shape[(fe_degree-i)*stride+q] =
781 0.5 * (shapes_1d(i,q) - shapes_1d(i,fe_degree-q));
783 if (fe_degree % 2 == 0)
784 for (
unsigned int q=0; q<(fe_degree+2)/2; ++q)
785 inverse_shape[fe_degree/2*stride+q] = shapes_1d(fe_degree/2,q);
790 template <
int dim,
int fe_degree,
int n_components,
typename Number>
797 Assert(inverse_jxw.size() > 0 &&
798 inverse_jxw.size() % dofs_per_cell == 0,
799 ExcMessage(
"Expected diagonal to be a multiple of scalar dof per cells"));
803 const unsigned int previous_size = inverse_jxw.size();
804 inverse_jxw.resize(dofs_per_cell);
805 fe_eval.fill_JxW_values(inverse_jxw);
808 inverse_jxw.resize_fast(previous_size);
809 for (
unsigned int q=0; q<dofs_per_cell; ++q)
810 inverse_jxw[q] = 1./inverse_jxw[q];
812 for (
unsigned int q=dofs_per_cell; q<previous_size; )
813 for (
unsigned int i=0; i<dofs_per_cell; ++i, ++q)
814 inverse_jxw[q] = inverse_jxw[i];
819 template <
int dim,
int fe_degree,
int n_components,
typename Number>
824 const unsigned int n_actual_components,
829 Assert(inverse_coefficients.size() > 0 &&
830 inverse_coefficients.size() % dofs_per_cell == 0,
831 ExcMessage(
"Expected diagonal to be a multiple of scalar dof per cells"));
832 if (inverse_coefficients.size() != dofs_per_cell)
833 AssertDimension(n_actual_components * dofs_per_cell, inverse_coefficients.size());
839 evaluator(inverse_shape, inverse_shape, inverse_shape);
841 const unsigned int shift_coefficient =
842 inverse_coefficients.size() > dofs_per_cell ? dofs_per_cell : 0;
845 for (
unsigned int d=0; d<n_actual_components; ++d)
851 evaluator.template hessians<0,false,false> (in, temp_data_field);
852 evaluator.template hessians<1,false,false> (temp_data_field, out);
856 evaluator.template hessians<2,false,false> (out, temp_data_field);
857 for (
unsigned int q=0; q<dofs_per_cell; ++q)
858 temp_data_field[q] *= inv_coefficient[q];
859 evaluator.template hessians<2,true,false> (temp_data_field, out);
862 for (
unsigned int q=0; q<dofs_per_cell; ++q)
863 out[q] *= inv_coefficient[q];
865 evaluator.template hessians<1,true,false>(out, temp_data_field);
866 evaluator.template hessians<0,true,false>(temp_data_field, out);
868 inv_coefficient += shift_coefficient;
873 template <
int dim,
typename VectorType>
880 template <
int dim,
typename VectorType>
884 have_interface_matrices(false)
893 template <
int dim,
typename VectorType>
897 Assert(data.get() != NULL,
900 for (
unsigned int i=0; i<selected_rows.size(); ++i)
901 total_size += data->get_vector_partitioner(selected_rows[i])->size();
907 template <
int dim,
typename VectorType>
911 Assert(data.get() != NULL,
914 for (
unsigned int i=0; i<selected_columns.size(); ++i)
915 total_size += data->get_vector_partitioner(selected_columns[i])->size();
921 template <
int dim,
typename VectorType>
926 inverse_diagonal_entries.reset();
931 template <
int dim,
typename VectorType>
934 const unsigned int col)
const 938 Assert (inverse_diagonal_entries.get() != NULL &&
940 return 1.0/(*inverse_diagonal_entries)(row,row);
945 template <
int dim,
typename VectorType>
949 Assert(data.get() != NULL,
952 for (
unsigned int i = 0; i < n_blocks(vec); ++i)
954 const unsigned int index = selected_rows[i];
955 if (!subblock(vec,index).partitioners_are_compatible
956 (*data->get_dof_info(index).vector_partitioner))
957 data->initialize_dof_vector(subblock(vec,index));
958 Assert(subblock(vec,index).partitioners_are_globally_compatible
959 (*data->get_dof_info(index).vector_partitioner),
967 template <
int dim,
typename VectorType>
971 const std::vector<unsigned int> &given_row_selection,
972 const std::vector<unsigned int> &given_column_selection)
976 selected_rows.clear();
977 selected_columns.clear();
978 if (given_row_selection.empty())
979 for (
unsigned int i=0; i<data_->n_components(); ++i)
980 selected_rows.push_back(i);
983 for (
unsigned int i=0; i<given_row_selection.size(); ++i)
986 for (
unsigned int j=0; j<given_row_selection.size(); ++j)
988 Assert(given_row_selection[j] != given_row_selection[i],
989 ExcMessage(
"Given row indices must be unique"));
991 selected_rows.push_back(given_row_selection[i]);
994 if (given_column_selection.size() == 0)
995 selected_columns = selected_rows;
998 for (
unsigned int i=0; i<given_column_selection.size(); ++i)
1001 for (
unsigned int j=0; j<given_column_selection.size(); ++j)
1003 Assert(given_column_selection[j] != given_column_selection[i],
1004 ExcMessage(
"Given column indices must be unique"));
1006 selected_columns.push_back(given_column_selection[i]);
1010 edge_constrained_indices.clear();
1011 edge_constrained_indices.resize(selected_rows.size());
1012 edge_constrained_values.clear();
1013 edge_constrained_values.resize(selected_rows.size());
1014 have_interface_matrices =
false;
1019 template <
int dim,
typename VectorType>
1024 const unsigned int level,
1025 const std::vector<unsigned int> &given_row_selection)
1027 std::vector<MGConstrainedDoFs> mg_constrained_dofs_vector(1, mg_constrained_dofs);
1028 initialize(data_, mg_constrained_dofs_vector, level, given_row_selection);
1033 template <
int dim,
typename VectorType>
1037 const std::vector<MGConstrainedDoFs> &mg_constrained_dofs,
1038 const unsigned int level,
1039 const std::vector<unsigned int> &given_row_selection)
1044 selected_rows.clear();
1045 selected_columns.clear();
1046 if (given_row_selection.empty())
1047 for (
unsigned int i=0; i<data_->n_components(); ++i)
1048 selected_rows.push_back(i);
1051 for (
unsigned int i=0; i<given_row_selection.size(); ++i)
1054 for (
unsigned int j=0; j<given_row_selection.size(); ++j)
1056 Assert(given_row_selection[j] != given_row_selection[i],
1057 ExcMessage(
"Given row indices must be unique"));
1059 selected_rows.push_back(given_row_selection[i]);
1062 selected_columns = selected_rows;
1065 edge_constrained_indices.clear();
1066 edge_constrained_indices.resize(selected_rows.size());
1067 edge_constrained_values.clear();
1068 edge_constrained_values.resize(selected_rows.size());
1072 for (
unsigned int j = 0; j < selected_rows.size(); ++j)
1074 if (data_->n_macro_cells() > 0)
1077 data_->get_cell_iterator(0,0,j)->level());
1081 std::vector<types::global_dof_index> interface_indices;
1082 mg_constrained_dofs[j].get_refinement_edge_indices(level)
1083 .fill_index_vector(interface_indices);
1084 edge_constrained_indices[j].clear();
1085 edge_constrained_indices[j].reserve(interface_indices.size());
1086 edge_constrained_values[j].resize(interface_indices.size());
1087 const IndexSet &locally_owned = data->get_dof_handler(selected_rows[j]).locally_owned_mg_dofs(level);
1088 for (
unsigned int i=0; i<interface_indices.size(); ++i)
1089 if (locally_owned.
is_element(interface_indices[i]))
1090 edge_constrained_indices[j].push_back
1092 have_interface_matrices |=
1094 data->get_vector_partitioner()->get_mpi_communicator()) > 0;
1100 template <
int dim,
typename VectorType>
1104 for (
unsigned int j=0; j<n_blocks(dst); ++j)
1106 const std::vector<unsigned int> &
1107 constrained_dofs = data->get_constrained_dofs(selected_rows[j]);
1108 for (
unsigned int i=0; i<constrained_dofs.size(); ++i)
1109 subblock(dst,j).local_element(constrained_dofs[i]) = 1.;
1110 for (
unsigned int i=0; i<edge_constrained_indices[j].size(); ++i)
1111 subblock(dst,j).local_element(edge_constrained_indices[j][i]) = 1.;
1117 template <
int dim,
typename VectorType>
1120 const VectorType &src)
const 1124 vmult_add (dst, src);
1129 template <
int dim,
typename VectorType>
1132 const VectorType &src)
const 1134 mult_add (dst, src,
false);
1139 template <
int dim,
typename VectorType>
1142 const VectorType &src)
const 1144 mult_add (dst, src,
true);
1149 template <
int dim,
typename VectorType>
1152 const bool is_row)
const 1155 for (
unsigned int i=0; i<n_blocks(src); ++i)
1157 const unsigned int mf_component = is_row ? selected_rows[i] : selected_columns[i];
1159 if (subblock(src,i).get_partitioner().get() ==
1160 data->get_dof_info(mf_component).vector_partitioner.get())
1165 Assert(subblock(src,i).get_partitioner()->local_size() ==
1166 data->get_dof_info(mf_component).vector_partitioner->local_size(),
1167 ExcMessage(
"The vector passed to the vmult() function does not have " 1168 "the correct size for compatibility with MatrixFree."));
1173 subblock(src,i).begin());
1175 subblock(const_cast<VectorType &>(src),i).
1176 reinit(data->get_dof_info(mf_component).vector_partitioner);
1178 subblock(src,i).begin());
1185 template <
int dim,
typename VectorType>
1188 const VectorType &src,
1189 const bool transpose)
const 1195 adjust_ghost_range_if_necessary(src,
false);
1196 adjust_ghost_range_if_necessary(dst,
true);
1200 for (
unsigned int j=0; j<n_blocks(dst); ++j)
1202 for (
unsigned int i=0; i<edge_constrained_indices[j].size(); ++i)
1204 edge_constrained_values[j][i] =
1205 std::pair<Number,Number>
1206 (subblock(src,j).local_element(edge_constrained_indices[j][i]),
1207 subblock(dst,j).local_element(edge_constrained_indices[j][i]));
1208 subblock(const_cast<VectorType &>(src),j).local_element(edge_constrained_indices[j][i]) = 0.;
1213 Tapply_add(dst,src);
1217 for (
unsigned int j=0; j<n_blocks(dst); ++j)
1219 const std::vector<unsigned int> &
1220 constrained_dofs = data->get_constrained_dofs(selected_rows[j]);
1221 for (
unsigned int i=0; i<constrained_dofs.size(); ++i)
1222 subblock(dst,j).local_element(constrained_dofs[i]) +=
1223 subblock(src,j).local_element(constrained_dofs[i]);
1228 for (
unsigned int j=0; j<n_blocks(dst); ++j)
1230 for (
unsigned int i=0; i<edge_constrained_indices[j].size(); ++i)
1232 subblock(const_cast<VectorType &>(src),j).local_element
1233 (edge_constrained_indices[j][i])
1234 = edge_constrained_values[j][i].first;
1235 subblock(dst,j).local_element(edge_constrained_indices[j][i])
1236 = edge_constrained_values[j][i].second +
1237 edge_constrained_values[j][i].first;
1244 template <
int dim,
typename VectorType>
1248 const VectorType &src)
const 1252 adjust_ghost_range_if_necessary(src,
false);
1253 adjust_ghost_range_if_necessary(dst,
true);
1257 if (!have_interface_matrices)
1262 for (
unsigned int j=0; j<n_blocks(dst); ++j)
1263 for (
unsigned int i=0; i<edge_constrained_indices[j].size(); ++i)
1265 edge_constrained_values[j][i] =
1266 std::pair<Number,Number>
1267 (subblock(src,j).local_element(edge_constrained_indices[j][i]),
1268 subblock(dst,j).local_element(edge_constrained_indices[j][i]));
1269 subblock(const_cast<VectorType &>(src),j)
1270 .local_element(edge_constrained_indices[j][i]) = 0.;
1275 for (
unsigned int j=0; j<n_blocks(dst); ++j)
1278 for (
unsigned int i=0; i<edge_constrained_indices[j].size(); ++i)
1280 for ( ; c<edge_constrained_indices[j][i]; ++c)
1281 subblock(dst,j).local_element(c) = 0.;
1285 subblock(const_cast<VectorType &>(src),j)
1286 .local_element(edge_constrained_indices[j][i])
1287 = edge_constrained_values[j][i].first;
1289 for ( ; c<dst.local_size(); ++c)
1290 subblock(dst,j).local_element(c) = 0.;
1296 template <
int dim,
typename VectorType>
1300 const VectorType &src)
const 1304 adjust_ghost_range_if_necessary(src,
false);
1305 adjust_ghost_range_if_necessary(dst,
true);
1309 if (!have_interface_matrices)
1312 VectorType src_cpy (src);
1313 for (
unsigned int j=0; j<n_blocks(dst); ++j)
1316 for (
unsigned int i=0; i<edge_constrained_indices[j].size(); ++i)
1318 for ( ; c<edge_constrained_indices[j][i]; ++c)
1319 subblock(src_cpy,j).local_element(c) = 0.;
1322 for ( ; c<subblock(src_cpy,j).local_size(); ++c)
1323 subblock(src_cpy,j).local_element(c) = 0.;
1326 apply_add(dst,src_cpy);
1328 for (
unsigned int j=0; j<n_blocks(dst); ++j)
1329 for (
unsigned int i=0; i<edge_constrained_indices[j].size(); ++i)
1330 subblock(dst,j).local_element(edge_constrained_indices[j][i]) = 0.;
1335 template <
int dim,
typename VectorType>
1338 const VectorType &src)
const 1342 Tvmult_add (dst,src);
1347 template <
int dim,
typename VectorType>
1351 return inverse_diagonal_entries.get() != NULL ? inverse_diagonal_entries->memory_consumption() :
sizeof(*this);
1356 template <
int dim,
typename VectorType>
1357 std_cxx11::shared_ptr<const MatrixFree<dim,typename Base<dim,VectorType>::value_type> >
1365 template <
int dim,
typename VectorType>
1366 const std_cxx11::shared_ptr<DiagonalMatrix<VectorType> > &
1369 Assert(inverse_diagonal_entries.get() != NULL &&
1371 return inverse_diagonal_entries;
1376 template <
int dim,
typename VectorType>
1379 const VectorType &src)
const 1386 template <
int dim,
typename VectorType>
1389 const VectorType &src,
1392 Assert(inverse_diagonal_entries.get() &&
1394 inverse_diagonal_entries->vmult(dst,src);
1402 template <
typename OperatorType>
1406 mf_base_operator(NULL)
1412 template <
typename OperatorType>
1416 mf_base_operator = NULL;
1421 template <
typename OperatorType>
1425 mf_base_operator = &operator_in;
1430 template <
typename OperatorType>
1431 template <
typename VectorType>
1434 const VectorType &src)
const 1436 #ifdef DEAL_II_WITH_CXX11 1437 #ifndef DEAL_II_MSVC 1438 static_assert (std::is_same<typename VectorType::value_type,value_type>::value,
1439 "The vector type must be based on the same value type as this" 1444 Assert(mf_base_operator != NULL,
1447 mf_base_operator->vmult_interface_down(dst, src);
1452 template <
typename OperatorType>
1453 template <
typename VectorType>
1456 const VectorType &src)
const 1458 #ifdef DEAL_II_WITH_CXX11 1459 #ifndef DEAL_II_MSVC 1460 static_assert (std::is_same<typename VectorType::value_type,value_type>::value,
1461 "The vector type must be based on the same value type as this" 1466 Assert(mf_base_operator != NULL,
1469 mf_base_operator->vmult_interface_up(dst, src);
1476 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
int n_components,
typename VectorType>
1480 Base<dim, VectorType>()
1485 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
int n_components,
typename VectorType>
1493 this->inverse_diagonal_entries.
1495 VectorType &inverse_diagonal_vector = this->inverse_diagonal_entries->get_vector();
1497 this->initialize_dof_vector(ones);
1498 this->initialize_dof_vector(inverse_diagonal_vector);
1500 apply_add(inverse_diagonal_vector, ones);
1502 this->set_constrained_entries_to_one(inverse_diagonal_vector);
1504 const unsigned int local_size = inverse_diagonal_vector.local_size();
1505 for (
unsigned int i=0; i<local_size; ++i)
1506 inverse_diagonal_vector.local_element(i)
1507 =1./inverse_diagonal_vector.local_element(i);
1509 inverse_diagonal_vector.update_ghost_values();
1514 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
int n_components,
typename VectorType>
1518 const VectorType &src)
const 1526 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
int n_components,
typename VectorType>
1531 const VectorType &src,
1532 const std::pair<unsigned int,unsigned int> &cell_range)
const 1536 for (
unsigned int cell=cell_range.first; cell<cell_range.second; ++cell)
1541 for (
unsigned int q=0; q<phi.
n_q_points; ++q)
1551 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
int n_components,
typename VectorType>
1555 Base<dim, VectorType>()
1561 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
int n_components,
typename VectorType>
1567 scalar_coefficient.reset();
1572 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
int n_components,
typename VectorType>
1577 scalar_coefficient = scalar_coefficient_;
1582 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
int n_components,
typename VectorType>
1583 std_cxx11::shared_ptr< Table<2, VectorizedArray< typename LaplaceOperator<dim, fe_degree, n_q_points_1d, n_components, VectorType>::value_type> > >
1587 Assert (scalar_coefficient.get(),
1589 return scalar_coefficient;
1594 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
int n_components,
typename VectorType>
1602 unsigned int dummy = 0;
1603 this->inverse_diagonal_entries.
1605 VectorType &inverse_diagonal_vector = this->inverse_diagonal_entries->get_vector();
1606 this->initialize_dof_vector(inverse_diagonal_vector);
1609 this, inverse_diagonal_vector, dummy);
1610 this->set_constrained_entries_to_one(inverse_diagonal_vector);
1612 for (
unsigned int i=0; i<inverse_diagonal_vector.local_size(); ++i)
1613 if (std::abs(inverse_diagonal_vector.local_element(i)) > std::sqrt(std::numeric_limits<Number>::epsilon()))
1614 inverse_diagonal_vector.local_element(i) = 1./inverse_diagonal_vector.local_element(i);
1616 inverse_diagonal_vector.local_element(i) = 1.;
1618 inverse_diagonal_vector.update_ghost_values();
1623 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
int n_components,
typename VectorType>
1627 const VectorType &src)
const 1635 template<
typename Number>
1639 for (
unsigned int v=0; v<VectorizedArray<Number>::n_array_elements; ++v)
1649 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
int n_components,
typename VectorType>
1653 const unsigned int cell)
const 1655 phi.evaluate (
false,
true,
false);
1656 if (scalar_coefficient.get())
1658 for (
unsigned int q=0; q<phi.n_q_points; ++q)
1660 Assert (non_negative((*scalar_coefficient)(cell,q)),
1661 ExcMessage(
"Coefficient must be non-negative"));
1662 phi.submit_gradient ((*scalar_coefficient)(cell,q)*phi.get_gradient(q), q);
1667 for (
unsigned int q=0; q<phi.n_q_points; ++q)
1669 phi.submit_gradient (phi.get_gradient(q), q);
1672 phi.integrate (
false,
true);
1678 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
int n_components,
typename VectorType>
1683 const VectorType &src,
1684 const std::pair<unsigned int,unsigned int> &cell_range)
const 1688 for (
unsigned int cell=cell_range.first; cell<cell_range.second; ++cell)
1692 do_operation_on_cell(phi,cell);
1698 template <
int dim,
int fe_degree,
int n_q_po
ints_1d,
int n_components,
typename VectorType>
1703 const unsigned int &,
1704 const std::pair<unsigned int,unsigned int> &cell_range)
const 1708 for (
unsigned int cell=cell_range.first; cell<cell_range.second; ++cell)
1717 do_operation_on_cell(phi,cell);
1720 for (
unsigned int i=0; i<phi.tensor_dofs_per_cell; ++i)
1730 DEAL_II_NAMESPACE_CLOSE
std_cxx11::shared_ptr< Table< 2, VectorizedArray< value_type > > > get_coefficient()
void initialize(std_cxx11::shared_ptr< const MatrixFree< dim, value_type > > 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 do_operation_on_cell(FEEvaluation< dim, fe_degree, n_q_points_1d, n_components, value_type > &phi, const unsigned int cell) const
const std_cxx11::shared_ptr< DiagonalMatrix< VectorType > > & get_matrix_diagonal_inverse() const
static const unsigned int invalid_unsigned_int
void vmult_interface_down(VectorType &dst, const VectorType &src) const
OperatorType::value_type value_type
std_cxx11::shared_ptr< DiagonalMatrix< VectorType > > inverse_diagonal_entries
std::vector< unsigned int > selected_rows
#define AssertDimension(dim1, dim2)
void reinit(const unsigned int cell)
virtual void apply_add(VectorType &dst, const VectorType &src) const
void local_apply_cell(const MatrixFree< dim, value_type > &data, VectorType &dst, const VectorType &src, const std::pair< unsigned int, unsigned int > &cell_range) const
void apply(const AlignedVector< VectorizedArray< Number > > &inverse_coefficient, const unsigned int n_actual_components, const VectorizedArray< Number > *in_array, VectorizedArray< Number > *out_array) const
const unsigned int dofs_per_cell
void evaluate(const bool evaluate_values, const bool evaluate_gradients, const bool evaluate_hessians=false)
void read_dof_values(const VectorType &src, const unsigned int first_index=0)
void integrate(const bool integrate_values, const bool integrate_gradients)
void precondition_Jacobi(VectorType &dst, const VectorType &src, const value_type omega) const
#define AssertIndexRange(index, range)
const FEEvaluationBase< dim, n_components, Number > & fe_eval
void mult_add(VectorType &dst, const VectorType &src, const bool transpose) const
static ::ExceptionBase & ExcNotInitialized()
void vmult_interface_up(VectorType &dst, const VectorType &src) const
#define AssertThrow(cond, exc)
void set_constrained_entries_to_one(VectorType &dst) const
std_cxx11::shared_ptr< const MatrixFree< dim, value_type > > get_matrix_free() const
const VectorizedArray< Number > * begin_dof_values() const
Base< dim, VectorType >::value_type value_type
virtual void compute_diagonal()
std::vector< unsigned int > selected_columns
value_type el(const unsigned int row, const unsigned int col) const
VectorType::value_type value_type
static ::ExceptionBase & ExcMessage(std::string arg1)
const unsigned int n_q_points
virtual void Tapply_add(VectorType &dst, const VectorType &src) const
void set_coefficient(const std_cxx11::shared_ptr< Table< 2, VectorizedArray< value_type > > > &scalar_coefficient)
void Tvmult(VectorType &dst, const VectorType &src) const
#define Assert(cond, exc)
void local_diagonal_cell(const MatrixFree< dim, value_type > &data, VectorType &dst, const unsigned int &, const std::pair< unsigned int, unsigned int > &cell_range) const
size_type index_within_set(const size_type global_index) const
virtual void compute_diagonal()=0
std_cxx11::shared_ptr< const MatrixFree< dim, value_type > > data
VectorType::size_type size_type
std::vector< std::vector< unsigned int > > edge_constrained_indices
AlignedVector< VectorizedArray< Number > > inverse_shape
void Tvmult_add(VectorType &dst, const VectorType &src) const
virtual void apply_add(VectorType &dst, const VectorType &src) const
Base< dim, VectorType >::size_type size_type
void vmult(VectorType &dst, const VectorType &src) const
void distribute_local_to_global(VectorType &dst, const unsigned int first_index=0) const
virtual void apply_add(VectorType &dst, const VectorType &src) const =0
void vmult(VectorType &dst, const VectorType &src) const
virtual void compute_diagonal()
void submit_value(const value_type val_in, const unsigned int q_point)
CellwiseInverseMassMatrix(const FEEvaluationBase< dim, n_components, Number > &fe_eval)
std::vector< std::vector< std::pair< value_type, value_type > > > edge_constrained_values
void local_apply_cell(const MatrixFree< dim, value_type > &data, VectorType &dst, const VectorType &src, const std::pair< unsigned int, unsigned int > &cell_range) const
virtual std::size_t memory_consumption() const
void fill_inverse_JxW_values(AlignedVector< VectorizedArray< Number > > &inverse_jxw) const
const internal::MatrixFreeFunctions::ShapeInfo< Number > & get_shape_info() const
static ::ExceptionBase & ExcNotImplemented()
std_cxx11::shared_ptr< Table< 2, VectorizedArray< value_type > > > scalar_coefficient
Base< dim, VectorType >::value_type value_type
bool is_element(const size_type index) const
void adjust_ghost_range_if_necessary(const VectorType &vec, const bool is_row) const
void vmult_add(VectorType &dst, const VectorType &src) const
void initialize(const OperatorType &operator_in)
SmartPointer< const OperatorType > mf_base_operator
value_type get_value(const unsigned int q_point) const
Base< dim, VectorType >::size_type size_type
void Tvmult(VectorType &dst, const VectorType &src) const
void initialize_dof_vector(VectorType &vec) const
T max(const T &t, const MPI_Comm &mpi_communicator)
bool have_interface_matrices
static ::ExceptionBase & ExcInternalError()