deal.II version 9.7.0
\(\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-89.h
Go to the documentation of this file.
1{},
349 /* interior_face_operation= */ face_function,
350 /* boundary_face_operation= */{},
351 dst, src);
352@endcode
353
354The code to do the same with FERemoteEvaluation is shown below.
355For brevity, we assume all boundary faces are somehow connected via
356non-conforming interfaces for FERemoteEvaluation.
357
358@code
359// Initialize FERemoteEvaluation: Note, that FERemoteEvaluation internally manages
360// the memory to store precomputed values. Therefore, FERemoteEvaluation
361// should be initialized only once to avoid frequent memory
362// allocation/deallocation. At this point, remote_communicator is assumed
363// to be initialized.
364FERemoteEvaluation<dim,Number> u_p_evaluator(remote_communicator);
365
366// Precompute the interpolated values of the finite element solution at all
367// the quadrature points outside the loop, invoking the vector entries and
368// respective basis function at possibly remote MPI processes before communication.
369u_p_evaluator.gather_evaluate(src, EvaluationFlags::values);
370
371const auto boundary_function =
372 [&](const MatrixFree &data, VectorType &dst, const VectorType &src,
373 const std::pair<unsigned int, unsigned int> face_range) {
374
375 FEFaceEvaluation phi_m(data, true);
376 // To access the values in a thread safe way each thread has
377 // to create a own accessor object. A small helper function
378 // provides the accessor.
379 internal::PrecomputedEvaluationDataAccessor u_p = u_p_evaluator.get_data_accessor();
380
381 for (unsigned int f = face_range.first; f < face_range.second; ++f)
382 {
383 phi_m.reinit(f);
384 u_p.reinit(f);
385
386 for (unsigned int q = 0; q < phi_m.n_q_points; ++q)
387 phi_m.submit_value(u_p.get_value(q), q); // access values with phi_p
388
389 phi_m.integrate_scatter(EvaluationFlags::values, dst);
390 }
391 };
392
393matrix_free.template loop<VectorType, VectorType>({}, {}, boundary_function, dst, src);
394@endcode
395The object @c remote_communicator is of type `FERemoteEvaluationCommunicator`
396and assumed to be correctly initialized prior to the above code snippet.
397`FERemoteEvaluationCommunicator` internally manages the update of ghost values
398over non-matching interfaces and keeps track of the mapping between quadrature
399point index and corresponding values/gradients. As mentioned above, the update
400of the values/gradients happens <i>before</i> the actual matrix-free loop.
402differently for the given template parameter @c value_type. If we want to access
403values at arbitrary points (e.g. in combination with @c FEPointEvaluation), then
404we need to choose @c value_type=Number. If the values are defined at quadrature
405points of a @c FEEvaluation object it is possible to get the values at the
406quadrature points of <i>batches</i> and we need to choose
407@c value_type=VectorizedArray<Number>.
408
409
410<a name="step_89-Overviewofthetestcase"></a><h3>Overview of the test case</h3>
411
412
413In this program, we implemented both the point-to-point interpolation and
414Nitsche-type mortaring mentioned in the introduction.
415
416At first we are considering the test case of a vibrating membrane, see e.g.
417@cite nguyen2011high. Standing waves of length @f$\lambda=2/M@f$ are oscillating
418with a time period of @f$T=2 / (M \sqrt{d} c)@f$ where @f$d@f$ is the dimension of the
419space in which our domain is located and @f$M@f$ is the number of modes per meter,
420i.e. the number of half-waves per meter. The corresponding analytical solution
421reads as
422
423@f{align*}{
424 p &=\cos(M \sqrt{d} \pi c t)\prod_{i=1}^{d} \sin(M \pi x_i),\\
425 u_i&=-\frac{\sin(M \sqrt{d} \pi c t)}{\sqrt{d}\rho c}
426 \cos(M \pi x_i)\prod_{j=1,j\neq i}^{d} \sin(M \pi x_j),
427@f}
428
429For simplicity, we are using homogeneous pressure Dirichlet boundary conditions
430within this tutorial. To be able to do so we have to tailor the domain size as
431well as the number of modes to conform with the homogeneous pressure Dirichlet
432boundary conditions. Within this tutorial we are using @f$M=10@f$ and a domain
433@f$\Omega=(0,1)^2@f$. The domain will be meshed so that the left and right parts of
434the domain consist of separate meshes that do not match at the interface.
435
436As will become clear from the results, the point-to-point interpolation will
437result in aliasing, which can be resolved using Nitsche-type mortaring.
438
439In a more realistic second example, we apply this implementation to a test case
440in which a wave is propagating from one fluid into another fluid. The speed of
441sound in the left part of the domain is @f$c=1@f$ and in the right part it is @f$c=3@f$.
442Since the wavelength is directly proportional to the speed of sound, three times
443larger elements can be used in the right part of the domain to resolve waves up
444to the same frequency. A test case like this has been simulated with a different
445domain and different initial conditions, e.g., in @cite bangerth2010adaptive.
446 *
447 *
448 * <a name="step_89-CommProg"></a>
449 * <h1> The commented program</h1>
450 *
451 *
452 * <a name="step_89-Includefiles"></a>
453 * <h3>Include files</h3>
454 *
455
456 *
457 * The program starts with including all the relevant header files.
458 *
459 * @code
460 *   #include <deal.II/base/conditional_ostream.h>
461 *   #include <deal.II/base/mpi.h>
462 *  
463 *   #include <deal.II/distributed/tria.h>
464 *  
465 *   #include <deal.II/fe/fe_dgq.h>
466 *   #include <deal.II/fe/fe_system.h>
467 *   #include <deal.II/fe/fe_values.h>
468 *   #include <deal.II/fe/mapping_q1.h>
469 *  
470 *   #include <deal.II/grid/grid_generator.h>
471 *   #include <deal.II/grid/grid_tools.h>
472 *   #include <deal.II/grid/grid_tools_cache.h>
473 *  
474 *   #include <deal.II/matrix_free/fe_evaluation.h>
475 *   #include <deal.II/matrix_free/matrix_free.h>
476 *   #include <deal.II/matrix_free/operators.h>
477 *  
478 *   #include <deal.II/non_matching/mapping_info.h>
479 *  
480 *   #include <deal.II/numerics/data_out.h>
481 *   #include <deal.II/numerics/vector_tools.h>
482 * @endcode
483 *
484 * The following header file provides the class FERemoteEvaluation, which allows
485 * to access values and/or gradients at remote triangulations similar to
486 * FEEvaluation.
487 *
488 * @code
489 *   #include <deal.II/matrix_free/fe_remote_evaluation.h>
490 *  
491 * @endcode
492 *
493 * We pack everything that is specific for this program into a namespace
494 * of its own.
495 *
496 * @code
497 *   namespace Step89
498 *   {
499 *   using namespace dealii;
500 *  
501 * @endcode
502 *
503 *
504 * <a name="step_89-Initialconditionsforvibratingmembrane"></a>
505 * <h3>Initial conditions for vibrating membrane</h3>
506 *
507
508 *
509 * In the following, let us first define a function that provides
510 * the initial condition for the vibrating membrane test case. It
511 * implements both the initial pressure (component 0) and velocity
512 * (components 1 to 1 + dim).
513 *
514
515 *
516 * There is also a function that computes the duration of one
517 * oscillation.
518 *
519 * @code
520 *   template <int dim>
521 *   class InitialConditionVibratingMembrane : public Function<dim>
522 *   {
523 *   public:
524 *   InitialConditionVibratingMembrane(const double modes);
525 *  
526 *   double value(const Point<dim> &p, const unsigned int comp) const final;
527 *  
528 *   double get_period_duration(const double speed_of_sound) const;
529 *  
530 *   private:
531 *   const double M;
532 *   };
533 *  
534 *   template <int dim>
535 *   InitialConditionVibratingMembrane<dim>::InitialConditionVibratingMembrane(
536 *   const double modes)
537 *   : Function<dim>(dim + 1, 0.0)
538 *   , M(modes)
539 *   {
540 *   static_assert(dim == 2, "Only implemented for dim==2");
541 *   }
542 *  
543 *   template <int dim>
544 *   double
545 *   InitialConditionVibratingMembrane<dim>::value(const Point<dim> &p,
546 *   const unsigned int comp) const
547 *   {
548 *   if (comp == 0)
549 *   return std::sin(M * numbers::PI * p[0]) *
550 *   std::sin(M * numbers::PI * p[1]);
551 *  
552 *   return 0.0;
553 *   }
554 *  
555 *   template <int dim>
556 *   double InitialConditionVibratingMembrane<dim>::get_period_duration(
557 *   const double speed_of_sound) const
558 *   {
559 *   return 2.0 / (M * std::sqrt(dim) * speed_of_sound);
560 *   }
561 *  
562 * @endcode
563 *
564 *
565 * <a name="step_89-Gausspulse"></a>
566 * <h3>Gauss pulse</h3>
567 *
568
569 *
570 * Next up is a function that provides the values of a pressure
571 * Gauss pulse. As with the previous function, it implements both
572 * the initial pressure (component 0) and velocity (components 1 to
573 * 1 + dim).
574 *
575 * @code
576 *   template <int dim>
577 *   class GaussPulse : public Function<dim>
578 *   {
579 *   public:
580 *   GaussPulse(const double shift_x, const double shift_y);
581 *  
582 *   double value(const Point<dim> &p, const unsigned int comp) const final;
583 *  
584 *   private:
585 *   const double shift_x;
586 *   const double shift_y;
587 *   };
588 *  
589 *   template <int dim>
590 *   GaussPulse<dim>::GaussPulse(const double shift_x, const double shift_y)
591 *   : Function<dim>(dim + 1, 0.0)
592 *   , shift_x(shift_x)
593 *   , shift_y(shift_y)
594 *   {
595 *   static_assert(dim == 2, "Only implemented for dim==2");
596 *   }
597 *  
598 *   template <int dim>
599 *   double GaussPulse<dim>::value(const Point<dim> &p,
600 *   const unsigned int comp) const
601 *   {
602 *   if (comp == 0)
603 *   return std::exp(-1000.0 * ((std::pow(p[0] - shift_x, 2)) +
604 *   (std::pow(p[1] - shift_y, 2))));
605 *  
606 *   return 0.0;
607 *   }
608 *  
609 * @endcode
610 *
611 *
612 * <a name="step_89-Helperfunctions"></a>
613 * <h3>Helper functions</h3>
614 *
615
616 *
617 * The following namespace contains free helper functions that are used in the
618 * tutorial.
619 *
620 * @code
621 *   namespace HelperFunctions
622 *   {
623 * @endcode
624 *
625 * Helper function to check if a boundary ID is related to a non-matching
626 * face. A @c std::set that contains all non-matching boundary IDs is
627 * handed over in addition to the face ID under question.
628 *
629 * @code
630 *   bool is_non_matching_face(
631 *   const std::set<types::boundary_id> &non_matching_face_ids,
632 *   const types::boundary_id face_id)
633 *   {
634 *   return non_matching_face_ids.find(face_id) != non_matching_face_ids.end();
635 *   }
636 *  
637 * @endcode
638 *
639 * Helper function to set the initial conditions for the vibrating membrane
640 * test case.
641 *
642 * @code
643 *   template <int dim, typename Number, typename VectorType>
644 *   void set_initial_condition(MatrixFree<dim, Number> matrix_free,
645 *   const Function<dim> &initial_solution,
646 *   VectorType &dst)
647 *   {
648 *   VectorTools::interpolate(*matrix_free.get_mapping_info().mapping,
649 *   matrix_free.get_dof_handler(),
650 *   initial_solution,
651 *   dst);
652 *   }
653 *  
654 * @endcode
655 *
656 * Helper function to compute the time step size according to the CFL
657 * condition.
658 *
659 * @code
660 *   double
661 *   compute_dt_cfl(const double hmin, const unsigned int degree, const double c)
662 *   {
663 *   return hmin / (std::pow(degree, 1.5) * c);
664 *   }
665 *  
666 * @endcode
667 *
668 * Helper function that writes vtu output.
669 *
670 * @code
671 *   template <typename VectorType, int dim>
672 *   void write_vtu(const VectorType &solution,
673 *   const DoFHandler<dim> &dof_handler,
674 *   const Mapping<dim> &mapping,
675 *   const unsigned int degree,
676 *   const std::string &name_prefix)
677 *   {
678 *   DataOut<dim> data_out;
679 *   DataOutBase::VtkFlags flags;
680 *   flags.write_higher_order_cells = true;
681 *   data_out.set_flags(flags);
682 *  
683 *   std::vector<DataComponentInterpretation::DataComponentInterpretation>
684 *   interpretation(
686 *   std::vector<std::string> names(dim + 1, "U");
687 *  
689 *   names[0] = "P";
690 *  
691 *   data_out.add_data_vector(dof_handler, solution, names, interpretation);
692 *  
693 *   data_out.build_patches(mapping, degree, DataOut<dim>::curved_inner_cells);
694 *   data_out.write_vtu_in_parallel(name_prefix + ".vtu",
695 *   dof_handler.get_mpi_communicator());
696 *   }
697 *   } // namespace HelperFunctions
698 *  
699 * @endcode
700 *
701 *
702 * <a name="step_89-Materialparameterdescription"></a>
703 * <h3>Material parameter description</h3>
704 *
705
706 *
707 * The following class stores the information if the fluid is
708 * homogeneous as well as the material properties at every cell.
709 * This class helps access the correct values without accessing a
710 * large vector of materials in the homogeneous case.
711 *
712
713 *
714 * The class is provided a map from material ids to material
715 * properties (given as a pair of values for the speed of sound and
716 * the density). If the map has only one entry, the material is
717 * homogeneous -- using the same values everywhere -- and we can
718 * remember that fact in the `homogeneous` member variable and use it
719 * to optimize some code paths below. If the material is not
720 * homogeneous, we will fill a vector with the correct materials for
721 * each batch of cells; this information can then be access via
723 *
724
725 *
726 * As is usual when working with the MatrixFree framework, we will
727 * not only access material parameters at a single site, but for
728 * whole "batches" of cells. As a consequence, the functions below
729 * returned VectorizedArray objects for a batch at a time.
730 *
731 * @code
732 *   template <typename Number>
733 *   class CellwiseMaterialData
734 *   {
735 *   public:
736 *   template <int dim>
737 *   CellwiseMaterialData(
738 *   const MatrixFree<dim, Number, VectorizedArray<Number>> &matrix_free,
739 *   const std::map<types::material_id, std::pair<double, double>>
740 *   &material_id_map)
741 *   : homogeneous(material_id_map.size() == 1)
742 *   {
743 *   Assert(material_id_map.size() > 0,
744 *   ExcMessage("No materials given to CellwiseMaterialData."));
745 *  
746 *   if (homogeneous)
747 *   {
748 *   speed_of_sound_homogeneous = material_id_map.begin()->second.first;
749 *   density_homogeneous = material_id_map.begin()->second.second;
750 *   }
751 *   else
752 *   {
753 *   const auto n_cell_batches =
754 *   matrix_free.n_cell_batches() + matrix_free.n_ghost_cell_batches();
755 *  
756 *   speed_of_sound.resize(n_cell_batches);
757 *   density.resize(n_cell_batches);
758 *  
759 *   for (unsigned int cell = 0; cell < n_cell_batches; ++cell)
760 *   {
761 *   speed_of_sound[cell] = 1.;
762 *   density[cell] = 1.;
763 *   for (unsigned int v = 0;
764 *   v < matrix_free.n_active_entries_per_cell_batch(cell);
765 *   ++v)
766 *   {
767 *   const auto material_id =
768 *   matrix_free.get_cell_iterator(cell, v)->material_id();
769 *  
770 *   speed_of_sound[cell][v] =
771 *   material_id_map.at(material_id).first;
772 *   density[cell][v] = material_id_map.at(material_id).second;
773 *   }
774 *   }
775 *   }
776 *   }
777 *  
778 *   bool is_homogeneous() const
779 *   {
780 *   return homogeneous;
781 *   }
782 *  
783 *   const AlignedVector<VectorizedArray<Number>> &get_speed_of_sound() const
784 *   {
785 *   Assert(!homogeneous, ExcMessage("Use get_homogeneous_speed_of_sound()."));
786 *   return speed_of_sound;
787 *   }
788 *  
789 *   const AlignedVector<VectorizedArray<Number>> &get_density() const
790 *   {
791 *   Assert(!homogeneous, ExcMessage("Use get_homogeneous_density()."));
792 *   return density;
793 *   }
794 *  
795 *   VectorizedArray<Number> get_homogeneous_speed_of_sound() const
796 *   {
797 *   Assert(homogeneous, ExcMessage("Use get_speed_of_sound()."));
798 *   return speed_of_sound_homogeneous;
799 *   }
800 *  
801 *   VectorizedArray<Number> get_homogeneous_density() const
802 *   {
803 *   Assert(homogeneous, ExcMessage("Use get_density()."));
804 *   return density_homogeneous;
805 *   }
806 *  
807 *   private:
808 *   const bool homogeneous;
809 *  
810 *   /* Materials in the inhomogeneous case. */
811 *   AlignedVector<VectorizedArray<Number>> speed_of_sound;
813 *  
814 *   /* Materials in the homogeneous case. */
815 *   VectorizedArray<Number> speed_of_sound_homogeneous;
816 *   VectorizedArray<Number> density_homogeneous;
817 *   };
818 *  
819 * @endcode
820 *
821 * To be able to access the material data in every cell in a thread
822 * safe way, the following class @c MaterialEvaluation is
823 * used. Similar to @c FEEvaluation, functions below will create
824 * their own instances of this class; thus, there can be no race
825 * conditions even if these functions run multiple times in
826 * parallel. For inhomogeneous materials, a @c reinit_cell() or @c
827 * reinit_face() function is used to set the correct material at the
828 * current cell batch. In the homogeneous case the @c _reinit()
829 * functions don't have to reset the materials.
830 *
831 * @code
832 *   template <int dim, typename Number>
833 *   class MaterialEvaluation
834 *   {
835 *   public:
836 *   MaterialEvaluation(
837 *   const MatrixFree<dim, Number, VectorizedArray<Number>> &matrix_free,
838 *   const CellwiseMaterialData<Number> &material_data)
839 *   : phi(matrix_free)
840 *   , phi_face(matrix_free, true)
841 *   , material_data(material_data)
842 *   {
843 *   if (material_data.is_homogeneous())
844 *   {
845 *   speed_of_sound = material_data.get_homogeneous_speed_of_sound();
846 *   density = material_data.get_homogeneous_density();
847 *   }
848 *   }
849 *  
850 *   bool is_homogeneous() const
851 *   {
852 *   return material_data.is_homogeneous();
853 *   }
854 *  
855 * @endcode
856 *
857 * The following two functions update the data for the current
858 * cell or face, given a cell batch index. If the material is
859 * homogeneous, there is nothing to do. Otherwise, we reinit the
860 * FEEvaluation object and store the data for the current object.
861 *
862 * @code
863 *   void reinit_cell(const unsigned int cell)
864 *   {
865 *   if (!material_data.is_homogeneous())
866 *   {
867 *   phi.reinit(cell);
868 *   speed_of_sound =
869 *   phi.read_cell_data(material_data.get_speed_of_sound());
870 *   density = phi.read_cell_data(material_data.get_density());
871 *   }
872 *   }
873 *  
874 *   void reinit_face(const unsigned int face)
875 *   {
876 *   if (!material_data.is_homogeneous())
877 *   {
878 *   phi_face.reinit(face);
879 *   speed_of_sound =
880 *   phi_face.read_cell_data(material_data.get_speed_of_sound());
881 *   density = phi_face.read_cell_data(material_data.get_density());
882 *   }
883 *   }
884 *  
885 * @endcode
886 *
887 * The following functions then return the speed of sound and
888 * density for the current cell batch.
889 *
890 * @code
891 *   VectorizedArray<Number> get_speed_of_sound() const
892 *   {
893 *   return speed_of_sound;
894 *   }
895 *  
896 *   VectorizedArray<Number> get_density() const
897 *   {
898 *   return density;
899 *   }
900 *  
901 *   private:
902 *   /* Members needed for the inhomogeneous case. */
903 *   FEEvaluation<dim, -1, 0, 1, Number> phi;
904 *   FEFaceEvaluation<dim, -1, 0, 1, Number> phi_face;
905 *  
906 *   /* Material defined at every cell. */
907 *   const CellwiseMaterialData<Number> &material_data;
908 *  
909 *   /* Materials at current cell. */
910 *   VectorizedArray<Number> speed_of_sound;
911 *   VectorizedArray<Number> density;
912 *   };
913 *  
914 *  
915 * @endcode
916 *
917 *
918 * <a name="step_89-Boundaryconditions"></a>
919 * <h3>Boundary conditions</h3>
920 *
921
922 *
923 * To be able to use the same kernel, for all face integrals we define
924 * a class that returns the needed values at boundaries. In this tutorial
925 * homogeneous pressure Dirichlet boundary conditions are applied via
926 * the mirror principle, i.e. @f$p_h^+=-p_h^- + 2g@f$ with @f$g=0@f$.
927 *
928 * @code
929 *   template <int dim, typename Number>
930 *   class BCEvaluationP
931 *   {
932 *   public:
933 *   BCEvaluationP(const FEFaceEvaluation<dim, -1, 0, 1, Number> &pressure_m)
934 *   : pressure_m(pressure_m)
935 *   {}
936 *  
937 *   typename FEFaceEvaluation<dim, -1, 0, 1, Number>::value_type
938 *   get_value(const unsigned int q) const
939 *   {
940 *   return -pressure_m.get_value(q);
941 *   }
942 *  
943 *   private:
944 *   const FEFaceEvaluation<dim, -1, 0, 1, Number> &pressure_m;
945 *   };
946 *  
947 * @endcode
948 *
949 * We don't have to apply boundary conditions for the velocity, i.e.
950 * @f$\mathbf{u}_h^+=\mathbf{u}_h^-@f$.
951 *
952 * @code
953 *   template <int dim, typename Number>
954 *   class BCEvaluationU
955 *   {
956 *   public:
957 *   BCEvaluationU(const FEFaceEvaluation<dim, -1, 0, dim, Number> &velocity_m)
958 *   : velocity_m(velocity_m)
959 *   {}
960 *  
961 *   typename FEFaceEvaluation<dim, -1, 0, dim, Number>::value_type
962 *   get_value(const unsigned int q) const
963 *   {
964 *   return velocity_m.get_value(q);
965 *   }
966 *  
967 *   private:
968 *   const FEFaceEvaluation<dim, -1, 0, dim, Number> &velocity_m;
969 *   };
970 *  
971 * @endcode
972 *
973 *
974 * <a name="step_89-Acousticoperator"></a>
975 * <h3>Acoustic operator</h3>
976 *
977
978 *
979 * The following class then defines the acoustic operator. The class is
980 * heavily based on matrix-free methods. For a better understanding in
981 * matrix-free methods please refer to @ref step_67 "step-67".
982 *
983
984 *
985 * At the top we define a flag that decides whether we want to use
986 * mortaring. If the remote evaluators are set up with a
987 * VectorizedArray we are using point-to-point interpolation;
988 * otherwise we make use of Nitsche-type mortaring. The decision is
989 * made using `std::is_floating_point_v`, which is a variable that
990 * is true if the template argument is a floating point type (such
991 * as `float` or `double`) and false otherwise (in particular if the
992 * template argument is a vectorized array type).
993 *
994 * @code
995 *   template <int dim, typename Number, typename remote_value_type>
996 *   class AcousticOperator
997 *   {
998 *   static constexpr bool use_mortaring =
999 *   std::is_floating_point_v<remote_value_type>;
1000 *  
1001 *   public:
1002 * @endcode
1003 *
1004 * In case of Nitsche-type mortaring, `NonMatching::MappingInfo` has to
1005 * be provided in the constructor.
1006 *
1007 * @code
1008 *   AcousticOperator(
1009 *   const MatrixFree<dim, Number> &matrix_free,
1010 *   std::shared_ptr<CellwiseMaterialData<Number>> material_data,
1011 *   const std::set<types::boundary_id> &remote_face_ids,
1012 *   std::shared_ptr<FERemoteEvaluation<dim, 1, remote_value_type>>
1013 *   pressure_r_eval,
1014 *   std::shared_ptr<FERemoteEvaluation<dim, dim, remote_value_type>>
1015 *   velocity_r_eval,
1016 *   std::shared_ptr<FERemoteEvaluation<dim, 1, remote_value_type>> c_r_eval,
1017 *   std::shared_ptr<FERemoteEvaluation<dim, 1, remote_value_type>> rho_r_eval,
1018 *   std::shared_ptr<NonMatching::MappingInfo<dim, dim, Number>> nm_info =
1019 *   nullptr)
1020 *   : matrix_free(matrix_free)
1021 *   , material_data(material_data)
1022 *   , remote_face_ids(remote_face_ids)
1023 *   , pressure_r_eval(pressure_r_eval)
1024 *   , velocity_r_eval(velocity_r_eval)
1025 *   , c_r_eval(c_r_eval)
1026 *   , rho_r_eval(rho_r_eval)
1027 *   , nm_mapping_info(nm_info)
1028 *   {
1029 *   if (use_mortaring)
1030 *   Assert(nm_info,
1031 *   ExcMessage(
1032 *   "In case of Nitsche-type mortaring NonMatching::MappingInfo \
1033 *   has to be provided."));
1034 *   }
1035 *  
1036 * @endcode
1037 *
1038 * The following function then evaluates the acoustic operator.
1039 * It first updates the precomputed values in corresponding the
1040 * FERemoteEvaluation objects. The material parameters do not change and
1041 * thus, we do not have to update precomputed values in @c c_r_eval and @c
1042 * rho_r_eval.
1043 *
1044
1045 *
1046 * It then either performs a matrix-free loop with Nitsche-type
1047 * mortaring at non-matching faces (if `use_mortaring` is true) or
1048 * with point-to-point interpolation at non-matching faces (in the
1049 * `else` branch). The difference is only in the third argument to
1050 * the loop function, denoting the function object that is
1051 * executed at boundary faces.
1052 *
1053 * @code
1054 *   template <typename VectorType>
1055 *   void evaluate(VectorType &dst, const VectorType &src) const
1056 *   {
1057 *   pressure_r_eval->gather_evaluate(src, EvaluationFlags::values);
1058 *   velocity_r_eval->gather_evaluate(src, EvaluationFlags::values);
1059 *  
1060 *   if constexpr (use_mortaring)
1061 *   {
1062 *   matrix_free.loop(
1063 *   &AcousticOperator::local_apply_cell<VectorType>,
1064 *   &AcousticOperator::local_apply_face<VectorType>,
1065 *   &AcousticOperator::local_apply_boundary_face_mortaring<VectorType>,
1066 *   this,
1067 *   dst,
1068 *   src,
1069 *   true,
1070 *   MatrixFree<dim, Number>::DataAccessOnFaces::values,
1071 *   MatrixFree<dim, Number>::DataAccessOnFaces::values);
1072 *   }
1073 *   else
1074 *   {
1075 *   matrix_free.loop(
1076 *   &AcousticOperator::local_apply_cell<VectorType>,
1077 *   &AcousticOperator::local_apply_face<VectorType>,
1078 *   &AcousticOperator::local_apply_boundary_face_point_to_point<
1079 *   VectorType>,
1080 *   this,
1081 *   dst,
1082 *   src,
1083 *   true,
1084 *   MatrixFree<dim, Number>::DataAccessOnFaces::values,
1085 *   MatrixFree<dim, Number>::DataAccessOnFaces::values);
1086 *   }
1087 *   }
1088 *  
1089 * @endcode
1090 *
1091 * In the `private` section of the class, we then define the
1092 * functions that evaluate volume, interior face, and boundary
1093 * face integrals. The concrete terms these functions evaluate are
1094 * stated in the documentation at the top of this tutorial
1095 * program. Each of these functions has its own object of type
1096 * `MaterialEvaluation` that provides access to the material at
1097 * each cell or face.
1098 *
1099 * @code
1100 *   private:
1101 *   template <typename VectorType>
1102 *   void local_apply_cell(
1103 *   const MatrixFree<dim, Number> &matrix_free,
1104 *   VectorType &dst,
1105 *   const VectorType &src,
1106 *   const std::pair<unsigned int, unsigned int> &cell_range) const
1107 *   {
1108 *   FEEvaluation<dim, -1, 0, 1, Number> pressure(matrix_free, 0, 0, 0);
1109 *   FEEvaluation<dim, -1, 0, dim, Number> velocity(matrix_free, 0, 0, 1);
1110 *  
1111 *   MaterialEvaluation material(matrix_free, *material_data);
1112 *  
1113 *   for (unsigned int cell = cell_range.first; cell < cell_range.second;
1114 *   ++cell)
1115 *   {
1116 *   velocity.reinit(cell);
1117 *   pressure.reinit(cell);
1118 *  
1119 *   pressure.gather_evaluate(src, EvaluationFlags::gradients);
1120 *   velocity.gather_evaluate(src, EvaluationFlags::gradients);
1121 *  
1122 * @endcode
1123 *
1124 * Get the materials on the corresponding cell. Since we
1125 * introduced @c MaterialEvaluation we can write the code
1126 * independent of whether the material is homogeneous or
1127 * inhomogeneous.
1128 *
1129 * @code
1130 *   material.reinit_cell(cell);
1131 *   const auto c = material.get_speed_of_sound();
1132 *   const auto rho = material.get_density();
1133 *  
1134 *   for (const unsigned int q : pressure.quadrature_point_indices())
1135 *   {
1136 *   pressure.submit_value(rho * c * c * velocity.get_divergence(q),
1137 *   q);
1138 *   velocity.submit_value(1.0 / rho * pressure.get_gradient(q), q);
1139 *   }
1140 *  
1141 *   pressure.integrate_scatter(EvaluationFlags::values, dst);
1142 *   velocity.integrate_scatter(EvaluationFlags::values, dst);
1143 *   }
1144 *   }
1145 *  
1146 * @endcode
1147 *
1148 * This next function evaluates the fluxes at faces between cells with the
1149 * same material. If boundary faces are under consideration fluxes into
1150 * neighboring faces do not have to be considered which is enforced via
1151 * `weight_neighbor=false`. For non-matching faces the fluxes into
1152 * neighboring faces are not considered as well. This is because we iterate
1153 * over each side of the non-matching face separately (similar to a cell
1154 * centric loop).
1155 *
1156
1157 *
1158 * In this and following functions, we also introduce the factors
1159 * @f$\tau@f$ and @f$\gamma@f$ that are computed from @f$\rho@f$ and @f$c@f$ along
1160 * interfaces and that appear in the bilinear forms. Their
1161 * definitions are provided in the introduction.
1162 *
1163 * @code
1164 *   template <bool weight_neighbor,
1165 *   typename InternalFaceIntegratorPressure,
1166 *   typename InternalFaceIntegratorVelocity,
1167 *   typename ExternalFaceIntegratorPressure,
1168 *   typename ExternalFaceIntegratorVelocity>
1169 *   inline DEAL_II_ALWAYS_INLINE void evaluate_face_kernel(
1170 *   InternalFaceIntegratorPressure &pressure_m,
1171 *   InternalFaceIntegratorVelocity &velocity_m,
1172 *   ExternalFaceIntegratorPressure &pressure_p,
1173 *   ExternalFaceIntegratorVelocity &velocity_p,
1174 *   const typename InternalFaceIntegratorPressure::value_type c,
1175 *   const typename InternalFaceIntegratorPressure::value_type rho) const
1176 *   {
1177 *   const auto tau = 0.5 * rho * c;
1178 *   const auto gamma = 0.5 / (rho * c);
1179 *  
1180 *   for (const unsigned int q : pressure_m.quadrature_point_indices())
1181 *   {
1182 *   const auto n = pressure_m.normal_vector(q);
1183 *   const auto pm = pressure_m.get_value(q);
1184 *   const auto um = velocity_m.get_value(q);
1185 *  
1186 *   const auto pp = pressure_p.get_value(q);
1187 *   const auto up = velocity_p.get_value(q);
1188 *  
1189 * @endcode
1190 *
1191 * Compute homogeneous local Lax-Friedrichs fluxes and submit the
1192 * corresponding values to the integrators.
1193 *
1194 * @code
1195 *   const auto momentum_flux =
1196 *   0.5 * (pm + pp) + 0.5 * tau * (um - up) * n;
1197 *   velocity_m.submit_value(1.0 / rho * (momentum_flux - pm) * n, q);
1198 *   if constexpr (weight_neighbor)
1199 *   velocity_p.submit_value(1.0 / rho * (momentum_flux - pp) * (-n), q);
1200 *  
1201 *   const auto mass_flux = 0.5 * (um + up) + 0.5 * gamma * (pm - pp) * n;
1202 *   pressure_m.submit_value(rho * c * c * (mass_flux - um) * n, q);
1203 *   if constexpr (weight_neighbor)
1204 *   pressure_p.submit_value(rho * c * c * (mass_flux - up) * (-n), q);
1205 *   }
1206 *   }
1207 *  
1208 * @endcode
1209 *
1210 * This function evaluates the fluxes at faces between cells with different
1211 * materials. This can only happen over non-matching interfaces. Therefore,
1212 * it is implicitly known that `weight_neighbor=false` and we can omit the
1213 * parameter.
1214 *
1215 * @code
1216 *   template <typename InternalFaceIntegratorPressure,
1217 *   typename InternalFaceIntegratorVelocity,
1218 *   typename ExternalFaceIntegratorPressure,
1219 *   typename ExternalFaceIntegratorVelocity,
1220 *   typename MaterialIntegrator>
1221 *   void evaluate_face_kernel_inhomogeneous(
1222 *   InternalFaceIntegratorPressure &pressure_m,
1223 *   InternalFaceIntegratorVelocity &velocity_m,
1224 *   const ExternalFaceIntegratorPressure &pressure_p,
1225 *   const ExternalFaceIntegratorVelocity &velocity_p,
1226 *   const typename InternalFaceIntegratorPressure::value_type c,
1227 *   const typename InternalFaceIntegratorPressure::value_type rho,
1228 *   const MaterialIntegrator &c_r,
1229 *   const MaterialIntegrator &rho_r) const
1230 *   {
1231 *   const auto tau_m = 0.5 * rho * c;
1232 *   const auto gamma_m = 0.5 / (rho * c);
1233 *  
1234 *   for (const unsigned int q : pressure_m.quadrature_point_indices())
1235 *   {
1236 *   const auto c_p = c_r.get_value(q);
1237 *   const auto rho_p = rho_r.get_value(q);
1238 *   const auto tau_p = 0.5 * rho_p * c_p;
1239 *   const auto gamma_p = 0.5 / (rho_p * c_p);
1240 *   const auto tau_sum_inv = 1.0 / (tau_m + tau_p);
1241 *   const auto gamma_sum_inv = 1.0 / (gamma_m + gamma_p);
1242 *  
1243 *   const auto n = pressure_m.normal_vector(q);
1244 *   const auto pm = pressure_m.get_value(q);
1245 *   const auto um = velocity_m.get_value(q);
1246 *  
1247 *   const auto pp = pressure_p.get_value(q);
1248 *   const auto up = velocity_p.get_value(q);
1249 *  
1250 *  
1251 * @endcode
1252 *
1253 * Compute inhomogeneous fluxes and submit the corresponding values
1254 * to the integrators.
1255 *
1256 * @code
1257 *   const auto momentum_flux =
1258 *   pm - tau_m * tau_sum_inv * (pm - pp) +
1259 *   tau_m * tau_p * tau_sum_inv * (um - up) * n;
1260 *   velocity_m.submit_value(1.0 / rho * (momentum_flux - pm) * n, q);
1261 *  
1262 *  
1263 *   const auto mass_flux =
1264 *   um - gamma_m * gamma_sum_inv * (um - up) +
1265 *   gamma_m * gamma_p * gamma_sum_inv * (pm - pp) * n;
1266 *  
1267 *   pressure_m.submit_value(rho * c * c * (mass_flux - um) * n, q);
1268 *   }
1269 *   }
1270 *  
1271 * @endcode
1272 *
1273 * This function evaluates the inner face integrals.
1274 *
1275 * @code
1276 *   template <typename VectorType>
1277 *   void local_apply_face(
1278 *   const MatrixFree<dim, Number> &matrix_free,
1279 *   VectorType &dst,
1280 *   const VectorType &src,
1281 *   const std::pair<unsigned int, unsigned int> &face_range) const
1282 *   {
1283 *   FEFaceEvaluation<dim, -1, 0, 1, Number> pressure_m(
1284 *   matrix_free, true, 0, 0, 0);
1285 *   FEFaceEvaluation<dim, -1, 0, 1, Number> pressure_p(
1286 *   matrix_free, false, 0, 0, 0);
1287 *   FEFaceEvaluation<dim, -1, 0, dim, Number> velocity_m(
1288 *   matrix_free, true, 0, 0, 1);
1289 *   FEFaceEvaluation<dim, -1, 0, dim, Number> velocity_p(
1290 *   matrix_free, false, 0, 0, 1);
1291 *  
1292 *   MaterialEvaluation material(matrix_free, *material_data);
1293 *  
1294 *   for (unsigned int face = face_range.first; face < face_range.second;
1295 *   ++face)
1296 *   {
1297 *   velocity_m.reinit(face);
1298 *   velocity_p.reinit(face);
1299 *  
1300 *   pressure_m.reinit(face);
1301 *   pressure_p.reinit(face);
1302 *  
1303 *   pressure_m.gather_evaluate(src, EvaluationFlags::values);
1304 *   pressure_p.gather_evaluate(src, EvaluationFlags::values);
1305 *  
1306 *   velocity_m.gather_evaluate(src, EvaluationFlags::values);
1307 *   velocity_p.gather_evaluate(src, EvaluationFlags::values);
1308 *  
1309 *   material.reinit_face(face);
1310 *   evaluate_face_kernel<true>(pressure_m,
1311 *   velocity_m,
1312 *   pressure_p,
1313 *   velocity_p,
1314 *   material.get_speed_of_sound(),
1315 *   material.get_density());
1316 *  
1317 *   pressure_m.integrate_scatter(EvaluationFlags::values, dst);
1318 *   pressure_p.integrate_scatter(EvaluationFlags::values, dst);
1319 *   velocity_m.integrate_scatter(EvaluationFlags::values, dst);
1320 *   velocity_p.integrate_scatter(EvaluationFlags::values, dst);
1321 *   }
1322 *   }
1323 *  
1324 *  
1325 * @endcode
1326 *
1327 *
1328 * <a name="step_89-Matrixfreeboundaryfunctionforpointtopointinterpolation"></a>
1329 * <h4>Matrix-free boundary function for point-to-point interpolation</h4>
1330 *
1331
1332 *
1333 * The following function then evaluates the boundary face integrals and the
1334 * non-matching face integrals using point-to-point interpolation.
1335 *
1336 * @code
1337 *   template <typename VectorType>
1338 *   void local_apply_boundary_face_point_to_point(
1339 *   const MatrixFree<dim, Number> &matrix_free,
1340 *   VectorType &dst,
1341 *   const VectorType &src,
1342 *   const std::pair<unsigned int, unsigned int> &face_range) const
1343 *   {
1344 *   FEFaceEvaluation<dim, -1, 0, 1, Number> pressure_m(
1345 *   matrix_free, true, 0, 0, 0);
1346 *   FEFaceEvaluation<dim, -1, 0, dim, Number> velocity_m(
1347 *   matrix_free, true, 0, 0, 1);
1348 *  
1349 *   BCEvaluationP pressure_bc(pressure_m);
1350 *   BCEvaluationU velocity_bc(velocity_m);
1351 *  
1352 *   MaterialEvaluation material(matrix_free, *material_data);
1353 *  
1354 * @endcode
1355 *
1356 * Remote evaluators.
1357 *
1358 * @code
1359 *   auto pressure_r = pressure_r_eval->get_data_accessor();
1360 *   auto velocity_r = velocity_r_eval->get_data_accessor();
1361 *   auto c_r = c_r_eval->get_data_accessor();
1362 *   auto rho_r = rho_r_eval->get_data_accessor();
1363 *  
1364 *   for (unsigned int face = face_range.first; face < face_range.second;
1365 *   ++face)
1366 *   {
1367 *   velocity_m.reinit(face);
1368 *   pressure_m.reinit(face);
1369 *  
1370 *   pressure_m.gather_evaluate(src, EvaluationFlags::values);
1371 *   velocity_m.gather_evaluate(src, EvaluationFlags::values);
1372 *  
1373 *   if (HelperFunctions::is_non_matching_face(
1374 *   remote_face_ids, matrix_free.get_boundary_id(face)))
1375 *   {
1376 * @endcode
1377 *
1378 * If @c face is non-matching we have to query values via the
1379 * FERemoteEvaluaton objects. This is done by passing the
1380 * corresponding FERemoteEvaluaton objects to the function that
1381 * evaluates the kernel. As mentioned above, each side of the
1382 * non-matching interface is traversed separately and we do not
1383 * have to consider the neighbor in the kernel. Note, that the
1384 * values in the FERemoteEvaluaton objects are already updated at
1385 * this point.
1386 *
1387
1388 *
1389 * For point-to-point interpolation we simply use the
1390 * corresponding FERemoteEvaluaton objects in combination with the
1391 * standard FEFaceEvaluation objects.
1392 *
1393 * @code
1394 *   velocity_r.reinit(face);
1395 *   pressure_r.reinit(face);
1396 *  
1397 *   material.reinit_face(face);
1398 *  
1399 * @endcode
1400 *
1401 * If we are considering a homogeneous material, do not use the
1402 * inhomogeneous fluxes. While it would be possible
1403 * to use the inhomogeneous fluxes they are more expensive to
1404 * compute.
1405 *
1406 * @code
1407 *   if (material.is_homogeneous())
1408 *   {
1409 *   evaluate_face_kernel<false>(pressure_m,
1410 *   velocity_m,
1411 *   pressure_r,
1412 *   velocity_r,
1413 *   material.get_speed_of_sound(),
1414 *   material.get_density());
1415 *   }
1416 *   else
1417 *   {
1418 *   c_r.reinit(face);
1419 *   rho_r.reinit(face);
1420 *   evaluate_face_kernel_inhomogeneous(
1421 *   pressure_m,
1422 *   velocity_m,
1423 *   pressure_r,
1424 *   velocity_r,
1425 *   material.get_speed_of_sound(),
1426 *   material.get_density(),
1427 *   c_r,
1428 *   rho_r);
1429 *   }
1430 *   }
1431 *   else
1432 *   {
1433 * @endcode
1434 *
1435 * If @c face is a standard boundary face, evaluate the integral
1436 * as usual in the matrix free context. To be able to use the same
1437 * kernel as for inner faces we pass the boundary condition
1438 * objects to the function that evaluates the kernel. As detailed
1439 * above `weight_neighbor=false`.
1440 *
1441 * @code
1442 *   material.reinit_face(face);
1443 *   evaluate_face_kernel<false>(pressure_m,
1444 *   velocity_m,
1445 *   pressure_bc,
1446 *   velocity_bc,
1447 *   material.get_speed_of_sound(),
1448 *   material.get_density());
1449 *   }
1450 *  
1451 *   pressure_m.integrate_scatter(EvaluationFlags::values, dst);
1452 *   velocity_m.integrate_scatter(EvaluationFlags::values, dst);
1453 *   }
1454 *   }
1455 *  
1456 * @endcode
1457 *
1458 *
1459 * <a name="step_89-MatrixfreeboundaryfunctionforNitschetypemortaring"></a>
1460 * <h4>Matrix-free boundary function for Nitsche-type mortaring</h4>
1461 *
1462
1463 *
1464 * This function evaluates the boundary face integrals and the
1465 * non-matching face integrals using Nitsche-type mortaring.
1466 *
1467 * @code
1468 *   template <typename VectorType>
1469 *   void local_apply_boundary_face_mortaring(
1470 *   const MatrixFree<dim, Number> &matrix_free,
1471 *   VectorType &dst,
1472 *   const VectorType &src,
1473 *   const std::pair<unsigned int, unsigned int> &face_range) const
1474 *   {
1475 *   FEFaceEvaluation<dim, -1, 0, 1, Number> pressure_m(
1476 *   matrix_free, true, 0, 0, 0);
1477 *   FEFaceEvaluation<dim, -1, 0, dim, Number> velocity_m(
1478 *   matrix_free, true, 0, 0, 1);
1479 *  
1480 * @endcode
1481 *
1482 * For Nitsche-type mortaring we are evaluating the integrals over
1483 * intersections. This is why, quadrature points are arbitrarily
1484 * distributed on every face. Thus, we can not make use of face batches
1485 * and FEFaceEvaluation but have to consider each face individually and
1486 * make use of @c FEFacePointEvaluation to evaluate the integrals in the
1487 * arbitrarily distributed quadrature points.
1488 * Since the setup of FEFacePointEvaluation is more expensive than that of
1489 * FEEvaluation we do the setup only once. For this we are using the
1490 * helper function @c get_thread_safe_fe_face_point_evaluation_object().
1491 *
1492 * @code
1493 *   FEFacePointEvaluation<1, dim, dim, Number> &pressure_m_mortar =
1494 *   get_thread_safe_fe_face_point_evaluation_object<1>(
1495 *   thread_local_pressure_m_mortar, 0);
1496 *   FEFacePointEvaluation<dim, dim, dim, Number> &velocity_m_mortar =
1497 *   get_thread_safe_fe_face_point_evaluation_object<dim>(
1498 *   thread_local_velocity_m_mortar, 1);
1499 *  
1500 *   BCEvaluationP pressure_bc(pressure_m);
1501 *   BCEvaluationU velocity_bc(velocity_m);
1502 *  
1503 *   MaterialEvaluation material(matrix_free, *material_data);
1504 *  
1505 *   auto pressure_r_mortar = pressure_r_eval->get_data_accessor();
1506 *   auto velocity_r_mortar = velocity_r_eval->get_data_accessor();
1507 *   auto c_r = c_r_eval->get_data_accessor();
1508 *   auto rho_r = rho_r_eval->get_data_accessor();
1509 *  
1510 *   for (unsigned int face = face_range.first; face < face_range.second;
1511 *   ++face)
1512 *   {
1513 *   if (HelperFunctions::is_non_matching_face(
1514 *   remote_face_ids, matrix_free.get_boundary_id(face)))
1515 *   {
1516 *   material.reinit_face(face);
1517 *  
1518 * @endcode
1519 *
1520 * First fetch the DoF values with standard FEFaceEvaluation
1521 * objects.
1522 *
1523 * @code
1524 *   pressure_m.reinit(face);
1525 *   velocity_m.reinit(face);
1526 *  
1527 *   pressure_m.read_dof_values(src);
1528 *   velocity_m.read_dof_values(src);
1529 *  
1530 * @endcode
1531 *
1532 * Project the internally stored values into the face DoFs
1533 * of the current face.
1534 *
1535 * @code
1536 *   pressure_m.project_to_face(EvaluationFlags::values);
1537 *   velocity_m.project_to_face(EvaluationFlags::values);
1538 *  
1539 * @endcode
1540 *
1541 * For mortaring, we have to consider every face from the face
1542 * batches separately and have to use the FEFacePointEvaluation
1543 * objects to be able to evaluate the integrals with the
1544 * arbitrarily distributed quadrature points.
1545 *
1546 * @code
1547 *   for (unsigned int v = 0;
1548 *   v < matrix_free.n_active_entries_per_face_batch(face);
1549 *   ++v)
1550 *   {
1551 *   constexpr unsigned int n_lanes =
1552 *   VectorizedArray<Number>::size();
1553 *   velocity_m_mortar.reinit(face * n_lanes + v);
1554 *   pressure_m_mortar.reinit(face * n_lanes + v);
1555 *  
1556 * @endcode
1557 *
1558 * Evaluate using FEFacePointEvaluation. As buffer,
1559 * simply use the internal buffers from the
1560 * FEFaceEvaluation objects.
1561 *
1562 * @code
1563 *   velocity_m_mortar.evaluate_in_face(
1564 *   &velocity_m.get_scratch_data().begin()[0][v],
1565 *   EvaluationFlags::values);
1566 *  
1567 *   pressure_m_mortar.evaluate_in_face(
1568 *   &pressure_m.get_scratch_data().begin()[0][v],
1569 *   EvaluationFlags::values);
1570 *  
1571 *   velocity_r_mortar.reinit(face * n_lanes + v);
1572 *   pressure_r_mortar.reinit(face * n_lanes + v);
1573 *  
1574 * @endcode
1575 *
1576 * As above, if we are considering a homogeneous
1577 * material, do not use the inhomogeneous
1578 * fluxes. Since we are operating on face @c v we
1579 * call @c material.get_density()[v].
1580 *
1581 * @code
1582 *   if (material.is_homogeneous())
1583 *   {
1584 *   evaluate_face_kernel<false>(
1585 *   pressure_m_mortar,
1586 *   velocity_m_mortar,
1587 *   pressure_r_mortar,
1588 *   velocity_r_mortar,
1589 *   material.get_speed_of_sound()[v],
1590 *   material.get_density()[v]);
1591 *   }
1592 *   else
1593 *   {
1594 *   c_r.reinit(face * n_lanes + v);
1595 *   rho_r.reinit(face * n_lanes + v);
1596 *  
1597 *   evaluate_face_kernel_inhomogeneous(
1598 *   pressure_m_mortar,
1599 *   velocity_m_mortar,
1600 *   pressure_r_mortar,
1601 *   velocity_r_mortar,
1602 *   material.get_speed_of_sound()[v],
1603 *   material.get_density()[v],
1604 *   c_r,
1605 *   rho_r);
1606 *   }
1607 *  
1608 * @endcode
1609 *
1610 * Integrate using FEFacePointEvaluation. As buffer,
1611 * simply use the internal buffers from the
1612 * FEFaceEvaluation objects.
1613 *
1614 * @code
1615 *   velocity_m_mortar.integrate_in_face(
1616 *   &velocity_m.get_scratch_data().begin()[0][v],
1617 *   EvaluationFlags::values);
1618 *  
1619 *   pressure_m_mortar.integrate_in_face(
1620 *   &pressure_m.get_scratch_data().begin()[0][v],
1621 *   EvaluationFlags::values);
1622 *   }
1623 *  
1624 * @endcode
1625 *
1626 * Collect the contributions from the face DoFs to
1627 * the internal cell DoFs to be able to use the
1628 * member function @c distribute_local_to_global().
1629 *
1630 * @code
1631 *   pressure_m.collect_from_face(EvaluationFlags::values);
1632 *   velocity_m.collect_from_face(EvaluationFlags::values);
1633 *  
1634 *   pressure_m.distribute_local_to_global(dst);
1635 *   velocity_m.distribute_local_to_global(dst);
1636 *   }
1637 *   else
1638 *   {
1639 *   velocity_m.reinit(face);
1640 *   pressure_m.reinit(face);
1641 *  
1642 *   pressure_m.gather_evaluate(src, EvaluationFlags::values);
1643 *   velocity_m.gather_evaluate(src, EvaluationFlags::values);
1644 *  
1645 *   material.reinit_face(face);
1646 *   evaluate_face_kernel<false>(pressure_m,
1647 *   velocity_m,
1648 *   pressure_bc,
1649 *   velocity_bc,
1650 *   material.get_speed_of_sound(),
1651 *   material.get_density());
1652 *  
1653 *   pressure_m.integrate_scatter(EvaluationFlags::values, dst);
1654 *   velocity_m.integrate_scatter(EvaluationFlags::values, dst);
1655 *   }
1656 *   }
1657 *   }
1658 *  
1659 *   const MatrixFree<dim, Number> &matrix_free;
1660 *  
1661 *   const std::shared_ptr<CellwiseMaterialData<Number>> material_data;
1662 *  
1663 *   const std::set<types::boundary_id> remote_face_ids;
1664 *  
1665 * @endcode
1666 *
1667 * FERemoteEvaluation objects are stored as shared pointers. This way,
1668 * they can also be used for other operators without precomputing the values
1669 * multiple times.
1670 *
1671 * @code
1672 *   const std::shared_ptr<FERemoteEvaluation<dim, 1, remote_value_type>>
1673 *   pressure_r_eval;
1674 *   const std::shared_ptr<FERemoteEvaluation<dim, dim, remote_value_type>>
1675 *   velocity_r_eval;
1676 *  
1677 *   const std::shared_ptr<FERemoteEvaluation<dim, 1, remote_value_type>>
1678 *   c_r_eval;
1679 *   const std::shared_ptr<FERemoteEvaluation<dim, 1, remote_value_type>>
1680 *   rho_r_eval;
1681 *  
1682 *   const std::shared_ptr<NonMatching::MappingInfo<dim, dim, Number>>
1683 *   nm_mapping_info;
1684 *  
1685 * @endcode
1686 *
1687 * We store FEFacePointEvaluation objects as members in a thread local
1688 * way, since its creation is more expensive compared to FEEvaluation
1689 * objects.
1690 *
1691 * @code
1692 *   mutable Threads::ThreadLocalStorage<
1693 *   std::unique_ptr<FEFacePointEvaluation<1, dim, dim, Number>>>
1694 *   thread_local_pressure_m_mortar;
1695 *  
1696 *   mutable Threads::ThreadLocalStorage<
1697 *   std::unique_ptr<FEFacePointEvaluation<dim, dim, dim, Number>>>
1698 *   thread_local_velocity_m_mortar;
1699 *  
1700 * @endcode
1701 *
1702 * Helper function to create and get FEFacePointEvaluation objects in a
1703 * thread safe way. On each thread, FEFacePointEvaluation is created if it
1704 * has not been created by now. After that, simply return the object
1705 * corresponding to the thread under consideration.
1706 *
1707 * @code
1708 *   template <int n_components>
1709 *   FEFacePointEvaluation<n_components, dim, dim, Number> &
1710 *   get_thread_safe_fe_face_point_evaluation_object(
1711 *   Threads::ThreadLocalStorage<
1712 *   std::unique_ptr<FEFacePointEvaluation<n_components, dim, dim, Number>>>
1713 *   &fe_face_point_eval_thread_local,
1714 *   unsigned int fist_selected_comp) const
1715 *   {
1716 *   if (fe_face_point_eval_thread_local.get() == nullptr)
1717 *   {
1718 *   fe_face_point_eval_thread_local = std::make_unique<
1719 *   FEFacePointEvaluation<n_components, dim, dim, Number>>(
1720 *   *nm_mapping_info,
1721 *   matrix_free.get_dof_handler().get_fe(),
1722 *   true,
1723 *   fist_selected_comp);
1724 *   }
1725 *   return *fe_face_point_eval_thread_local.get();
1726 *   }
1727 *   };
1728 *  
1729 * @endcode
1730 *
1731 *
1732 * <a name="step_89-Inversemassoperator"></a>
1733 * <h3>Inverse mass operator</h3>
1734 *
1735
1736 *
1737 * For the time stepping methods below, we need the inverse mass
1738 * operator. We apply it via a loop over all (batches of) cells as
1739 * always, where the contribution of each cell is computed in a
1740 * matrix-free way:
1741 *
1742 * @code
1743 *   template <int dim, typename Number>
1744 *   class InverseMassOperator
1745 *   {
1746 *   public:
1747 *   InverseMassOperator(const MatrixFree<dim, Number> &matrix_free)
1748 *   : matrix_free(matrix_free)
1749 *   {}
1750 *  
1751 *   template <typename VectorType>
1752 *   void apply(VectorType &dst, const VectorType &src) const
1753 *   {
1754 *   dst.zero_out_ghost_values();
1755 *   matrix_free.cell_loop(&InverseMassOperator::local_apply_cell<VectorType>,
1756 *   this,
1757 *   dst,
1758 *   src);
1759 *   }
1760 *  
1761 *   private:
1762 *   template <typename VectorType>
1763 *   void local_apply_cell(
1764 *   const MatrixFree<dim, Number> &mf,
1765 *   VectorType &dst,
1766 *   const VectorType &src,
1767 *   const std::pair<unsigned int, unsigned int> &cell_range) const
1768 *   {
1769 *   FEEvaluation<dim, -1, 0, dim + 1, Number> phi(mf);
1770 *   MatrixFreeOperators::CellwiseInverseMassMatrix<dim, -1, dim + 1, Number>
1771 *   minv(phi);
1772 *  
1773 *   for (unsigned int cell = cell_range.first; cell < cell_range.second;
1774 *   ++cell)
1775 *   {
1776 *   phi.reinit(cell);
1777 *   phi.read_dof_values(src);
1778 *   minv.apply(phi.begin_dof_values(), phi.begin_dof_values());
1779 *   phi.set_dof_values(dst);
1780 *   }
1781 *   }
1782 *  
1783 *   const MatrixFree<dim, Number> &matrix_free;
1784 *   };
1785 *  
1786 * @endcode
1787 *
1788 *
1789 * <a name="step_89-RungeKuttatimestepping"></a>
1790 * <h3>Runge-Kutta time-stepping</h3>
1791 *
1792
1793 *
1794 * This class implements a Runge-Kutta scheme of order 2.
1795 *
1796 * @code
1797 *   template <int dim, typename Number, typename remote_value_type>
1798 *   class RungeKutta2
1799 *   {
1800 *   using VectorType = LinearAlgebra::distributed::Vector<Number>;
1801 *  
1802 *   public:
1803 *   RungeKutta2(
1804 *   const std::shared_ptr<InverseMassOperator<dim, Number>>
1805 *   inverse_mass_operator,
1806 *   const std::shared_ptr<AcousticOperator<dim, Number, remote_value_type>>
1807 *   acoustic_operator)
1808 *   : inverse_mass_operator(inverse_mass_operator)
1809 *   , acoustic_operator(acoustic_operator)
1810 *   {}
1811 *  
1812 * @endcode
1813 *
1814 * The `run()` function of this class contains the time loop. It
1815 * starts by initializing some member variables (such as
1816 * short-cuts to objects that describe the finite element, its
1817 * properties, and the mapping; an by initializing vectors). It
1818 * then computes the time step size via minimum element edge
1819 * length. Assuming non-distorted elements, we can compute the
1820 * edge length as the distance between two vertices. From this,
1821 * we can then obtain the time step size via the CFL condition.
1822 *
1823 * @code
1824 *   void run(const MatrixFree<dim, Number> &matrix_free,
1825 *   const double cr,
1826 *   const double end_time,
1827 *   const double speed_of_sound,
1828 *   const Function<dim> &initial_condition,
1829 *   const std::string &vtk_prefix)
1830 *   {
1831 *   const auto &dof_handler = matrix_free.get_dof_handler();
1832 *   const auto &mapping = *matrix_free.get_mapping_info().mapping;
1833 *   const auto degree = dof_handler.get_fe().degree;
1834 *  
1835 *   VectorType solution;
1836 *   matrix_free.initialize_dof_vector(solution);
1837 *   VectorType solution_temp;
1838 *   matrix_free.initialize_dof_vector(solution_temp);
1839 *  
1840 *   HelperFunctions::set_initial_condition(matrix_free,
1841 *   initial_condition,
1842 *   solution);
1843 *  
1844 *   double h_local_min = std::numeric_limits<double>::max();
1845 *   for (const auto &cell : dof_handler.active_cell_iterators())
1846 *   h_local_min =
1847 *   std::min(h_local_min,
1848 *   (cell->vertex(1) - cell->vertex(0)).norm_square());
1849 *   h_local_min = std::sqrt(h_local_min);
1850 *   const double h_min =
1851 *   Utilities::MPI::min(h_local_min, dof_handler.get_mpi_communicator());
1852 *  
1853 *   const double dt =
1854 *   cr * HelperFunctions::compute_dt_cfl(h_min, degree, speed_of_sound);
1855 *  
1856 * @endcode
1857 *
1858 * We can then perform the time integration loop:
1859 *
1860 * @code
1861 *   double time = 0.0;
1862 *   unsigned int timestep = 0;
1863 *   while (time < end_time)
1864 *   {
1865 *   HelperFunctions::write_vtu(solution,
1866 *   matrix_free.get_dof_handler(),
1867 *   mapping,
1868 *   degree,
1869 *   "step_89-" + vtk_prefix +
1870 *   std::to_string(timestep));
1871 *  
1872 *   std::swap(solution, solution_temp);
1873 *   time += dt;
1874 *   ++timestep;
1875 *   perform_time_step(dt, solution, solution_temp);
1876 *   }
1877 *   }
1878 *  
1879 * @endcode
1880 *
1881 * The main work of this class is done by a `private` member
1882 * function that performs one Runge-Kutta 2 time step. Recall that
1883 * this method requires two sub-steps ("stages") computing
1884 * intermediate values `k1` and `k2`, after which the intermediate
1885 * values are summed with weights to obtain the new solution at
1886 * the end of the time step. The RK2 method allows for the
1887 * elimination of the intermediate vector `k2` by utilizing the
1888 * `dst` vector as temporary storage.
1889 *
1890 * @code
1891 *   private:
1892 *   void
1893 *   perform_time_step(const double dt, VectorType &dst, const VectorType &src)
1894 *   {
1895 *   VectorType k1 = src;
1896 *  
1897 *   /* First stage. */
1898 *   evaluate_stage(k1, src);
1899 *  
1900 *   /* Second stage. */
1901 *   k1.sadd(0.5 * dt, 1.0, src);
1902 *   evaluate_stage(dst, k1);
1903 *  
1904 *   /* Summing things into the output vector. */
1905 *   dst.sadd(dt, 1.0, src);
1906 *   }
1907 *  
1908 * @endcode
1909 *
1910 * Evaluating a single Runge-Kutta stage is a straightforward step
1911 * that really only requires applying the operator once, and then
1912 * applying the inverse of the mass matrix.
1913 *
1914 * @code
1915 *   void evaluate_stage(VectorType &dst, const VectorType &src)
1916 *   {
1917 *   acoustic_operator->evaluate(dst, src);
1918 *   dst *= -1.0;
1919 *   inverse_mass_operator->apply(dst, dst);
1920 *   }
1921 *  
1922 *   const std::shared_ptr<InverseMassOperator<dim, Number>>
1923 *   inverse_mass_operator;
1924 *   const std::shared_ptr<AcousticOperator<dim, Number, remote_value_type>>
1925 *   acoustic_operator;
1926 *   };
1927 *  
1928 *  
1929 * @endcode
1930 *
1931 *
1932 * <a name="step_89-Constructionofnonmatchingtriangulations"></a>
1933 * <h3>Construction of non-matching triangulations</h3>
1934 *
1935
1936 *
1937 * Let us now make our way to the higher-level functions of this program.
1938 *
1939
1940 *
1941 * The first of these functions creates a two dimensional square
1942 * triangulation that spans from @f$(0,0)@f$ to @f$(1,1)@f$. It consists of
1943 * two sub-domains. The left sub-domain spans from @f$(0,0)@f$ to
1944 * @f$(0.525,1)@f$. The right sub-domain spans from @f$(0.525,0)@f$ to
1945 * @f$(1,1)@f$. The left sub-domain has elements that are three times
1946 * smaller compared to the ones for the right sub-domain.
1947 *
1948
1949 *
1950 * At non-matching interfaces, we need to provide different boundary
1951 * IDs for the cells that make up the two parts (because, while they
1952 * may be physically adjacent, they are not logically neighbors
1953 * given that the faces of cells on both sides do not match, and the
1954 * Triangulation class will therefore treat the interface between
1955 * the two parts as a "boundary"). These boundary IDs have to differ
1956 * because later on RemotePointEvaluation has to search for remote
1957 * points for each face, that are defined in the same mesh (since we
1958 * merge the mesh) but not on the same side of the non-matching
1959 * interface. As a consequence, we declare at the top symbolic names
1960 * for these boundary indicators, and ensure that we return a set
1961 * with these values to the caller for later use.
1962 *
1963
1964 *
1965 * The actual mesh is then constructed by first constructing the
1966 * left and right parts separately (setting material ids to zero and
1967 * one, respectively), and using the appropriate boundary ids for
1968 * all parts of the mesh. We then use
1969 * GridGenerator::merge_triangulations() to combine them into one
1970 * (non-matching) mesh. We have to pay attention that should the two
1971 * sub-triangulations have vertices at the same locations, that they
1972 * are not merged (connecting the two triangulations logically)
1973 * since we want the interface to be an actual boundary. We achieve
1974 * this by providing a tolerance of zero for the merge, see the
1975 * documentation of the function
1976 * GridGenerator::merge_triangulations().
1977 *
1978 * @code
1979 *   template <int dim>
1980 *   void build_non_matching_triangulation(
1981 *   Triangulation<dim> &tria,
1982 *   std::set<types::boundary_id> &non_matching_faces,
1983 *   const unsigned int refinements)
1984 *   {
1985 *   const double length = 1.0;
1986 *  
1987 *   const types::boundary_id non_matching_id_left = 98;
1988 *   const types::boundary_id non_matching_id_right = 99;
1989 *  
1990 *   non_matching_faces = {non_matching_id_left, non_matching_id_right};
1991 *  
1992 *   /* Construct left part of mesh. */
1993 *   Triangulation<dim> tria_left;
1994 *   const unsigned int subdiv_left = 11;
1995 *   GridGenerator::subdivided_hyper_rectangle(tria_left,
1996 *   {subdiv_left, 2 * subdiv_left},
1997 *   {0.0, 0.0},
1998 *   {0.525 * length, length});
1999 *  
2000 *   for (const auto &cell : tria_left.active_cell_iterators())
2001 *   cell->set_material_id(0);
2002 *   for (const auto &face : tria_left.active_face_iterators())
2003 *   if (face->at_boundary())
2004 *   {
2005 *   face->set_boundary_id(0);
2006 *   if (face->center()[0] > 0.525 * length - 1e-6)
2007 *   face->set_boundary_id(non_matching_id_left);
2008 *   }
2009 *  
2010 *   /* Construct right part of mesh. */
2011 *   Triangulation<dim> tria_right;
2012 *   const unsigned int subdiv_right = 4;
2013 *   GridGenerator::subdivided_hyper_rectangle(tria_right,
2014 *   {subdiv_right, 2 * subdiv_right},
2015 *   {0.525 * length, 0.0},
2016 *   {length, length});
2017 *  
2018 *   for (const auto &cell : tria_right.active_cell_iterators())
2019 *   cell->set_material_id(1);
2020 *   for (const auto &face : tria_right.active_face_iterators())
2021 *   if (face->at_boundary())
2022 *   {
2023 *   face->set_boundary_id(0);
2024 *   if (face->center()[0] < 0.525 * length + 1e-6)
2025 *   face->set_boundary_id(non_matching_id_right);
2026 *   }
2027 *  
2028 *   /* Merge triangulations. */
2029 *   GridGenerator::merge_triangulations(tria_left,
2030 *   tria_right,
2031 *   tria,
2032 *   /*tolerance*/ 0.,
2033 *   /*copy_manifold_ids*/ false,
2034 *   /*copy_boundary_ids*/ true);
2035 *  
2036 *   /* Refine the result. */
2037 *   tria.refine_global(refinements);
2038 *   }
2039 *  
2040 * @endcode
2041 *
2042 *
2043 * <a name="step_89-Setupandrunningofthetwoschemes"></a>
2044 * <h3>Set up and running of the two schemes</h3>
2045 *
2046
2047 *
2048 *
2049 * <a name="step_89-Setupandrunningofthepointtopointinterpolationscheme"></a>
2050 * <h4>Setup and running of the point-to-point interpolation scheme</h4>
2051 *
2052
2053 *
2054 * We are now at the two functions that run the overall schemes (the
2055 * point-to-point and the mortaring schemes). The first of these
2056 * functions fills a FERemoteEvaluationCommunicator object that is
2057 * needed for point-to-point interpolation. Additionally, the
2058 * corresponding remote evaluators are set up using this remote
2059 * communicator. Eventually, the operators are handed to the time
2060 * integrator that runs the simulation.
2061 *
2062 * @code
2063 *   template <int dim, typename Number>
2064 *   void run_with_point_to_point_interpolation(
2065 *   const MatrixFree<dim, Number> &matrix_free,
2066 *   const std::set<types::boundary_id> &non_matching_faces,
2067 *   const std::map<types::material_id, std::pair<double, double>> &materials,
2068 *   const double end_time,
2069 *   const Function<dim> &initial_condition,
2070 *   const std::string &vtk_prefix)
2071 *   {
2072 *   const auto &dof_handler = matrix_free.get_dof_handler();
2073 *   const auto &tria = dof_handler.get_triangulation();
2074 *  
2075 * @endcode
2076 *
2077 * Communication objects know about the communication pattern. That is,
2078 * they know about the cells and quadrature points that have to be
2079 * evaluated at remote faces. This information is given via
2080 * RemotePointEvaluation. Additionally, the communication objects
2081 * have to be able to match the quadrature points of the remote
2082 * points (that provide exterior information) to the quadrature points
2083 * defined at the interior cell. In case of point-to-point interpolation
2084 * a vector of pairs with face batch Ids and the number of faces in the
2085 * batch is needed. @c FERemoteCommunicationObjectEntityBatches
2086 * is a container to store this information.
2087 *
2088
2089 *
2090 * The information is filled outside of the actual class since in some cases
2091 * the information is available from some heuristic and
2092 * it is possible to skip some expensive operations. This is for example
2093 * the case for sliding rotating interfaces with equally spaced elements on
2094 * both sides of the non-matching interface @cite duerrwaechter2021an.
2095 *
2096
2097 *
2098 * For the standard case of point to point-to-point interpolation without
2099 * any heuristic we make use of the utility function
2100 * Utilities::compute_remote_communicator_faces_point_to_point_interpolation().
2101 * Please refer to the documentation of this function to see how to manually
2102 * set up the remote communicator from outside.
2103 *
2104
2105 *
2106 * Among the inputs for the remote communicator we need a list of
2107 * boundary ids for the non-matching faces, along with a function
2108 * object for each boundary id that returns a vector of true/false
2109 * flags in which exactly the vertices of cells of the
2110 * triangulation are marked that have a face at the boundary id in
2111 * question.
2112 *
2113 * @code
2114 *   std::vector<
2115 *   std::pair<types::boundary_id, std::function<std::vector<bool>()>>>
2116 *   non_matching_faces_marked_vertices;
2117 *   for (const auto &nm_face : non_matching_faces)
2118 *   {
2119 *   auto marked_vertices = [&nm_face, &tria]() -> std::vector<bool> {
2120 *   std::vector<bool> mask(tria.n_vertices(), true);
2121 *  
2122 *   for (const auto &cell : tria.active_cell_iterators())
2123 *   for (auto const &f : cell->face_indices())
2124 *   if (cell->face(f)->at_boundary() &&
2125 *   cell->face(f)->boundary_id() == nm_face)
2126 *   for (const auto v : cell->vertex_indices())
2127 *   mask[cell->vertex_index(v)] = false;
2128 *  
2129 *   return mask;
2130 *   };
2131 *  
2132 *   non_matching_faces_marked_vertices.emplace_back(
2133 *   std::make_pair(nm_face, marked_vertices));
2134 *   }
2135 *  
2136 *   auto remote_communicator =
2137 *   Utilities::compute_remote_communicator_faces_point_to_point_interpolation(
2138 *   matrix_free, non_matching_faces_marked_vertices);
2139 *  
2140 * @endcode
2141 *
2142 * We are using point-to-point interpolation and can therefore
2143 * easily access the corresponding data at face batches. This
2144 * is why we use a @c VectorizedArray as @c remote_value_type:
2145 *
2146 * @code
2147 *   using remote_value_type = VectorizedArray<Number>;
2148 *  
2149 * @endcode
2150 *
2151 * We then set up FERemoteEvaluation objects that access the
2152 * pressure and velocity at remote faces, along with an object to
2153 * describe cell-wise material data.
2154 *
2155 * @code
2156 *   const auto pressure_r =
2157 *   std::make_shared<FERemoteEvaluation<dim, 1, remote_value_type>>(
2158 *   remote_communicator, dof_handler, /*first_selected_component*/ 0);
2159 *  
2160 *   const auto velocity_r =
2161 *   std::make_shared<FERemoteEvaluation<dim, dim, remote_value_type>>(
2162 *   remote_communicator, dof_handler, /*first_selected_component*/ 1);
2163 *  
2164 *   const auto material_data =
2165 *   std::make_shared<CellwiseMaterialData<Number>>(matrix_free, materials);
2166 *  
2167 * @endcode
2168 *
2169 * If we have an inhomogeneous problem, we have to set up the
2170 * material handler that accesses the materials at remote faces.
2171 *
2172 * @code
2173 *   const auto c_r =
2174 *   std::make_shared<FERemoteEvaluation<dim, 1, remote_value_type>>(
2175 *   remote_communicator,
2176 *   matrix_free.get_dof_handler().get_triangulation(),
2177 *   /*first_selected_component*/ 0);
2178 *   const auto rho_r =
2179 *   std::make_shared<FERemoteEvaluation<dim, 1, remote_value_type>>(
2180 *   remote_communicator,
2181 *   matrix_free.get_dof_handler().get_triangulation(),
2182 *   /*first_selected_component*/ 0);
2183 *  
2184 * @endcode
2185 *
2186 * If the domain is not homogeneous, i.e., if material parameters
2187 * change from cell to cell, we initialize and fill DoF vectors
2188 * that contain the material properties. Materials do not change
2189 * during the simulation, therefore there is no need to ever
2190 * compute the values after the first @c gather_evaluate() (at the
2191 * end of the following block) again.
2192 *
2193 * @code
2194 *   if (!material_data->is_homogeneous())
2195 *   {
2196 *   Vector<Number> c(
2197 *   matrix_free.get_dof_handler().get_triangulation().n_active_cells());
2198 *   Vector<Number> rho(
2199 *   matrix_free.get_dof_handler().get_triangulation().n_active_cells());
2200 *  
2201 *   for (const auto &cell : matrix_free.get_dof_handler()
2202 *   .get_triangulation()
2203 *   .active_cell_iterators())
2204 *   {
2205 *   c[cell->active_cell_index()] =
2206 *   materials.at(cell->material_id()).first;
2207 *   rho[cell->active_cell_index()] =
2208 *   materials.at(cell->material_id()).second;
2209 *   }
2210 *  
2211 *   c_r->gather_evaluate(c, EvaluationFlags::values);
2212 *   rho_r->gather_evaluate(rho, EvaluationFlags::values);
2213 *   }
2214 *  
2215 *  
2216 * @endcode
2217 *
2218 * Next, we set up the inverse mass operator and the acoustic
2219 * operator. Using `remote_value_type=VectorizedArray<Number>`
2220 * makes the operator use point-to-point interpolation. These two
2221 * objects are then used to create a `RungeKutta2` object to
2222 * perform the time integration.
2223 *
2224
2225 *
2226 * We also compute the maximum speed of sound, needed for the
2227 * computation of the time-step size, and then run the time integrator.
2228 * For the examples considered here, we found a limiting Courant number of
2229 * @f$\mathrm{Cr}\approx 0.36@f$ to maintain stability. To ensure, the
2230 * error of the temporal discretization is small, we use a considerably
2231 * smaller Courant number of @f$0.2@f$.
2232 *
2233 * @code
2234 *   const auto inverse_mass_operator =
2235 *   std::make_shared<InverseMassOperator<dim, Number>>(matrix_free);
2236 *  
2237 *   const auto acoustic_operator =
2238 *   std::make_shared<AcousticOperator<dim, Number, remote_value_type>>(
2239 *   matrix_free,
2240 *   material_data,
2241 *   non_matching_faces,
2242 *   pressure_r,
2243 *   velocity_r,
2244 *   c_r,
2245 *   rho_r);
2246 *  
2247 *   RungeKutta2<dim, Number, remote_value_type> time_integrator(
2248 *   inverse_mass_operator, acoustic_operator);
2249 *  
2250 *   double speed_of_sound_max = 0.0;
2251 *   for (const auto &mat : materials)
2252 *   speed_of_sound_max = std::max(speed_of_sound_max, mat.second.first);
2253 *  
2254 *  
2255 *   time_integrator.run(matrix_free,
2256 *   /*Cr*/ 0.2,
2257 *   end_time,
2258 *   speed_of_sound_max,
2259 *   initial_condition,
2260 *   vtk_prefix);
2261 *   }
2262 *  
2263 * @endcode
2264 *
2265 *
2266 * <a name="step_89-SetupandrunningoftheNitschetypemortaringscheme"></a>
2267 * <h4>Setup and running of the Nitsche-type mortaring scheme</h4>
2268 *
2269
2270 *
2271 * The alternative to the previous function is to use the mortaring
2272 * scheme -- implemented in the following function. This function
2273 * can only be run when deal.II is configured using CGAL (but the
2274 * previous function can still be used without CGAL), and so errors
2275 * out if CGAL is not available.
2276 *
2277 * @code
2278 *   template <int dim, typename Number>
2279 *   void run_with_nitsche_type_mortaring(
2280 *   const MatrixFree<dim, Number> &matrix_free,
2281 *   const std::set<types::boundary_id> &non_matching_faces,
2282 *   const std::map<types::material_id, std::pair<double, double>> &materials,
2283 *   const double end_time,
2284 *   const Function<dim> &initial_condition,
2285 *   const std::string &vtk_prefix)
2286 *   {
2287 *   #ifndef DEAL_II_WITH_CGAL
2288 *   (void)matrix_free;
2289 *   (void)non_matching_faces;
2290 *   (void)materials;
2291 *   (void)end_time;
2292 *   (void)initial_condition;
2293 *   (void)vtk_prefix;
2294 *  
2295 *   ConditionalOStream pcout(
2296 *   std::cout, (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0));
2297 *  
2298 *   pcout << "In this function, mortars are computed using CGAL. "
2299 *   "Configure deal.II with DEAL_II_WITH_CGAL to run this function.\n";
2300 *  
2301 *   return;
2302 *   #else
2303 *  
2304 *   const auto &dof_handler = matrix_free.get_dof_handler();
2305 *   const auto &tria = dof_handler.get_triangulation();
2306 *   const auto &mapping = *matrix_free.get_mapping_info().mapping;
2307 *   const auto n_quadrature_pnts = matrix_free.get_quadrature().size();
2308 *  
2309 *   std::vector<
2310 *   std::pair<types::boundary_id, std::function<std::vector<bool>()>>>
2311 *   non_matching_faces_marked_vertices;
2312 *  
2313 *   for (const auto &nm_face : non_matching_faces)
2314 *   {
2315 *   auto marked_vertices = [&]() {
2316 *   std::vector<bool> mask(tria.n_vertices(), true);
2317 *  
2318 *   for (const auto &cell : tria.active_cell_iterators())
2319 *   for (auto const &f : cell->face_indices())
2320 *   if (cell->face(f)->at_boundary() &&
2321 *   cell->face(f)->boundary_id() == nm_face)
2322 *   for (const auto v : cell->vertex_indices())
2323 *   mask[cell->vertex_index(v)] = false;
2324 *  
2325 *   return mask;
2326 *   };
2327 *  
2328 *   non_matching_faces_marked_vertices.emplace_back(
2329 *   std::make_pair(nm_face, marked_vertices));
2330 *   }
2331 *  
2332 * @endcode
2333 *
2334 * The only parts in this function that are functionally different
2335 * from the previous one follow here.
2336 *
2337
2338 *
2339 * First, quadrature points are arbitrarily distributed on each non-matching
2340 * face. Therefore, we have to make use of FEFacePointEvaluation.
2341 * FEFacePointEvaluation needs NonMatching::MappingInfo to work at the
2342 * correct quadrature points that are in sync with used FERemoteEvaluation
2343 * object. Using
2344 * `compute_remote_communicator_faces_nitsche_type_mortaring()` to reinit
2345 * NonMatching::MappingInfo ensures this. In the case of mortaring, we have
2346 * to use the weights provided by the quadrature rules that are used to set
2347 * up NonMatching::MappingInfo. Therefore we set the flag @c
2348 * use_global_weights.
2349 *
2350 * @code
2351 *   typename NonMatching::MappingInfo<dim, dim, Number>::AdditionalData
2352 *   additional_data;
2353 *   additional_data.use_global_weights = true;
2354 *  
2355 * @endcode
2356 *
2357 * Set up NonMatching::MappingInfo with needed update flags and
2358 * @c additional_data.
2359 *
2360 * @code
2361 *   auto nm_mapping_info =
2362 *   std::make_shared<NonMatching::MappingInfo<dim, dim, Number>>(
2363 *   mapping,
2364 *   update_values | update_JxW_values | update_normal_vectors |
2365 *   update_quadrature_points,
2366 *   additional_data);
2367 *  
2368 *   auto remote_communicator =
2369 *   Utilities::compute_remote_communicator_faces_nitsche_type_mortaring(
2370 *   matrix_free,
2371 *   non_matching_faces_marked_vertices,
2372 *   n_quadrature_pnts,
2373 *   0,
2374 *   nm_mapping_info.get());
2375 *  
2376 * @endcode
2377 *
2378 * Second, since quadrature points are arbitrarily distributed we
2379 * have to consider each face in a batch separately and can not
2380 * make use of @c VecorizedArray.
2381 *
2382 * @code
2383 *   using remote_value_type = Number;
2384 *  
2385 * @endcode
2386 *
2387 * The rest of the code is then identical to what we had in the
2388 * previous function (though it functions differently because of
2389 * the difference in `remote_value_type`).
2390 *
2391 * @code
2392 *   const auto pressure_r =
2393 *   std::make_shared<FERemoteEvaluation<dim, 1, remote_value_type>>(
2394 *   remote_communicator, dof_handler, /*first_selected_component*/ 0);
2395 *  
2396 *   const auto velocity_r =
2397 *   std::make_shared<FERemoteEvaluation<dim, dim, remote_value_type>>(
2398 *   remote_communicator, dof_handler, /*first_selected_component*/ 1);
2399 *  
2400 *   const auto material_data =
2401 *   std::make_shared<CellwiseMaterialData<Number>>(matrix_free, materials);
2402 *  
2403 *   const auto c_r =
2404 *   std::make_shared<FERemoteEvaluation<dim, 1, remote_value_type>>(
2405 *   remote_communicator,
2406 *   matrix_free.get_dof_handler().get_triangulation(),
2407 *   /*first_selected_component*/ 0);
2408 *   const auto rho_r =
2409 *   std::make_shared<FERemoteEvaluation<dim, 1, remote_value_type>>(
2410 *   remote_communicator,
2411 *   matrix_free.get_dof_handler().get_triangulation(),
2412 *   /*first_selected_component*/ 0);
2413 *  
2414 *   if (!material_data->is_homogeneous())
2415 *   {
2416 *   Vector<Number> c(
2417 *   matrix_free.get_dof_handler().get_triangulation().n_active_cells());
2418 *   Vector<Number> rho(
2419 *   matrix_free.get_dof_handler().get_triangulation().n_active_cells());
2420 *  
2421 *   for (const auto &cell : matrix_free.get_dof_handler()
2422 *   .get_triangulation()
2423 *   .active_cell_iterators())
2424 *   {
2425 *   c[cell->active_cell_index()] =
2426 *   materials.at(cell->material_id()).first;
2427 *   rho[cell->active_cell_index()] =
2428 *   materials.at(cell->material_id()).second;
2429 *   }
2430 *  
2431 *   c_r->gather_evaluate(c, EvaluationFlags::values);
2432 *   rho_r->gather_evaluate(rho, EvaluationFlags::values);
2433 *   }
2434 *  
2435 *   const auto inverse_mass_operator =
2436 *   std::make_shared<InverseMassOperator<dim, Number>>(matrix_free);
2437 *  
2438 *   const auto acoustic_operator =
2439 *   std::make_shared<AcousticOperator<dim, Number, remote_value_type>>(
2440 *   matrix_free,
2441 *   material_data,
2442 *   non_matching_faces,
2443 *   pressure_r,
2444 *   velocity_r,
2445 *   c_r,
2446 *   rho_r,
2447 *   nm_mapping_info);
2448 *  
2449 *   RungeKutta2<dim, Number, remote_value_type> time_integrator(
2450 *   inverse_mass_operator, acoustic_operator);
2451 *  
2452 *   double speed_of_sound_max = 0.0;
2453 *   for (const auto &mat : materials)
2454 *   speed_of_sound_max = std::max(speed_of_sound_max, mat.second.first);
2455 *  
2456 *   time_integrator.run(matrix_free,
2457 *   /*Cr*/ 0.2,
2458 *   end_time,
2459 *   speed_of_sound_max,
2460 *   initial_condition,
2461 *   vtk_prefix);
2462 *   #endif
2463 *   }
2464 *   } // namespace Step89
2465 *  
2466 *  
2467 * @endcode
2468 *
2469 *
2470 * <a name="step_89-main"></a>
2471 * <h3>main()</h3>
2472 *
2473
2474 *
2475 * Finally, the `main()` function executes the different versions of handling
2476 * non-matching interfaces.
2477 *
2478
2479 *
2480 * Similar to @ref step_87 "step-87", the minimum requirement of this tutorial is MPI.
2481 * The parallel::distributed::Triangulation class is used if deal.II is
2482 * configured with p4est. Otherwise parallel::shared::Triangulation is used.
2483 *
2484 * @code
2485 *   int main(int argc, char *argv[])
2486 *   {
2487 *   using namespace dealii;
2488 *   constexpr int dim = 2;
2489 *   using Number = double;
2490 *  
2491 *   Utilities::MPI::MPI_InitFinalize mpi(argc, argv);
2492 *   std::cout.precision(5);
2493 *   ConditionalOStream pcout(std::cout,
2494 *   (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) ==
2495 *   0));
2496 *  
2497 *   const unsigned int refinements = 1;
2498 *   const unsigned int degree = 3;
2499 *  
2500 *   #ifdef DEAL_II_WITH_P4EST
2501 *   parallel::distributed::Triangulation<dim> tria(MPI_COMM_WORLD);
2502 *   #else
2503 *   parallel::shared::Triangulation<dim> tria(MPI_COMM_WORLD);
2504 *   #endif
2505 *  
2506 *   pcout << "Create non-matching grid..." << std::endl;
2507 *  
2508 *   std::set<types::boundary_id> non_matching_faces;
2509 *   Step89::build_non_matching_triangulation(tria,
2510 *   non_matching_faces,
2511 *   refinements);
2512 *  
2513 *   pcout << " - Refinement level: " << refinements << std::endl;
2514 *   pcout << " - Number of cells: " << tria.n_cells() << std::endl;
2515 *  
2516 *   pcout << "Create DoFHandler..." << std::endl;
2517 *   DoFHandler<dim> dof_handler(tria);
2518 *   dof_handler.distribute_dofs(FESystem<dim>(FE_DGQ<dim>(degree) ^ (dim + 1)));
2519 *   pcout << " - Number of DoFs: " << dof_handler.n_dofs() << std::endl;
2520 *  
2521 *   AffineConstraints<Number> constraints;
2522 *   constraints.close();
2523 *  
2524 *   pcout << "Set up MatrixFree..." << std::endl;
2525 *   typename MatrixFree<dim, Number>::AdditionalData data;
2526 *   data.mapping_update_flags = update_gradients | update_values;
2527 *   data.mapping_update_flags_inner_faces = update_values;
2528 *   data.mapping_update_flags_boundary_faces =
2529 *   update_quadrature_points | update_values;
2530 *  
2531 *   MatrixFree<dim, Number> matrix_free;
2532 *   matrix_free.reinit(
2533 *   MappingQ1<dim>(), dof_handler, constraints, QGauss<dim>(degree + 1), data);
2534 *  
2535 *  
2536 * @endcode
2537 *
2538 *
2539 * <a name="step_89-RunvibratingmembranetestcaseHomogeneouspressure"></a>
2540 * <h4>Run vibrating membrane test case} Homogeneous pressure</h4>
2541 * Dirichlet boundary conditions are applied for
2542 * simplicity. Therefore, modes can not be chosen arbitrarily.
2543 *
2544 * @code
2545 *   pcout << "Run vibrating membrane test case..." << std::endl;
2546 *   const double modes = 10.0;
2547 *   std::map<types::material_id, std::pair<double, double>> homogeneous_material;
2548 *   homogeneous_material[numbers::invalid_material_id] = std::make_pair(1.0, 1.0);
2549 *   const auto initial_solution_membrane =
2550 *   Step89::InitialConditionVibratingMembrane<dim>(modes);
2551 *  
2552 *   /* Run vibrating membrane test case using point-to-point interpolation: */
2553 *   pcout << " - Point-to-point interpolation: " << std::endl;
2554 *   Step89::run_with_point_to_point_interpolation(
2555 *   matrix_free,
2556 *   non_matching_faces,
2557 *   homogeneous_material,
2558 *   8.0 * initial_solution_membrane.get_period_duration(
2559 *   homogeneous_material.begin()->second.first),
2560 *   initial_solution_membrane,
2561 *   "vm-p2p");
2562 *  
2563 *   /* Run vibrating membrane test case using Nitsche-type mortaring: */
2564 *   pcout << " - Nitsche-type mortaring: " << std::endl;
2565 *   Step89::run_with_nitsche_type_mortaring(
2566 *   matrix_free,
2567 *   non_matching_faces,
2568 *   homogeneous_material,
2569 *   8.0 * initial_solution_membrane.get_period_duration(
2570 *   homogeneous_material.begin()->second.first),
2571 *   initial_solution_membrane,
2572 *   "vm-nitsche");
2573 *  
2574 * @endcode
2575 *
2576 *
2577 * <a name="step_89-Runtestcasewithinhomogeneousmaterial"></a>
2578 * <h4>Run test case with inhomogeneous material</h4>
2579 * Run simple test case with inhomogeneous material and Nitsche-type
2580 * mortaring:
2581 *
2582 * @code
2583 *   pcout << "Run test case with inhomogeneous material..." << std::endl;
2584 *   std::map<types::material_id, std::pair<double, double>>
2585 *   inhomogeneous_material;
2586 *   inhomogeneous_material[0] = std::make_pair(1.0, 1.0);
2587 *   inhomogeneous_material[1] = std::make_pair(3.0, 1.0);
2588 *   Step89::run_with_nitsche_type_mortaring(matrix_free,
2589 *   non_matching_faces,
2590 *   inhomogeneous_material,
2591 *   /*runtime*/ 0.3,
2592 *   Step89::GaussPulse<dim>(0.3, 0.5),
2593 *   "inhomogeneous");
2594 *  
2595 *  
2596 *   return 0;
2597 *   }
2598 * @endcode
2599<a name="step_89-Results"></a><h1>Results</h1>
2600
2601
2602<a name="step_89-VibratingmembranePointtopointinterpolationvsNitschetypemortaring"></a><h3>Vibrating membrane: Point-to-point interpolation vs. Nitsche-type mortaring</h3>
2603
2604
2605We compare the results of the simulations after the last time step, i.e. at @f$t=8T@f$.
2606The @f$y@f$-component of the velocity field using Nitsche-type mortaring is depicted on the left.
2607The same field using point-to-point interpolation is depicted on the right.
2608
2609<table align="center" class="doxtable">
2610 <tr>
2611 <td>
2612 @image html https://www.dealii.org/images/steps/developer/step_89_membrane_test_case_mortaring_velocity_Y.png "" width=60%
2613 </td>
2614 <td>
2615 @image html https://www.dealii.org/images/steps/developer/step_89_membrane_test_case_point_to_point_velocity_Y.png "" width=60%
2616 </td>
2617 </tr>
2618</table>
2619
2620Besides this, the results for the pressure and the velocity in @f$y@f$ direction
2621are plotted along the horizontal line that spans from (0,0.69) to (1,0.69).
2622
2623<table align="center" class="doxtable">
2624 <tr>
2625 <td>
2626 @image html https://www.dealii.org/images/steps/developer/step_89_membrane_test_case_mortaring_vs_point_to_point_pressure.svg "" width=150%
2627 </td>
2628 <td>
2629 @image html https://www.dealii.org/images/steps/developer/step_89_membrane_test_case_mortaring_vs_point_to_point_velocity_Y.svg "" width=150%
2630 </td>
2631 </tr>
2632</table>
2633
2634While the results of the pressure are similar, @f$u_y@f$ clearly differs. At certain
2635positions we can see aliasing errors for the point-to-point interpolation.
2636For different mesh configurations and/or longer run times, the aliasing effects
2637of the point-to-point simulation accumulate and the simulation becomes instable.
2638To keep the tutorial short we have chosen one mesh that can be used for all
2639examples. For a configuration that yields instable results for a wide range of
2640polynomial degrees, see @cite heinz2022high.
2641
2642<a name="step_89-Wavepropagationthroughinhomogeneousfluid"></a><h3>Wave propagation through in-homogeneous fluid</h3>
2643
2644
2645The example that follows is just one example in which non-matching discretizations can be efficiently
2646used to reduce the number of DoFs. The example is nice, since results for a similar
2647test case are shown in multiple publications. As before, we slightly adapted the
2648test case to be able to use the same mesh for all simulations. The pressure field
2649at @f$t=0.3@f$ is depicted below.
2650
2651@image html https://www.dealii.org/images/steps/developer/step_89_inhomogenous_test_case_pressure.png "" width=30%
2652
2653As expected, we can easily see that the wave length in the right domain is roughly
2654three times times the wave length in the left domain. Hence, the wave can be
2655resolved with a coarser discretization.
2656
2657Using the same element size in the whole domain, we can compute a reference result.
2658The displayed reference result is obtained by choosing the same subdivision level
2659for both sub-domains, i.e. @c subdiv_right = 11. In this particular example the
2660reference result uses @f$92,928@f$ DoFs, while the non-matching result uses @f$52,608@f$ DoFs.
2661The pressure result is plotted along the horizontal line that spans from (0,0.5) to (1,0.5).
2662
2663@image html https://www.dealii.org/images/steps/developer/step_89_inhomogenous_test_case_conforming_vs_nonmatching.svg "" width=60%
2664
2665The results computed using the non-matching discretization are clearly in good agreement with
2666the reference result.
2667
2668<a name="step_89-Possibilitiesforextensions"></a><h3>Possibilities for extensions</h3>
2669
2670
2671All the implementations are done with overlapping triangulations in mind. In particular the
2672intersections in the mortaring case are constructed such that they are computed correctly
2673for overlapping triangulations. For this the intersection requests are of dimension @f$dim-1@f$.
2674The cells which are intersected with the intersection requests are of dimension @f$dim@f$. For the
2675simple case of non-conforming interfaces it would be sufficient to compute the intersections
2676between @f$dim-1@f$ and @f$dim-1@f$ entities. Furthermore, the lambda could be adapted, such that cells are
2677only marked if they are connected to a certain boundary ID (in this case, e.g., 99) instead of
2678marking every cell that is <i>not</i> connected to the opposite boundary ID (in this case, e.g., 98).
2679Marking fewer cells can reduce the setup costs significantly.
2680
2681Note that the use of inhomogeneous materials in this procedure is questionable, since it is not clear which
2682material is present in the overlapping part of the mesh.
2683 *
2684 *
2685<a name="step_89-PlainProg"></a>
2686<h1> The plain program</h1>
2687@include "step-89.cc"
2688*/
@ curved_inner_cells
Definition data_out.h:211
T read_cell_data(const AlignedVector< T > &array) const
Abstract base class for mapping classes.
Definition mapping.h:320
Definition point.h:113
const value_type get_value(const unsigned int q) const
const unsigned int DoFAccessor< structdim, dim, spacedim, level_dof_access >::dimension
#define Assert(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
const bool IsBlockVector< VectorType >::value
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
Definition mapping_q1.h:104
const Event initial
Definition event.cc:68
void write_vtu(const std::vector< Patch< dim, spacedim > > &patches, const std::vector< std::string > &data_names, const std::vector< std::tuple< unsigned int, unsigned int, std::string, DataComponentInterpretation::DataComponentInterpretation > > &nonscalar_data_ranges, const VtkFlags &flags, std::ostream &out)
@ matrix
Contents is actually a matrix.
constexpr char A
constexpr types::blas_int one
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition utilities.cc:193
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)
void interpolate(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Function< spacedim, typename VectorType::value_type > &function, VectorType &vec, const ComponentMask &component_mask={}, const unsigned int level=numbers::invalid_unsigned_int)
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)
bool check(const ConstraintKinds kind_in, const unsigned int dim)
int(& functions)(const void *v1, const void *v2)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
constexpr double PI
Definition numbers.h:239
::VectorizedArray< Number, width > exp(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sin(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
unsigned int material_id
Definition types.h:184
unsigned int boundary_id
Definition types.h:161