39template <
int dim,
typename RtreeType>
52 return (
pow >=
sizeof(
unsigned int) * 8) ? 0 :
65 std::vector<types::global_cell_index>
nodes;
66 std::vector<std::vector<types::global_cell_index>>
adjacency;
71 std::cout <<
"Graph information" << std::endl;
72 for (
const auto &node :
nodes)
74 std::cout <<
"Neighbors for node " << node << std::endl;
75 for (
const auto &neigh_node :
adjacency[node])
76 std::cout << neigh_node << std::endl;
95 template <
int dim,
int spacedim,
typename MatrixType>
102 using NumberType =
typename MatrixType::value_type;
104 static constexpr bool is_trilinos_matrix =
105 std::is_same_v<TrilinosWrappers::SparseMatrix, MatrixType>;
106 static constexpr bool is_serial_matrix =
107 std::is_same_v<SparseMatrix<NumberType>, MatrixType>;
108 static constexpr bool is_supported_matrix =
109 is_trilinos_matrix || is_serial_matrix;
110 static_assert(is_supported_matrix);
111 if constexpr (is_trilinos_matrix)
112 Assert(matrix.m() == 0, ExcInternalError());
113 if constexpr (is_serial_matrix)
116 "The destination matrix and its sparsity pattern must the empty "
117 "upon calling this function."));
132 const IndexSet &locally_owned_dofs_fine =
134 const IndexSet locally_relevant_dofs_fine =
137 const IndexSet &locally_owned_dofs_coarse =
140 std::conditional_t<is_trilinos_matrix,
145 if constexpr (is_trilinos_matrix)
146 dsp.
reinit(locally_owned_dofs_fine,
147 locally_owned_dofs_coarse,
150 dsp.reinit(fine_agglo_dh.
n_dofs(),
152 locally_relevant_dofs_fine);
155 std::vector<types::global_dof_index> agglo_dof_indices(dofs_per_cell);
156 std::vector<types::global_dof_index> standard_dof_indices(dofs_per_cell);
157 std::vector<types::global_dof_index> output_dof_indices(dofs_per_cell);
159 const std::vector<Point<dim>> &unit_support_points =
167 std::vector<types::global_dof_index> local_dof_indices_coarse(
169 std::vector<types::global_dof_index> local_dof_indices_child(dofs_per_cell);
172 if (polytope->is_locally_owned())
174 polytope->get_dof_indices(local_dof_indices_coarse);
177 const auto &children_polytopes = polytope->children();
182 child_dh->get_dof_indices(local_dof_indices_child);
183 for (
const auto row : local_dof_indices_child)
185 local_dof_indices_coarse.begin(),
186 local_dof_indices_coarse.end());
190 const auto assemble_injection_matrix = [&]() {
192 std::vector<Point<dim>> reference_q_points(dofs_per_cell);
199 if (polytope->is_locally_owned())
201 polytope->get_dof_indices(local_dof_indices_coarse);
203 coarse_bboxes[polytope->index()];
206 const auto &children_polytopes = polytope->children();
212 child_dh->get_dof_indices(local_dof_indices_child);
217 std::vector<Point<dim>> real_qpoints;
218 real_qpoints.reserve(unit_support_points.size());
219 for (
const Point<dim> &p : unit_support_points)
222 for (
unsigned int i = 0; i < local_dof_indices_coarse.size();
225 const auto &p = coarse_bbox.
real_to_unit(real_qpoints[i]);
226 for (
unsigned int j = 0; j < local_dof_indices_child.size();
234 local_dof_indices_child,
235 local_dof_indices_coarse,
242 if constexpr (is_trilinos_matrix)
246 assemble_injection_matrix();
248 else if constexpr (is_serial_matrix)
252 assemble_injection_matrix();
263 "This injection does not support PETSc types."));
286 template <
typename VectorType,
typename MatrixType,
typename SolverType>
295 static constexpr bool is_serial_matrix =
296 std::is_same_v<SolverType, SparseDirectUMFPACK>;
297 static constexpr bool is_trilinos_matrix =
298 std::is_same_v<SolverType, TrilinosWrappers::SolverDirect>;
299 [[maybe_unused]]
static constexpr bool is_matrix_type_supported =
300 is_serial_matrix || is_trilinos_matrix;
301 Assert(is_matrix_type_supported, ExcNotImplemented());
305 if constexpr (is_serial_matrix)
308 std::is_same_v<VectorType,
309 dealii::Vector<typename MatrixType::value_type>>;
312 else if constexpr (is_trilinos_matrix)
314 [[maybe_unused]]
static constexpr bool is_supported_parallel_vector =
315 std::is_same_v<VectorType,
317 typename MatrixType::value_type>> ||
318 std::is_same_v<VectorType, TrilinosWrappers::MPI::Vector>;
319 Assert(is_supported_parallel_vector, ExcNotImplemented());
337 const VectorType &src)
const override
379 typename number =
double>
442 *
data.get_dof_info(0).vector_partitioner))
445 src_copy.
reinit(
data.get_dof_info().vector_partitioner);
451 *
data.get_dof_info(0).vector_partitioner))
454 dst_copy.
reinit(
data.get_dof_info().vector_partitioner);
470 return data.get_vector_partitioner()->size();
476 return data.get_vector_partitioner()->size();
480 el(
const unsigned int row,
const unsigned int col)
const
485 ExcMessage(
"Matrix-free does not allow for entry access"));
493 data.initialize_dof_vector(vector);
500 return data.get_dof_handler();
507 return data.get_dof_handler().get_triangulation();
530 const auto cell_operation = [&](
auto &phi) {
532 for (
unsigned int q = 0; q < phi.n_q_points; ++q)
533 phi.submit_gradient(phi.get_gradient(q), q);
537 const auto face_operation = [&](
auto &phi_m,
auto &phi_p) {
543 (phi_m.normal_vector(0) * phi_m.inverse_jacobian(0))[dim - 1]) +
545 (phi_m.normal_vector(0) * phi_p.inverse_jacobian(0))[dim - 1])) *
548 for (
unsigned int q = 0; q < phi_m.n_q_points; ++q)
551 (phi_m.get_value(q) - phi_p.get_value(q)) * 0.5;
553 phi_m.get_normal_derivative(q) + phi_p.get_normal_derivative(q);
555 jump_value * 2. * sigmaF - average_valgrad * 0.5;
556 phi_m.submit_normal_derivative(-jump_value, q);
557 phi_p.submit_normal_derivative(-jump_value, q);
558 phi_m.submit_value(average_valgrad, q);
559 phi_p.submit_value(-average_valgrad, q);
565 const auto boundary_operation = [&](
auto &phi_m) {
569 (phi_m.normal_vector(0) * phi_m.inverse_jacobian(0))[dim - 1]) *
572 for (
unsigned int q = 0; q < phi_m.n_q_points; ++q)
576 -phi_m.get_normal_derivative(q);
577 average_valgrad += jump_value * sigmaF * 2.0;
578 phi_m.submit_normal_derivative(-jump_value, q);
579 phi_m.submit_value(average_valgrad, q);
593 const auto &dof_handler =
data.get_dof_handler();
597 dof_handler.locally_owned_mg_dofs(
data.get_mg_level()) :
598 dof_handler.locally_owned_dofs(),
599 data.get_task_info().communicator);
640 const auto cell_operation = [&](
auto &phi) {
642 for (
unsigned int q = 0; q < phi.n_q_points; ++q)
643 phi.submit_gradient(phi.get_gradient(q), q);
647 const auto face_operation = [&](
auto &phi_m,
auto &phi_p) {
653 (phi_m.normal_vector(0) * phi_m.inverse_jacobian(0))[dim - 1]) +
655 (phi_m.normal_vector(0) * phi_p.inverse_jacobian(0))[dim - 1])) *
658 for (
unsigned int q = 0; q < phi_m.n_q_points; ++q)
661 (phi_m.get_value(q) - phi_p.get_value(q)) * 0.5;
663 phi_m.get_normal_derivative(q) + phi_p.get_normal_derivative(q);
665 jump_value * 2. * sigmaF - average_valgrad * 0.5;
666 phi_m.submit_normal_derivative(-jump_value, q);
667 phi_p.submit_normal_derivative(-jump_value, q);
668 phi_m.submit_value(average_valgrad, q);
669 phi_p.submit_value(-average_valgrad, q);
675 const auto boundary_operation = [&](
auto &phi_m) {
679 (phi_m.normal_vector(0) * phi_m.inverse_jacobian(0))[dim - 1]) *
682 for (
unsigned int q = 0; q < phi_m.n_q_points; ++q)
686 -phi_m.get_normal_derivative(q);
687 average_valgrad += jump_value * sigmaF * 2.0;
688 phi_m.submit_normal_derivative(-jump_value, q);
689 phi_m.submit_value(average_valgrad, q);
706 const IndexSet system_relevant_set =
712 system_relevant_set);
718 data.get_task_info().communicator,
719 system_relevant_set);
720 mg_matrix.
reinit(system_partitioning,
722 data.get_task_info().communicator);
748 .template cell_loop<LinearAlgebra::distributed::Vector<number>,
int>(
749 [](
const auto &matrix_free,
753 FEEvaluation<dim, -1, 0, n_components, number> phi(matrix_free,
755 for (
unsigned int cell = cells.first; cell < cells.second; ++cell)
758 for (
unsigned int q = 0; q < phi.n_q_points; ++q)
760 if constexpr (n_components == 1)
762 phi.submit_value(1.0, q);
767 for (
unsigned int v = 0;
768 v < VectorizedArray<number>::size();
771 for (
unsigned int i = 0; i < n_components; i++)
774 phi.submit_value(temp, q);
792 unsigned int dummy = 0;
800 for (
unsigned int i = 0;
814 const std::pair<unsigned int, unsigned int> &cell_range)
const
818 for (
unsigned int cell = cell_range.first; cell < cell_range.second;
822 phi.read_dof_values(src);
824 for (
unsigned int q = 0; q < phi.n_q_points; ++q)
825 phi.submit_gradient(phi.get_gradient(q), q);
827 phi.distribute_local_to_global(dst);
836 const std::pair<unsigned int, unsigned int> &face_range)
const
841 for (
unsigned int face = face_range.first; face < face_range.second;
844 fe_eval.reinit(face);
845 fe_eval_neighbor.reinit(face);
847 fe_eval.read_dof_values(src);
850 fe_eval_neighbor.read_dof_values(src);
854 (
std::abs((fe_eval.normal_vector(0) *
855 fe_eval.inverse_jacobian(0))[dim - 1]) +
856 std::abs((fe_eval.normal_vector(0) *
857 fe_eval_neighbor.inverse_jacobian(0))[dim - 1])) *
860 for (
unsigned int q = 0; q < fe_eval.n_q_points; ++q)
863 (fe_eval.get_value(q) - fe_eval_neighbor.get_value(q)) * 0.5;
865 fe_eval.get_normal_derivative(q) +
866 fe_eval_neighbor.get_normal_derivative(q);
868 jump_value * 2. * sigmaF - average_valgrad * 0.5;
869 fe_eval.submit_normal_derivative(-jump_value, q);
870 fe_eval_neighbor.submit_normal_derivative(-jump_value, q);
871 fe_eval.submit_value(average_valgrad, q);
872 fe_eval_neighbor.submit_value(-average_valgrad, q);
876 fe_eval.distribute_local_to_global(dst);
879 fe_eval_neighbor.distribute_local_to_global(dst);
888 const std::pair<unsigned int, unsigned int> &face_range)
const
891 for (
unsigned int face = face_range.first; face < face_range.second;
894 fe_eval.reinit(face);
895 fe_eval.read_dof_values(src);
899 std::abs((fe_eval.normal_vector(0) *
900 fe_eval.inverse_jacobian(0))[dim - 1]) *
903 for (
unsigned int q = 0; q < fe_eval.n_q_points; ++q)
907 -fe_eval.get_normal_derivative(q);
908 average_valgrad += jump_value * sigmaF * 2.0;
909 fe_eval.submit_normal_derivative(-jump_value, q);
910 fe_eval.submit_value(average_valgrad, q);
915 fe_eval.distribute_local_to_global(dst);
924 const unsigned int &,
925 const std::pair<unsigned int, unsigned int> &cell_range)
const
931 for (
unsigned int cell = cell_range.first; cell < cell_range.second;
936 for (
unsigned int i = 0; i < phi.dofs_per_cell; ++i)
938 for (
unsigned int j = 0; j < phi.dofs_per_cell; ++j)
940 phi.begin_dof_values()[i] = 1.;
942 for (
unsigned int q = 0; q < phi.n_q_points; ++q)
943 phi.submit_gradient(phi.get_gradient(q), q);
945 local_diagonal_vector[i] = phi.begin_dof_values()[i];
947 for (
unsigned int i = 0; i < phi.dofs_per_cell; ++i)
948 phi.begin_dof_values()[i] = local_diagonal_vector[i];
949 phi.distribute_local_to_global(dst);
957 const unsigned int &,
958 const std::pair<unsigned int, unsigned int> &face_range)
const
965 for (
unsigned int face = face_range.first; face < face_range.second;
969 phi_outer.reinit(face);
973 (phi.normal_vector(0) * phi.inverse_jacobian(0))[dim - 1]) +
975 phi_outer.inverse_jacobian(0))[dim - 1])) *
979 for (
unsigned int j = 0; j < phi.dofs_per_cell; ++j)
983 for (
unsigned int i = 0; i < phi.dofs_per_cell; ++i)
985 for (
unsigned int j = 0; j < phi.dofs_per_cell; ++j)
987 phi.begin_dof_values()[i] = 1.;
991 for (
unsigned int q = 0; q < phi.n_q_points; ++q)
994 (phi.get_value(q) - phi_outer.get_value(q)) * 0.5;
996 phi.get_normal_derivative(q) +
997 phi_outer.get_normal_derivative(q);
999 jump_value * 2. * sigmaF - average_valgrad * 0.5;
1000 phi.submit_normal_derivative(-jump_value, q);
1001 phi.submit_value(average_valgrad, q);
1005 local_diagonal_vector[i] = phi.begin_dof_values()[i];
1007 for (
unsigned int i = 0; i < phi.dofs_per_cell; ++i)
1008 phi.begin_dof_values()[i] = local_diagonal_vector[i];
1009 phi.distribute_local_to_global(dst);
1012 for (
unsigned int j = 0; j < phi.dofs_per_cell; ++j)
1015 for (
unsigned int i = 0; i < phi.dofs_per_cell; ++i)
1017 for (
unsigned int j = 0; j < phi.dofs_per_cell; ++j)
1019 phi_outer.begin_dof_values()[i] = 1.;
1023 for (
unsigned int q = 0; q < phi.n_q_points; ++q)
1026 (phi.get_value(q) - phi_outer.get_value(q)) * 0.5;
1028 phi.get_normal_derivative(q) +
1029 phi_outer.get_normal_derivative(q);
1031 jump_value * 2. * sigmaF - average_valgrad * 0.5;
1032 phi_outer.submit_normal_derivative(-jump_value, q);
1033 phi_outer.submit_value(-average_valgrad, q);
1037 local_diagonal_vector[i] = phi_outer.begin_dof_values()[i];
1039 for (
unsigned int i = 0; i < phi.dofs_per_cell; ++i)
1040 phi_outer.begin_dof_values()[i] = local_diagonal_vector[i];
1041 phi_outer.distribute_local_to_global(dst);
1049 const unsigned int &,
1050 const std::pair<unsigned int, unsigned int> &face_range)
const
1056 for (
unsigned int face = face_range.first; face < face_range.second;
1063 (phi.normal_vector(0) * phi.inverse_jacobian(0))[dim - 1]) *
1066 for (
unsigned int i = 0; i < phi.dofs_per_cell; ++i)
1068 for (
unsigned int j = 0; j < phi.dofs_per_cell; ++j)
1070 phi.begin_dof_values()[i] = 1.;
1074 for (
unsigned int q = 0; q < phi.n_q_points; ++q)
1078 -phi.get_normal_derivative(q);
1079 average_valgrad += jump_value * sigmaF * 2.0;
1080 phi.submit_normal_derivative(-jump_value, q);
1081 phi.submit_value(average_valgrad, q);
1086 local_diagonal_vector[i] = phi.begin_dof_values()[i];
1088 for (
unsigned int i = 0; i < phi.dofs_per_cell; ++i)
1089 phi.begin_dof_values()[i] = local_diagonal_vector[i];
1090 phi.distribute_local_to_global(dst);
1112 typename number =
double>
1133 "Degree of the operator must match the degree of the DoFHandler"));
1183 *
data.get_dof_info(0).vector_partitioner))
1186 src_copy.
reinit(
data.get_dof_info().vector_partitioner);
1192 *
data.get_dof_info(0).vector_partitioner))
1195 dst_copy.
reinit(
data.get_dof_info().vector_partitioner);
1214 return data.get_vector_partitioner()->size();
1220 return data.get_vector_partitioner()->size();
1224 el(
const unsigned int row,
const unsigned int col)
const
1229 ExcMessage(
"Matrix-free does not allow for entry access"));
1237 data.initialize_dof_vector(vector);
1244 return data.get_dof_handler();
1251 return data.get_dof_handler().get_triangulation();
1276 const auto cell_operation = [&](
auto &phi) {
1278 for (
unsigned int q = 0; q < phi.n_q_points; ++q)
1280 phi.submit_value(phi.get_value(q) * factor, q);
1281 phi.submit_gradient(phi.get_gradient(q) *
parameters.sigma, q);
1286 const auto face_operation = [&](
auto &phi_m,
auto &phi_p) {
1292 (phi_m.normal_vector(0) * phi_m.inverse_jacobian(0))[dim - 1]) +
1294 (phi_m.normal_vector(0) * phi_p.inverse_jacobian(0))[dim - 1])) *
1297 for (
unsigned int q = 0; q < phi_m.n_q_points; ++q)
1300 (phi_m.get_value(q) - phi_p.get_value(q)) * 0.5;
1302 phi_m.get_normal_derivative(q) + phi_p.get_normal_derivative(q);
1304 jump_value * 2. * sigmaF - average_valgrad * 0.5;
1305 phi_m.submit_normal_derivative(-
parameters.sigma * jump_value, q);
1306 phi_p.submit_normal_derivative(-
parameters.sigma * jump_value, q);
1307 phi_m.submit_value(
parameters.sigma * average_valgrad, q);
1308 phi_p.submit_value(-
parameters.sigma * average_valgrad, q);
1331 const auto &dof_handler =
data.get_dof_handler();
1335 dof_handler.locally_owned_mg_dofs(
data.get_mg_level()) :
1336 dof_handler.locally_owned_dofs(),
1337 data.get_task_info().communicator);
1342 data.get_mg_level(),
1380 const auto cell_operation = [&](
auto &phi) {
1382 for (
unsigned int q = 0; q < phi.n_q_points; ++q)
1384 phi.submit_value(phi.get_value(q) * factor, q);
1385 phi.submit_gradient(phi.get_gradient(q) *
parameters.sigma, q);
1390 const auto face_operation = [&](
auto &phi_m,
auto &phi_p) {
1396 (phi_m.normal_vector(0) * phi_m.inverse_jacobian(0))[dim - 1]) +
1398 (phi_m.normal_vector(0) * phi_p.inverse_jacobian(0))[dim - 1])) *
1401 for (
unsigned int q = 0; q < phi_m.n_q_points; ++q)
1404 (phi_m.get_value(q) - phi_p.get_value(q)) * 0.5;
1406 phi_m.get_normal_derivative(q) + phi_p.get_normal_derivative(q);
1408 jump_value * 2. * sigmaF - average_valgrad * 0.5;
1409 phi_m.submit_normal_derivative(-
parameters.sigma * jump_value, q);
1410 phi_p.submit_normal_derivative(-
parameters.sigma * jump_value, q);
1411 phi_m.submit_value(
parameters.sigma * average_valgrad, q);
1412 phi_p.submit_value(-
parameters.sigma * average_valgrad, q);
1423 ExcInternalError());
1429 const IndexSet system_relevant_set =
1435 system_relevant_set);
1440 system_partitioning,
1441 data.get_task_info().communicator,
1442 system_relevant_set);
1443 mg_matrix.
reinit(system_partitioning,
1445 data.get_task_info().communicator);
1467 unsigned int dummy = 0;
1475 for (
unsigned int i = 0;
1493 for (
unsigned int cell = 0; cell <
data.n_cell_batches(); ++cell)
1498 for (
unsigned int q = 0; q < phi.
n_q_points; ++q)
1515 for (
unsigned int cell = 0; cell <
data.n_cell_batches(); ++cell)
1520 for (
unsigned int q = 0; q < phi.
n_q_points; ++q)
1527 for (
unsigned int v = 0; v < VectorizedArray<number>::size();
1531 for (
unsigned int d = 0;
d < dim; ++
d)
1532 p[
d] = p_vect[
d][v];
1533 applied_current_value[v] = Iext.
value(p);
1539 applied_current_value,
1555 const std::pair<unsigned int, unsigned int> &cell_range)
const
1560 for (
unsigned int cell = cell_range.first; cell < cell_range.second;
1566 for (
unsigned int q = 0; q < phi.
n_q_points; ++q)
1581 const std::pair<unsigned int, unsigned int> &face_range)
const
1586 fe_eval_neighbor(
data,
false);
1588 for (
unsigned int face = face_range.first; face < face_range.second;
1592 fe_eval_neighbor.
reinit(face);
1607 for (
unsigned int q = 0; q < fe_eval.
n_q_points; ++q)
1615 jump_value * 2. * sigmaF - average_valgrad * 0.5;
1640 const std::pair<unsigned int, unsigned int> &face_range)
const
1654 const unsigned int &,
1655 const std::pair<unsigned int, unsigned int> &cell_range)
const
1662 for (
unsigned int cell = cell_range.first; cell < cell_range.second;
1674 for (
unsigned int q = 0; q < phi.
n_q_points; ++q)
1694 const unsigned int &,
1695 const std::pair<unsigned int, unsigned int> &face_range)
const
1700 phi_outer(
data,
false);
1704 for (
unsigned int face = face_range.first; face < face_range.second;
1730 for (
unsigned int q = 0; q < phi.
n_q_points; ++q)
1738 jump_value * 2. * sigmaF - average_valgrad * 0.5;
1763 for (
unsigned int q = 0; q < phi.
n_q_points; ++q)
1771 jump_value * 2. * sigmaF - average_valgrad * 0.5;
1792 const unsigned int &,
1793 const std::pair<unsigned int, unsigned int> &face_range)
const
1824 return std::distance(v.begin(), std::find(v.begin(), v.end(), index));
1835 Assert(cells.size() > 0, ExcMessage(
"No cells to be agglomerated."));
1836 const unsigned int n_faces = cells[0]->n_faces();
1838 std::vector<types::global_cell_index> vec_cells(cells.size());
1839 for (
size_t i = 0; i < cells.size(); i++)
1840 vec_cells[i] = cells[i]->active_cell_index();
1843 for (
const auto &cell : cells)
1849 g.
nodes.push_back(
get_index(vec_cells, cell->active_cell_index()));
1850 for (
unsigned int f = 0; f < n_faces; ++f)
1852 const auto &neigh = cell->neighbor(f);
1853 if (neigh.state() == IteratorState::IteratorStates::valid &&
1854 std::find(cells.begin(), cells.end(), neigh) != std::end(cells))
1856 .push_back(
get_index(vec_cells, neigh->active_cell_index()));
1864 dfs(std::vector<types::global_cell_index> &comp,
1865 std::vector<bool> &visited,
1874 dfs(comp, visited, g, u);
1883 std::vector<std::vector<types::global_cell_index>> &connected_components)
1885 Assert(g.
nodes.size() > 0, ExcMessage(
"No nodes in this graph."));
1887 connected_components.size() == 0,
1889 "Connected components have to be computed by the present function."));
1891 std::vector<bool> visited(g.
nodes.size());
1892 std::fill(visited.begin(), visited.end(), 0);
1899 connected_components.emplace_back();
1900 dfs(connected_components.back(), visited, g, v);
void distribute_local_to_global(const InVector &local_vector, const std::vector< size_type > &local_dof_indices, OutVector &global_vector) const
Point< spacedim, Number > real_to_unit(const Point< spacedim, Number > &point) const
Point< spacedim, Number > unit_to_real(const Point< spacedim, Number > &point) const
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
const IndexSet & locally_owned_dofs() const
types::global_dof_index n_dofs() const
void reinit(const size_type m, const size_type n, const IndexSet &rowset=IndexSet())
gradient_type get_gradient(const unsigned int q_point) const
void submit_gradient(const gradient_type grad_in, const unsigned int q_point)
void submit_value(const value_type val_in, const unsigned int q_point)
void distribute_local_to_global(VectorType &dst, const unsigned int first_index=0, const std::bitset< n_lanes > &mask=std::bitset< n_lanes >().flip()) const
void submit_normal_derivative(const value_type grad_in, const unsigned int q_point)
value_type get_normal_derivative(const unsigned int q_point) const
value_type get_value(const unsigned int q_point) const
void read_dof_values(const VectorType &src, const unsigned int first_index=0, const std::bitset< n_lanes > &mask=std::bitset< n_lanes >().flip())
Tensor< 2, dim, VectorizedArrayType > inverse_jacobian(const unsigned int q_point) const
Tensor< 1, dim, VectorizedArrayType > normal_vector(const unsigned int q_point) const
Point< dim, VectorizedArrayType > quadrature_point(const unsigned int q) const
const VectorizedArrayType * begin_dof_values() const
const unsigned int n_q_points
void reinit(const unsigned int cell_batch_index)
void evaluate(const EvaluationFlags::EvaluationFlags evaluation_flag)
const unsigned int dofs_per_cell
void integrate(const EvaluationFlags::EvaluationFlags integration_flag)
void reinit(const unsigned int face_batch_number)
const unsigned int n_q_points
const unsigned int dofs_per_cell
void evaluate(const EvaluationFlags::EvaluationFlags evaluation_flag)
void integrate(const EvaluationFlags::EvaluationFlags integration_flag, const bool sum_into_values=false)
const unsigned int degree
const unsigned int dofs_per_cell
const std::vector< Point< dim > > & get_unit_support_points() const
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
virtual RangeNumberType value(const Point< dim > &p, const unsigned int component=0) const
void zero_out_ghost_values() const
void swap(Vector< Number, MemorySpace > &v) noexcept
void compress(VectorOperation::values operation)
bool partitioners_are_globally_compatible(const Utilities::MPI::Partitioner &part) const
void reinit(const size_type size, const bool omit_zeroing_entries=false)
void copy_from(const size_type n_rows, const size_type n_cols, const ForwardIterator begin, const ForwardIterator end)
virtual MPI_Comm get_mpi_communicator() const
void reinit(const SparsityPatternType &sparsity_pattern)
const DoFHandler< dim, spacedim >::active_cell_iterator polytope_to_dh_iterator(const types::global_cell_index polytope_index) const
types::global_dof_index n_dofs() const noexcept
const Triangulation< dim, spacedim > & get_triangulation() const
const Mapping< dim > & get_mapping() const
const std::vector< BoundingBox< dim > > & get_local_bboxes() const
const FiniteElement< dim, spacedim > & get_fe() const
IteratorRange< agglomeration_iterator > polytope_iterators() const
DoFHandler< dim, spacedim > agglo_dh
~MGCoarseDirect()=default
const MatrixType & coarse_matrix
MGCoarseDirect(const MatrixType &matrix)
void operator()(const unsigned int, VectorType &dst, const VectorType &src) const override
void Tvmult(LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src) const
LinearAlgebra::distributed::Vector< number > inverse_diagonal_entries
AffineConstraints< number > constraints
void rhs(LinearAlgebra::distributed::Vector< number > &b) const
void vmult_add(LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src) const
void local_diagonal_face(const MatrixFree< dim, number > &data, LinearAlgebra::distributed::Vector< number > &dst, const unsigned int &, const std::pair< unsigned int, unsigned int > &face_range) const
void local_diagonal_boundary(const MatrixFree< dim, number > &data, LinearAlgebra::distributed::Vector< number > &dst, const unsigned int &, const std::pair< unsigned int, unsigned int > &face_range) const
TrilinosWrappers::SparseMatrix system_matrix
number el(const unsigned int row, const unsigned int col) const
void initialize_dof_vector(LinearAlgebra::distributed::Vector< number > &vector) const
void reinit(const Mapping< dim > &mapping, const DoFHandler< dim > &dof_handler, const unsigned int level=numbers::invalid_unsigned_int)
VectorizedArray< number > VectorizedArrayType
const Triangulation< dim > & get_triangulation() const
MatrixFree< dim, number > data
const LinearAlgebra::distributed::Vector< number > & get_matrix_diagonal_inverse() const
const DoFHandler< dim > & get_dof_handler() const
const MatrixFree< dim, number > * get_matrix_free() const
void local_diagonal_cell(const MatrixFree< dim, number > &data, LinearAlgebra::distributed::Vector< number > &dst, const unsigned int &, const std::pair< unsigned int, unsigned int > &cell_range) const
void vmult(LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src) const
void get_system_matrix(TrilinosWrappers::SparseMatrix &mg_matrix) const
void compute_inverse_diagonal()
void local_apply_boundary(const MatrixFree< dim, number > &data, LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src, const std::pair< unsigned int, unsigned int > &face_range) const
types::global_dof_index n() const
const TrilinosWrappers::SparseMatrix & get_system_matrix() const
void Tvmult_add(LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src) const
LinearAlgebra::distributed::Vector< number > VectorType
types::global_dof_index m() const
void local_apply(const MatrixFree< dim, number > &data, LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src, const std::pair< unsigned int, unsigned int > &cell_range) const
void local_apply_face(const MatrixFree< dim, number > &data, LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src, const std::pair< unsigned int, unsigned int > &face_range) const
Physics::BilinearFormParameters parameters
void compute_inverse_diagonal()
void rhs(LinearAlgebra::distributed::Vector< number > &rhs, const LinearAlgebra::distributed::Vector< number > &solution_minus_ion, const Function< dim > &Iext) const
MonodomainOperatorDG(const Physics::BilinearFormParameters ¶meters_)
const DoFHandler< dim > & get_dof_handler() const
void local_apply_face(const MatrixFree< dim, number > &data, LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src, const std::pair< unsigned int, unsigned int > &face_range) const
void local_apply(const MatrixFree< dim, number > &data, LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src, const std::pair< unsigned int, unsigned int > &cell_range) const
const LinearAlgebra::distributed::Vector< number > & get_matrix_diagonal_inverse() const
void local_diagonal_face(const MatrixFree< dim, number > &data, LinearAlgebra::distributed::Vector< number > &dst, const unsigned int &, const std::pair< unsigned int, unsigned int > &face_range) const
const MatrixFree< dim, number > * get_matrix_free() const
MatrixFree< dim, number > data
void reinit(const Mapping< dim > &mapping, const DoFHandler< dim > &dof_handler, const unsigned int level=numbers::invalid_unsigned_int)
TrilinosWrappers::SparseMatrix system_matrix
LinearAlgebra::distributed::Vector< number > inverse_diagonal_entries
void local_diagonal_cell(const MatrixFree< dim, number > &data, LinearAlgebra::distributed::Vector< number > &dst, const unsigned int &, const std::pair< unsigned int, unsigned int > &cell_range) const
VectorizedArray< number > VectorizedArrayType
AffineConstraints< number > constraints
void Tvmult_add(LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src) const
void apply_mass_term(LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src) const
const TrilinosWrappers::SparseMatrix & get_system_matrix() const
number el(const unsigned int row, const unsigned int col) const
LinearAlgebra::distributed::Vector< number > VectorType
types::global_dof_index m() const
void local_apply_boundary(const MatrixFree< dim, number > &data, LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src, const std::pair< unsigned int, unsigned int > &face_range) const
void get_system_matrix(TrilinosWrappers::SparseMatrix &mg_matrix) const
void local_diagonal_boundary(const MatrixFree< dim, number > &data, LinearAlgebra::distributed::Vector< number > &dst, const unsigned int &, const std::pair< unsigned int, unsigned int > &face_range) const
void vmult_add(LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src) const
void vmult(LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src) const
types::global_dof_index n() const
void initialize_dof_vector(LinearAlgebra::distributed::Vector< number > &vector) const
void Tvmult(LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src) const
const Triangulation< dim > & get_triangulation() const
EnableObserverPointer Subscriptor
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertThrow(cond, exc)
void make_flux_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity_pattern)
Expression pow(const Expression &base, const Expression &exponent)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void create_graph_from_agglomerate(const std::vector< typename Triangulation< dim >::active_cell_iterator > &cells, Graph &g)
void dfs(std::vector< types::global_cell_index > &comp, std::vector< bool > &visited, const Graph &g, const types::global_cell_index v)
constexpr T constexpr_pow(T num, unsigned int pow)
void fill_injection_matrix(const AgglomerationHandler< dim, spacedim > &coarse_ah, const AgglomerationHandler< dim, spacedim > &fine_ah, SparsityPattern &sp, MatrixType &matrix)
void compute_connected_components(Graph &g, std::vector< std::vector< types::global_cell_index > > &connected_components)
types::global_cell_index get_index(const std::vector< types::global_cell_index > &v, const types::global_cell_index index)
constexpr unsigned int invalid_unsigned_int
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
unsigned int global_dof_index
unsigned int global_cell_index
void swap(ObserverPointer< T, P > &t1, ObserverPointer< T, Q > &t2)
TasksParallelScheme tasks_parallel_scheme
UpdateFlags mapping_update_flags_inner_faces
UpdateFlags mapping_update_flags_boundary_faces
UpdateFlags mapping_update_flags
unsigned int tasks_block_size
std::vector< std::vector< types::global_cell_index > > adjacency
std::vector< types::global_cell_index > nodes