PolyDEAL
 
Loading...
Searching...
No Matches
utils.h
Go to the documentation of this file.
1// -----------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception OR LGPL-2.1-or-later
4// Copyright (C) XXXX - YYYY by the polyDEAL authors
5//
6// This file is part of the polyDEAL library.
7//
8// Detailed license information governing the source code
9// can be found in LICENSE.md at the top level directory.
10//
11// -----------------------------------------------------------------------------
12
13
14#ifndef utils_h
15#define utils_h
16
17#include <deal.II/base/point.h>
18
20
28
32
34
35#include <multigrid_amg.h>
36
37
38
39template <int dim, typename RtreeType>
41template <int, int>
43
44
45
46namespace Utils
47{
48 template <typename T>
49 inline constexpr T
50 constexpr_pow(T num, unsigned int pow)
51 {
52 return (pow >= sizeof(unsigned int) * 8) ? 0 :
53 pow == 0 ? 1 :
54 num * constexpr_pow(num, pow - 1);
55 }
56
57
58
63 struct Graph
64 {
65 std::vector<types::global_cell_index> nodes;
66 std::vector<std::vector<types::global_cell_index>> adjacency;
67
68 void
70 {
71 std::cout << "Graph information" << std::endl;
72 for (const auto &node : nodes)
73 {
74 std::cout << "Neighbors for node " << node << std::endl;
75 for (const auto &neigh_node : adjacency[node])
76 std::cout << neigh_node << std::endl;
77 }
78 }
79 };
80
81
82
95 template <int dim, int spacedim, typename MatrixType>
96 void
100 MatrixType &matrix)
101 {
102 using NumberType = typename MatrixType::value_type;
103 // First, check that we support the matrix types
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)
114 Assert(sp.empty() && matrix.empty(),
115 ExcMessage(
116 "The destination matrix and its sparsity pattern must the empty "
117 "upon calling this function."));
118
119 Assert(coarse_ah.n_dofs() < fine_ah.n_dofs(), ExcInternalError());
120 AssertDimension(dim, spacedim);
121
122 // Get information from the handlers
123 const DoFHandler<dim, spacedim> &coarse_agglo_dh = coarse_ah.agglo_dh;
124 const DoFHandler<dim, spacedim> &fine_agglo_dh = fine_ah.agglo_dh;
125
126 const Mapping<dim, spacedim> &mapping = fine_ah.get_mapping();
127 const FiniteElement<dim, spacedim> &fe = coarse_ah.get_fe();
128 const Triangulation<dim, spacedim> &tria = coarse_ah.get_triangulation();
129 const auto &fine_bboxes = fine_ah.get_local_bboxes();
130 const auto &coarse_bboxes = coarse_ah.get_local_bboxes();
131
132 const IndexSet &locally_owned_dofs_fine =
133 fine_agglo_dh.locally_owned_dofs();
134 const IndexSet locally_relevant_dofs_fine =
136
137 const IndexSet &locally_owned_dofs_coarse =
138 coarse_agglo_dh.locally_owned_dofs();
139
140 std::conditional_t<is_trilinos_matrix,
143 dsp;
144
145 if constexpr (is_trilinos_matrix)
146 dsp.reinit(locally_owned_dofs_fine,
147 locally_owned_dofs_coarse,
148 tria.get_mpi_communicator());
149 else
150 dsp.reinit(fine_agglo_dh.n_dofs(),
151 coarse_agglo_dh.n_dofs(),
152 locally_relevant_dofs_fine);
153
154 const unsigned int dofs_per_cell = fe.dofs_per_cell;
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);
158
159 const std::vector<Point<dim>> &unit_support_points =
161 Quadrature<dim> quad(unit_support_points);
162 FEValues<dim, spacedim> output_fe_values(mapping,
163 fe,
164 quad,
166
167 std::vector<types::global_dof_index> local_dof_indices_coarse(
168 dofs_per_cell);
169 std::vector<types::global_dof_index> local_dof_indices_child(dofs_per_cell);
170
171 for (const auto &polytope : coarse_ah.polytope_iterators())
172 if (polytope->is_locally_owned())
173 {
174 polytope->get_dof_indices(local_dof_indices_coarse);
175
176 // Get local children and their DoFs
177 const auto &children_polytopes = polytope->children();
178 for (const types::global_cell_index child_idx : children_polytopes)
179 {
180 const typename DoFHandler<dim>::active_cell_iterator &child_dh =
181 fine_ah.polytope_to_dh_iterator(child_idx);
182 child_dh->get_dof_indices(local_dof_indices_child);
183 for (const auto row : local_dof_indices_child)
184 dsp.add_entries(row,
185 local_dof_indices_coarse.begin(),
186 local_dof_indices_coarse.end());
187 }
188 }
189
190 const auto assemble_injection_matrix = [&]() {
191 FullMatrix<NumberType> local_matrix(dofs_per_cell, dofs_per_cell);
192 std::vector<Point<dim>> reference_q_points(dofs_per_cell);
193
194 // Dummy AffineConstraints, only needed for loc2glb
196 c.close();
197
198 for (const auto &polytope : coarse_ah.polytope_iterators())
199 if (polytope->is_locally_owned())
200 {
201 polytope->get_dof_indices(local_dof_indices_coarse);
202 const BoundingBox<dim> &coarse_bbox =
203 coarse_bboxes[polytope->index()];
204
205 // Get local children of the present polytope
206 const auto &children_polytopes = polytope->children();
207 for (const types::global_cell_index child_idx : children_polytopes)
208 {
209 const BoundingBox<dim> &fine_bbox = fine_bboxes[child_idx];
210 const typename DoFHandler<dim>::active_cell_iterator &child_dh =
211 fine_ah.polytope_to_dh_iterator(child_idx);
212 child_dh->get_dof_indices(local_dof_indices_child);
213
214 local_matrix = 0.;
215
216 // compute real location of support points
217 std::vector<Point<dim>> real_qpoints;
218 real_qpoints.reserve(unit_support_points.size());
219 for (const Point<dim> &p : unit_support_points)
220 real_qpoints.push_back(fine_bbox.unit_to_real(p));
221
222 for (unsigned int i = 0; i < local_dof_indices_coarse.size();
223 ++i)
224 {
225 const auto &p = coarse_bbox.real_to_unit(real_qpoints[i]);
226 for (unsigned int j = 0; j < local_dof_indices_child.size();
227 ++j)
228 {
229 local_matrix(i, j) = fe.shape_value(j, p);
230 }
231 }
232
233 c.distribute_local_to_global(local_matrix,
234 local_dof_indices_child,
235 local_dof_indices_coarse,
236 matrix);
237 }
238 }
239 };
240
241
242 if constexpr (is_trilinos_matrix)
243 {
244 dsp.compress();
245 matrix.reinit(dsp);
246 assemble_injection_matrix();
247 }
248 else if constexpr (is_serial_matrix)
249 {
250 sp.copy_from(dsp);
251 matrix.reinit(sp);
252 assemble_injection_matrix();
253 }
254 else
255 {
256 // PETSc types not implemented.
257 (void)coarse_ah;
258 (void)fine_ah;
259 (void)sp;
260 (void)matrix;
261 AssertThrow(false,
262 ExcNotImplemented(
263 "This injection does not support PETSc types."));
264 }
265
266 // If tria is distributed
267 if (dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
268 &tria) != nullptr)
269 matrix.compress(VectorOperation::add);
270 }
271
272
273
286 template <typename VectorType, typename MatrixType, typename SolverType>
287 class MGCoarseDirect : public MGCoarseGridBase<VectorType>
288 {
289 public:
290 explicit MGCoarseDirect(const MatrixType &matrix)
291 : direct_solver(TrilinosWrappers::SolverDirect::AdditionalData())
292 , coarse_matrix(matrix)
293 {
294 // Check if matrix types are supported
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());
302
303 // Check on the vector types: standard deal.II vectors, LA::d::V, or
304 // Trilinos vectors.
305 if constexpr (is_serial_matrix)
306 {
307 [[maybe_unused]] static constexpr bool is_serial_vector =
308 std::is_same_v<VectorType,
309 dealii::Vector<typename MatrixType::value_type>>;
310 Assert(is_serial_vector, ExcNotImplemented());
311 }
312 else if constexpr (is_trilinos_matrix)
313 {
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());
320 }
321 else
322 {
323 AssertThrow(false, ExcNotImplemented());
324 // DEAL_II_NOT_IMPLEMENTED(); //not available in older releases
325 }
326
327
328 // regardless of UMFPACK or Trilinos, both direct solvers need a call
329 // to `initialize()`.
330 direct_solver.initialize(coarse_matrix);
331 }
332
333
334 void
335 operator()(const unsigned int,
336 VectorType &dst,
337 const VectorType &src) const override
338 {
339 AssertDimension(coarse_matrix.n(), src.size());
340 AssertDimension(coarse_matrix.m(), dst.size());
341 const_cast<SolverType &>(direct_solver).solve(dst, src);
342 }
343
344
345 ~MGCoarseDirect() = default;
346
347 private:
348 SolverType direct_solver;
349 const MatrixType &coarse_matrix;
350 };
351
352
353
354 namespace Physics
355 {
357 {
358 double dt = 1e-4;
359 double penalty_constant = 10.;
360 double chi = 1e5;
361 double Cm = 1e-2;
362 double sigma = 12e-2;
363 };
364
365 } // namespace Physics
366
368 {
375 template <int dim,
376 int degree,
377 int n_q_points,
378 int n_components,
379 typename number = double>
381 {
382 public:
383 using value_type = number;
386
387
389
390 void
391 reinit(const Mapping<dim> &mapping,
392 const DoFHandler<dim> &dof_handler,
393 const unsigned int level = numbers::invalid_unsigned_int)
394 {
395 fe_degree = dof_handler.get_fe().degree;
396
397 const QGauss<1> quad(n_q_points);
398 typename MatrixFree<dim, number>::AdditionalData addit_data;
399 addit_data.tasks_parallel_scheme =
401 addit_data.tasks_block_size = 3;
402 addit_data.mg_level = level;
407 constraints.close();
408
409 data.reinit(mapping, dof_handler, constraints, quad, addit_data);
410
412 }
413
414 void
417 {
418 dst = 0;
419 vmult_add(dst, src);
420 }
421
422 void
425 {
426 dst = 0;
427 vmult_add(dst, src);
428 }
429
430 void
436
437 void
440 {
442 *data.get_dof_info(0).vector_partitioner))
443 {
445 src_copy.reinit(data.get_dof_info().vector_partitioner);
446 src_copy = src;
448 src_copy);
449 }
451 *data.get_dof_info(0).vector_partitioner))
452 {
454 dst_copy.reinit(data.get_dof_info().vector_partitioner);
455 dst_copy = dst;
456 dst.swap(dst_copy);
457 }
462 this,
463 dst,
464 src);
465 }
466
468 m() const
469 {
470 return data.get_vector_partitioner()->size();
471 }
472
474 n() const
475 {
476 return data.get_vector_partitioner()->size();
477 }
478
479 number
480 el(const unsigned int row, const unsigned int col) const
481 {
482 (void)row;
483 (void)col;
484 AssertThrow(false,
485 ExcMessage("Matrix-free does not allow for entry access"));
486 return number();
487 }
488
489 void
492 {
493 data.initialize_dof_vector(vector);
494 }
495
496
497 const DoFHandler<dim> &
499 {
500 return data.get_dof_handler();
501 }
502
503
504 const Triangulation<dim> &
506 {
507 return data.get_dof_handler().get_triangulation();
508 }
509
515
516
519 {
520 return &data;
521 }
522
523
524
527 {
528 // Boilerplate for SIP-DG form. TODO: unify interface.
530 const auto cell_operation = [&](auto &phi) {
531 phi.evaluate(EvaluationFlags::gradients);
532 for (unsigned int q = 0; q < phi.n_q_points; ++q)
533 phi.submit_gradient(phi.get_gradient(q), q);
534 phi.integrate(EvaluationFlags::gradients);
535 };
536
537 const auto face_operation = [&](auto &phi_m, auto &phi_p) {
540
541 VectorizedArrayType sigmaF =
542 (std::abs(
543 (phi_m.normal_vector(0) * phi_m.inverse_jacobian(0))[dim - 1]) +
544 std::abs(
545 (phi_m.normal_vector(0) * phi_p.inverse_jacobian(0))[dim - 1])) *
546 (number)(std::max(fe_degree, 1) * (fe_degree + 1.0));
547
548 for (unsigned int q = 0; q < phi_m.n_q_points; ++q)
549 {
550 VectorizedArrayType jump_value =
551 (phi_m.get_value(q) - phi_p.get_value(q)) * 0.5;
552 VectorizedArrayType average_valgrad =
553 phi_m.get_normal_derivative(q) + phi_p.get_normal_derivative(q);
554 average_valgrad =
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);
560 }
563 };
564
565 const auto boundary_operation = [&](auto &phi_m) {
567 VectorizedArrayType sigmaF =
568 std::abs(
569 (phi_m.normal_vector(0) * phi_m.inverse_jacobian(0))[dim - 1]) *
570 number(std::max(fe_degree, 1) * (fe_degree + 1.0)) * 2.0;
571
572 for (unsigned int q = 0; q < phi_m.n_q_points; ++q)
573 {
574 VectorizedArrayType jump_value = phi_m.get_value(q);
575 VectorizedArrayType average_valgrad =
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);
580 }
581
583 };
584
585
587
588
589 // Check if matrix has already been set up.
590 if (system_matrix.m() == 0 && system_matrix.n() == 0)
591 {
592 // Set up sparsity pattern of system matrix.
593 const auto &dof_handler = data.get_dof_handler();
594
596 data.get_mg_level() != numbers::invalid_unsigned_int ?
597 dof_handler.locally_owned_mg_dofs(data.get_mg_level()) :
598 dof_handler.locally_owned_dofs(),
599 data.get_task_info().communicator);
600
601 if (data.get_mg_level() != numbers::invalid_unsigned_int)
603 dsp,
604 data.get_mg_level(),
606 else
608 dsp,
610
611 dsp.compress();
612 system_matrix.reinit(dsp);
613
614 // Assemble system matrix. Notice that degree 1 has been hardcoded.
615
617 degree,
618 n_q_points,
619 n_components,
620 number,
622 data,
625 cell_operation,
626 face_operation,
627 boundary_operation);
628 }
629
630 return system_matrix;
631 }
632
633
634
635 void
637 {
638 // Boilerplate for SIP-DG form. TODO: unify interface.
640 const auto cell_operation = [&](auto &phi) {
641 phi.evaluate(EvaluationFlags::gradients);
642 for (unsigned int q = 0; q < phi.n_q_points; ++q)
643 phi.submit_gradient(phi.get_gradient(q), q);
644 phi.integrate(EvaluationFlags::gradients);
645 };
646
647 const auto face_operation = [&](auto &phi_m, auto &phi_p) {
650
651 VectorizedArrayType sigmaF =
652 (std::abs(
653 (phi_m.normal_vector(0) * phi_m.inverse_jacobian(0))[dim - 1]) +
654 std::abs(
655 (phi_m.normal_vector(0) * phi_p.inverse_jacobian(0))[dim - 1])) *
656 (number)(std::max(fe_degree, 1) * (fe_degree + 1.0));
657
658 for (unsigned int q = 0; q < phi_m.n_q_points; ++q)
659 {
660 VectorizedArrayType jump_value =
661 (phi_m.get_value(q) - phi_p.get_value(q)) * 0.5;
662 VectorizedArrayType average_valgrad =
663 phi_m.get_normal_derivative(q) + phi_p.get_normal_derivative(q);
664 average_valgrad =
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);
670 }
673 };
674
675 const auto boundary_operation = [&](auto &phi_m) {
677 VectorizedArrayType sigmaF =
678 std::abs(
679 (phi_m.normal_vector(0) * phi_m.inverse_jacobian(0))[dim - 1]) *
680 number(std::max(fe_degree, 1) * (fe_degree + 1.0)) * 2.0;
681
682 for (unsigned int q = 0; q < phi_m.n_q_points; ++q)
683 {
684 VectorizedArrayType jump_value = phi_m.get_value(q);
685 VectorizedArrayType average_valgrad =
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);
690 }
691
693 };
694
695
697
698 // Check if matrix has already been set up.
699 AssertThrow((mg_matrix.m() == 0 && mg_matrix.n() == 0),
700 ExcInternalError());
701
702 // Set up sparsity pattern of system matrix.
703 const DoFHandler<dim> &dof_handler = data.get_dof_handler();
704
705 const IndexSet &system_partitioning = dof_handler.locally_owned_dofs();
706 const IndexSet system_relevant_set =
708
709
710 DynamicSparsityPattern dsp(dof_handler.n_dofs(),
711 dof_handler.n_dofs(),
712 system_relevant_set);
713 DoFTools::make_flux_sparsity_pattern(dof_handler, dsp);
714
716 dsp,
717 system_partitioning,
718 data.get_task_info().communicator,
719 system_relevant_set);
720 mg_matrix.reinit(system_partitioning,
721 dsp,
722 data.get_task_info().communicator);
723
724 // Assemble system matrix.
726 degree,
727 n_q_points,
728 n_components,
729 number,
731 data,
733 mg_matrix,
734 cell_operation,
735 face_operation,
736 boundary_operation);
737 }
738
739
740
741 void
743 {
744 const int dummy = 0;
745
746
747 data
748 .template cell_loop<LinearAlgebra::distributed::Vector<number>, int>(
749 [](const auto &matrix_free,
750 auto &dst,
751 const auto &,
752 const auto cells) {
753 FEEvaluation<dim, -1, 0, n_components, number> phi(matrix_free,
754 cells);
755 for (unsigned int cell = cells.first; cell < cells.second; ++cell)
756 {
757 phi.reinit(cell);
758 for (unsigned int q = 0; q < phi.n_q_points; ++q)
759 {
760 if constexpr (n_components == 1)
761 {
762 phi.submit_value(1.0, q);
763 }
764 else
765 {
767 for (unsigned int v = 0;
768 v < VectorizedArray<number>::size();
769 ++v)
770 {
771 for (unsigned int i = 0; i < n_components; i++)
772 temp[i][v] = 1.;
773 }
774 phi.submit_value(temp, q);
775 }
776 }
777
778 phi.integrate_scatter(EvaluationFlags::values, dst);
779 }
780 },
781 b,
782 dummy,
783 true);
784 }
785
786
787
788 void
790 {
791 data.initialize_dof_vector(inverse_diagonal_entries);
792 unsigned int dummy = 0;
796 this,
798 dummy);
799
800 for (unsigned int i = 0;
801 i < inverse_diagonal_entries.locally_owned_size();
802 ++i)
803 if (std::abs(inverse_diagonal_entries.local_element(i)) > 1e-10)
804 inverse_diagonal_entries.local_element(i) =
805 1. / inverse_diagonal_entries.local_element(i);
806 }
807
808
809 private:
810 void
814 const std::pair<unsigned int, unsigned int> &cell_range) const
815 {
816 FEEvaluation<dim, -1, 0, 1, number> phi(data);
817
818 for (unsigned int cell = cell_range.first; cell < cell_range.second;
819 ++cell)
820 {
821 phi.reinit(cell);
822 phi.read_dof_values(src);
823 phi.evaluate(EvaluationFlags::gradients);
824 for (unsigned int q = 0; q < phi.n_q_points; ++q)
825 phi.submit_gradient(phi.get_gradient(q), q);
826 phi.integrate(EvaluationFlags::gradients);
827 phi.distribute_local_to_global(dst);
828 }
829 }
830
831 void
836 const std::pair<unsigned int, unsigned int> &face_range) const
837 {
838 FEFaceEvaluation<dim, -1, 0, 1, number> fe_eval(data, true);
839 FEFaceEvaluation<dim, -1, 0, 1, number> fe_eval_neighbor(data, false);
840
841 for (unsigned int face = face_range.first; face < face_range.second;
842 ++face)
843 {
844 fe_eval.reinit(face);
845 fe_eval_neighbor.reinit(face);
846
847 fe_eval.read_dof_values(src);
848 fe_eval.evaluate(EvaluationFlags::values |
850 fe_eval_neighbor.read_dof_values(src);
851 fe_eval_neighbor.evaluate(EvaluationFlags::values |
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])) *
858 (number)(std::max(fe_degree, 1) * (fe_degree + 1.0));
859
860 for (unsigned int q = 0; q < fe_eval.n_q_points; ++q)
861 {
862 VectorizedArray<number> jump_value =
863 (fe_eval.get_value(q) - fe_eval_neighbor.get_value(q)) * 0.5;
864 VectorizedArray<number> average_valgrad =
865 fe_eval.get_normal_derivative(q) +
866 fe_eval_neighbor.get_normal_derivative(q);
867 average_valgrad =
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);
873 }
874 fe_eval.integrate(EvaluationFlags::values |
876 fe_eval.distribute_local_to_global(dst);
877 fe_eval_neighbor.integrate(EvaluationFlags::values |
879 fe_eval_neighbor.distribute_local_to_global(dst);
880 }
881 }
882
883 void
888 const std::pair<unsigned int, unsigned int> &face_range) const
889 {
890 FEFaceEvaluation<dim, -1, 0, 1, number> fe_eval(data, true);
891 for (unsigned int face = face_range.first; face < face_range.second;
892 ++face)
893 {
894 fe_eval.reinit(face);
895 fe_eval.read_dof_values(src);
896 fe_eval.evaluate(EvaluationFlags::values |
899 std::abs((fe_eval.normal_vector(0) *
900 fe_eval.inverse_jacobian(0))[dim - 1]) *
901 number(std::max(fe_degree, 1) * (fe_degree + 1.0)) * 2.0;
902
903 for (unsigned int q = 0; q < fe_eval.n_q_points; ++q)
904 {
905 VectorizedArray<number> jump_value = fe_eval.get_value(q);
906 VectorizedArray<number> average_valgrad =
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);
911 }
912
913 fe_eval.integrate(EvaluationFlags::values |
915 fe_eval.distribute_local_to_global(dst);
916 }
917 }
918
919
920 void
924 const unsigned int &,
925 const std::pair<unsigned int, unsigned int> &cell_range) const
926 {
927 FEEvaluation<dim, -1, 0, 1, number> phi(data);
928 AlignedVector<VectorizedArray<number>> local_diagonal_vector(
929 phi.dofs_per_cell);
930
931 for (unsigned int cell = cell_range.first; cell < cell_range.second;
932 ++cell)
933 {
934 phi.reinit(cell);
935
936 for (unsigned int i = 0; i < phi.dofs_per_cell; ++i)
937 {
938 for (unsigned int j = 0; j < phi.dofs_per_cell; ++j)
939 phi.begin_dof_values()[j] = VectorizedArray<number>();
940 phi.begin_dof_values()[i] = 1.;
941 phi.evaluate(EvaluationFlags::gradients);
942 for (unsigned int q = 0; q < phi.n_q_points; ++q)
943 phi.submit_gradient(phi.get_gradient(q), q);
944 phi.integrate(EvaluationFlags::gradients);
945 local_diagonal_vector[i] = phi.begin_dof_values()[i];
946 }
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);
950 }
951 }
952
953 void
957 const unsigned int &,
958 const std::pair<unsigned int, unsigned int> &face_range) const
959 {
960 FEFaceEvaluation<dim, -1, 0, 1, number> phi(data, true);
961 FEFaceEvaluation<dim, -1, 0, 1, number> phi_outer(data, false);
962 AlignedVector<VectorizedArray<number>> local_diagonal_vector(
963 phi.dofs_per_cell);
964
965 for (unsigned int face = face_range.first; face < face_range.second;
966 ++face)
967 {
968 phi.reinit(face);
969 phi_outer.reinit(face);
970
972 (std::abs(
973 (phi.normal_vector(0) * phi.inverse_jacobian(0))[dim - 1]) +
974 std::abs((phi.normal_vector(0) *
975 phi_outer.inverse_jacobian(0))[dim - 1])) *
976 (number)(std::max(fe_degree, 1) * (fe_degree + 1.0));
977
978 // Compute phi part
979 for (unsigned int j = 0; j < phi.dofs_per_cell; ++j)
980 phi_outer.begin_dof_values()[j] = VectorizedArray<number>();
981 phi_outer.evaluate(EvaluationFlags::values |
983 for (unsigned int i = 0; i < phi.dofs_per_cell; ++i)
984 {
985 for (unsigned int j = 0; j < phi.dofs_per_cell; ++j)
986 phi.begin_dof_values()[j] = VectorizedArray<number>();
987 phi.begin_dof_values()[i] = 1.;
988 phi.evaluate(EvaluationFlags::values |
990
991 for (unsigned int q = 0; q < phi.n_q_points; ++q)
992 {
993 VectorizedArray<number> jump_value =
994 (phi.get_value(q) - phi_outer.get_value(q)) * 0.5;
995 VectorizedArray<number> average_valgrad =
996 phi.get_normal_derivative(q) +
997 phi_outer.get_normal_derivative(q);
998 average_valgrad =
999 jump_value * 2. * sigmaF - average_valgrad * 0.5;
1000 phi.submit_normal_derivative(-jump_value, q);
1001 phi.submit_value(average_valgrad, q);
1002 }
1003 phi.integrate(EvaluationFlags::values |
1005 local_diagonal_vector[i] = phi.begin_dof_values()[i];
1006 }
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);
1010
1011 // Compute phi_outer part
1012 for (unsigned int j = 0; j < phi.dofs_per_cell; ++j)
1013 phi.begin_dof_values()[j] = VectorizedArray<number>();
1015 for (unsigned int i = 0; i < phi.dofs_per_cell; ++i)
1016 {
1017 for (unsigned int j = 0; j < phi.dofs_per_cell; ++j)
1018 phi_outer.begin_dof_values()[j] = VectorizedArray<number>();
1019 phi_outer.begin_dof_values()[i] = 1.;
1020 phi_outer.evaluate(EvaluationFlags::values |
1022
1023 for (unsigned int q = 0; q < phi.n_q_points; ++q)
1024 {
1025 VectorizedArray<number> jump_value =
1026 (phi.get_value(q) - phi_outer.get_value(q)) * 0.5;
1027 VectorizedArray<number> average_valgrad =
1028 phi.get_normal_derivative(q) +
1029 phi_outer.get_normal_derivative(q);
1030 average_valgrad =
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);
1034 }
1035 phi_outer.integrate(EvaluationFlags::values |
1037 local_diagonal_vector[i] = phi_outer.begin_dof_values()[i];
1038 }
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);
1042 }
1043 }
1044
1045 void
1049 const unsigned int &,
1050 const std::pair<unsigned int, unsigned int> &face_range) const
1051 {
1052 FEFaceEvaluation<dim, -1, 0, 1, number> phi(data);
1053 AlignedVector<VectorizedArray<number>> local_diagonal_vector(
1054 phi.dofs_per_cell);
1055
1056 for (unsigned int face = face_range.first; face < face_range.second;
1057 ++face)
1058 {
1059 phi.reinit(face);
1060
1062 std::abs(
1063 (phi.normal_vector(0) * phi.inverse_jacobian(0))[dim - 1]) *
1064 number(std::max(fe_degree, 1) * (fe_degree + 1.0)) * 2.0;
1065
1066 for (unsigned int i = 0; i < phi.dofs_per_cell; ++i)
1067 {
1068 for (unsigned int j = 0; j < phi.dofs_per_cell; ++j)
1069 phi.begin_dof_values()[j] = VectorizedArray<number>();
1070 phi.begin_dof_values()[i] = 1.;
1071 phi.evaluate(EvaluationFlags::values |
1073
1074 for (unsigned int q = 0; q < phi.n_q_points; ++q)
1075 {
1076 VectorizedArray<number> jump_value = phi.get_value(q);
1077 VectorizedArray<number> average_valgrad =
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);
1082 }
1083
1084 phi.integrate(EvaluationFlags::values |
1086 local_diagonal_vector[i] = phi.begin_dof_values()[i];
1087 }
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);
1091 }
1092 }
1093
1094
1095
1101 };
1102
1103
1104
1108 template <int dim,
1109 int degree,
1110 int n_q_points,
1111 int n_components,
1112 typename number = double>
1114 {
1115 public:
1116 using value_type = number;
1119
1121 {
1122 parameters = parameters_;
1123 };
1124
1125 void
1126 reinit(const Mapping<dim> &mapping,
1127 const DoFHandler<dim> &dof_handler,
1128 const unsigned int level = numbers::invalid_unsigned_int)
1129 {
1130 Assert(
1131 degree == dof_handler.get_fe().degree,
1132 ExcMessage(
1133 "Degree of the operator must match the degree of the DoFHandler"));
1134 fe_degree = dof_handler.get_fe().degree;
1135
1136 const QGauss<1> quad(n_q_points);
1137 typename MatrixFree<dim, number>::AdditionalData addit_data;
1138 addit_data.tasks_parallel_scheme =
1140 addit_data.tasks_block_size = 3;
1141 addit_data.mg_level = level;
1142 addit_data.mapping_update_flags =
1148 constraints.close();
1149
1150 data.reinit(mapping, dof_handler, constraints, quad, addit_data);
1151
1153 }
1154
1155 void
1158 {
1159 dst = 0;
1160 vmult_add(dst, src);
1161 }
1162
1163 void
1166 {
1167 dst = 0;
1168 vmult_add(dst, src);
1169 }
1170
1171 void
1177
1178 void
1181 {
1183 *data.get_dof_info(0).vector_partitioner))
1184 {
1186 src_copy.reinit(data.get_dof_info().vector_partitioner);
1187 src_copy = src;
1189 src_copy);
1190 }
1192 *data.get_dof_info(0).vector_partitioner))
1193 {
1195 dst_copy.reinit(data.get_dof_info().vector_partitioner);
1196 dst_copy = dst;
1197 dst.swap(dst_copy);
1198 }
1203 this,
1204 dst,
1205 src,
1206 /*zero_dst =*/false,
1209 }
1210
1212 m() const
1213 {
1214 return data.get_vector_partitioner()->size();
1215 }
1216
1218 n() const
1219 {
1220 return data.get_vector_partitioner()->size();
1221 }
1222
1223 number
1224 el(const unsigned int row, const unsigned int col) const
1225 {
1226 (void)row;
1227 (void)col;
1228 AssertThrow(false,
1229 ExcMessage("Matrix-free does not allow for entry access"));
1230 return number();
1231 }
1232
1233 void
1236 {
1237 data.initialize_dof_vector(vector);
1238 }
1239
1240
1241 const DoFHandler<dim> &
1243 {
1244 return data.get_dof_handler();
1245 }
1246
1247
1248 const Triangulation<dim> &
1250 {
1251 return data.get_dof_handler().get_triangulation();
1252 }
1253
1256 {
1258 }
1259
1260
1263 {
1264 return &data;
1265 }
1266
1267
1268
1271 {
1272 // Boilerplate for SIP-DG form. TODO: unify interface.
1274 const double factor = (parameters.chi * parameters.Cm) / parameters.dt;
1275
1276 const auto cell_operation = [&](auto &phi) {
1278 for (unsigned int q = 0; q < phi.n_q_points; ++q)
1279 {
1280 phi.submit_value(phi.get_value(q) * factor, q);
1281 phi.submit_gradient(phi.get_gradient(q) * parameters.sigma, q);
1282 }
1284 };
1285
1286 const auto face_operation = [&](auto &phi_m, auto &phi_p) {
1289
1290 VectorizedArrayType sigmaF =
1291 (std::abs(
1292 (phi_m.normal_vector(0) * phi_m.inverse_jacobian(0))[dim - 1]) +
1293 std::abs(
1294 (phi_m.normal_vector(0) * phi_p.inverse_jacobian(0))[dim - 1])) *
1295 (number)(std::max(fe_degree, 1) * (fe_degree + 1.0));
1296
1297 for (unsigned int q = 0; q < phi_m.n_q_points; ++q)
1298 {
1299 VectorizedArrayType jump_value =
1300 (phi_m.get_value(q) - phi_p.get_value(q)) * 0.5;
1301 VectorizedArrayType average_valgrad =
1302 phi_m.get_normal_derivative(q) + phi_p.get_normal_derivative(q);
1303 average_valgrad =
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);
1309 }
1312 };
1313
1314 // const auto boundary_operation = [&](auto &phi_m) {
1315 // (void)phi_m;
1316 // AssertThrow(
1317 // false,
1318 // ExcMessage(
1319 // "We have Neumann BCs only, so this should never be called."));
1320 // // Do nothing since we have homogeneous Neumann BCs
1321 // };
1322
1323
1325
1326
1327 // Check if matrix has already been set up.
1328 if (system_matrix.m() == 0 && system_matrix.n() == 0)
1329 {
1330 // Set up sparsity pattern of system matrix.
1331 const auto &dof_handler = data.get_dof_handler();
1332
1334 data.get_mg_level() != numbers::invalid_unsigned_int ?
1335 dof_handler.locally_owned_mg_dofs(data.get_mg_level()) :
1336 dof_handler.locally_owned_dofs(),
1337 data.get_task_info().communicator);
1338
1339 if (data.get_mg_level() != numbers::invalid_unsigned_int)
1341 dsp,
1342 data.get_mg_level(),
1343 constraints);
1344 else
1346 dsp,
1347 constraints);
1348
1349 dsp.compress();
1350 system_matrix.reinit(dsp);
1351
1352 // Assemble system matrix. Notice that degree 1 has been hardcoded.
1353
1355 degree,
1356 n_q_points,
1357 n_components,
1358 number,
1360 data,
1363 cell_operation,
1364 face_operation,
1365 nullptr /*boundary_operation*/);
1366 }
1367
1368 return system_matrix;
1369 }
1370
1371
1372
1373 void
1375 {
1376 // Boilerplate for SIP-DG form. TODO: unify interface.
1378 const double factor = (parameters.chi * parameters.Cm) / parameters.dt;
1379
1380 const auto cell_operation = [&](auto &phi) {
1382 for (unsigned int q = 0; q < phi.n_q_points; ++q)
1383 {
1384 phi.submit_value(phi.get_value(q) * factor, q);
1385 phi.submit_gradient(phi.get_gradient(q) * parameters.sigma, q);
1386 }
1388 };
1389
1390 const auto face_operation = [&](auto &phi_m, auto &phi_p) {
1393
1394 VectorizedArrayType sigmaF =
1395 (std::abs(
1396 (phi_m.normal_vector(0) * phi_m.inverse_jacobian(0))[dim - 1]) +
1397 std::abs(
1398 (phi_m.normal_vector(0) * phi_p.inverse_jacobian(0))[dim - 1])) *
1399 (number)(std::max(fe_degree, 1) * (fe_degree + 1.0));
1400
1401 for (unsigned int q = 0; q < phi_m.n_q_points; ++q)
1402 {
1403 VectorizedArrayType jump_value =
1404 (phi_m.get_value(q) - phi_p.get_value(q)) * 0.5;
1405 VectorizedArrayType average_valgrad =
1406 phi_m.get_normal_derivative(q) + phi_p.get_normal_derivative(q);
1407 average_valgrad =
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);
1413 }
1416 };
1417
1418
1420
1421 // Check if matrix has already been set up.
1422 AssertThrow((mg_matrix.m() == 0 && mg_matrix.n() == 0),
1423 ExcInternalError());
1424
1425 // Set up sparsity pattern of system matrix.
1426 const DoFHandler<dim> &dof_handler = data.get_dof_handler();
1427
1428 const IndexSet &system_partitioning = dof_handler.locally_owned_dofs();
1429 const IndexSet system_relevant_set =
1431
1432
1433 DynamicSparsityPattern dsp(dof_handler.n_dofs(),
1434 dof_handler.n_dofs(),
1435 system_relevant_set);
1436 DoFTools::make_flux_sparsity_pattern(dof_handler, dsp);
1437
1439 dsp,
1440 system_partitioning,
1441 data.get_task_info().communicator,
1442 system_relevant_set);
1443 mg_matrix.reinit(system_partitioning,
1444 dsp,
1445 data.get_task_info().communicator);
1446
1447 // Assemble system matrix.
1449 degree,
1450 n_q_points,
1451 n_components,
1452 number,
1454 data,
1456 mg_matrix,
1457 cell_operation,
1458 face_operation,
1459 nullptr /*boundary_operation*/);
1460 }
1461
1462
1463 void
1465 {
1466 data.initialize_dof_vector(inverse_diagonal_entries);
1467 unsigned int dummy = 0;
1471 this,
1473 dummy);
1474
1475 for (unsigned int i = 0;
1476 i < inverse_diagonal_entries.locally_owned_size();
1477 ++i)
1478 if (std::abs(inverse_diagonal_entries.local_element(i)) > 1e-10)
1479 inverse_diagonal_entries.local_element(i) =
1480 1. / inverse_diagonal_entries.local_element(i);
1481 }
1482
1483
1484
1485 void
1489 {
1491 const double factor = (parameters.chi * parameters.Cm) / parameters.dt;
1492
1493 for (unsigned int cell = 0; cell < data.n_cell_batches(); ++cell)
1494 {
1495 phi.reinit(cell);
1496 phi.read_dof_values(src);
1498 for (unsigned int q = 0; q < phi.n_q_points; ++q)
1499 phi.submit_value(phi.get_value(q) * factor, q);
1502 }
1504 }
1505
1506
1507
1508 void
1510 const LinearAlgebra::distributed::Vector<number> &solution_minus_ion,
1511 const Function<dim> &Iext) const
1512 {
1514
1515 for (unsigned int cell = 0; cell < data.n_cell_batches(); ++cell)
1516 {
1517 phi.reinit(cell);
1518 phi.read_dof_values(solution_minus_ion);
1520 for (unsigned int q = 0; q < phi.n_q_points; ++q)
1521 {
1523 phi.quadrature_point(q);
1524 // evaluate the external current for each component in
1525 // VectorizedArray
1526 VectorizedArray<number> applied_current_value = 0.0;
1527 for (unsigned int v = 0; v < VectorizedArray<number>::size();
1528 ++v)
1529 {
1530 Point<dim> p;
1531 for (unsigned int d = 0; d < dim; ++d)
1532 p[d] = p_vect[d][v];
1533 applied_current_value[v] = Iext.value(p);
1534 }
1535
1536 // reaction_term = (parameters.chi * phi.get_value(q)) +
1537 // applied_current_value;
1538 phi.submit_value((parameters.chi * phi.get_value(q)) +
1539 applied_current_value,
1540 q);
1541 }
1544 }
1545 rhs.compress(VectorOperation::add);
1546 }
1547
1548
1549
1550 private:
1551 void
1555 const std::pair<unsigned int, unsigned int> &cell_range) const
1556 {
1558 const double factor = (parameters.chi * parameters.Cm) / parameters.dt;
1559
1560 for (unsigned int cell = cell_range.first; cell < cell_range.second;
1561 ++cell)
1562 {
1563 phi.reinit(cell);
1564 phi.read_dof_values(src);
1566 for (unsigned int q = 0; q < phi.n_q_points; ++q)
1567 {
1568 phi.submit_value(phi.get_value(q) * factor, q);
1569 phi.submit_gradient(phi.get_gradient(q) * parameters.sigma, q);
1570 }
1573 }
1574 }
1575
1576 void
1581 const std::pair<unsigned int, unsigned int> &face_range) const
1582 {
1584 data, true);
1586 fe_eval_neighbor(data, false);
1587
1588 for (unsigned int face = face_range.first; face < face_range.second;
1589 ++face)
1590 {
1591 fe_eval.reinit(face);
1592 fe_eval_neighbor.reinit(face);
1593
1594 fe_eval.read_dof_values(src);
1597 fe_eval_neighbor.read_dof_values(src);
1598 fe_eval_neighbor.evaluate(EvaluationFlags::values |
1601 (std::abs((fe_eval.normal_vector(0) *
1602 fe_eval.inverse_jacobian(0))[dim - 1]) +
1603 std::abs((fe_eval.normal_vector(0) *
1604 fe_eval_neighbor.inverse_jacobian(0))[dim - 1])) *
1605 (number)(std::max(fe_degree, 1) * (fe_degree + 1.0));
1606
1607 for (unsigned int q = 0; q < fe_eval.n_q_points; ++q)
1608 {
1609 VectorizedArray<number> jump_value =
1610 (fe_eval.get_value(q) - fe_eval_neighbor.get_value(q)) * 0.5;
1611 VectorizedArray<number> average_valgrad =
1612 fe_eval.get_normal_derivative(q) +
1613 fe_eval_neighbor.get_normal_derivative(q);
1614 average_valgrad =
1615 jump_value * 2. * sigmaF - average_valgrad * 0.5;
1616 fe_eval.submit_normal_derivative(-parameters.sigma * jump_value,
1617 q);
1618 fe_eval_neighbor.submit_normal_derivative(-parameters.sigma *
1619 jump_value,
1620 q);
1621 fe_eval.submit_value(parameters.sigma * average_valgrad, q);
1622 fe_eval_neighbor.submit_value(-parameters.sigma *
1623 average_valgrad,
1624 q);
1625 }
1628 fe_eval.distribute_local_to_global(dst);
1629 fe_eval_neighbor.integrate(EvaluationFlags::values |
1631 fe_eval_neighbor.distribute_local_to_global(dst);
1632 }
1633 }
1634
1635 void
1640 const std::pair<unsigned int, unsigned int> &face_range) const
1641 {
1642 (void)data;
1643 (void)dst;
1644 (void)src;
1645 (void)face_range;
1646 // Do nothing since we have homogeneous Neumann BCs
1647 }
1648
1649
1650 void
1654 const unsigned int &,
1655 const std::pair<unsigned int, unsigned int> &cell_range) const
1656 {
1658 AlignedVector<VectorizedArray<number>> local_diagonal_vector(
1659 phi.dofs_per_cell);
1660 const double factor = (parameters.chi * parameters.Cm) / parameters.dt;
1661
1662 for (unsigned int cell = cell_range.first; cell < cell_range.second;
1663 ++cell)
1664 {
1665 phi.reinit(cell);
1666
1667 for (unsigned int i = 0; i < phi.dofs_per_cell; ++i)
1668 {
1669 for (unsigned int j = 0; j < phi.dofs_per_cell; ++j)
1671 phi.begin_dof_values()[i] = 1.;
1674 for (unsigned int q = 0; q < phi.n_q_points; ++q)
1675 {
1676 phi.submit_value(phi.get_value(q) * factor, q);
1677 phi.submit_gradient(phi.get_gradient(q) * parameters.sigma,
1678 q);
1679 }
1682 local_diagonal_vector[i] = phi.begin_dof_values()[i];
1683 }
1684 for (unsigned int i = 0; i < phi.dofs_per_cell; ++i)
1685 phi.begin_dof_values()[i] = local_diagonal_vector[i];
1687 }
1688 }
1689
1690 void
1694 const unsigned int &,
1695 const std::pair<unsigned int, unsigned int> &face_range) const
1696 {
1698 data, true);
1700 phi_outer(data, false);
1701 AlignedVector<VectorizedArray<number>> local_diagonal_vector(
1702 phi.dofs_per_cell);
1703
1704 for (unsigned int face = face_range.first; face < face_range.second;
1705 ++face)
1706 {
1707 phi.reinit(face);
1708 phi_outer.reinit(face);
1709
1711 (std::abs(
1712 (phi.normal_vector(0) * phi.inverse_jacobian(0))[dim - 1]) +
1713 std::abs((phi.normal_vector(0) *
1714 phi_outer.inverse_jacobian(0))[dim - 1])) *
1715 (number)(std::max(fe_degree, 1) * (fe_degree + 1.0));
1716
1717 // Compute phi part
1718 for (unsigned int j = 0; j < phi.dofs_per_cell; ++j)
1719 phi_outer.begin_dof_values()[j] = VectorizedArray<number>();
1722 for (unsigned int i = 0; i < phi.dofs_per_cell; ++i)
1723 {
1724 for (unsigned int j = 0; j < phi.dofs_per_cell; ++j)
1726 phi.begin_dof_values()[i] = 1.;
1729
1730 for (unsigned int q = 0; q < phi.n_q_points; ++q)
1731 {
1732 VectorizedArray<number> jump_value =
1733 (phi.get_value(q) - phi_outer.get_value(q)) * 0.5;
1734 VectorizedArray<number> average_valgrad =
1735 phi.get_normal_derivative(q) +
1736 phi_outer.get_normal_derivative(q);
1737 average_valgrad =
1738 jump_value * 2. * sigmaF - average_valgrad * 0.5;
1739 phi.submit_normal_derivative(-parameters.sigma * jump_value,
1740 q);
1741 phi.submit_value(parameters.sigma * average_valgrad, q);
1742 }
1745 local_diagonal_vector[i] = phi.begin_dof_values()[i];
1746 }
1747 for (unsigned int i = 0; i < phi.dofs_per_cell; ++i)
1748 phi.begin_dof_values()[i] = local_diagonal_vector[i];
1750
1751 // Compute phi_outer part
1752 for (unsigned int j = 0; j < phi.dofs_per_cell; ++j)
1755 for (unsigned int i = 0; i < phi.dofs_per_cell; ++i)
1756 {
1757 for (unsigned int j = 0; j < phi.dofs_per_cell; ++j)
1758 phi_outer.begin_dof_values()[j] = VectorizedArray<number>();
1759 phi_outer.begin_dof_values()[i] = 1.;
1762
1763 for (unsigned int q = 0; q < phi.n_q_points; ++q)
1764 {
1765 VectorizedArray<number> jump_value =
1766 (phi.get_value(q) - phi_outer.get_value(q)) * 0.5;
1767 VectorizedArray<number> average_valgrad =
1768 phi.get_normal_derivative(q) +
1769 phi_outer.get_normal_derivative(q);
1770 average_valgrad =
1771 jump_value * 2. * sigmaF - average_valgrad * 0.5;
1772 phi_outer.submit_normal_derivative(-parameters.sigma *
1773 jump_value,
1774 q);
1775 phi_outer.submit_value(-parameters.sigma * average_valgrad,
1776 q);
1777 }
1780 local_diagonal_vector[i] = phi_outer.begin_dof_values()[i];
1781 }
1782 for (unsigned int i = 0; i < phi.dofs_per_cell; ++i)
1783 phi_outer.begin_dof_values()[i] = local_diagonal_vector[i];
1784 phi_outer.distribute_local_to_global(dst);
1785 }
1786 }
1787
1788 void
1792 const unsigned int &,
1793 const std::pair<unsigned int, unsigned int> &face_range) const
1794 {
1795 (void)data;
1796 (void)dst;
1797 (void)face_range;
1798 // Do nothing since we have homogeneous Neumann BCs
1799 }
1800
1801
1802
1809 };
1810
1811
1812
1813 } // namespace MatrixFreeOperators
1814
1815
1816
1821 get_index(const std::vector<types::global_cell_index> &v,
1822 const types::global_cell_index index)
1823 {
1824 return std::distance(v.begin(), std::find(v.begin(), v.end(), index));
1825 }
1826
1827
1828
1829 template <int dim>
1830 void
1832 const std::vector<typename Triangulation<dim>::active_cell_iterator> &cells,
1833 Graph &g)
1834 {
1835 Assert(cells.size() > 0, ExcMessage("No cells to be agglomerated."));
1836 const unsigned int n_faces = cells[0]->n_faces();
1837
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();
1841
1842 g.adjacency.resize(cells.size());
1843 for (const auto &cell : cells)
1844 {
1845 // std::cout << "Cell idx: " << cell->active_cell_index() <<
1846 // std::endl; std::cout << "new idx: "
1847 // << get_index(vec_cells, cell->active_cell_index())
1848 // << std::endl;
1849 g.nodes.push_back(get_index(vec_cells, cell->active_cell_index()));
1850 for (unsigned int f = 0; f < n_faces; ++f)
1851 {
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))
1855 g.adjacency[get_index(vec_cells, cell->active_cell_index())]
1856 .push_back(get_index(vec_cells, neigh->active_cell_index()));
1857 }
1858 }
1859 }
1860
1861
1862
1863 inline void
1864 dfs(std::vector<types::global_cell_index> &comp,
1865 std::vector<bool> &visited,
1866 const Graph &g,
1868 {
1869 visited[v] = true;
1870 comp.push_back(v);
1871 for (const types::global_cell_index u : g.adjacency[v])
1872 {
1873 if (!visited[u])
1874 dfs(comp, visited, g, u);
1875 }
1876 }
1877
1878
1879
1880 void
1882 Graph &g,
1883 std::vector<std::vector<types::global_cell_index>> &connected_components)
1884 {
1885 Assert(g.nodes.size() > 0, ExcMessage("No nodes in this graph."));
1886 Assert(
1887 connected_components.size() == 0,
1888 ExcMessage(
1889 "Connected components have to be computed by the present function."));
1890
1891 std::vector<bool> visited(g.nodes.size()); // register visited node
1892 std::fill(visited.begin(), visited.end(), 0);
1893
1894
1896 {
1897 if (!visited[v])
1898 {
1899 connected_components.emplace_back();
1900 dfs(connected_components.back(), visited, g, v);
1901 }
1902 }
1903 }
1904
1905
1906} // namespace Utils
1907
1908
1909
1910#endif
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 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)
bool empty() const
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
const MatrixType & coarse_matrix
Definition utils.h:349
MGCoarseDirect(const MatrixType &matrix)
Definition utils.h:290
void operator()(const unsigned int, VectorType &dst, const VectorType &src) const override
Definition utils.h:335
SolverType direct_solver
Definition utils.h:348
void Tvmult(LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src) const
Definition utils.h:423
LinearAlgebra::distributed::Vector< number > inverse_diagonal_entries
Definition utils.h:1097
AffineConstraints< number > constraints
Definition utils.h:1100
void rhs(LinearAlgebra::distributed::Vector< number > &b) const
Definition utils.h:742
void vmult_add(LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src) const
Definition utils.h:438
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
Definition utils.h:954
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
Definition utils.h:1046
TrilinosWrappers::SparseMatrix system_matrix
Definition utils.h:1099
number el(const unsigned int row, const unsigned int col) const
Definition utils.h:480
void initialize_dof_vector(LinearAlgebra::distributed::Vector< number > &vector) const
Definition utils.h:490
void reinit(const Mapping< dim > &mapping, const DoFHandler< dim > &dof_handler, const unsigned int level=numbers::invalid_unsigned_int)
Definition utils.h:391
VectorizedArray< number > VectorizedArrayType
Definition utils.h:384
const Triangulation< dim > & get_triangulation() const
Definition utils.h:505
const LinearAlgebra::distributed::Vector< number > & get_matrix_diagonal_inverse() const
Definition utils.h:511
const DoFHandler< dim > & get_dof_handler() const
Definition utils.h:498
const MatrixFree< dim, number > * get_matrix_free() const
Definition utils.h:518
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
Definition utils.h:921
void vmult(LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src) const
Definition utils.h:415
void get_system_matrix(TrilinosWrappers::SparseMatrix &mg_matrix) const
Definition utils.h:636
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
Definition utils.h:884
types::global_dof_index n() const
Definition utils.h:474
const TrilinosWrappers::SparseMatrix & get_system_matrix() const
Definition utils.h:526
void Tvmult_add(LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src) const
Definition utils.h:431
LinearAlgebra::distributed::Vector< number > VectorType
Definition utils.h:385
types::global_dof_index m() const
Definition utils.h:468
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
Definition utils.h:811
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
Definition utils.h:832
Physics::BilinearFormParameters parameters
Definition utils.h:1808
void rhs(LinearAlgebra::distributed::Vector< number > &rhs, const LinearAlgebra::distributed::Vector< number > &solution_minus_ion, const Function< dim > &Iext) const
Definition utils.h:1509
MonodomainOperatorDG(const Physics::BilinearFormParameters &parameters_)
Definition utils.h:1120
const DoFHandler< dim > & get_dof_handler() const
Definition utils.h:1242
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
Definition utils.h:1577
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
Definition utils.h:1552
const LinearAlgebra::distributed::Vector< number > & get_matrix_diagonal_inverse() const
Definition utils.h:1255
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
Definition utils.h:1691
const MatrixFree< dim, number > * get_matrix_free() const
Definition utils.h:1262
void reinit(const Mapping< dim > &mapping, const DoFHandler< dim > &dof_handler, const unsigned int level=numbers::invalid_unsigned_int)
Definition utils.h:1126
TrilinosWrappers::SparseMatrix system_matrix
Definition utils.h:1806
LinearAlgebra::distributed::Vector< number > inverse_diagonal_entries
Definition utils.h:1804
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
Definition utils.h:1651
VectorizedArray< number > VectorizedArrayType
Definition utils.h:1117
AffineConstraints< number > constraints
Definition utils.h:1807
void Tvmult_add(LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src) const
Definition utils.h:1172
void apply_mass_term(LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src) const
Definition utils.h:1486
const TrilinosWrappers::SparseMatrix & get_system_matrix() const
Definition utils.h:1270
number el(const unsigned int row, const unsigned int col) const
Definition utils.h:1224
LinearAlgebra::distributed::Vector< number > VectorType
Definition utils.h:1118
types::global_dof_index m() const
Definition utils.h:1212
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
Definition utils.h:1636
void get_system_matrix(TrilinosWrappers::SparseMatrix &mg_matrix) const
Definition utils.h:1374
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
Definition utils.h:1789
void vmult_add(LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src) const
Definition utils.h:1179
void vmult(LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src) const
Definition utils.h:1156
types::global_dof_index n() const
Definition utils.h:1218
void initialize_dof_vector(LinearAlgebra::distributed::Vector< number > &vector) const
Definition utils.h:1234
void Tvmult(LinearAlgebra::distributed::Vector< number > &dst, const LinearAlgebra::distributed::Vector< number > &src) const
Definition utils.h:1164
const Triangulation< dim > & get_triangulation() const
Definition utils.h:1249
EnableObserverPointer Subscriptor
unsigned int level
#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)
update_values
update_JxW_values
update_gradients
update_quadrature_points
Expression pow(const Expression &base, const Expression &exponent)
IndexSet extract_locally_relevant_dofs(const DoFHandler< dim, spacedim > &dof_handler)
void make_flux_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity, const unsigned int level, const AffineConstraints< number > &constraints={}, const bool keep_constrained_dofs=true)
void compute_matrix(const MatrixFree< dim, Number, VectorizedArrayType > &matrix_free, const AffineConstraints< Number > &constraints, MatrixType &matrix, const std::function< void(FEEvaluation< dim, fe_degree, n_q_points_1d, n_components, Number, VectorizedArrayType > &)> &cell_operation, const unsigned int dof_handler_index=0, const unsigned int quadrature_index=0, const unsigned int first_selected_component=0)
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 distribute_sparsity_pattern(DynamicSparsityPattern &dsp, const IndexSet &locally_owned_rows, const MPI_Comm mpi_comm, const IndexSet &locally_relevant_rows)
Definition utils.h:47
void create_graph_from_agglomerate(const std::vector< typename Triangulation< dim >::active_cell_iterator > &cells, Graph &g)
Definition utils.h:1831
void dfs(std::vector< types::global_cell_index > &comp, std::vector< bool > &visited, const Graph &g, const types::global_cell_index v)
Definition utils.h:1864
constexpr T constexpr_pow(T num, unsigned int pow)
Definition utils.h:50
void fill_injection_matrix(const AgglomerationHandler< dim, spacedim > &coarse_ah, const AgglomerationHandler< dim, spacedim > &fine_ah, SparsityPattern &sp, MatrixType &matrix)
Definition utils.h:97
void compute_connected_components(Graph &g, std::vector< std::vector< types::global_cell_index > > &connected_components)
Definition utils.h:1881
types::global_cell_index get_index(const std::vector< types::global_cell_index > &v, const types::global_cell_index index)
Definition utils.h:1821
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
void print_graph() const
Definition utils.h:69
std::vector< std::vector< types::global_cell_index > > adjacency
Definition utils.h:66
std::vector< types::global_cell_index > nodes
Definition utils.h:65