deal.II version GIT relicensing-1721-g8100761196 2024-08-31 12:30:00+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
step-76.h
Go to the documentation of this file.
1true);
280 FEFaceEvaluation<dim, degree, degree + 1, 1, Number> phi_p(data, /*is_interior_face=*/false);
281
282 for (unsigned int face = range.first; face < range.second; ++face)
283 {
284 phi_m.reinit(face);
285 phi_m.gather_evaluate(src, face_evaluation_flags);
286 phi_p.reinit(face);
287 phi_p.gather_evaluate(src, face_evaluation_flags);
288
289 // some operations on the face quadrature points
290
291 phi_m.integrate_scatter(face_evaluation_flags, dst);
292 phi_p.integrate_scatter(face_evaluation_flags, dst);
293 }
294 },
295 [&](const auto &data, auto &dst, const auto &src, const auto range) {
296 // operation performed boundary faces
297
298 FEFaceEvaluation<dim, degree, degree + 1, 1, Number> phi_m(data, /*is_interior_face=*/true);
299
300 for (unsigned int face = range.first; face < range.second; ++face)
301 {
302 phi_m.reinit(face);
303 phi_m.gather_evaluate(src, face_evaluation_flags);
304
305 // some operations on the face quadrature points
306
307 phi_m.integrate_scatter(face_evaluation_flags, dst);
308 }
309 },
310 dst,
311 src);
312@endcode
313
314in the following way:
315
316@code
317matrix_free.template loop_cell_centric<VectorType, VectorType>(
318 [&](const auto &data, auto &dst, const auto &src, const auto range) {
320 FEFaceEvaluation<dim, degree, degree + 1, 1, Number> phi_m(data, /*is_interior_face=*/true);
321 FEFaceEvaluation<dim, degree, degree + 1, 1, Number> phi_p(data, /*is_interior_face=*/false);
322
323 for (unsigned int cell = range.first; cell < range.second; ++cell)
324 {
325 phi.reinit(cell);
326 phi.gather_evaluate(src, cell_evaluation_flags);
327
328 // some operations on the cell quadrature points
329
330 phi.integrate_scatter(cell_evaluation_flags, dst);
331
332 // loop over all faces of cell
333 for (unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell; ++face)
334 {
335 if (data.get_faces_by_cells_boundary_id(cell, face)[0] ==
337 {
338 // internal face
339 phi_m.reinit(cell, face);
340 phi_m.gather_evaluate(src, face_evaluation_flags);
341 phi_p.reinit(cell, face);
342 phi_p.gather_evaluate(src, face_evaluation_flags);
343
344 // some operations on the face quadrature points
345
346 phi_m.integrate_scatter(face_evaluation_flags, dst);
347 }
348 else
349 {
350 // boundary face
351 phi_m.reinit(cell, face);
352 phi_m.gather_evaluate(src, face_evaluation_flags);
353
354 // some operations on the face quadrature points
355
356 phi_m.integrate_scatter(face_evaluation_flags, dst);
357 }
358 }
359 }
360 },
361 dst,
362 src);
363@endcode
364
365It should be noted that FEFaceEvaluation is initialized now with two numbers,
366the cell number and the local face number. The given example only
367highlights how to transform face-centric loops into cell-centric loops and
368is by no means efficient, since data is read and written multiple times
369from and to the global vector as well as computations are performed
370redundantly. Below, we will discuss advanced techniques that target these issues.
371
372To be able to use MatrixFree::loop_cell_centric(), following flags of MatrixFree::AdditionalData
373have to be enabled:
374
375@code
376typename MatrixFree<dim, Number>::AdditionalData additional_data;
377
378// set flags as usual (not shown)
379
380additional_data.hold_all_faces_to_owned_cells = true;
381additional_data.mapping_update_flags_faces_by_cells =
382 additional_data.mapping_update_flags_inner_faces |
383 additional_data.mapping_update_flags_boundary_faces;
384
385data.reinit(mapping, dof_handler, constraint, quadrature, additional_data);
386@endcode
387
388In particular, these flags enable that the internal data structures are set up
389for all faces of the cells.
390
391Currently, cell-centric loops in deal.II only work for uniformly refined meshes
392and if no constraints are applied (which is the standard case DG is normally
393used).
394
395
396<a name="step_76-ProvidinglambdastoMatrixFreeloops"></a><h3>Providing lambdas to MatrixFree loops</h3>
397
398
399The examples given above have already used lambdas, which have been provided to
400matrix-free loops. The following short examples present how to transform functions between
401a version where a class and a pointer to one of its methods are used and a
402variant where lambdas are utilized.
403
404In the following code, a class and a pointer to one of its methods, which should
405be interpreted as cell integral, are passed to MatrixFree::loop():
406
407@code
408void
409local_apply_cell(const MatrixFree<dim, Number> & data,
410 VectorType & dst,
411 const VectorType & src,
412 const std::pair<unsigned int, unsigned int> &range) const
413{
415 for (unsigned int cell = range.first; cell < range.second; ++cell)
416 {
417 phi.reinit(cell);
418 phi.gather_evaluate(src, cell_evaluation_flags);
419
420 // some operations on the quadrature points
421
422 phi.integrate_scatter(cell_evaluation_flags, dst);
423 }
424}
425@endcode
426
427@code
428matrix_free.cell_loop(&Operator::local_apply_cell, this, dst, src);
429@endcode
430
431However, it is also possible to pass an anonymous function via a lambda function
432with the same result:
433
434@code
435matrix_free.template cell_loop<VectorType, VectorType>(
436 [&](const auto &data, auto &dst, const auto &src, const auto range) {
438 for (unsigned int cell = range.first; cell < range.second; ++cell)
439 {
440 phi.reinit(cell);
441 phi.gather_evaluate(src, cell_evaluation_flags);
442
443 // some operations on the quadrature points
444
445 phi.integrate_scatter(cell_evaluation_flags, dst);
446 }
447 },
448 dst,
449 src);
450@endcode
451
452<a name="step_76-VectorizedArrayType"></a><h3>VectorizedArrayType</h3>
453
454
455The class VectorizedArray<Number> is a key component to achieve the high
456node-level performance of the matrix-free algorithms in deal.II.
457It is a wrapper class around a short vector of @f$n@f$ entries of type Number and
458maps arithmetic operations to appropriate single-instruction/multiple-data
459(SIMD) concepts by intrinsic functions. The length of the vector can be
460queried by VectorizedArray::size() and its underlying number type by
462
463In the default case (<code>VectorizedArray<Number></code>), the vector length is
464set at compile time of the library to
465match the highest value supported by the given processor architecture.
466However, also a second optional template argument can be
467specified as <code>VectorizedArray<Number, size></code>, where <code>size</code> explicitly
468controls the vector length within the capabilities of a particular instruction
469set. A full list of supported vector lengths is presented in the following table:
470
471<table align="center" class="doxtable">
472 <tr>
473 <th>double</th>
474 <th>float</th>
475 <th>ISA</th>
476 </tr>
477 <tr>
478 <td><code>VectorizedArray<double, 1></code></td>
479 <td><code>VectorizedArray<float, 1></code></td>
480 <td>(auto-vectorization)</td>
481 </tr>
482 <tr>
483 <td><code>VectorizedArray<double, 2></code></td>
484 <td><code>VectorizedArray<float, 4></code></td>
485 <td>SSE2/AltiVec</td>
486 </tr>
487 <tr>
488 <td><code>VectorizedArray<double, 4></code></td>
489 <td><code>VectorizedArray<float, 8></code></td>
490 <td>AVX/AVX2</td>
491 </tr>
492 <tr>
493 <td><code>VectorizedArray<double, 8></code></td>
494 <td><code>VectorizedArray<float, 16></code></td>
495 <td>AVX-512</td>
496 </tr>
497</table>
498
499This allows users to select the vector length/ISA and, as a consequence, the
500number of cells to be processed at once in matrix-free operator evaluations,
501possibly reducing the pressure on the caches, an severe issue for very high
502degrees (and dimensions).
503
504A possible further reason to reduce the number of filled lanes
505is to simplify debugging: instead of having to look at, e.g., 8
506cells, one can concentrate on a single cell.
507
508The interface of VectorizedArray also enables the replacement by any type with
509a matching interface. Specifically, this prepares deal.II for the <code>std::simd</code>
510class that is planned to become part of the C++23 standard. The following table
511compares the deal.II-specific SIMD classes and the equivalent C++23 classes:
512
513
514<table align="center" class="doxtable">
515 <tr>
516 <th>VectorizedArray (deal.II)</th>
517 <th>std::simd (C++23)</th>
518 </tr>
519 <tr>
520 <td><code>VectorizedArray<Number></code></td>
521 <td><code>std::experimental::native_simd<Number></code></td>
522 </tr>
523 <tr>
524 <td><code>VectorizedArray<Number, size></code></td>
525 <td><code>std::experimental::fixed_size_simd<Number, size></code></td>
526 </tr>
527</table>
528 *
529 *
530 * <a name="step_76-CommProg"></a>
531 * <h1> The commented program</h1>
532 *
533 *
534 * <a name="step_76-Parametersandutilityfunctions"></a>
535 * <h3>Parameters and utility functions</h3>
536 *
537
538 *
539 * The same includes as in @ref step_67 "step-67":
540 *
541 * @code
542 *   #include <deal.II/base/conditional_ostream.h>
543 *   #include <deal.II/base/function.h>
544 *   #include <deal.II/base/time_stepping.h>
545 *   #include <deal.II/base/timer.h>
546 *   #include <deal.II/base/utilities.h>
547 *   #include <deal.II/base/vectorization.h>
548 *  
549 *   #include <deal.II/distributed/tria.h>
550 *  
551 *   #include <deal.II/dofs/dof_handler.h>
552 *  
553 *   #include <deal.II/fe/fe_dgq.h>
554 *   #include <deal.II/fe/fe_system.h>
555 *  
556 *   #include <deal.II/grid/grid_generator.h>
557 *   #include <deal.II/grid/tria.h>
558 *   #include <deal.II/grid/tria_accessor.h>
559 *   #include <deal.II/grid/tria_iterator.h>
560 *  
561 *   #include <deal.II/lac/affine_constraints.h>
562 *   #include <deal.II/lac/la_parallel_vector.h>
563 *  
564 *   #include <deal.II/matrix_free/fe_evaluation.h>
565 *   #include <deal.II/matrix_free/matrix_free.h>
566 *   #include <deal.II/matrix_free/operators.h>
567 *  
568 *   #include <deal.II/numerics/data_out.h>
569 *  
570 *   #include <fstream>
571 *   #include <iomanip>
572 *   #include <iostream>
573 *  
574 * @endcode
575 *
576 * A new include for categorizing of cells according to their boundary IDs:
577 *
578 * @code
579 *   #include <deal.II/matrix_free/tools.h>
580 *  
581 *  
582 *  
583 *   namespace Euler_DG
584 *   {
585 *   using namespace dealii;
586 *  
587 * @endcode
588 *
589 * The same input parameters as in @ref step_67 "step-67":
590 *
591 * @code
592 *   constexpr unsigned int testcase = 1;
593 *   constexpr unsigned int dimension = 2;
594 *   constexpr unsigned int n_global_refinements = 2;
595 *   constexpr unsigned int fe_degree = 5;
596 *   constexpr unsigned int n_q_points_1d = fe_degree + 2;
597 *  
598 * @endcode
599 *
600 * This parameter specifies the size of the shared-memory group. Currently,
601 * only the values 1 and numbers::invalid_unsigned_int is possible, leading
602 * to the options that the memory features can be turned off or all processes
603 * having access to the same shared-memory domain are grouped together.
604 *
605 * @code
606 *   constexpr unsigned int group_size = numbers::invalid_unsigned_int;
607 *  
608 *   using Number = double;
609 *  
610 * @endcode
611 *
612 * Here, the type of the data structure is chosen for vectorization. In the
613 * default case, VectorizedArray<Number> is used, i.e., the highest
614 * instruction-set-architecture extension available on the given hardware with
615 * the maximum number of vector lanes is used. However, one might reduce
616 * the number of filled lanes, e.g., by writing
617 * <code>using VectorizedArrayType = VectorizedArray<Number, 4></code> to only
618 * process 4 cells.
619 *
620 * @code
621 *   using VectorizedArrayType = VectorizedArray<Number>;
622 *  
623 * @endcode
624 *
625 * The following parameters have not changed:
626 *
627 * @code
628 *   constexpr double gamma = 1.4;
629 *   constexpr double final_time = testcase == 0 ? 10 : 2.0;
630 *   constexpr double output_tick = testcase == 0 ? 1 : 0.05;
631 *  
632 *   const double courant_number = 0.15 / std::pow(fe_degree, 1.5);
633 *  
634 * @endcode
635 *
636 * Specify max number of time steps useful for performance studies.
637 *
638 * @code
639 *   constexpr unsigned int max_time_steps = numbers::invalid_unsigned_int;
640 *  
641 * @endcode
642 *
643 * Runge-Kutta-related functions copied from @ref step_67 "step-67" and slightly modified
644 * with the purpose to minimize global vector access:
645 *
646 * @code
647 *   enum LowStorageRungeKuttaScheme
648 *   {
649 *   stage_3_order_3,
650 *   stage_5_order_4,
651 *   stage_7_order_4,
652 *   stage_9_order_5,
653 *   };
654 *   constexpr LowStorageRungeKuttaScheme lsrk_scheme = stage_5_order_4;
655 *  
656 *  
657 *  
658 *   class LowStorageRungeKuttaIntegrator
659 *   {
660 *   public:
661 *   LowStorageRungeKuttaIntegrator(const LowStorageRungeKuttaScheme scheme)
662 *   {
664 *   switch (scheme)
665 *   {
666 *   case stage_3_order_3:
668 *   break;
669 *   case stage_5_order_4:
671 *   break;
672 *   case stage_7_order_4:
674 *   break;
675 *   case stage_9_order_5:
677 *   break;
678 *  
679 *   default:
680 *   AssertThrow(false, ExcNotImplemented());
681 *   }
684 *   rk_integrator(lsrk);
685 *   std::vector<double> ci; // not used
686 *   rk_integrator.get_coefficients(ai, bi, ci);
687 *   }
688 *  
689 *   unsigned int n_stages() const
690 *   {
691 *   return bi.size();
692 *   }
693 *  
694 *   template <typename VectorType, typename Operator>
695 *   void perform_time_step(const Operator &pde_operator,
696 *   const double current_time,
697 *   const double time_step,
698 *   VectorType &solution,
699 *   VectorType &vec_ri,
700 *   VectorType &vec_ki) const
701 *   {
702 *   AssertDimension(ai.size() + 1, bi.size());
703 *  
704 *   vec_ki.swap(solution);
705 *  
706 *   double sum_previous_bi = 0;
707 *   for (unsigned int stage = 0; stage < bi.size(); ++stage)
708 *   {
709 *   const double c_i = stage == 0 ? 0 : sum_previous_bi + ai[stage - 1];
710 *  
711 *   pde_operator.perform_stage(stage,
712 *   current_time + c_i * time_step,
713 *   bi[stage] * time_step,
714 *   (stage == bi.size() - 1 ?
715 *   0 :
716 *   ai[stage] * time_step),
717 *   (stage % 2 == 0 ? vec_ki : vec_ri),
718 *   (stage % 2 == 0 ? vec_ri : vec_ki),
719 *   solution);
720 *  
721 *   if (stage > 0)
722 *   sum_previous_bi += bi[stage - 1];
723 *   }
724 *   }
725 *  
726 *   private:
727 *   std::vector<double> bi;
728 *   std::vector<double> ai;
729 *   };
730 *  
731 *  
732 * @endcode
733 *
734 * Euler-specific utility functions from @ref step_67 "step-67":
735 *
736 * @code
737 *   enum EulerNumericalFlux
738 *   {
739 *   lax_friedrichs_modified,
740 *   harten_lax_vanleer,
741 *   };
742 *   constexpr EulerNumericalFlux numerical_flux_type = lax_friedrichs_modified;
743 *  
744 *  
745 *  
746 *   template <int dim>
747 *   class ExactSolution : public Function<dim>
748 *   {
749 *   public:
750 *   ExactSolution(const double time)
751 *   : Function<dim>(dim + 2, time)
752 *   {}
753 *  
754 *   virtual double value(const Point<dim> &p,
755 *   const unsigned int component = 0) const override;
756 *   };
757 *  
758 *  
759 *  
760 *   template <int dim>
761 *   double ExactSolution<dim>::value(const Point<dim> &x,
762 *   const unsigned int component) const
763 *   {
764 *   const double t = this->get_time();
765 *  
766 *   switch (testcase)
767 *   {
768 *   case 0:
769 *   {
770 *   Assert(dim == 2, ExcNotImplemented());
771 *   const double beta = 5;
772 *  
773 *   Point<dim> x0;
774 *   x0[0] = 5.;
775 *   const double radius_sqr =
776 *   (x - x0).norm_square() - 2. * (x[0] - x0[0]) * t + t * t;
777 *   const double factor =
778 *   beta / (numbers::PI * 2) * std::exp(1. - radius_sqr);
779 *   const double density_log = std::log2(
780 *   std::abs(1. - (gamma - 1.) / gamma * 0.25 * factor * factor));
781 *   const double density = std::exp2(density_log * (1. / (gamma - 1.)));
782 *   const double u = 1. - factor * (x[1] - x0[1]);
783 *   const double v = factor * (x[0] - t - x0[0]);
784 *  
785 *   if (component == 0)
786 *   return density;
787 *   else if (component == 1)
788 *   return density * u;
789 *   else if (component == 2)
790 *   return density * v;
791 *   else
792 *   {
793 *   const double pressure =
794 *   std::exp2(density_log * (gamma / (gamma - 1.)));
795 *   return pressure / (gamma - 1.) +
796 *   0.5 * (density * u * u + density * v * v);
797 *   }
798 *   }
799 *  
800 *   case 1:
801 *   {
802 *   if (component == 0)
803 *   return 1.;
804 *   else if (component == 1)
805 *   return 0.4;
806 *   else if (component == dim + 1)
807 *   return 3.097857142857143;
808 *   else
809 *   return 0.;
810 *   }
811 *  
812 *   default:
814 *   return 0.;
815 *   }
816 *   }
817 *  
818 *  
819 *  
820 *   template <int dim, typename Number>
821 *   inline DEAL_II_ALWAYS_INLINE
823 *   euler_velocity(const Tensor<1, dim + 2, Number> &conserved_variables)
824 *   {
825 *   const Number inverse_density = Number(1.) / conserved_variables[0];
826 *  
827 *   Tensor<1, dim, Number> velocity;
828 *   for (unsigned int d = 0; d < dim; ++d)
829 *   velocity[d] = conserved_variables[1 + d] * inverse_density;
830 *  
831 *   return velocity;
832 *   }
833 *  
834 *   template <int dim, typename Number>
835 *   inline DEAL_II_ALWAYS_INLINE
836 *   Number
837 *   euler_pressure(const Tensor<1, dim + 2, Number> &conserved_variables)
838 *   {
839 *   const Tensor<1, dim, Number> velocity =
840 *   euler_velocity<dim>(conserved_variables);
841 *  
842 *   Number rho_u_dot_u = conserved_variables[1] * velocity[0];
843 *   for (unsigned int d = 1; d < dim; ++d)
844 *   rho_u_dot_u += conserved_variables[1 + d] * velocity[d];
845 *  
846 *   return (gamma - 1.) * (conserved_variables[dim + 1] - 0.5 * rho_u_dot_u);
847 *   }
848 *  
849 *   template <int dim, typename Number>
850 *   inline DEAL_II_ALWAYS_INLINE
852 *   euler_flux(const Tensor<1, dim + 2, Number> &conserved_variables)
853 *   {
854 *   const Tensor<1, dim, Number> velocity =
855 *   euler_velocity<dim>(conserved_variables);
856 *   const Number pressure = euler_pressure<dim>(conserved_variables);
857 *  
859 *   for (unsigned int d = 0; d < dim; ++d)
860 *   {
861 *   flux[0][d] = conserved_variables[1 + d];
862 *   for (unsigned int e = 0; e < dim; ++e)
863 *   flux[e + 1][d] = conserved_variables[e + 1] * velocity[d];
864 *   flux[d + 1][d] += pressure;
865 *   flux[dim + 1][d] =
866 *   velocity[d] * (conserved_variables[dim + 1] + pressure);
867 *   }
868 *  
869 *   return flux;
870 *   }
871 *  
872 *   template <int n_components, int dim, typename Number>
873 *   inline DEAL_II_ALWAYS_INLINE
875 *   operator*(const Tensor<1, n_components, Tensor<1, dim, Number>> &matrix,
876 *   const Tensor<1, dim, Number> &vector)
877 *   {
879 *   for (unsigned int d = 0; d < n_components; ++d)
880 *   result[d] = matrix[d] * vector;
881 *   return result;
882 *   }
883 *  
884 *   template <int dim, typename Number>
885 *   inline DEAL_II_ALWAYS_INLINE
887 *   euler_numerical_flux(const Tensor<1, dim + 2, Number> &u_m,
888 *   const Tensor<1, dim + 2, Number> &u_p,
889 *   const Tensor<1, dim, Number> &normal)
890 *   {
891 *   const auto velocity_m = euler_velocity<dim>(u_m);
892 *   const auto velocity_p = euler_velocity<dim>(u_p);
893 *  
894 *   const auto pressure_m = euler_pressure<dim>(u_m);
895 *   const auto pressure_p = euler_pressure<dim>(u_p);
896 *  
897 *   const auto flux_m = euler_flux<dim>(u_m);
898 *   const auto flux_p = euler_flux<dim>(u_p);
899 *  
900 *   switch (numerical_flux_type)
901 *   {
902 *   case lax_friedrichs_modified:
903 *   {
904 *   const auto lambda =
905 *   0.5 * std::sqrt(std::max(velocity_p.norm_square() +
906 *   gamma * pressure_p * (1. / u_p[0]),
907 *   velocity_m.norm_square() +
908 *   gamma * pressure_m * (1. / u_m[0])));
909 *  
910 *   return 0.5 * (flux_m * normal + flux_p * normal) +
911 *   0.5 * lambda * (u_m - u_p);
912 *   }
913 *  
914 *   case harten_lax_vanleer:
915 *   {
916 *   const auto avg_velocity_normal =
917 *   0.5 * ((velocity_m + velocity_p) * normal);
918 *   const auto avg_c = std::sqrt(std::abs(
919 *   0.5 * gamma *
920 *   (pressure_p * (1. / u_p[0]) + pressure_m * (1. / u_m[0]))));
921 *   const Number s_pos =
922 *   std::max(Number(), avg_velocity_normal + avg_c);
923 *   const Number s_neg =
924 *   std::min(Number(), avg_velocity_normal - avg_c);
925 *   const Number inverse_s = Number(1.) / (s_pos - s_neg);
926 *  
927 *   return inverse_s *
928 *   ((s_pos * (flux_m * normal) - s_neg * (flux_p * normal)) -
929 *   s_pos * s_neg * (u_m - u_p));
930 *   }
931 *  
932 *   default:
933 *   {
935 *   return {};
936 *   }
937 *   }
938 *   }
939 *  
940 *  
941 *  
942 * @endcode
943 *
944 * General-purpose utility functions from @ref step_67 "step-67":
945 *
946 * @code
947 *   template <int dim, typename VectorizedArrayType>
948 *   VectorizedArrayType
949 *   evaluate_function(const Function<dim> &function,
950 *   const Point<dim, VectorizedArrayType> &p_vectorized,
951 *   const unsigned int component)
952 *   {
953 *   VectorizedArrayType result;
954 *   for (unsigned int v = 0; v < VectorizedArrayType::size(); ++v)
955 *   {
956 *   Point<dim> p;
957 *   for (unsigned int d = 0; d < dim; ++d)
958 *   p[d] = p_vectorized[d][v];
959 *   result[v] = function.value(p, component);
960 *   }
961 *   return result;
962 *   }
963 *  
964 *  
965 *   template <int dim, typename VectorizedArrayType, int n_components = dim + 2>
967 *   evaluate_function(const Function<dim> &function,
968 *   const Point<dim, VectorizedArrayType> &p_vectorized)
969 *   {
970 *   AssertDimension(function.n_components, n_components);
972 *   for (unsigned int v = 0; v < VectorizedArrayType::size(); ++v)
973 *   {
974 *   Point<dim> p;
975 *   for (unsigned int d = 0; d < dim; ++d)
976 *   p[d] = p_vectorized[d][v];
977 *   for (unsigned int d = 0; d < n_components; ++d)
978 *   result[d][v] = function.value(p, d);
979 *   }
980 *   return result;
981 *   }
982 *  
983 *  
984 * @endcode
985 *
986 *
987 * <a name="step_76-EuleroperatorusingacellcentricloopandMPI30sharedmemory"></a>
988 * <h3>Euler operator using a cell-centric loop and MPI-3.0 shared memory</h3>
989 *
990
991 *
992 * Euler operator from @ref step_67 "step-67" with some changes as detailed below:
993 *
994 * @code
995 *   template <int dim, int degree, int n_points_1d>
996 *   class EulerOperator
997 *   {
998 *   public:
999 *   static constexpr unsigned int n_quadrature_points_1d = n_points_1d;
1000 *  
1001 *   EulerOperator(TimerOutput &timer_output);
1002 *  
1003 *   ~EulerOperator();
1004 *  
1005 *   void reinit(const Mapping<dim> &mapping,
1006 *   const DoFHandler<dim> &dof_handler);
1007 *  
1008 *   void set_inflow_boundary(const types::boundary_id boundary_id,
1009 *   std::unique_ptr<Function<dim>> inflow_function);
1010 *  
1011 *   void set_subsonic_outflow_boundary(
1012 *   const types::boundary_id boundary_id,
1013 *   std::unique_ptr<Function<dim>> outflow_energy);
1014 *  
1015 *   void set_wall_boundary(const types::boundary_id boundary_id);
1016 *  
1017 *   void set_body_force(std::unique_ptr<Function<dim>> body_force);
1018 *  
1019 *   void
1020 *   perform_stage(const unsigned int stage,
1021 *   const Number cur_time,
1022 *   const Number bi,
1023 *   const Number ai,
1024 *   const LinearAlgebra::distributed::Vector<Number> &current_ri,
1026 *   LinearAlgebra::distributed::Vector<Number> &solution) const;
1027 *  
1028 *   void project(const Function<dim> &function,
1029 *   LinearAlgebra::distributed::Vector<Number> &solution) const;
1030 *  
1031 *   std::array<double, 3> compute_errors(
1032 *   const Function<dim> &function,
1033 *   const LinearAlgebra::distributed::Vector<Number> &solution) const;
1034 *  
1035 *   double compute_cell_transport_speed(
1036 *   const LinearAlgebra::distributed::Vector<Number> &solution) const;
1037 *  
1038 *   void
1039 *   initialize_vector(LinearAlgebra::distributed::Vector<Number> &vector) const;
1040 *  
1041 *   private:
1042 * @endcode
1043 *
1044 * Instance of SubCommunicatorWrapper containing the sub-communicator, which
1045 * we need to pass to MatrixFree::reinit() to be able to exploit MPI-3.0
1046 * shared-memory capabilities:
1047 *
1048 * @code
1049 *   MPI_Comm subcommunicator;
1050 *  
1051 *   MatrixFree<dim, Number, VectorizedArrayType> data;
1052 *  
1053 *   TimerOutput &timer;
1054 *  
1055 *   std::map<types::boundary_id, std::unique_ptr<Function<dim>>>
1056 *   inflow_boundaries;
1057 *   std::map<types::boundary_id, std::unique_ptr<Function<dim>>>
1058 *   subsonic_outflow_boundaries;
1059 *   std::set<types::boundary_id> wall_boundaries;
1060 *   std::unique_ptr<Function<dim>> body_force;
1061 *   };
1062 *  
1063 *  
1064 *  
1065 * @endcode
1066 *
1067 * New constructor, which creates a sub-communicator. The user can specify
1068 * the size of the sub-communicator via the global parameter group_size. If
1069 * the size is set to -1, all MPI processes of a
1070 * shared-memory domain are combined to a group. The specified size is
1071 * decisive for the benefit of the shared-memory capabilities of MatrixFree
1072 * and, therefore, setting the <code>size</code> to <code>-1</code> is a
1073 * reasonable choice. By setting, the size to <code>1</code> users explicitly
1074 * disable the MPI-3.0 shared-memory features of MatrixFree and rely
1075 * completely on MPI-2.0 features, like <code>MPI_Isend</code> and
1076 * <code>MPI_Irecv</code>.
1077 *
1078 * @code
1079 *   template <int dim, int degree, int n_points_1d>
1080 *   EulerOperator<dim, degree, n_points_1d>::EulerOperator(TimerOutput &timer)
1081 *   : timer(timer)
1082 *   {
1083 *   #ifdef DEAL_II_WITH_MPI
1084 *   if (group_size == 1)
1085 *   {
1086 *   this->subcommunicator = MPI_COMM_SELF;
1087 *   }
1088 *   else if (group_size == numbers::invalid_unsigned_int)
1089 *   {
1090 *   const auto rank = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
1091 *  
1092 *   MPI_Comm_split_type(MPI_COMM_WORLD,
1093 *   MPI_COMM_TYPE_SHARED,
1094 *   rank,
1095 *   MPI_INFO_NULL,
1096 *   &subcommunicator);
1097 *   }
1098 *   else
1099 *   {
1101 *   }
1102 *   #else
1103 *   (void)subcommunicator;
1104 *   (void)group_size;
1105 *   this->subcommunicator = MPI_COMM_SELF;
1106 *   #endif
1107 *   }
1108 *  
1109 *  
1110 * @endcode
1111 *
1112 * New destructor responsible for freeing of the sub-communicator.
1113 *
1114 * @code
1115 *   template <int dim, int degree, int n_points_1d>
1116 *   EulerOperator<dim, degree, n_points_1d>::~EulerOperator()
1117 *   {
1118 *   #ifdef DEAL_II_WITH_MPI
1119 *   if (this->subcommunicator != MPI_COMM_SELF)
1120 *   MPI_Comm_free(&subcommunicator);
1121 *   #endif
1122 *   }
1123 *  
1124 *  
1125 * @endcode
1126 *
1127 * Modified reinit() function to set up the internal data structures in
1128 * MatrixFree in a way that it is usable by the cell-centric loops and
1129 * the MPI-3.0 shared-memory capabilities are used:
1130 *
1131 * @code
1132 *   template <int dim, int degree, int n_points_1d>
1133 *   void EulerOperator<dim, degree, n_points_1d>::reinit(
1134 *   const Mapping<dim> &mapping,
1135 *   const DoFHandler<dim> &dof_handler)
1136 *   {
1137 *   const std::vector<const DoFHandler<dim> *> dof_handlers = {&dof_handler};
1138 *   const AffineConstraints<double> dummy;
1139 *   const std::vector<const AffineConstraints<double> *> constraints = {&dummy};
1140 *   const std::vector<Quadrature<1>> quadratures = {QGauss<1>(n_q_points_1d),
1141 *   QGauss<1>(fe_degree + 1)};
1142 *  
1144 *   additional_data;
1145 *   additional_data.mapping_update_flags =
1147 *   update_values);
1148 *   additional_data.mapping_update_flags_inner_faces =
1150 *   update_values);
1151 *   additional_data.mapping_update_flags_boundary_faces =
1153 *   update_values);
1154 *   additional_data.tasks_parallel_scheme =
1156 *  
1157 * @endcode
1158 *
1159 * Categorize cells so that all lanes have the same boundary IDs for each
1160 * face. This is strictly not necessary, however, allows to write simpler
1161 * code in EulerOperator::perform_stage() without masking, since it is
1162 * guaranteed that all cells grouped together (in a VectorizedArray)
1163 * have to perform exactly the same operation also on the faces.
1164 *
1165 * @code
1166 *   MatrixFreeTools::categorize_by_boundary_ids(dof_handler.get_triangulation(),
1167 *   additional_data);
1168 *  
1169 * @endcode
1170 *
1171 * Enable MPI-3.0 shared-memory capabilities within MatrixFree by providing
1172 * the sub-communicator:
1173 *
1174 * @code
1175 *   additional_data.communicator_sm = subcommunicator;
1176 *  
1177 *   data.reinit(
1178 *   mapping, dof_handlers, constraints, quadratures, additional_data);
1179 *   }
1180 *  
1181 *  
1182 * @endcode
1183 *
1184 * The following function does an entire stage of a Runge--Kutta update
1185 * and is, alongside the slightly modified setup, the heart of this tutorial
1186 * compared to @ref step_67 "step-67".
1187 *
1188
1189 *
1190 * In contrast to @ref step_67 "step-67", we are not executing the advection step
1191 * (using MatrixFree::loop()) and the inverse mass-matrix step
1192 * (using MatrixFree::cell_loop()) in sequence, but evaluate everything in
1193 * one go inside of MatrixFree::loop_cell_centric(). This function expects
1194 * a single function that is executed on each locally-owned (macro) cell as
1195 * parameter so that we need to loop over all faces of that cell and perform
1196 * needed integration steps on our own.
1197 *
1198
1199 *
1200 * The following function contains to a large extent copies of the following
1201 * functions from @ref step_67 "step-67" so that comments related the evaluation of the weak
1202 * form are skipped here:
1203 * - <code>EulerDG::EulerOperator::local_apply_cell</code>
1204 * - <code>EulerDG::EulerOperator::local_apply_face</code>
1205 * - <code>EulerDG::EulerOperator::local_apply_boundary_face</code>
1206 * - <code>EulerDG::EulerOperator::local_apply_inverse_mass_matrix</code>
1207 *
1208 * @code
1209 *   template <int dim, int degree, int n_points_1d>
1210 *   void EulerOperator<dim, degree, n_points_1d>::perform_stage(
1211 *   const unsigned int stage,
1212 *   const Number current_time,
1213 *   const Number bi,
1214 *   const Number ai,
1215 *   const LinearAlgebra::distributed::Vector<Number> &current_ri,
1216 *   LinearAlgebra::distributed::Vector<Number> &vec_ki,
1217 *   LinearAlgebra::distributed::Vector<Number> &solution) const
1218 *   {
1219 *   for (auto &i : inflow_boundaries)
1220 *   i.second->set_time(current_time);
1221 *   for (auto &i : subsonic_outflow_boundaries)
1222 *   i.second->set_time(current_time);
1223 *  
1224 * @endcode
1225 *
1226 * Run a cell-centric loop by calling MatrixFree::loop_cell_centric() and
1227 * providing a lambda containing the effects of the cell, face and
1228 * boundary-face integrals:
1229 *
1230 * @code
1231 *   data.template loop_cell_centric<LinearAlgebra::distributed::Vector<Number>,
1232 *   LinearAlgebra::distributed::Vector<Number>>(
1233 *   [&](const auto &data, auto &dst, const auto &src, const auto cell_range) {
1234 *   using FECellIntegral = FEEvaluation<dim,
1235 *   degree,
1236 *   n_points_1d,
1237 *   dim + 2,
1238 *   Number,
1239 *   VectorizedArrayType>;
1240 *   using FEFaceIntegral = FEFaceEvaluation<dim,
1241 *   degree,
1242 *   n_points_1d,
1243 *   dim + 2,
1244 *   Number,
1245 *   VectorizedArrayType>;
1246 *  
1247 *   FECellIntegral phi(data);
1248 *   FECellIntegral phi_temp(data);
1249 *   FEFaceIntegral phi_m(data, true);
1250 *   FEFaceIntegral phi_p(data, false);
1251 *  
1252 *   Tensor<1, dim, VectorizedArrayType> constant_body_force;
1253 *   const Functions::ConstantFunction<dim> *constant_function =
1254 *   dynamic_cast<Functions::ConstantFunction<dim> *>(body_force.get());
1255 *  
1256 *   if (constant_function)
1257 *   constant_body_force =
1258 *   evaluate_function<dim, VectorizedArrayType, dim>(
1259 *   *constant_function, Point<dim, VectorizedArrayType>());
1260 *  
1263 *   dim,
1264 *   n_points_1d,
1265 *   n_points_1d,
1266 *   VectorizedArrayType,
1267 *   Number>
1268 *   eval({},
1269 *   data.get_shape_info().data[0].shape_gradients_collocation_eo,
1270 *   {});
1271 *  
1272 *   AlignedVector<VectorizedArrayType> buffer(phi.static_n_q_points *
1273 *   phi.n_components);
1274 *  
1275 * @endcode
1276 *
1277 * Loop over all cell batches:
1278 *
1279 * @code
1280 *   for (unsigned int cell = cell_range.first; cell < cell_range.second;
1281 *   ++cell)
1282 *   {
1283 *   phi.reinit(cell);
1284 *  
1285 *   if (ai != Number())
1286 *   phi_temp.reinit(cell);
1287 *  
1288 * @endcode
1289 *
1290 * Read values from global vector and compute the values at the
1291 * quadrature points:
1292 *
1293 * @code
1294 *   if (ai != Number() && stage == 0)
1295 *   {
1296 *   phi.read_dof_values(src);
1297 *  
1298 *   for (unsigned int i = 0;
1299 *   i < phi.static_dofs_per_component * (dim + 2);
1300 *   ++i)
1301 *   phi_temp.begin_dof_values()[i] = phi.begin_dof_values()[i];
1302 *  
1303 *   phi.evaluate(EvaluationFlags::values);
1304 *   }
1305 *   else
1306 *   {
1307 *   phi.gather_evaluate(src, EvaluationFlags::values);
1308 *   }
1309 *  
1310 * @endcode
1311 *
1312 * Buffer the computed values at the quadrature points, since
1313 * these are overridden by FEEvaluation::submit_value() in the next
1314 * step, however, are needed later on for the face integrals:
1315 *
1316 * @code
1317 *   for (unsigned int i = 0; i < phi.static_n_q_points * (dim + 2); ++i)
1318 *   buffer[i] = phi.begin_values()[i];
1319 *  
1320 * @endcode
1321 *
1322 * Apply the cell integral at the cell quadrature points. See also
1323 * the function <code>EulerOperator::local_apply_cell()</code> from
1324 * @ref step_67 "step-67":
1325 *
1326 * @code
1327 *   for (const unsigned int q : phi.quadrature_point_indices())
1328 *   {
1329 *   const auto w_q = phi.get_value(q);
1330 *   phi.submit_gradient(euler_flux<dim>(w_q), q);
1331 *   if (body_force.get() != nullptr)
1332 *   {
1333 *   const Tensor<1, dim, VectorizedArrayType> force =
1334 *   constant_function ?
1335 *   constant_body_force :
1336 *   evaluate_function<dim, VectorizedArrayType, dim>(
1337 *   *body_force, phi.quadrature_point(q));
1338 *  
1340 *   for (unsigned int d = 0; d < dim; ++d)
1341 *   forcing[d + 1] = w_q[0] * force[d];
1342 *   for (unsigned int d = 0; d < dim; ++d)
1343 *   forcing[dim + 1] += force[d] * w_q[d + 1];
1344 *  
1345 *   phi.submit_value(forcing, q);
1346 *   }
1347 *   }
1348 *  
1349 * @endcode
1350 *
1351 * Test with the gradient of the test functions in the quadrature
1352 * points. We skip the interpolation back to the support points
1353 * of the element, since we first collect all contributions in the
1354 * cell quadrature points and only perform the interpolation back
1355 * as the final step.
1356 *
1357 * @code
1358 *   {
1359 *   auto *values_ptr = phi.begin_values();
1360 *   auto *gradient_ptr = phi.begin_gradients();
1361 *  
1362 *   for (unsigned int c = 0; c < dim + 2; ++c)
1363 *   {
1364 *   if (dim >= 1 && body_force.get() == nullptr)
1365 *   eval.template gradients<0, false, false, dim>(gradient_ptr,
1366 *   values_ptr);
1367 *   else if (dim >= 1)
1368 *   eval.template gradients<0, false, true, dim>(gradient_ptr,
1369 *   values_ptr);
1370 *   if (dim >= 2)
1371 *   eval.template gradients<1, false, true, dim>(gradient_ptr +
1372 *   1,
1373 *   values_ptr);
1374 *   if (dim >= 3)
1375 *   eval.template gradients<2, false, true, dim>(gradient_ptr +
1376 *   2,
1377 *   values_ptr);
1378 *  
1379 *   values_ptr += phi.static_n_q_points;
1380 *   gradient_ptr += phi.static_n_q_points * dim;
1381 *   }
1382 *   }
1383 *  
1384 * @endcode
1385 *
1386 * Loop over all faces of the current cell:
1387 *
1388 * @code
1389 *   for (unsigned int face = 0;
1390 *   face < GeometryInfo<dim>::faces_per_cell;
1391 *   ++face)
1392 *   {
1393 * @endcode
1394 *
1395 * Determine the boundary ID of the current face. Since we have
1396 * set up MatrixFree in a way that all filled lanes have
1397 * guaranteed the same boundary ID, we can select the
1398 * boundary ID of the first lane.
1399 *
1400 * @code
1401 *   const auto boundary_ids =
1402 *   data.get_faces_by_cells_boundary_id(cell, face);
1403 *  
1404 *   Assert(std::equal(boundary_ids.begin(),
1405 *   boundary_ids.begin() +
1406 *   data.n_active_entries_per_cell_batch(cell),
1407 *   boundary_ids.begin()),
1408 *   ExcMessage("Boundary IDs of lanes differ."));
1409 *  
1410 *   const auto boundary_id = boundary_ids[0];
1411 *  
1412 *   phi_m.reinit(cell, face);
1413 *  
1414 * @endcode
1415 *
1416 * Interpolate the values from the cell quadrature points to the
1417 * quadrature points of the current face via a simple 1d
1418 * interpolation:
1419 *
1420 * @code
1422 *   n_points_1d - 1,
1423 *   VectorizedArrayType>::
1424 *   template interpolate_quadrature<true, false>(
1425 *   dim + 2,
1427 *   data.get_shape_info(),
1428 *   buffer.data(),
1429 *   phi_m.begin_values(),
1430 *   face);
1431 *  
1432 * @endcode
1433 *
1434 * Check if the face is an internal or a boundary face and
1435 * select a different code path based on this information:
1436 *
1437 * @code
1438 *   if (boundary_id == numbers::internal_face_boundary_id)
1439 *   {
1440 * @endcode
1441 *
1442 * Process and internal face. The following lines of code
1443 * are a copy of the function
1444 * <code>EulerDG::EulerOperator::local_apply_face</code>
1445 * from @ref step_67 "step-67":
1446 *
1447 * @code
1448 *   phi_p.reinit(cell, face);
1449 *   phi_p.gather_evaluate(src, EvaluationFlags::values);
1450 *  
1451 *   for (const unsigned int q :
1452 *   phi_m.quadrature_point_indices())
1453 *   {
1454 *   const auto numerical_flux =
1455 *   euler_numerical_flux<dim>(phi_m.get_value(q),
1456 *   phi_p.get_value(q),
1457 *   phi_m.normal_vector(q));
1458 *   phi_m.submit_value(-numerical_flux, q);
1459 *   }
1460 *   }
1461 *   else
1462 *   {
1463 * @endcode
1464 *
1465 * Process a boundary face. These following lines of code
1466 * are a copy of the function
1467 * <code>EulerDG::EulerOperator::local_apply_boundary_face</code>
1468 * from @ref step_67 "step-67":
1469 *
1470 * @code
1471 *   for (const unsigned int q :
1472 *   phi_m.quadrature_point_indices())
1473 *   {
1474 *   const auto w_m = phi_m.get_value(q);
1475 *   const auto normal = phi_m.normal_vector(q);
1476 *  
1477 *   auto rho_u_dot_n = w_m[1] * normal[0];
1478 *   for (unsigned int d = 1; d < dim; ++d)
1479 *   rho_u_dot_n += w_m[1 + d] * normal[d];
1480 *  
1481 *   bool at_outflow = false;
1482 *  
1484 *  
1485 *   if (wall_boundaries.find(boundary_id) !=
1486 *   wall_boundaries.end())
1487 *   {
1488 *   w_p[0] = w_m[0];
1489 *   for (unsigned int d = 0; d < dim; ++d)
1490 *   w_p[d + 1] =
1491 *   w_m[d + 1] - 2. * rho_u_dot_n * normal[d];
1492 *   w_p[dim + 1] = w_m[dim + 1];
1493 *   }
1494 *   else if (inflow_boundaries.find(boundary_id) !=
1495 *   inflow_boundaries.end())
1496 *   w_p = evaluate_function(
1497 *   *inflow_boundaries.find(boundary_id)->second,
1498 *   phi_m.quadrature_point(q));
1499 *   else if (subsonic_outflow_boundaries.find(
1500 *   boundary_id) !=
1501 *   subsonic_outflow_boundaries.end())
1502 *   {
1503 *   w_p = w_m;
1504 *   w_p[dim + 1] =
1505 *   evaluate_function(*subsonic_outflow_boundaries
1506 *   .find(boundary_id)
1507 *   ->second,
1508 *   phi_m.quadrature_point(q),
1509 *   dim + 1);
1510 *   at_outflow = true;
1511 *   }
1512 *   else
1513 *   AssertThrow(false,
1514 *   ExcMessage(
1515 *   "Unknown boundary id, did "
1516 *   "you set a boundary condition for "
1517 *   "this part of the domain boundary?"));
1518 *  
1519 *   auto flux = euler_numerical_flux<dim>(w_m, w_p, normal);
1520 *  
1521 *   if (at_outflow)
1522 *   for (unsigned int v = 0;
1523 *   v < VectorizedArrayType::size();
1524 *   ++v)
1525 *   {
1526 *   if (rho_u_dot_n[v] < -1e-12)
1527 *   for (unsigned int d = 0; d < dim; ++d)
1528 *   flux[d + 1][v] = 0.;
1529 *   }
1530 *  
1531 *   phi_m.submit_value(-flux, q);
1532 *   }
1533 *   }
1534 *  
1535 * @endcode
1536 *
1537 * Evaluate local integrals related to cell by quadrature and
1538 * add into cell contribution via a simple 1d interpolation:
1539 *
1540 * @code
1542 *   n_points_1d - 1,
1543 *   VectorizedArrayType>::
1544 *   template interpolate_quadrature<false, true>(
1545 *   dim + 2,
1547 *   data.get_shape_info(),
1548 *   phi_m.begin_values(),
1549 *   phi.begin_values(),
1550 *   face);
1551 *   }
1552 *  
1553 * @endcode
1554 *
1555 * Apply inverse mass matrix in the cell quadrature points. See
1556 * also the function
1557 * <code>EulerDG::EulerOperator::local_apply_inverse_mass_matrix()</code>
1558 * from @ref step_67 "step-67":
1559 *
1560 * @code
1561 *   for (unsigned int q = 0; q < phi.static_n_q_points; ++q)
1562 *   {
1563 *   const auto factor = VectorizedArrayType(1.0) / phi.JxW(q);
1564 *   for (unsigned int c = 0; c < dim + 2; ++c)
1565 *   phi.begin_values()[c * phi.static_n_q_points + q] =
1566 *   phi.begin_values()[c * phi.static_n_q_points + q] * factor;
1567 *   }
1568 *  
1569 * @endcode
1570 *
1571 * Transform values from collocation space to the original
1572 * Gauss-Lobatto space:
1573 *
1574 * @code
1578 *   dim,
1579 *   degree + 1,
1580 *   n_points_1d>::do_backward(dim + 2,
1581 *   data.get_shape_info()
1582 *   .data[0]
1583 *   .inverse_shape_values_eo,
1584 *   false,
1585 *   phi.begin_values(),
1586 *   phi.begin_dof_values());
1587 *  
1588 * @endcode
1589 *
1590 * Perform Runge-Kutta update and write results back to global
1591 * vectors:
1592 *
1593 * @code
1594 *   if (ai == Number())
1595 *   {
1596 *   for (unsigned int q = 0; q < phi.static_dofs_per_cell; ++q)
1597 *   phi.begin_dof_values()[q] = bi * phi.begin_dof_values()[q];
1598 *   phi.distribute_local_to_global(solution);
1599 *   }
1600 *   else
1601 *   {
1602 *   if (stage != 0)
1603 *   phi_temp.read_dof_values(solution);
1604 *  
1605 *   for (unsigned int q = 0; q < phi.static_dofs_per_cell; ++q)
1606 *   {
1607 *   const auto K_i = phi.begin_dof_values()[q];
1608 *  
1609 *   phi.begin_dof_values()[q] =
1610 *   phi_temp.begin_dof_values()[q] + (ai * K_i);
1611 *  
1612 *   phi_temp.begin_dof_values()[q] += bi * K_i;
1613 *   }
1614 *   phi.set_dof_values(dst);
1615 *   phi_temp.set_dof_values(solution);
1616 *   }
1617 *   }
1618 *   },
1619 *   vec_ki,
1620 *   current_ri,
1621 *   true,
1623 *   }
1624 *  
1625 *  
1626 *  
1627 * @endcode
1628 *
1629 * From here, the code of @ref step_67 "step-67" has not changed.
1630 *
1631 * @code
1632 *   template <int dim, int degree, int n_points_1d>
1633 *   void EulerOperator<dim, degree, n_points_1d>::initialize_vector(
1635 *   {
1636 *   data.initialize_dof_vector(vector);
1637 *   }
1638 *  
1639 *  
1640 *  
1641 *   template <int dim, int degree, int n_points_1d>
1642 *   void EulerOperator<dim, degree, n_points_1d>::set_inflow_boundary(
1643 *   const types::boundary_id boundary_id,
1644 *   std::unique_ptr<Function<dim>> inflow_function)
1645 *   {
1646 *   AssertThrow(subsonic_outflow_boundaries.find(boundary_id) ==
1647 *   subsonic_outflow_boundaries.end() &&
1648 *   wall_boundaries.find(boundary_id) == wall_boundaries.end(),
1649 *   ExcMessage("You already set the boundary with id " +
1650 *   std::to_string(static_cast<int>(boundary_id)) +
1651 *   " to another type of boundary before now setting " +
1652 *   "it as inflow"));
1653 *   AssertThrow(inflow_function->n_components == dim + 2,
1654 *   ExcMessage("Expected function with dim+2 components"));
1655 *  
1656 *   inflow_boundaries[boundary_id] = std::move(inflow_function);
1657 *   }
1658 *  
1659 *  
1660 *  
1661 *   template <int dim, int degree, int n_points_1d>
1662 *   void EulerOperator<dim, degree, n_points_1d>::set_subsonic_outflow_boundary(
1663 *   const types::boundary_id boundary_id,
1664 *   std::unique_ptr<Function<dim>> outflow_function)
1665 *   {
1666 *   AssertThrow(inflow_boundaries.find(boundary_id) ==
1667 *   inflow_boundaries.end() &&
1668 *   wall_boundaries.find(boundary_id) == wall_boundaries.end(),
1669 *   ExcMessage("You already set the boundary with id " +
1670 *   std::to_string(static_cast<int>(boundary_id)) +
1671 *   " to another type of boundary before now setting " +
1672 *   "it as subsonic outflow"));
1673 *   AssertThrow(outflow_function->n_components == dim + 2,
1674 *   ExcMessage("Expected function with dim+2 components"));
1675 *  
1676 *   subsonic_outflow_boundaries[boundary_id] = std::move(outflow_function);
1677 *   }
1678 *  
1679 *  
1680 *  
1681 *   template <int dim, int degree, int n_points_1d>
1682 *   void EulerOperator<dim, degree, n_points_1d>::set_wall_boundary(
1683 *   const types::boundary_id boundary_id)
1684 *   {
1685 *   AssertThrow(inflow_boundaries.find(boundary_id) ==
1686 *   inflow_boundaries.end() &&
1687 *   subsonic_outflow_boundaries.find(boundary_id) ==
1688 *   subsonic_outflow_boundaries.end(),
1689 *   ExcMessage("You already set the boundary with id " +
1690 *   std::to_string(static_cast<int>(boundary_id)) +
1691 *   " to another type of boundary before now setting " +
1692 *   "it as wall boundary"));
1693 *  
1694 *   wall_boundaries.insert(boundary_id);
1695 *   }
1696 *  
1697 *  
1698 *  
1699 *   template <int dim, int degree, int n_points_1d>
1700 *   void EulerOperator<dim, degree, n_points_1d>::set_body_force(
1701 *   std::unique_ptr<Function<dim>> body_force)
1702 *   {
1703 *   AssertDimension(body_force->n_components, dim);
1704 *  
1705 *   this->body_force = std::move(body_force);
1706 *   }
1707 *  
1708 *  
1709 *  
1710 *   template <int dim, int degree, int n_points_1d>
1711 *   void EulerOperator<dim, degree, n_points_1d>::project(
1712 *   const Function<dim> &function,
1713 *   LinearAlgebra::distributed::Vector<Number> &solution) const
1714 *   {
1716 *   phi(data, 0, 1);
1718 *   degree,
1719 *   dim + 2,
1720 *   Number,
1721 *   VectorizedArrayType>
1722 *   inverse(phi);
1723 *   solution.zero_out_ghost_values();
1724 *   for (unsigned int cell = 0; cell < data.n_cell_batches(); ++cell)
1725 *   {
1726 *   phi.reinit(cell);
1727 *   for (const unsigned int q : phi.quadrature_point_indices())
1728 *   phi.submit_dof_value(evaluate_function(function,
1729 *   phi.quadrature_point(q)),
1730 *   q);
1731 *   inverse.transform_from_q_points_to_basis(dim + 2,
1732 *   phi.begin_dof_values(),
1733 *   phi.begin_dof_values());
1734 *   phi.set_dof_values(solution);
1735 *   }
1736 *   }
1737 *  
1738 *  
1739 *  
1740 *   template <int dim, int degree, int n_points_1d>
1741 *   std::array<double, 3> EulerOperator<dim, degree, n_points_1d>::compute_errors(
1742 *   const Function<dim> &function,
1743 *   const LinearAlgebra::distributed::Vector<Number> &solution) const
1744 *   {
1745 *   TimerOutput::Scope t(timer, "compute errors");
1746 *   double errors_squared[3] = {};
1748 *   phi(data, 0, 0);
1749 *  
1750 *   for (unsigned int cell = 0; cell < data.n_cell_batches(); ++cell)
1751 *   {
1752 *   phi.reinit(cell);
1753 *   phi.gather_evaluate(solution, EvaluationFlags::values);
1754 *   VectorizedArrayType local_errors_squared[3] = {};
1755 *   for (const unsigned int q : phi.quadrature_point_indices())
1756 *   {
1757 *   const auto error =
1758 *   evaluate_function(function, phi.quadrature_point(q)) -
1759 *   phi.get_value(q);
1760 *   const auto JxW = phi.JxW(q);
1761 *  
1762 *   local_errors_squared[0] += error[0] * error[0] * JxW;
1763 *   for (unsigned int d = 0; d < dim; ++d)
1764 *   local_errors_squared[1] += (error[d + 1] * error[d + 1]) * JxW;
1765 *   local_errors_squared[2] += (error[dim + 1] * error[dim + 1]) * JxW;
1766 *   }
1767 *   for (unsigned int v = 0; v < data.n_active_entries_per_cell_batch(cell);
1768 *   ++v)
1769 *   for (unsigned int d = 0; d < 3; ++d)
1770 *   errors_squared[d] += local_errors_squared[d][v];
1771 *   }
1772 *  
1773 *   Utilities::MPI::sum(errors_squared, MPI_COMM_WORLD, errors_squared);
1774 *  
1775 *   std::array<double, 3> errors;
1776 *   for (unsigned int d = 0; d < 3; ++d)
1777 *   errors[d] = std::sqrt(errors_squared[d]);
1778 *  
1779 *   return errors;
1780 *   }
1781 *  
1782 *  
1783 *  
1784 *   template <int dim, int degree, int n_points_1d>
1785 *   double EulerOperator<dim, degree, n_points_1d>::compute_cell_transport_speed(
1786 *   const LinearAlgebra::distributed::Vector<Number> &solution) const
1787 *   {
1788 *   TimerOutput::Scope t(timer, "compute transport speed");
1789 *   Number max_transport = 0;
1791 *   phi(data, 0, 1);
1792 *  
1793 *   for (unsigned int cell = 0; cell < data.n_cell_batches(); ++cell)
1794 *   {
1795 *   phi.reinit(cell);
1796 *   phi.gather_evaluate(solution, EvaluationFlags::values);
1797 *   VectorizedArrayType local_max = 0.;
1798 *   for (const unsigned int q : phi.quadrature_point_indices())
1799 *   {
1800 *   const auto solution = phi.get_value(q);
1801 *   const auto velocity = euler_velocity<dim>(solution);
1802 *   const auto pressure = euler_pressure<dim>(solution);
1803 *  
1804 *   const auto inverse_jacobian = phi.inverse_jacobian(q);
1805 *   const auto convective_speed = inverse_jacobian * velocity;
1806 *   VectorizedArrayType convective_limit = 0.;
1807 *   for (unsigned int d = 0; d < dim; ++d)
1808 *   convective_limit =
1809 *   std::max(convective_limit, std::abs(convective_speed[d]));
1810 *  
1811 *   const auto speed_of_sound =
1812 *   std::sqrt(gamma * pressure * (1. / solution[0]));
1813 *  
1814 *   Tensor<1, dim, VectorizedArrayType> eigenvector;
1815 *   for (unsigned int d = 0; d < dim; ++d)
1816 *   eigenvector[d] = 1.;
1817 *   for (unsigned int i = 0; i < 5; ++i)
1818 *   {
1819 *   eigenvector = transpose(inverse_jacobian) *
1820 *   (inverse_jacobian * eigenvector);
1821 *   VectorizedArrayType eigenvector_norm = 0.;
1822 *   for (unsigned int d = 0; d < dim; ++d)
1823 *   eigenvector_norm =
1824 *   std::max(eigenvector_norm, std::abs(eigenvector[d]));
1825 *   eigenvector /= eigenvector_norm;
1826 *   }
1827 *   const auto jac_times_ev = inverse_jacobian * eigenvector;
1828 *   const auto max_eigenvalue = std::sqrt(
1829 *   (jac_times_ev * jac_times_ev) / (eigenvector * eigenvector));
1830 *   local_max =
1831 *   std::max(local_max,
1832 *   max_eigenvalue * speed_of_sound + convective_limit);
1833 *   }
1834 *  
1835 *   for (unsigned int v = 0; v < data.n_active_entries_per_cell_batch(cell);
1836 *   ++v)
1837 *   max_transport = std::max(max_transport, local_max[v]);
1838 *   }
1839 *  
1840 *   max_transport = Utilities::MPI::max(max_transport, MPI_COMM_WORLD);
1841 *  
1842 *   return max_transport;
1843 *   }
1844 *  
1845 *  
1846 *  
1847 *   template <int dim>
1848 *   class EulerProblem
1849 *   {
1850 *   public:
1851 *   EulerProblem();
1852 *  
1853 *   void run();
1854 *  
1855 *   private:
1856 *   void make_grid_and_dofs();
1857 *  
1858 *   void output_results(const unsigned int result_number);
1859 *  
1861 *  
1862 *   ConditionalOStream pcout;
1863 *  
1864 *   #ifdef DEAL_II_WITH_P4EST
1866 *   #else
1868 *   #endif
1869 *  
1870 *   const FESystem<dim> fe;
1871 *   const MappingQ<dim> mapping;
1872 *   DoFHandler<dim> dof_handler;
1873 *  
1874 *   TimerOutput timer;
1875 *  
1876 *   EulerOperator<dim, fe_degree, n_q_points_1d> euler_operator;
1877 *  
1878 *   double time, time_step;
1879 *  
1880 *   class Postprocessor : public DataPostprocessor<dim>
1881 *   {
1882 *   public:
1883 *   Postprocessor();
1884 *  
1885 *   virtual void evaluate_vector_field(
1886 *   const DataPostprocessorInputs::Vector<dim> &inputs,
1887 *   std::vector<Vector<double>> &computed_quantities) const override;
1888 *  
1889 *   virtual std::vector<std::string> get_names() const override;
1890 *  
1891 *   virtual std::vector<
1893 *   get_data_component_interpretation() const override;
1894 *  
1895 *   virtual UpdateFlags get_needed_update_flags() const override;
1896 *  
1897 *   private:
1898 *   const bool do_schlieren_plot;
1899 *   };
1900 *   };
1901 *  
1902 *  
1903 *  
1904 *   template <int dim>
1905 *   EulerProblem<dim>::Postprocessor::Postprocessor()
1906 *   : do_schlieren_plot(dim == 2)
1907 *   {}
1908 *  
1909 *  
1910 *  
1911 *   template <int dim>
1912 *   void EulerProblem<dim>::Postprocessor::evaluate_vector_field(
1913 *   const DataPostprocessorInputs::Vector<dim> &inputs,
1914 *   std::vector<Vector<double>> &computed_quantities) const
1915 *   {
1916 *   const unsigned int n_evaluation_points = inputs.solution_values.size();
1917 *  
1918 *   if (do_schlieren_plot == true)
1919 *   Assert(inputs.solution_gradients.size() == n_evaluation_points,
1920 *   ExcInternalError());
1921 *  
1922 *   Assert(computed_quantities.size() == n_evaluation_points,
1923 *   ExcInternalError());
1924 *   Assert(inputs.solution_values[0].size() == dim + 2, ExcInternalError());
1925 *   Assert(computed_quantities[0].size() ==
1926 *   dim + 2 + (do_schlieren_plot == true ? 1 : 0),
1927 *   ExcInternalError());
1928 *  
1929 *   for (unsigned int p = 0; p < n_evaluation_points; ++p)
1930 *   {
1931 *   Tensor<1, dim + 2> solution;
1932 *   for (unsigned int d = 0; d < dim + 2; ++d)
1933 *   solution[d] = inputs.solution_values[p](d);
1934 *  
1935 *   const double density = solution[0];
1936 *   const Tensor<1, dim> velocity = euler_velocity<dim>(solution);
1937 *   const double pressure = euler_pressure<dim>(solution);
1938 *  
1939 *   for (unsigned int d = 0; d < dim; ++d)
1940 *   computed_quantities[p](d) = velocity[d];
1941 *   computed_quantities[p](dim) = pressure;
1942 *   computed_quantities[p](dim + 1) = std::sqrt(gamma * pressure / density);
1943 *  
1944 *   if (do_schlieren_plot == true)
1945 *   computed_quantities[p](dim + 2) =
1946 *   inputs.solution_gradients[p][0] * inputs.solution_gradients[p][0];
1947 *   }
1948 *   }
1949 *  
1950 *  
1951 *  
1952 *   template <int dim>
1953 *   std::vector<std::string> EulerProblem<dim>::Postprocessor::get_names() const
1954 *   {
1955 *   std::vector<std::string> names;
1956 *   for (unsigned int d = 0; d < dim; ++d)
1957 *   names.emplace_back("velocity");
1958 *   names.emplace_back("pressure");
1959 *   names.emplace_back("speed_of_sound");
1960 *  
1961 *   if (do_schlieren_plot == true)
1962 *   names.emplace_back("schlieren_plot");
1963 *  
1964 *   return names;
1965 *   }
1966 *  
1967 *  
1968 *  
1969 *   template <int dim>
1970 *   std::vector<DataComponentInterpretation::DataComponentInterpretation>
1971 *   EulerProblem<dim>::Postprocessor::get_data_component_interpretation() const
1972 *   {
1973 *   std::vector<DataComponentInterpretation::DataComponentInterpretation>
1974 *   interpretation;
1975 *   for (unsigned int d = 0; d < dim; ++d)
1976 *   interpretation.push_back(
1978 *   interpretation.push_back(DataComponentInterpretation::component_is_scalar);
1979 *   interpretation.push_back(DataComponentInterpretation::component_is_scalar);
1980 *  
1981 *   if (do_schlieren_plot == true)
1982 *   interpretation.push_back(
1984 *  
1985 *   return interpretation;
1986 *   }
1987 *  
1988 *  
1989 *  
1990 *   template <int dim>
1991 *   UpdateFlags EulerProblem<dim>::Postprocessor::get_needed_update_flags() const
1992 *   {
1993 *   if (do_schlieren_plot == true)
1994 *   return update_values | update_gradients;
1995 *   else
1996 *   return update_values;
1997 *   }
1998 *  
1999 *  
2000 *  
2001 *   template <int dim>
2002 *   EulerProblem<dim>::EulerProblem()
2003 *   : pcout(std::cout, Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
2004 *   #ifdef DEAL_II_WITH_P4EST
2005 *   , triangulation(MPI_COMM_WORLD)
2006 *   #endif
2007 *   , fe(FE_DGQ<dim>(fe_degree), dim + 2)
2008 *   , mapping(fe_degree)
2009 *   , dof_handler(triangulation)
2010 *   , timer(pcout, TimerOutput::never, TimerOutput::wall_times)
2011 *   , euler_operator(timer)
2012 *   , time(0)
2013 *   , time_step(0)
2014 *   {}
2015 *  
2016 *  
2017 *  
2018 *   template <int dim>
2019 *   void EulerProblem<dim>::make_grid_and_dofs()
2020 *   {
2021 *   switch (testcase)
2022 *   {
2023 *   case 0:
2024 *   {
2025 *   Point<dim> lower_left;
2026 *   for (unsigned int d = 1; d < dim; ++d)
2027 *   lower_left[d] = -5;
2028 *  
2029 *   Point<dim> upper_right;
2030 *   upper_right[0] = 10;
2031 *   for (unsigned int d = 1; d < dim; ++d)
2032 *   upper_right[d] = 5;
2033 *  
2035 *   lower_left,
2036 *   upper_right);
2037 *   triangulation.refine_global(2);
2038 *  
2039 *   euler_operator.set_inflow_boundary(
2040 *   0, std::make_unique<ExactSolution<dim>>(0));
2041 *  
2042 *   break;
2043 *   }
2044 *  
2045 *   case 1:
2046 *   {
2048 *   triangulation, 0.03, 1, 0, true);
2049 *  
2050 *   euler_operator.set_inflow_boundary(
2051 *   0, std::make_unique<ExactSolution<dim>>(0));
2052 *   euler_operator.set_subsonic_outflow_boundary(
2053 *   1, std::make_unique<ExactSolution<dim>>(0));
2054 *  
2055 *   euler_operator.set_wall_boundary(2);
2056 *   euler_operator.set_wall_boundary(3);
2057 *  
2058 *   if (dim == 3)
2059 *   euler_operator.set_body_force(
2060 *   std::make_unique<Functions::ConstantFunction<dim>>(
2061 *   std::vector<double>({0., 0., -0.2})));
2062 *  
2063 *   break;
2064 *   }
2065 *  
2066 *   default:
2068 *   }
2069 *  
2070 *   triangulation.refine_global(n_global_refinements);
2071 *  
2072 *   dof_handler.distribute_dofs(fe);
2073 *  
2074 *   euler_operator.reinit(mapping, dof_handler);
2075 *   euler_operator.initialize_vector(solution);
2076 *  
2077 *   std::locale s = pcout.get_stream().getloc();
2078 *   pcout.get_stream().imbue(std::locale(""));
2079 *   pcout << "Number of degrees of freedom: " << dof_handler.n_dofs()
2080 *   << " ( = " << (dim + 2) << " [vars] x "
2081 *   << triangulation.n_global_active_cells() << " [cells] x "
2082 *   << Utilities::pow(fe_degree + 1, dim) << " [dofs/cell/var] )"
2083 *   << std::endl;
2084 *   pcout.get_stream().imbue(s);
2085 *   }
2086 *  
2087 *  
2088 *  
2089 *   template <int dim>
2090 *   void EulerProblem<dim>::output_results(const unsigned int result_number)
2091 *   {
2092 *   const std::array<double, 3> errors =
2093 *   euler_operator.compute_errors(ExactSolution<dim>(time), solution);
2094 *   const std::string quantity_name = testcase == 0 ? "error" : "norm";
2095 *  
2096 *   pcout << "Time:" << std::setw(8) << std::setprecision(3) << time
2097 *   << ", dt: " << std::setw(8) << std::setprecision(2) << time_step
2098 *   << ", " << quantity_name << " rho: " << std::setprecision(4)
2099 *   << std::setw(10) << errors[0] << ", rho * u: " << std::setprecision(4)
2100 *   << std::setw(10) << errors[1] << ", energy:" << std::setprecision(4)
2101 *   << std::setw(10) << errors[2] << std::endl;
2102 *  
2103 *   {
2104 *   TimerOutput::Scope t(timer, "output");
2105 *  
2106 *   Postprocessor postprocessor;
2107 *   DataOut<dim> data_out;
2108 *  
2109 *   DataOutBase::VtkFlags flags;
2110 *   flags.write_higher_order_cells = true;
2111 *   data_out.set_flags(flags);
2112 *  
2113 *   data_out.attach_dof_handler(dof_handler);
2114 *   {
2115 *   std::vector<std::string> names;
2116 *   names.emplace_back("density");
2117 *   for (unsigned int d = 0; d < dim; ++d)
2118 *   names.emplace_back("momentum");
2119 *   names.emplace_back("energy");
2120 *  
2121 *   std::vector<DataComponentInterpretation::DataComponentInterpretation>
2122 *   interpretation;
2123 *   interpretation.push_back(
2125 *   for (unsigned int d = 0; d < dim; ++d)
2126 *   interpretation.push_back(
2128 *   interpretation.push_back(
2130 *  
2131 *   data_out.add_data_vector(dof_handler, solution, names, interpretation);
2132 *   }
2133 *   data_out.add_data_vector(solution, postprocessor);
2134 *  
2136 *   if (testcase == 0 && dim == 2)
2137 *   {
2138 *   reference.reinit(solution);
2139 *   euler_operator.project(ExactSolution<dim>(time), reference);
2140 *   reference.sadd(-1., 1, solution);
2141 *   std::vector<std::string> names;
2142 *   names.emplace_back("error_density");
2143 *   for (unsigned int d = 0; d < dim; ++d)
2144 *   names.emplace_back("error_momentum");
2145 *   names.emplace_back("error_energy");
2146 *  
2147 *   std::vector<DataComponentInterpretation::DataComponentInterpretation>
2148 *   interpretation;
2149 *   interpretation.push_back(
2151 *   for (unsigned int d = 0; d < dim; ++d)
2152 *   interpretation.push_back(
2154 *   interpretation.push_back(
2156 *  
2157 *   data_out.add_data_vector(dof_handler,
2158 *   reference,
2159 *   names,
2160 *   interpretation);
2161 *   }
2162 *  
2163 *   Vector<double> mpi_owner(triangulation.n_active_cells());
2164 *   mpi_owner = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
2165 *   data_out.add_data_vector(mpi_owner, "owner");
2166 *  
2167 *   data_out.build_patches(mapping,
2168 *   fe.degree,
2170 *  
2171 *   const std::string filename =
2172 *   "solution_" + Utilities::int_to_string(result_number, 3) + ".vtu";
2173 *   data_out.write_vtu_in_parallel(filename, MPI_COMM_WORLD);
2174 *   }
2175 *   }
2176 *  
2177 *  
2178 *  
2179 *   template <int dim>
2180 *   void EulerProblem<dim>::run()
2181 *   {
2182 *   {
2183 *   const unsigned int n_vect_number = VectorizedArrayType::size();
2184 *   const unsigned int n_vect_bits = 8 * sizeof(Number) * n_vect_number;
2185 *  
2186 *   pcout << "Running with "
2187 *   << Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD)
2188 *   << " MPI processes" << std::endl;
2189 *   pcout << "Vectorization over " << n_vect_number << ' '
2190 *   << (std::is_same_v<Number, double> ? "doubles" : "floats") << " = "
2191 *   << n_vect_bits << " bits ("
2193 *   << std::endl;
2194 *   }
2195 *  
2196 *   make_grid_and_dofs();
2197 *  
2198 *   const LowStorageRungeKuttaIntegrator integrator(lsrk_scheme);
2199 *  
2202 *   rk_register_1.reinit(solution);
2203 *   rk_register_2.reinit(solution);
2204 *  
2205 *   euler_operator.project(ExactSolution<dim>(time), solution);
2206 *  
2207 *  
2208 *   double min_vertex_distance = std::numeric_limits<double>::max();
2209 *   for (const auto &cell : triangulation.active_cell_iterators())
2210 *   if (cell->is_locally_owned())
2211 *   min_vertex_distance =
2212 *   std::min(min_vertex_distance, cell->minimum_vertex_distance());
2213 *   min_vertex_distance =
2214 *   Utilities::MPI::min(min_vertex_distance, MPI_COMM_WORLD);
2215 *  
2216 *   time_step = courant_number * integrator.n_stages() /
2217 *   euler_operator.compute_cell_transport_speed(solution);
2218 *   pcout << "Time step size: " << time_step
2219 *   << ", minimal h: " << min_vertex_distance
2220 *   << ", initial transport scaling: "
2221 *   << 1. / euler_operator.compute_cell_transport_speed(solution)
2222 *   << std::endl
2223 *   << std::endl;
2224 *  
2225 *   output_results(0);
2226 *  
2227 *   unsigned int timestep_number = 0;
2228 *  
2229 *   while (time < final_time - 1e-12 && timestep_number < max_time_steps)
2230 *   {
2231 *   ++timestep_number;
2232 *   if (timestep_number % 5 == 0)
2233 *   time_step =
2234 *   courant_number * integrator.n_stages() /
2236 *   euler_operator.compute_cell_transport_speed(solution), 3);
2237 *  
2238 *   {
2239 *   TimerOutput::Scope t(timer, "rk time stepping total");
2240 *   integrator.perform_time_step(euler_operator,
2241 *   time,
2242 *   time_step,
2243 *   solution,
2244 *   rk_register_1,
2245 *   rk_register_2);
2246 *   }
2247 *  
2248 *   time += time_step;
2249 *  
2250 *   if (static_cast<int>(time / output_tick) !=
2251 *   static_cast<int>((time - time_step) / output_tick) ||
2252 *   time >= final_time - 1e-12)
2253 *   output_results(
2254 *   static_cast<unsigned int>(std::round(time / output_tick)));
2255 *   }
2256 *  
2257 *   timer.print_wall_time_statistics(MPI_COMM_WORLD);
2258 *   pcout << std::endl;
2259 *   }
2260 *  
2261 *   } // namespace Euler_DG
2262 *  
2263 *  
2264 *   int main(int argc, char **argv)
2265 *   {
2266 *   using namespace Euler_DG;
2267 *   using namespace dealii;
2268 *  
2269 *   Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
2270 *  
2271 *   try
2272 *   {
2273 *   EulerProblem<dimension> euler_problem;
2274 *   euler_problem.run();
2275 *   }
2276 *   catch (std::exception &exc)
2277 *   {
2278 *   std::cerr << std::endl
2279 *   << std::endl
2280 *   << "----------------------------------------------------"
2281 *   << std::endl;
2282 *   std::cerr << "Exception on processing: " << std::endl
2283 *   << exc.what() << std::endl
2284 *   << "Aborting!" << std::endl
2285 *   << "----------------------------------------------------"
2286 *   << std::endl;
2287 *  
2288 *   return 1;
2289 *   }
2290 *   catch (...)
2291 *   {
2292 *   std::cerr << std::endl
2293 *   << std::endl
2294 *   << "----------------------------------------------------"
2295 *   << std::endl;
2296 *   std::cerr << "Unknown exception!" << std::endl
2297 *   << "Aborting!" << std::endl
2298 *   << "----------------------------------------------------"
2299 *   << std::endl;
2300 *   return 1;
2301 *   }
2302 *  
2303 *   return 0;
2304 *   }
2305 * @endcode
2306<a name="step_76-Results"></a><h1>Results</h1>
2307
2308
2309Running the program with the default settings on a machine with 40 processes
2310produces the following output:
2311
2312@code
2313Running with 40 MPI processes
2314Vectorization over 8 doubles = 512 bits (AVX512)
2315Number of degrees of freedom: 27.648.000 ( = 5 [vars] x 25.600 [cells] x 216 [dofs/cell/var] )
2316Time step size: 0.000295952, minimal h: 0.0075, initial transport scaling: 0.00441179
2317Time: 0, dt: 0.0003, norm rho: 5.385e-16, rho * u: 1.916e-16, energy: 1.547e-15
2318+--------------------------------------+------------------+------------+------------------+
2319| Total wallclock time elapsed | 17.52s 10 | 17.52s | 17.52s 11 |
2320| | | |
2321| Section | no. calls | min time rank | avg time | max time rank |
2322+--------------------------------------+------------------+------------+------------------+
2323| compute errors | 1 | 0.009594s 16 | 0.009705s | 0.009819s 8 |
2324| compute transport speed | 22 | 0.1366s 0 | 0.1367s | 0.1368s 18 |
2325| output | 1 | 1.233s 0 | 1.233s | 1.233s 32 |
2326| rk time stepping total | 100 | 8.746s 35 | 8.746s | 8.746s 0 |
2327| rk_stage - integrals L_h | 500 | 8.742s 36 | 8.742s | 8.743s 2 |
2328+--------------------------------------+------------------+------------+------------------+
2329@endcode
2330
2331and the following visual output:
2332
2333<table align="center" class="doxtable" style="width:85%">
2334 <tr>
2335 <td>
2336 <img src="https://www.dealii.org/images/steps/developer/step-67.pressure_010.png" alt="" width="100%">
2337 </td>
2338 <td>
2339 <img src="https://www.dealii.org/images/steps/developer/step-67.pressure_025.png" alt="" width="100%">
2340 </td>
2341 </tr>
2342 <tr>
2343 <td>
2344 <img src="https://www.dealii.org/images/steps/developer/step-67.pressure_050.png" alt="" width="100%">
2345 </td>
2346 <td>
2347 <img src="https://www.dealii.org/images/steps/developer/step-67.pressure_100.png" alt="" width="100%">
2348 </td>
2349 </tr>
2350</table>
2351
2352As a reference, the results of @ref step_67 "step-67" using FCL are:
2353
2354@code
2355Running with 40 MPI processes
2356Vectorization over 8 doubles = 512 bits (AVX512)
2357Number of degrees of freedom: 27.648.000 ( = 5 [vars] x 25.600 [cells] x 216 [dofs/cell/var] )
2358Time step size: 0.000295952, minimal h: 0.0075, initial transport scaling: 0.00441179
2359Time: 0, dt: 0.0003, norm rho: 5.385e-16, rho * u: 1.916e-16, energy: 1.547e-15
2360+-------------------------------------------+------------------+------------+------------------+
2361| Total wallclock time elapsed | 13.33s 0 | 13.34s | 13.35s 34 |
2362| | | |
2363| Section | no. calls | min time rank | avg time | max time rank |
2364+-------------------------------------------+------------------+------------+------------------+
2365| compute errors | 1 | 0.007977s 10 | 0.008053s | 0.008161s 30 |
2366| compute transport speed | 22 | 0.1228s 34 | 0.2227s | 0.3845s 0 |
2367| output | 1 | 1.255s 3 | 1.257s | 1.259s 27 |
2368| rk time stepping total | 100 | 11.15s 0 | 11.32s | 11.42s 34 |
2369| rk_stage - integrals L_h | 500 | 8.719s 10 | 8.932s | 9.196s 0 |
2370| rk_stage - inv mass + vec upd | 500 | 1.944s 0 | 2.377s | 2.55s 10 |
2371+-------------------------------------------+------------------+------------+------------------+
2372@endcode
2373
2374By the modifications shown in this tutorial, we were able to achieve a speedup of
237527% for the Runge-Kutta stages.
2376
2377<a name="step_76-Possibilitiesforextensions"></a><h3>Possibilities for extensions</h3>
2378
2379
2380The algorithms are easily extendable to higher dimensions: a high-dimensional
2381<a href="https://github.com/hyperdeal/hyperdeal/blob/a9e67b4e625ff1dde2fed93ad91cdfacfaa3acdf/include/hyper.deal/operators/advection/advection_operation.h#L219-L569">advection operator based on cell-centric loops</a>
2382is part of the hyper.deal library. An extension of cell-centric loops
2383to locally-refined meshes is more involved.
2384
2385<a name="step_76-ExtensiontothecompressibleNavierStokesequations"></a><h4>Extension to the compressible Navier-Stokes equations</h4>
2386
2387
2388The solver presented in this tutorial program can also be extended to the
2389compressible Navier–Stokes equations by adding viscous terms, as also
2390suggested in @ref step_67 "step-67". To keep as much of the performance obtained here despite
2391the additional cost of elliptic terms, e.g. via an interior penalty method, that
2392tutorial has proposed to switch the basis from FE_DGQ to FE_DGQHermite like in
2393the @ref step_59 "step-59" tutorial program. The reasoning behind this switch is that in the
2394case of FE_DGQ all values of neighboring cells (i.e., @f$k+1@f$ layers) are needed,
2395whilst in the case of FE_DGQHermite only 2 layers, making the latter
2396significantly more suitable for higher degrees. The additional layers have to be,
2397on the one hand, loaded from main memory during flux computation and, one the
2398other hand, have to be communicated. Using the shared-memory capabilities
2399introduced in this tutorial, the second point can be eliminated on a single
2400compute node or its influence can be reduced in a hybrid context.
2401
2402<a name="step_76-BlockGaussSeidellikepreconditioners"></a><h4>Block Gauss-Seidel-like preconditioners</h4>
2403
2404
2405Cell-centric loops could be used to create block Gauss-Seidel preconditioners
2406that are multiplicative within one process and additive over processes. These
2407type of preconditioners use during flux computation, in contrast to Jacobi-type
2408preconditioners, already updated values from neighboring cells. The following
2409pseudo-code visualizes how this could in principal be achieved:
2410
2411@code
2412// vector monitor if cells have been updated or not
2413Vector<Number> visit_flags(data.n_cell_batches () + data.n_ghost_cell_batches ());
2414
2415// element centric loop with a modified kernel
2416data.template loop_cell_centric<VectorType, VectorType>(
2417 [&](const auto &data, auto &dst, const auto &src, const auto cell_range) {
2418
2419 for (unsigned int cell = cell_range.first; cell < cell_range.second; ++cell)
2420 {
2421 // cell integral as usual (not shown)
2422
2423 for (unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell; ++face)
2424 {
2425 const auto boundary_id = data.get_faces_by_cells_boundary_id(cell, face)[0];
2426
2427 if (boundary_id == numbers::internal_face_boundary_id)
2428 {
2429 phi_p.reinit(cell, face);
2430
2431 const auto flags = phi_p.read_cell_data(visit_flags);
2432 const auto all_neighbors_have_been_updated =
2433 std::min(flags.begin(),
2434 flags().begin() + data.n_active_entries_per_cell_batch(cell) == 1;
2435
2436 if(all_neighbors_have_been_updated)
2437 phi_p.gather_evaluate(dst, EvaluationFlags::values);
2438 else
2439 phi_p.gather_evaluate(src, EvaluationFlags::values);
2440
2441 // continue as usual (not shown)
2442 }
2443 else
2444 {
2445 // boundary integral as usual (not shown)
2446 }
2447 }
2448
2449 // continue as above and apply your favorite algorithm to invert
2450 // the cell-local operator (not shown)
2451
2452 // make cells as updated
2453 phi.set_cell_data(visit_flags, VectorizedArrayType(1.0));
2454 }
2455 },
2456 dst,
2457 src,
2458 true,
2460@endcode
2461
2462For this purpose, one can exploit the cell-data vector capabilities of
2463MatrixFree and the range-based iteration capabilities of VectorizedArray.
2464
2465Please note that in the given example we process <code>VectorizedArrayType::size()</code>
2466number of blocks, since each lane corresponds to one block. We consider blocks
2467as updated if all blocks processed by a vector register have been updated. In
2468the case of Cartesian meshes this is a reasonable approach, however, for
2469general unstructured meshes this conservative approach might lead to a decrease in the
2470efficiency of the preconditioner. A reduction of cells processed in parallel
2471by explicitly reducing the number of lanes used by <code>VectorizedArrayType</code>
2472might increase the quality of the preconditioner, but with the cost that each
2473iteration might be more expensive. This dilemma leads us to a further
2474"possibility for extension": vectorization within an element.
2475 *
2476 *
2477<a name="step_76-PlainProg"></a>
2478<h1> The plain program</h1>
2479@include "step-76.cc"
2480*/
virtual UpdateFlags get_needed_update_flags() const =0
virtual void evaluate_vector_field(const DataPostprocessorInputs::Vector< dim > &input_data, std::vector< Vector< double > > &computed_quantities) const
virtual std::vector< std::string > get_names() const =0
virtual std::vector< DataComponentInterpretation::DataComponentInterpretation > get_data_component_interpretation() const
void submit_value(const value_type val_in, const unsigned int q_point)
void reinit(const size_type size, const bool omit_zeroing_entries=false)
Abstract base class for mapping classes.
Definition mapping.h:318
void loop_cell_centric(void(CLASS::*cell_operation)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const, const CLASS *owning_class, OutVector &dst, const InVector &src, const bool zero_dst_vector=false, const DataAccessOnFaces src_vector_face_access=DataAccessOnFaces::unspecified) const
void reinit(const MappingType &mapping, const DoFHandler< dim > &dof_handler, const AffineConstraints< number2 > &constraint, const QuadratureType &quad, const AdditionalData &additional_data=AdditionalData())
Definition point.h:111
constexpr numbers::NumberTraits< Number >::real_type norm_square() const
@ wall_times
Definition timer.h:651
std::enable_if_t< std::is_floating_point_v< T > &&std::is_floating_point_v< U >, typename ProductType< std::complex< T >, std::complex< U > >::type > operator*(const std::complex< T > &left, const std::complex< U > &right)
#define DEAL_II_ALWAYS_INLINE
Definition config.h:108
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)
Point< 2 > second
Definition grid_out.cc:4624
Point< 2 > first
Definition grid_out.cc:4623
unsigned int level
Definition grid_out.cc:4626
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertThrow(cond, exc)
void loop(IteratorType begin, std_cxx20::type_identity_t< IteratorType > end, DOFINFO &dinfo, INFOBOX &info, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &)> &cell_worker, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &)> &boundary_worker, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, AssemblerType &assembler, const LoopControl &lctrl=LoopControl())
Definition loop.h:564
UpdateFlags
@ update_values
Shape function values.
@ update_normal_vectors
Normal vectors.
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
#define DEAL_II_NOT_IMPLEMENTED()
const Event initial
Definition event.cc:64
void hyper_rectangle(Triangulation< dim, spacedim > &tria, const Point< dim > &p1, const Point< dim > &p2, const bool colorize=false)
void channel_with_cylinder(Triangulation< dim > &tria, const double shell_region_width=0.03, const unsigned int n_shells=2, const double skewness=2.0, const bool colorize=false)
@ matrix
Contents is actually a matrix.
@ general
No special properties.
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
Definition divergence.h:471
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition utilities.cc:191
SymmetricTensor< 2, dim, Number > C(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
@ LOW_STORAGE_RK_STAGE9_ORDER5
@ LOW_STORAGE_RK_STAGE3_ORDER3
@ LOW_STORAGE_RK_STAGE7_ORDER4
@ LOW_STORAGE_RK_STAGE5_ORDER4
VectorType::value_type * begin(VectorType &V)
T sum(const T &t, const MPI_Comm mpi_communicator)
unsigned int n_mpi_processes(const MPI_Comm mpi_communicator)
Definition mpi.cc:92
T max(const T &t, const MPI_Comm mpi_communicator)
T min(const T &t, const MPI_Comm mpi_communicator)
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
Definition mpi.cc:107
T reduce(const T &local_value, const MPI_Comm comm, const std::function< T(const T &, const T &)> &combiner, const unsigned int root_process=0)
std::string get_time()
std::string get_current_vectorization_level()
Definition utilities.cc:938
Number truncate_to_n_digits(const Number number, const unsigned int n_digits)
Definition utilities.cc:578
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition utilities.cc:470
constexpr T pow(const T base, const int iexp)
Definition utilities.h:966
void project(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const AffineConstraints< typename VectorType::value_type > &constraints, const Quadrature< dim > &quadrature, const Function< spacedim, typename VectorType::value_type > &function, VectorType &vec, const bool enforce_zero_boundary=false, const Quadrature< dim - 1 > &q_boundary=(dim > 1 ? QGauss< dim - 1 >(2) :Quadrature< dim - 1 >()), const bool project_to_boundary_first=false)
void run(const Iterator &begin, const std_cxx20::type_identity_t< Iterator > &end, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length, const unsigned int chunk_size)
long double gamma(const unsigned int n)
void copy(const T *begin, const T *end, U *dest)
int(&) functions(const void *v1, const void *v2)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
static constexpr double PI
Definition numbers.h:254
const types::boundary_id internal_face_boundary_id
Definition types.h:312
static const unsigned int invalid_unsigned_int
Definition types.h:220
void transform(const InputIterator &begin_in, const InputIterator &end_in, OutputIterator out, const Function &function, const unsigned int grainsize)
Definition parallel.h:165
STL namespace.
::VectorizedArray< Number, width > exp(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
Definition types.h:32
unsigned int boundary_id
Definition types.h:144
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation