deal.II version GIT relicensing-2173-gae8fc9d14b 2024-11-24 06:40: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-43.h
Go to the documentation of this file.
1) const
740 *   {
741 *   return 1 - p[0];
742 *   }
743 *  
744 *  
745 *   template <int dim>
746 *   class SaturationBoundaryValues : public Function<dim>
747 *   {
748 *   public:
749 *   SaturationBoundaryValues()
750 *   : Function<dim>(1)
751 *   {}
752 *  
753 *   virtual double value(const Point<dim> &p,
754 *   const unsigned int component = 0) const override;
755 *   };
756 *  
757 *  
758 *  
759 *   template <int dim>
760 *   double
761 *   SaturationBoundaryValues<dim>::value(const Point<dim> &p,
762 *   const unsigned int /*component*/) const
763 *   {
764 *   if (p[0] == 0)
765 *   return 1;
766 *   else
767 *   return 0;
768 *   }
769 *  
770 *  
771 *   template <int dim>
772 *   class SaturationInitialValues : public Function<dim>
773 *   {
774 *   public:
775 *   SaturationInitialValues()
776 *   : Function<dim>(1)
777 *   {}
778 *  
779 *   virtual double value(const Point<dim> &p,
780 *   const unsigned int component = 0) const override;
781 *  
782 *   virtual void vector_value(const Point<dim> &p,
783 *   Vector<double> &value) const override;
784 *   };
785 *  
786 *  
787 *   template <int dim>
788 *   double
789 *   SaturationInitialValues<dim>::value(const Point<dim> & /*p*/,
790 *   const unsigned int /*component*/) const
791 *   {
792 *   return 0.2;
793 *   }
794 *  
795 *  
796 *   template <int dim>
797 *   void SaturationInitialValues<dim>::vector_value(const Point<dim> &p,
798 *   Vector<double> &values) const
799 *   {
800 *   for (unsigned int c = 0; c < this->n_components; ++c)
801 *   values(c) = SaturationInitialValues<dim>::value(p, c);
802 *   }
803 *  
804 *  
805 * @endcode
806 *
807 *
808 * <a name="step_43-Permeabilitymodels"></a>
809 * <h3>Permeability models</h3>
810 *
811
812 *
813 * In this tutorial, we still use the two permeability models previously
814 * used in @ref step_21 "step-21" so we again refrain from commenting in detail about them.
815 *
816 * @code
817 *   namespace SingleCurvingCrack
818 *   {
819 *   template <int dim>
820 *   class KInverse : public TensorFunction<2, dim>
821 *   {
822 *   public:
823 *   KInverse()
825 *   {}
826 *  
827 *   virtual void
828 *   value_list(const std::vector<Point<dim>> &points,
829 *   std::vector<Tensor<2, dim>> &values) const override;
830 *   };
831 *  
832 *  
833 *   template <int dim>
834 *   void KInverse<dim>::value_list(const std::vector<Point<dim>> &points,
835 *   std::vector<Tensor<2, dim>> &values) const
836 *   {
837 *   AssertDimension(points.size(), values.size());
838 *  
839 *   for (unsigned int p = 0; p < points.size(); ++p)
840 *   {
841 *   values[p].clear();
842 *  
843 *   const double distance_to_flowline =
844 *   std::fabs(points[p][1] - 0.5 - 0.1 * std::sin(10 * points[p][0]));
845 *  
846 *   const double permeability =
847 *   std::max(std::exp(-(distance_to_flowline * distance_to_flowline) /
848 *   (0.1 * 0.1)),
849 *   0.01);
850 *  
851 *   for (unsigned int d = 0; d < dim; ++d)
852 *   values[p][d][d] = 1. / permeability;
853 *   }
854 *   }
855 *   } // namespace SingleCurvingCrack
856 *  
857 *  
858 *   namespace RandomMedium
859 *   {
860 *   template <int dim>
861 *   class KInverse : public TensorFunction<2, dim>
862 *   {
863 *   public:
864 *   KInverse()
866 *   {}
867 *  
868 *   virtual void
869 *   value_list(const std::vector<Point<dim>> &points,
870 *   std::vector<Tensor<2, dim>> &values) const override;
871 *  
872 *   private:
873 *   static std::vector<Point<dim>> centers;
874 *   };
875 *  
876 *  
877 *  
878 *   template <int dim>
879 *   std::vector<Point<dim>> KInverse<dim>::centers = []() {
880 *   const unsigned int N =
881 *   (dim == 2 ? 40 : (dim == 3 ? 100 : throw ExcNotImplemented()));
882 *  
883 *   std::vector<Point<dim>> centers_list(N);
884 *   for (unsigned int i = 0; i < N; ++i)
885 *   for (unsigned int d = 0; d < dim; ++d)
886 *   centers_list[i][d] = static_cast<double>(rand()) / RAND_MAX;
887 *  
888 *   return centers_list;
889 *   }();
890 *  
891 *  
892 *  
893 *   template <int dim>
894 *   void KInverse<dim>::value_list(const std::vector<Point<dim>> &points,
895 *   std::vector<Tensor<2, dim>> &values) const
896 *   {
897 *   AssertDimension(points.size(), values.size());
898 *  
899 *   for (unsigned int p = 0; p < points.size(); ++p)
900 *   {
901 *   values[p].clear();
902 *  
903 *   double permeability = 0;
904 *   for (unsigned int i = 0; i < centers.size(); ++i)
905 *   permeability +=
906 *   std::exp(-(points[p] - centers[i]).norm_square() / (0.05 * 0.05));
907 *  
908 *   const double normalized_permeability =
909 *   std::min(std::max(permeability, 0.01), 4.);
910 *  
911 *   for (unsigned int d = 0; d < dim; ++d)
912 *   values[p][d][d] = 1. / normalized_permeability;
913 *   }
914 *   }
915 *   } // namespace RandomMedium
916 *  
917 *  
918 * @endcode
919 *
920 *
921 * <a name="step_43-Physicalquantities"></a>
922 * <h3>Physical quantities</h3>
923 *
924
925 *
926 * The implementations of all the physical quantities such as total mobility
927 * @f$\lambda_t@f$ and fractional flow of water @f$F@f$ are taken from @ref step_21 "step-21" so
928 * again we don't have do any comment about them. Compared to @ref step_21 "step-21" we
929 * have added checks that the saturation passed to these functions is in
930 * fact within the physically valid range. Furthermore, given that the
931 * wetting phase moves at speed @f$\mathbf u F'(S)@f$ it is clear that @f$F'(S)@f$
932 * must be greater or equal to zero, so we assert that as well to make sure
933 * that our calculations to get at the formula for the derivative made
934 * sense.
935 *
936 * @code
937 *   double mobility_inverse(const double S, const double viscosity)
938 *   {
939 *   return 1.0 / (1.0 / viscosity * S * S + (1 - S) * (1 - S));
940 *   }
941 *  
942 *  
943 *   double fractional_flow(const double S, const double viscosity)
944 *   {
945 *   Assert((S >= 0) && (S <= 1),
946 *   ExcMessage("Saturation is outside its physically valid range."));
947 *  
948 *   return S * S / (S * S + viscosity * (1 - S) * (1 - S));
949 *   }
950 *  
951 *  
952 *   double fractional_flow_derivative(const double S, const double viscosity)
953 *   {
954 *   Assert((S >= 0) && (S <= 1),
955 *   ExcMessage("Saturation is outside its physically valid range."));
956 *  
957 *   const double temp = (S * S + viscosity * (1 - S) * (1 - S));
958 *  
959 *   const double numerator =
960 *   2.0 * S * temp - S * S * (2.0 * S - 2.0 * viscosity * (1 - S));
961 *   const double denominator = Utilities::fixed_power<2>(temp);
962 *  
963 *   const double F_prime = numerator / denominator;
964 *  
965 *   Assert(F_prime >= 0, ExcInternalError());
966 *  
967 *   return F_prime;
968 *   }
969 *  
970 *  
971 * @endcode
972 *
973 *
974 * <a name="step_43-Helperclassesforsolversandpreconditioners"></a>
975 * <h3>Helper classes for solvers and preconditioners</h3>
976 *
977
978 *
979 * In this first part we define a number of classes that we need in the
980 * construction of linear solvers and preconditioners. This part is
981 * essentially the same as that used in @ref step_31 "step-31". The only difference is that
982 * the original variable name stokes_matrix is replaced by another name
983 * darcy_matrix to match our problem.
984 *
985 * @code
986 *   namespace LinearSolvers
987 *   {
988 *   template <class MatrixType, class PreconditionerType>
989 *   class InverseMatrix : public EnableObserverPointer
990 *   {
991 *   public:
992 *   InverseMatrix(const MatrixType &m,
993 *   const PreconditionerType &preconditioner);
994 *  
995 *  
996 *   template <typename VectorType>
997 *   void vmult(VectorType &dst, const VectorType &src) const;
998 *  
999 *   private:
1000 *   const ObserverPointer<const MatrixType> matrix;
1001 *   const PreconditionerType &preconditioner;
1002 *   };
1003 *  
1004 *  
1005 *   template <class MatrixType, class PreconditionerType>
1006 *   InverseMatrix<MatrixType, PreconditionerType>::InverseMatrix(
1007 *   const MatrixType &m,
1008 *   const PreconditionerType &preconditioner)
1009 *   : matrix(&m)
1010 *   , preconditioner(preconditioner)
1011 *   {}
1012 *  
1013 *  
1014 *  
1015 *   template <class MatrixType, class PreconditionerType>
1016 *   template <typename VectorType>
1017 *   void InverseMatrix<MatrixType, PreconditionerType>::vmult(
1018 *   VectorType &dst,
1019 *   const VectorType &src) const
1020 *   {
1021 *   SolverControl solver_control(src.size(), 1e-7 * src.l2_norm());
1022 *   SolverCG<VectorType> cg(solver_control);
1023 *  
1024 *   dst = 0;
1025 *  
1026 *   try
1027 *   {
1028 *   cg.solve(*matrix, dst, src, preconditioner);
1029 *   }
1030 *   catch (std::exception &e)
1031 *   {
1032 *   Assert(false, ExcMessage(e.what()));
1033 *   }
1034 *   }
1035 *  
1036 *   template <class PreconditionerTypeA, class PreconditionerTypeMp>
1037 *   class BlockSchurPreconditioner : public EnableObserverPointer
1038 *   {
1039 *   public:
1040 *   BlockSchurPreconditioner(
1041 *   const TrilinosWrappers::BlockSparseMatrix &S,
1042 *   const InverseMatrix<TrilinosWrappers::SparseMatrix,
1043 *   PreconditionerTypeMp> &Mpinv,
1044 *   const PreconditionerTypeA &Apreconditioner);
1045 *  
1046 *   void vmult(TrilinosWrappers::MPI::BlockVector &dst,
1047 *   const TrilinosWrappers::MPI::BlockVector &src) const;
1048 *  
1049 *   private:
1050 *   const ObserverPointer<const TrilinosWrappers::BlockSparseMatrix>
1051 *   darcy_matrix;
1052 *   const ObserverPointer<const InverseMatrix<TrilinosWrappers::SparseMatrix,
1053 *   PreconditionerTypeMp>>
1054 *   m_inverse;
1055 *   const PreconditionerTypeA &a_preconditioner;
1056 *  
1057 *   mutable TrilinosWrappers::MPI::Vector tmp;
1058 *   };
1059 *  
1060 *  
1061 *  
1062 *   template <class PreconditionerTypeA, class PreconditionerTypeMp>
1063 *   BlockSchurPreconditioner<PreconditionerTypeA, PreconditionerTypeMp>::
1064 *   BlockSchurPreconditioner(
1065 *   const TrilinosWrappers::BlockSparseMatrix &S,
1066 *   const InverseMatrix<TrilinosWrappers::SparseMatrix,
1067 *   PreconditionerTypeMp> &Mpinv,
1068 *   const PreconditionerTypeA &Apreconditioner)
1069 *   : darcy_matrix(&S)
1070 *   , m_inverse(&Mpinv)
1071 *   , a_preconditioner(Apreconditioner)
1072 *   , tmp(complete_index_set(darcy_matrix->block(1, 1).m()))
1073 *   {}
1074 *  
1075 *  
1076 *   template <class PreconditionerTypeA, class PreconditionerTypeMp>
1077 *   void
1078 *   BlockSchurPreconditioner<PreconditionerTypeA, PreconditionerTypeMp>::vmult(
1079 *   TrilinosWrappers::MPI::BlockVector &dst,
1080 *   const TrilinosWrappers::MPI::BlockVector &src) const
1081 *   {
1082 *   a_preconditioner.vmult(dst.block(0), src.block(0));
1083 *   darcy_matrix->block(1, 0).residual(tmp, dst.block(0), src.block(1));
1084 *   tmp *= -1;
1085 *   m_inverse->vmult(dst.block(1), tmp);
1086 *   }
1087 *   } // namespace LinearSolvers
1088 *  
1089 *  
1090 * @endcode
1091 *
1092 *
1093 * <a name="step_43-TheTwoPhaseFlowProblemclass"></a>
1094 * <h3>The TwoPhaseFlowProblem class</h3>
1095 *
1096
1097 *
1098 * The definition of the class that defines the top-level logic of solving
1099 * the time-dependent advection-dominated two-phase flow problem (or
1100 * Buckley-Leverett problem @cite Buckley1942) is mainly based on tutorial
1101 * programs @ref step_21 "step-21" and @ref step_33 "step-33", and in particular on @ref step_31 "step-31" where we have
1102 * used basically the same general structure as done here. As in @ref step_31 "step-31",
1103 * the key routines to look for in the implementation below are the
1104 * <code>run()</code> and <code>solve()</code> functions.
1105 *
1106
1107 *
1108 * The main difference to @ref step_31 "step-31" is that, since adaptive operator splitting
1109 * is considered, we need a couple more member variables to hold the last
1110 * two computed Darcy (velocity/pressure) solutions in addition to the
1111 * current one (which is either computed directly, or extrapolated from the
1112 * previous two), and we need to remember the last two times we computed the
1113 * Darcy solution. We also need a helper function that figures out whether
1114 * we do indeed need to recompute the Darcy solution.
1115 *
1116
1117 *
1118 * Unlike @ref step_31 "step-31", this step uses one more AffineConstraints object called
1119 * darcy_preconditioner_constraints. This constraint object is used only for
1120 * assembling the matrix for the Darcy preconditioner and includes hanging
1121 * node constraints as well as Dirichlet boundary value constraints for the
1122 * pressure variable. We need this because we are building a Laplace matrix
1123 * for the pressure as an approximation of the Schur complement) which is
1124 * only positive definite if boundary conditions are applied.
1125 *
1126
1127 *
1128 * The collection of member functions and variables thus declared in this
1129 * class is then rather similar to those in @ref step_31 "step-31":
1130 *
1131 * @code
1132 *   template <int dim>
1133 *   class TwoPhaseFlowProblem
1134 *   {
1135 *   public:
1136 *   TwoPhaseFlowProblem(const unsigned int degree);
1137 *   void run();
1138 *  
1139 *   private:
1140 *   void setup_dofs();
1141 *   void assemble_darcy_preconditioner();
1142 *   void build_darcy_preconditioner();
1143 *   void assemble_darcy_system();
1144 *   void assemble_saturation_system();
1145 *   void assemble_saturation_matrix();
1146 *   void assemble_saturation_rhs();
1147 *   void assemble_saturation_rhs_cell_term(
1148 *   const FEValues<dim> &saturation_fe_values,
1149 *   const FEValues<dim> &darcy_fe_values,
1150 *   const double global_max_u_F_prime,
1151 *   const double global_S_variation,
1152 *   const std::vector<types::global_dof_index> &local_dof_indices);
1153 *   void assemble_saturation_rhs_boundary_term(
1154 *   const FEFaceValues<dim> &saturation_fe_face_values,
1155 *   const FEFaceValues<dim> &darcy_fe_face_values,
1156 *   const std::vector<types::global_dof_index> &local_dof_indices);
1157 *   void solve();
1158 *   void refine_mesh(const unsigned int min_grid_level,
1159 *   const unsigned int max_grid_level);
1160 *   void output_results() const;
1161 *  
1162 * @endcode
1163 *
1164 * We follow with a number of helper functions that are used in a variety
1165 * of places throughout the program:
1166 *
1167 * @code
1168 *   double get_max_u_F_prime() const;
1169 *   std::pair<double, double> get_extrapolated_saturation_range() const;
1170 *   bool determine_whether_to_solve_for_pressure_and_velocity() const;
1171 *   void project_back_saturation();
1172 *   double compute_viscosity(
1173 *   const std::vector<double> &old_saturation,
1174 *   const std::vector<double> &old_old_saturation,
1175 *   const std::vector<Tensor<1, dim>> &old_saturation_grads,
1176 *   const std::vector<Tensor<1, dim>> &old_old_saturation_grads,
1177 *   const std::vector<Vector<double>> &present_darcy_values,
1178 *   const double global_max_u_F_prime,
1179 *   const double global_S_variation,
1180 *   const double cell_diameter) const;
1181 *  
1182 *  
1183 * @endcode
1184 *
1185 * This all is followed by the member variables, most of which are similar
1186 * to the ones in @ref step_31 "step-31", with the exception of the ones that pertain to
1187 * the macro time stepping for the velocity/pressure system:
1188 *
1189 * @code
1190 *   Triangulation<dim> triangulation;
1191 *   double global_Omega_diameter;
1192 *  
1193 *   const unsigned int degree;
1194 *  
1195 *   const unsigned int darcy_degree;
1196 *   const FESystem<dim> darcy_fe;
1197 *   DoFHandler<dim> darcy_dof_handler;
1198 *   AffineConstraints<double> darcy_constraints;
1199 *  
1200 *   AffineConstraints<double> darcy_preconditioner_constraints;
1201 *  
1202 *   TrilinosWrappers::BlockSparseMatrix darcy_matrix;
1203 *   TrilinosWrappers::BlockSparseMatrix darcy_preconditioner_matrix;
1204 *  
1205 *   TrilinosWrappers::MPI::BlockVector darcy_solution;
1206 *   TrilinosWrappers::MPI::BlockVector darcy_rhs;
1207 *  
1208 *   TrilinosWrappers::MPI::BlockVector last_computed_darcy_solution;
1209 *   TrilinosWrappers::MPI::BlockVector second_last_computed_darcy_solution;
1210 *  
1211 *  
1212 *   const unsigned int saturation_degree;
1213 *   const FE_Q<dim> saturation_fe;
1214 *   DoFHandler<dim> saturation_dof_handler;
1215 *   AffineConstraints<double> saturation_constraints;
1216 *  
1217 *   TrilinosWrappers::SparseMatrix saturation_matrix;
1218 *  
1219 *  
1220 *   TrilinosWrappers::MPI::Vector saturation_solution;
1221 *   TrilinosWrappers::MPI::Vector old_saturation_solution;
1222 *   TrilinosWrappers::MPI::Vector old_old_saturation_solution;
1223 *   TrilinosWrappers::MPI::Vector saturation_rhs;
1224 *  
1225 *   TrilinosWrappers::MPI::Vector
1226 *   saturation_matching_last_computed_darcy_solution;
1227 *  
1228 *   const double saturation_refinement_threshold;
1229 *  
1230 *   double time;
1231 *   const double end_time;
1232 *  
1233 *   double current_macro_time_step;
1234 *   double old_macro_time_step;
1235 *  
1236 *   double time_step;
1237 *   double old_time_step;
1238 *   unsigned int timestep_number;
1239 *  
1240 *   const double viscosity;
1241 *   const double porosity;
1242 *   const double AOS_threshold;
1243 *  
1244 *   std::shared_ptr<TrilinosWrappers::PreconditionIC> top_left_preconditioner;
1245 *   std::shared_ptr<TrilinosWrappers::PreconditionIC>
1246 *   bottom_right_preconditioner;
1247 *  
1248 *   bool rebuild_saturation_matrix;
1249 *  
1250 * @endcode
1251 *
1252 * At the very end we declare a variable that denotes the material
1253 * model. Compared to @ref step_21 "step-21", we do this here as a member variable since
1254 * we will want to use it in a variety of places and so having a central
1255 * place where such a variable is declared will make it simpler to replace
1256 * one class by another (e.g. replace RandomMedium::KInverse by
1257 * SingleCurvingCrack::KInverse).
1258 *
1259 * @code
1260 *   const RandomMedium::KInverse<dim> k_inverse;
1261 *   };
1262 *  
1263 *  
1264 * @endcode
1265 *
1266 *
1267 * <a name="step_43-TwoPhaseFlowProblemdimTwoPhaseFlowProblem"></a>
1268 * <h3>TwoPhaseFlowProblem<dim>::TwoPhaseFlowProblem</h3>
1269 *
1270
1271 *
1272 * The constructor of this class is an extension of the constructors in
1273 * @ref step_21 "step-21" and @ref step_31 "step-31". We need to add the various variables that concern
1274 * the saturation. As discussed in the introduction, we are going to use
1275 * @f$Q_2 \times Q_1@f$ (Taylor-Hood) elements again for the Darcy system, an
1276 * element combination that fulfills the Ladyzhenskaya-Babuska-Brezzi (LBB)
1277 * conditions [Brezzi and Fortin 1991, Chen 2005], and @f$Q_1@f$ elements for
1278 * the saturation. However, by using variables that store the polynomial
1279 * degree of the Darcy and temperature finite elements, it is easy to
1280 * consistently modify the degree of the elements as well as all quadrature
1281 * formulas used on them downstream. Moreover, we initialize the time
1282 * stepping variables related to operator splitting as well as the option
1283 * for matrix assembly and preconditioning:
1284 *
1285 * @code
1286 *   template <int dim>
1287 *   TwoPhaseFlowProblem<dim>::TwoPhaseFlowProblem(const unsigned int degree)
1288 *   : triangulation(Triangulation<dim>::maximum_smoothing)
1289 *   , global_Omega_diameter(std::numeric_limits<double>::quiet_NaN())
1290 *   , degree(degree)
1291 *   , darcy_degree(degree)
1292 *   , darcy_fe(FE_Q<dim>(darcy_degree + 1) ^ dim, FE_Q<dim>(darcy_degree))
1293 *   , darcy_dof_handler(triangulation)
1294 *   ,
1295 *  
1296 *   saturation_degree(degree + 1)
1297 *   , saturation_fe(saturation_degree)
1298 *   , saturation_dof_handler(triangulation)
1299 *   ,
1300 *  
1301 *   saturation_refinement_threshold(0.5)
1302 *   ,
1303 *  
1304 *   time(0)
1305 *   , end_time(10)
1306 *   ,
1307 *  
1308 *   current_macro_time_step(0)
1309 *   , old_macro_time_step(0)
1310 *   ,
1311 *  
1312 *   time_step(0)
1313 *   , old_time_step(0)
1314 *   , timestep_number(0)
1315 *   , viscosity(0.2)
1316 *   , porosity(1.0)
1317 *   , AOS_threshold(3.0)
1318 *   ,
1319 *  
1320 *   rebuild_saturation_matrix(true)
1321 *   {}
1322 *  
1323 *  
1324 * @endcode
1325 *
1326 *
1327 * <a name="step_43-TwoPhaseFlowProblemdimsetup_dofs"></a>
1328 * <h3>TwoPhaseFlowProblem<dim>::setup_dofs</h3>
1329 *
1330
1331 *
1332 * This is the function that sets up the DoFHandler objects we have here
1333 * (one for the Darcy part and one for the saturation part) as well as set
1334 * to the right sizes the various objects required for the linear algebra in
1335 * this program. Its basic operations are similar to what @ref step_31 "step-31" did.
1336 *
1337
1338 *
1339 * The body of the function first enumerates all degrees of freedom for the
1340 * Darcy and saturation systems. For the Darcy part, degrees of freedom are
1341 * then sorted to ensure that velocities precede pressure DoFs so that we
1342 * can partition the Darcy matrix into a @f$2 \times 2@f$ matrix.
1343 *
1344
1345 *
1346 * Then, we need to incorporate hanging node constraints and Dirichlet
1347 * boundary value constraints into darcy_preconditioner_constraints. The
1348 * boundary condition constraints are only set on the pressure component
1349 * since the Schur complement preconditioner that corresponds to the porous
1350 * media flow operator in non-mixed form, @f$-\nabla \cdot [\mathbf K
1351 * \lambda_t(S)]\nabla@f$, acts only on the pressure variable. Therefore, we
1352 * use a component_mask that filters out the velocity component, so that the
1353 * condensation is performed on pressure degrees of freedom only.
1354 *
1355
1356 *
1357 * After having done so, we count the number of degrees of freedom in the
1358 * various blocks. This information is then used to create the sparsity
1359 * pattern for the Darcy and saturation system matrices as well as the
1360 * preconditioner matrix from which we build the Darcy preconditioner. As in
1361 * @ref step_31 "step-31", we choose to create the pattern using the blocked version of
1362 * DynamicSparsityPattern. So, for this, we follow the same way as @ref step_31 "step-31"
1363 * did and we don't have to repeat descriptions again for the rest of the
1364 * member function.
1365 *
1366 * @code
1367 *   template <int dim>
1368 *   void TwoPhaseFlowProblem<dim>::setup_dofs()
1369 *   {
1370 *   std::vector<unsigned int> darcy_block_component(dim + 1, 0);
1371 *   darcy_block_component[dim] = 1;
1372 *   {
1373 *   darcy_dof_handler.distribute_dofs(darcy_fe);
1374 *   DoFRenumbering::Cuthill_McKee(darcy_dof_handler);
1375 *   DoFRenumbering::component_wise(darcy_dof_handler, darcy_block_component);
1376 *  
1377 *   darcy_constraints.clear();
1378 *   DoFTools::make_hanging_node_constraints(darcy_dof_handler,
1379 *   darcy_constraints);
1380 *   darcy_constraints.close();
1381 *   }
1382 *   {
1383 *   saturation_dof_handler.distribute_dofs(saturation_fe);
1384 *  
1385 *   saturation_constraints.clear();
1386 *   DoFTools::make_hanging_node_constraints(saturation_dof_handler,
1387 *   saturation_constraints);
1388 *   saturation_constraints.close();
1389 *   }
1390 *   {
1391 *   darcy_preconditioner_constraints.clear();
1392 *  
1393 *   const FEValuesExtractors::Scalar pressure(dim);
1394 *  
1395 *   DoFTools::make_hanging_node_constraints(darcy_dof_handler,
1396 *   darcy_preconditioner_constraints);
1397 *   DoFTools::make_zero_boundary_constraints(darcy_dof_handler,
1398 *   darcy_preconditioner_constraints,
1399 *   darcy_fe.component_mask(
1400 *   pressure));
1401 *  
1402 *   darcy_preconditioner_constraints.close();
1403 *   }
1404 *  
1405 *  
1406 *   const std::vector<types::global_dof_index> darcy_dofs_per_block =
1407 *   DoFTools::count_dofs_per_fe_block(darcy_dof_handler,
1408 *   darcy_block_component);
1409 *   const types::global_dof_index n_u = darcy_dofs_per_block[0],
1410 *   n_p = darcy_dofs_per_block[1],
1411 *   n_s = saturation_dof_handler.n_dofs();
1412 *  
1413 *   std::cout << "Number of active cells: " << triangulation.n_active_cells()
1414 *   << " (on " << triangulation.n_levels() << " levels)" << std::endl
1415 *   << "Number of degrees of freedom: " << n_u + n_p + n_s << " ("
1416 *   << n_u << '+' << n_p << '+' << n_s << ')' << std::endl
1417 *   << std::endl;
1418 *  
1419 *   {
1420 *   darcy_matrix.clear();
1421 *  
1422 *   BlockDynamicSparsityPattern dsp(darcy_dofs_per_block,
1423 *   darcy_dofs_per_block);
1424 *  
1425 *   Table<2, DoFTools::Coupling> coupling(dim + 1, dim + 1);
1426 *   for (unsigned int c = 0; c < dim + 1; ++c)
1427 *   for (unsigned int d = 0; d < dim + 1; ++d)
1428 *   if (!((c == dim) && (d == dim)))
1429 *   coupling[c][d] = DoFTools::always;
1430 *   else
1431 *   coupling[c][d] = DoFTools::none;
1432 *  
1433 *  
1435 *   darcy_dof_handler, coupling, dsp, darcy_constraints, false);
1436 *  
1437 *   darcy_matrix.reinit(dsp);
1438 *   }
1439 *  
1440 *   {
1441 *   top_left_preconditioner.reset();
1442 *   bottom_right_preconditioner.reset();
1443 *   darcy_preconditioner_matrix.clear();
1444 *  
1445 *   BlockDynamicSparsityPattern dsp(darcy_dofs_per_block,
1446 *   darcy_dofs_per_block);
1447 *  
1448 *   Table<2, DoFTools::Coupling> coupling(dim + 1, dim + 1);
1449 *   for (unsigned int c = 0; c < dim + 1; ++c)
1450 *   for (unsigned int d = 0; d < dim + 1; ++d)
1451 *   if (c == d)
1452 *   coupling[c][d] = DoFTools::always;
1453 *   else
1454 *   coupling[c][d] = DoFTools::none;
1455 *  
1457 *   darcy_dof_handler, coupling, dsp, darcy_constraints, false);
1458 *  
1459 *   darcy_preconditioner_matrix.reinit(dsp);
1460 *   }
1461 *  
1462 *  
1463 *   {
1464 *   saturation_matrix.clear();
1465 *  
1466 *   DynamicSparsityPattern dsp(n_s, n_s);
1467 *  
1468 *   DoFTools::make_sparsity_pattern(saturation_dof_handler,
1469 *   dsp,
1470 *   saturation_constraints,
1471 *   false);
1472 *  
1473 *  
1474 *   saturation_matrix.reinit(dsp);
1475 *   }
1476 *  
1477 *   const std::vector<IndexSet> darcy_partitioning = {complete_index_set(n_u),
1478 *   complete_index_set(n_p)};
1479 *  
1480 *   darcy_solution.reinit(darcy_partitioning, MPI_COMM_WORLD);
1481 *  
1482 *   last_computed_darcy_solution.reinit(darcy_partitioning, MPI_COMM_WORLD);
1483 *  
1484 *   second_last_computed_darcy_solution.reinit(darcy_partitioning,
1485 *   MPI_COMM_WORLD);
1486 *  
1487 *   darcy_rhs.reinit(darcy_partitioning, MPI_COMM_WORLD);
1488 *  
1489 *   const IndexSet saturation_partitioning = complete_index_set(n_s);
1490 *   saturation_solution.reinit(saturation_partitioning, MPI_COMM_WORLD);
1491 *   old_saturation_solution.reinit(saturation_partitioning, MPI_COMM_WORLD);
1492 *   old_old_saturation_solution.reinit(saturation_partitioning, MPI_COMM_WORLD);
1493 *  
1494 *   saturation_matching_last_computed_darcy_solution.reinit(
1495 *   saturation_partitioning, MPI_COMM_WORLD);
1496 *  
1497 *   saturation_rhs.reinit(saturation_partitioning, MPI_COMM_WORLD);
1498 *   }
1499 *  
1500 *  
1501 * @endcode
1502 *
1503 *
1504 * <a name="step_43-Assemblingmatricesandpreconditioners"></a>
1505 * <h3>Assembling matrices and preconditioners</h3>
1506 *
1507
1508 *
1509 * The next few functions are devoted to setting up the various system and
1510 * preconditioner matrices and right hand sides that we have to deal with in
1511 * this program.
1512 *
1513
1514 *
1515 *
1516 * <a name="step_43-TwoPhaseFlowProblemdimassemble_darcy_preconditioner"></a>
1517 * <h4>TwoPhaseFlowProblem<dim>::assemble_darcy_preconditioner</h4>
1518 *
1519
1520 *
1521 * This function assembles the matrix we use for preconditioning the Darcy
1522 * system. What we need are a vector @ref GlossMassMatrix "mass matrix" weighted by
1523 * @f$\left(\mathbf{K} \lambda_t\right)^{-1}@f$ on the velocity components and a
1524 * mass matrix weighted by @f$\left(\mathbf{K} \lambda_t\right)@f$ on the
1525 * pressure component. We start by generating a quadrature object of
1526 * appropriate order, the FEValues object that can give values and gradients
1527 * at the quadrature points (together with quadrature weights). Next we
1528 * create data structures for the cell matrix and the relation between local
1529 * and global DoFs. The vectors phi_u and grad_phi_p are going to hold the
1530 * values of the basis functions in order to faster build up the local
1531 * matrices, as was already done in @ref step_22 "step-22". Before we start the loop over
1532 * all active cells, we have to specify which components are pressure and
1533 * which are velocity.
1534 *
1535
1536 *
1537 * The creation of the local matrix is rather simple. There are only a term
1538 * weighted by @f$\left(\mathbf{K} \lambda_t\right)^{-1}@f$ (on the velocity)
1539 * and a Laplace matrix weighted by @f$\left(\mathbf{K} \lambda_t\right)@f$ to
1540 * be generated, so the creation of the local matrix is done in essentially
1541 * two lines. Since the material model functions at the top of this file
1542 * only provide the inverses of the permeability and mobility, we have to
1543 * compute @f$\mathbf K@f$ and @f$\lambda_t@f$ by hand from the given values, once
1544 * per quadrature point.
1545 *
1546
1547 *
1548 * Once the local matrix is ready (loop over rows and columns in the local
1549 * matrix on each quadrature point), we get the local DoF indices and write
1550 * the local information into the global matrix. We do this by directly
1551 * applying the constraints (i.e. darcy_preconditioner_constraints) that
1552 * takes care of hanging node and zero Dirichlet boundary condition
1553 * constraints. By doing so, we don't have to do that afterwards, and we
1554 * later don't have to use AffineConstraints::condense and
1555 * MatrixTools::apply_boundary_values, both functions that would need to
1556 * modify matrix and vector entries and so are difficult to write for the
1557 * Trilinos classes where we don't immediately have access to individual
1558 * memory locations.
1559 *
1560 * @code
1561 *   template <int dim>
1562 *   void TwoPhaseFlowProblem<dim>::assemble_darcy_preconditioner()
1563 *   {
1564 *   std::cout << " Rebuilding darcy preconditioner..." << std::endl;
1565 *  
1566 *   darcy_preconditioner_matrix = 0;
1567 *  
1568 *   const QGauss<dim> quadrature_formula(darcy_degree + 2);
1569 *   FEValues<dim> darcy_fe_values(darcy_fe,
1570 *   quadrature_formula,
1571 *   update_JxW_values | update_values |
1572 *   update_gradients |
1573 *   update_quadrature_points);
1574 *   FEValues<dim> saturation_fe_values(saturation_fe,
1575 *   quadrature_formula,
1576 *   update_values);
1577 *  
1578 *   const unsigned int dofs_per_cell = darcy_fe.n_dofs_per_cell();
1579 *   const unsigned int n_q_points = quadrature_formula.size();
1580 *  
1581 *   std::vector<Tensor<2, dim>> k_inverse_values(n_q_points);
1582 *  
1583 *   std::vector<double> old_saturation_values(n_q_points);
1584 *  
1585 *   FullMatrix<double> local_matrix(dofs_per_cell, dofs_per_cell);
1586 *   std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
1587 *  
1588 *   std::vector<Tensor<1, dim>> phi_u(dofs_per_cell);
1589 *   std::vector<Tensor<1, dim>> grad_phi_p(dofs_per_cell);
1590 *  
1591 *   const FEValuesExtractors::Vector velocities(0);
1592 *   const FEValuesExtractors::Scalar pressure(dim);
1593 *  
1594 *   auto cell = darcy_dof_handler.begin_active();
1595 *   const auto endc = darcy_dof_handler.end();
1596 *   auto saturation_cell = saturation_dof_handler.begin_active();
1597 *  
1598 *   for (; cell != endc; ++cell, ++saturation_cell)
1599 *   {
1600 *   darcy_fe_values.reinit(cell);
1601 *   saturation_fe_values.reinit(saturation_cell);
1602 *  
1603 *   local_matrix = 0;
1604 *  
1605 *   saturation_fe_values.get_function_values(old_saturation_solution,
1606 *   old_saturation_values);
1607 *  
1608 *   k_inverse.value_list(darcy_fe_values.get_quadrature_points(),
1609 *   k_inverse_values);
1610 *  
1611 *   for (unsigned int q = 0; q < n_q_points; ++q)
1612 *   {
1613 *   const double old_s = old_saturation_values[q];
1614 *  
1615 *   const double inverse_mobility = mobility_inverse(old_s, viscosity);
1616 *   const double mobility = 1.0 / inverse_mobility;
1617 *   const Tensor<2, dim> permeability = invert(k_inverse_values[q]);
1618 *  
1619 *   for (unsigned int k = 0; k < dofs_per_cell; ++k)
1620 *   {
1621 *   phi_u[k] = darcy_fe_values[velocities].value(k, q);
1622 *   grad_phi_p[k] = darcy_fe_values[pressure].gradient(k, q);
1623 *   }
1624 *  
1625 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
1626 *   for (unsigned int j = 0; j < dofs_per_cell; ++j)
1627 *   {
1628 *   local_matrix(i, j) +=
1629 *   (k_inverse_values[q] * inverse_mobility * phi_u[i] *
1630 *   phi_u[j] +
1631 *   permeability * mobility * grad_phi_p[i] * grad_phi_p[j]) *
1632 *   darcy_fe_values.JxW(q);
1633 *   }
1634 *   }
1635 *  
1636 *   cell->get_dof_indices(local_dof_indices);
1637 *   darcy_preconditioner_constraints.distribute_local_to_global(
1638 *   local_matrix, local_dof_indices, darcy_preconditioner_matrix);
1639 *   }
1640 *   }
1641 *  
1642 *  
1643 * @endcode
1644 *
1645 *
1646 * <a name="step_43-TwoPhaseFlowProblemdimbuild_darcy_preconditioner"></a>
1647 * <h4>TwoPhaseFlowProblem<dim>::build_darcy_preconditioner</h4>
1648 *
1649
1650 *
1651 * After calling the above functions to assemble the preconditioner matrix,
1652 * this function generates the inner preconditioners that are going to be
1653 * used for the Schur complement block preconditioner. The preconditioners
1654 * need to be regenerated at every saturation time step since they depend on
1655 * the saturation @f$S@f$ that varies with time.
1656 *
1657
1658 *
1659 * In here, we set up the preconditioner for the velocity-velocity matrix
1660 * @f$\mathbf{M}^{\mathbf{u}}@f$ and the Schur complement @f$\mathbf{S}@f$. As
1661 * explained in the introduction, we are going to use an IC preconditioner
1662 * based on the vector matrix @f$\mathbf{M}^{\mathbf{u}}@f$ and another based on
1663 * the scalar Laplace matrix @f$\tilde{\mathbf{S}}^p@f$ (which is spectrally
1664 * close to the Schur complement of the Darcy matrix). Usually, the
1665 * TrilinosWrappers::PreconditionIC class can be seen as a good black-box
1666 * preconditioner which does not need any special knowledge of the matrix
1667 * structure and/or the operator that's behind it.
1668 *
1669 * @code
1670 *   template <int dim>
1671 *   void TwoPhaseFlowProblem<dim>::build_darcy_preconditioner()
1672 *   {
1673 *   assemble_darcy_preconditioner();
1674 *  
1675 *   top_left_preconditioner =
1676 *   std::make_shared<TrilinosWrappers::PreconditionIC>();
1677 *   top_left_preconditioner->initialize(
1678 *   darcy_preconditioner_matrix.block(0, 0));
1679 *  
1680 *   bottom_right_preconditioner =
1681 *   std::make_shared<TrilinosWrappers::PreconditionIC>();
1682 *   bottom_right_preconditioner->initialize(
1683 *   darcy_preconditioner_matrix.block(1, 1));
1684 *   }
1685 *  
1686 *  
1687 * @endcode
1688 *
1689 *
1690 * <a name="step_43-TwoPhaseFlowProblemdimassemble_darcy_system"></a>
1691 * <h4>TwoPhaseFlowProblem<dim>::assemble_darcy_system</h4>
1692 *
1693
1694 *
1695 * This is the function that assembles the linear system for the Darcy
1696 * system.
1697 *
1698
1699 *
1700 * Regarding the technical details of implementation, the procedures are
1701 * similar to those in @ref step_22 "step-22" and @ref step_31 "step-31". We reset matrix and vector,
1702 * create a quadrature formula on the cells, and then create the respective
1703 * FEValues object.
1704 *
1705
1706 *
1707 * There is one thing that needs to be commented: since we have a separate
1708 * finite element and DoFHandler for the saturation, we need to generate a
1709 * second FEValues object for the proper evaluation of the saturation
1710 * solution. This isn't too complicated to realize here: just use the
1711 * saturation structures and set an update flag for the basis function
1712 * values which we need for evaluation of the saturation solution. The only
1713 * important part to remember here is that the same quadrature formula is
1714 * used for both FEValues objects to ensure that we get matching information
1715 * when we loop over the quadrature points of the two objects.
1716 *
1717
1718 *
1719 * The declarations proceed with some shortcuts for array sizes, the
1720 * creation of the local matrix, right hand side as well as the vector for
1721 * the indices of the local dofs compared to the global system.
1722 *
1723 * @code
1724 *   template <int dim>
1725 *   void TwoPhaseFlowProblem<dim>::assemble_darcy_system()
1726 *   {
1727 *   darcy_matrix = 0;
1728 *   darcy_rhs = 0;
1729 *  
1730 *   const QGauss<dim> quadrature_formula(darcy_degree + 2);
1731 *   const QGauss<dim - 1> face_quadrature_formula(darcy_degree + 2);
1732 *  
1733 *   FEValues<dim> darcy_fe_values(darcy_fe,
1734 *   quadrature_formula,
1735 *   update_values | update_gradients |
1736 *   update_quadrature_points |
1737 *   update_JxW_values);
1738 *  
1739 *   FEValues<dim> saturation_fe_values(saturation_fe,
1740 *   quadrature_formula,
1741 *   update_values);
1742 *  
1743 *   FEFaceValues<dim> darcy_fe_face_values(darcy_fe,
1744 *   face_quadrature_formula,
1745 *   update_values |
1746 *   update_normal_vectors |
1747 *   update_quadrature_points |
1748 *   update_JxW_values);
1749 *  
1750 *   const unsigned int dofs_per_cell = darcy_fe.n_dofs_per_cell();
1751 *  
1752 *   const unsigned int n_q_points = quadrature_formula.size();
1753 *   const unsigned int n_face_q_points = face_quadrature_formula.size();
1754 *  
1755 *   FullMatrix<double> local_matrix(dofs_per_cell, dofs_per_cell);
1756 *   Vector<double> local_rhs(dofs_per_cell);
1757 *  
1758 *   std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
1759 *  
1760 *   const Functions::ZeroFunction<dim> pressure_right_hand_side;
1761 *   const PressureBoundaryValues<dim> pressure_boundary_values;
1762 *  
1763 *   std::vector<double> pressure_rhs_values(n_q_points);
1764 *   std::vector<double> boundary_values(n_face_q_points);
1765 *   std::vector<Tensor<2, dim>> k_inverse_values(n_q_points);
1766 *  
1767 * @endcode
1768 *
1769 * Next we need a vector that will contain the values of the saturation
1770 * solution at the previous time level at the quadrature points to
1771 * assemble the saturation dependent coefficients in the Darcy equations.
1772 *
1773
1774 *
1775 * The set of vectors we create next hold the evaluations of the basis
1776 * functions as well as their gradients that will be used for creating the
1777 * matrices. Putting these into their own arrays rather than asking the
1778 * FEValues object for this information each time it is needed is an
1779 * optimization to accelerate the assembly process, see @ref step_22 "step-22" for
1780 * details.
1781 *
1782
1783 *
1784 * The last two declarations are used to extract the individual blocks
1785 * (velocity, pressure, saturation) from the total FE system.
1786 *
1787 * @code
1788 *   std::vector<double> old_saturation_values(n_q_points);
1789 *  
1790 *   std::vector<Tensor<1, dim>> phi_u(dofs_per_cell);
1791 *   std::vector<double> div_phi_u(dofs_per_cell);
1792 *   std::vector<double> phi_p(dofs_per_cell);
1793 *  
1794 *   const FEValuesExtractors::Vector velocities(0);
1795 *   const FEValuesExtractors::Scalar pressure(dim);
1796 *  
1797 * @endcode
1798 *
1799 * Now start the loop over all cells in the problem. We are working on two
1800 * different DoFHandlers for this assembly routine, so we must have two
1801 * different cell iterators for the two objects in use. This might seem a
1802 * bit peculiar, but since both the Darcy system and the saturation system
1803 * use the same grid we can assume that the two iterators run in sync over
1804 * the cells of the two DoFHandler objects.
1805 *
1806
1807 *
1808 * The first statements within the loop are again all very familiar, doing
1809 * the update of the finite element data as specified by the update flags,
1810 * zeroing out the local arrays and getting the values of the old solution
1811 * at the quadrature points. At this point we also have to get the values
1812 * of the saturation function of the previous time step at the quadrature
1813 * points. To this end, we can use the FEValues::get_function_values
1814 * (previously already used in @ref step_9 "step-9", @ref step_14 "step-14" and @ref step_15 "step-15"), a function
1815 * that takes a solution vector and returns a list of function values at
1816 * the quadrature points of the present cell. In fact, it returns the
1817 * complete vector-valued solution at each quadrature point, i.e. not only
1818 * the saturation but also the velocities and pressure.
1819 *
1820
1821 *
1822 * Then we are ready to loop over the quadrature points on the cell to do
1823 * the integration. The formula for this follows in a straightforward way
1824 * from what has been discussed in the introduction.
1825 *
1826
1827 *
1828 * Once this is done, we start the loop over the rows and columns of the
1829 * local matrix and feed the matrix with the relevant products.
1830 *
1831
1832 *
1833 * The last step in the loop over all cells is to enter the local
1834 * contributions into the global matrix and vector structures to the
1835 * positions specified in local_dof_indices. Again, we let the
1836 * AffineConstraints class do the insertion of the cell matrix
1837 * elements to the global matrix, which already condenses the hanging node
1838 * constraints.
1839 *
1840 * @code
1841 *   auto cell = darcy_dof_handler.begin_active();
1842 *   const auto endc = darcy_dof_handler.end();
1843 *   auto saturation_cell = saturation_dof_handler.begin_active();
1844 *  
1845 *   for (; cell != endc; ++cell, ++saturation_cell)
1846 *   {
1847 *   darcy_fe_values.reinit(cell);
1848 *   saturation_fe_values.reinit(saturation_cell);
1849 *  
1850 *   local_matrix = 0;
1851 *   local_rhs = 0;
1852 *  
1853 *   saturation_fe_values.get_function_values(old_saturation_solution,
1854 *   old_saturation_values);
1855 *  
1856 *   pressure_right_hand_side.value_list(
1857 *   darcy_fe_values.get_quadrature_points(), pressure_rhs_values);
1858 *   k_inverse.value_list(darcy_fe_values.get_quadrature_points(),
1859 *   k_inverse_values);
1860 *  
1861 *   for (unsigned int q = 0; q < n_q_points; ++q)
1862 *   {
1863 *   for (unsigned int k = 0; k < dofs_per_cell; ++k)
1864 *   {
1865 *   phi_u[k] = darcy_fe_values[velocities].value(k, q);
1866 *   div_phi_u[k] = darcy_fe_values[velocities].divergence(k, q);
1867 *   phi_p[k] = darcy_fe_values[pressure].value(k, q);
1868 *   }
1869 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
1870 *   {
1871 *   const double old_s = old_saturation_values[q];
1872 *   for (unsigned int j = 0; j <= i; ++j)
1873 *   {
1874 *   local_matrix(i, j) +=
1875 *   (phi_u[i] * k_inverse_values[q] *
1876 *   mobility_inverse(old_s, viscosity) * phi_u[j] -
1877 *   div_phi_u[i] * phi_p[j] - phi_p[i] * div_phi_u[j]) *
1878 *   darcy_fe_values.JxW(q);
1879 *   }
1880 *  
1881 *   local_rhs(i) +=
1882 *   (-phi_p[i] * pressure_rhs_values[q]) * darcy_fe_values.JxW(q);
1883 *   }
1884 *   }
1885 *  
1886 *   for (const auto &face : cell->face_iterators())
1887 *   if (face->at_boundary())
1888 *   {
1889 *   darcy_fe_face_values.reinit(cell, face);
1890 *  
1891 *   pressure_boundary_values.value_list(
1892 *   darcy_fe_face_values.get_quadrature_points(), boundary_values);
1893 *  
1894 *   for (unsigned int q = 0; q < n_face_q_points; ++q)
1895 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
1896 *   {
1897 *   const Tensor<1, dim> phi_i_u =
1898 *   darcy_fe_face_values[velocities].value(i, q);
1899 *  
1900 *   local_rhs(i) +=
1901 *   -(phi_i_u * darcy_fe_face_values.normal_vector(q) *
1902 *   boundary_values[q] * darcy_fe_face_values.JxW(q));
1903 *   }
1904 *   }
1905 *  
1906 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
1907 *   for (unsigned int j = i + 1; j < dofs_per_cell; ++j)
1908 *   local_matrix(i, j) = local_matrix(j, i);
1909 *  
1910 *   cell->get_dof_indices(local_dof_indices);
1911 *  
1912 *   darcy_constraints.distribute_local_to_global(
1913 *   local_matrix, local_rhs, local_dof_indices, darcy_matrix, darcy_rhs);
1914 *   }
1915 *   }
1916 *  
1917 *  
1918 * @endcode
1919 *
1920 *
1921 * <a name="step_43-TwoPhaseFlowProblemdimassemble_saturation_system"></a>
1922 * <h4>TwoPhaseFlowProblem<dim>::assemble_saturation_system</h4>
1923 *
1924
1925 *
1926 * This function is to assemble the linear system for the saturation
1927 * transport equation. It calls, if necessary, two other member functions:
1928 * assemble_saturation_matrix() and assemble_saturation_rhs(). The former
1929 * function then assembles the saturation matrix that only needs to be
1930 * changed occasionally. On the other hand, the latter function that
1931 * assembles the right hand side must be called at every saturation time
1932 * step.
1933 *
1934 * @code
1935 *   template <int dim>
1936 *   void TwoPhaseFlowProblem<dim>::assemble_saturation_system()
1937 *   {
1938 *   if (rebuild_saturation_matrix == true)
1939 *   {
1940 *   saturation_matrix = 0;
1941 *   assemble_saturation_matrix();
1942 *   }
1943 *  
1944 *   saturation_rhs = 0;
1945 *   assemble_saturation_rhs();
1946 *   }
1947 *  
1948 *  
1949 *  
1950 * @endcode
1951 *
1952 *
1953 * <a name="step_43-TwoPhaseFlowProblemdimassemble_saturation_matrix"></a>
1954 * <h4>TwoPhaseFlowProblem<dim>::assemble_saturation_matrix</h4>
1955 *
1956
1957 *
1958 * This function is easily understood since it only forms a simple mass
1959 * matrix for the left hand side of the saturation linear system by basis
1960 * functions phi_i_s and phi_j_s only. Finally, as usual, we enter the local
1961 * contribution into the global matrix by specifying the position in
1962 * local_dof_indices. This is done by letting the AffineConstraints class do
1963 * the insertion of the cell matrix elements to the global matrix, which
1964 * already condenses the hanging node constraints.
1965 *
1966 * @code
1967 *   template <int dim>
1968 *   void TwoPhaseFlowProblem<dim>::assemble_saturation_matrix()
1969 *   {
1970 *   const QGauss<dim> quadrature_formula(saturation_degree + 2);
1971 *  
1972 *   FEValues<dim> saturation_fe_values(saturation_fe,
1973 *   quadrature_formula,
1974 *   update_values | update_JxW_values);
1975 *  
1976 *   const unsigned int dofs_per_cell = saturation_fe.n_dofs_per_cell();
1977 *  
1978 *   const unsigned int n_q_points = quadrature_formula.size();
1979 *  
1980 *   FullMatrix<double> local_matrix(dofs_per_cell, dofs_per_cell);
1981 *   Vector<double> local_rhs(dofs_per_cell);
1982 *  
1983 *   std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
1984 *  
1985 *   for (const auto &cell : saturation_dof_handler.active_cell_iterators())
1986 *   {
1987 *   saturation_fe_values.reinit(cell);
1988 *   local_matrix = 0;
1989 *   local_rhs = 0;
1990 *  
1991 *   for (unsigned int q = 0; q < n_q_points; ++q)
1992 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
1993 *   {
1994 *   const double phi_i_s = saturation_fe_values.shape_value(i, q);
1995 *   for (unsigned int j = 0; j < dofs_per_cell; ++j)
1996 *   {
1997 *   const double phi_j_s = saturation_fe_values.shape_value(j, q);
1998 *   local_matrix(i, j) +=
1999 *   porosity * phi_i_s * phi_j_s * saturation_fe_values.JxW(q);
2000 *   }
2001 *   }
2002 *   cell->get_dof_indices(local_dof_indices);
2003 *  
2004 *   saturation_constraints.distribute_local_to_global(local_matrix,
2005 *   local_dof_indices,
2006 *   saturation_matrix);
2007 *   }
2008 *   }
2009 *  
2010 *  
2011 *  
2012 * @endcode
2013 *
2014 *
2015 * <a name="step_43-TwoPhaseFlowProblemdimassemble_saturation_rhs"></a>
2016 * <h4>TwoPhaseFlowProblem<dim>::assemble_saturation_rhs</h4>
2017 *
2018
2019 *
2020 * This function is to assemble the right hand side of the saturation
2021 * transport equation. Before going about it, we have to create two FEValues
2022 * objects for the Darcy and saturation systems respectively and, in
2023 * addition, two FEFaceValues objects for the two systems because we have a
2024 * boundary integral term in the weak form of saturation equation. For the
2025 * FEFaceValues object of the saturation system, we also require normal
2026 * vectors, which we request using the update_normal_vectors flag.
2027 *
2028
2029 *
2030 * Next, before looping over all the cells, we have to compute some
2031 * parameters (e.g. global_u_infty, global_S_variation, and
2032 * global_Omega_diameter) that the artificial viscosity @f$\nu@f$ needs. This is
2033 * largely the same as was done in @ref step_31 "step-31", so you may see there for more
2034 * information.
2035 *
2036
2037 *
2038 * The real works starts with the loop over all the saturation and Darcy
2039 * cells to put the local contributions into the global vector. In this
2040 * loop, in order to simplify the implementation, we split some of the work
2041 * into two helper functions: assemble_saturation_rhs_cell_term and
2042 * assemble_saturation_rhs_boundary_term. We note that we insert cell or
2043 * boundary contributions into the global vector in the two functions rather
2044 * than in this present function.
2045 *
2046 * @code
2047 *   template <int dim>
2048 *   void TwoPhaseFlowProblem<dim>::assemble_saturation_rhs()
2049 *   {
2050 *   const QGauss<dim> quadrature_formula(saturation_degree + 2);
2051 *   const QGauss<dim - 1> face_quadrature_formula(saturation_degree + 2);
2052 *  
2053 *   FEValues<dim> saturation_fe_values(saturation_fe,
2054 *   quadrature_formula,
2055 *   update_values | update_gradients |
2056 *   update_quadrature_points |
2057 *   update_JxW_values);
2058 *   FEValues<dim> darcy_fe_values(darcy_fe, quadrature_formula, update_values);
2059 *   FEFaceValues<dim> saturation_fe_face_values(saturation_fe,
2060 *   face_quadrature_formula,
2061 *   update_values |
2062 *   update_normal_vectors |
2063 *   update_quadrature_points |
2064 *   update_JxW_values);
2065 *   FEFaceValues<dim> darcy_fe_face_values(darcy_fe,
2066 *   face_quadrature_formula,
2067 *   update_values);
2068 *   FEFaceValues<dim> saturation_fe_face_values_neighbor(
2069 *   saturation_fe, face_quadrature_formula, update_values);
2070 *  
2071 *   const unsigned int dofs_per_cell =
2072 *   saturation_dof_handler.get_fe().n_dofs_per_cell();
2073 *   std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
2074 *  
2075 *   const double global_max_u_F_prime = get_max_u_F_prime();
2076 *   const std::pair<double, double> global_S_range =
2077 *   get_extrapolated_saturation_range();
2078 *   const double global_S_variation =
2079 *   global_S_range.second - global_S_range.first;
2080 *  
2081 *   auto cell = saturation_dof_handler.begin_active();
2082 *   const auto endc = saturation_dof_handler.end();
2083 *   auto darcy_cell = darcy_dof_handler.begin_active();
2084 *   for (; cell != endc; ++cell, ++darcy_cell)
2085 *   {
2086 *   saturation_fe_values.reinit(cell);
2087 *   darcy_fe_values.reinit(darcy_cell);
2088 *  
2089 *   cell->get_dof_indices(local_dof_indices);
2090 *  
2091 *   assemble_saturation_rhs_cell_term(saturation_fe_values,
2092 *   darcy_fe_values,
2093 *   global_max_u_F_prime,
2094 *   global_S_variation,
2095 *   local_dof_indices);
2096 *  
2097 *   for (const auto &face : cell->face_iterators())
2098 *   if (face->at_boundary())
2099 *   {
2100 *   darcy_fe_face_values.reinit(darcy_cell, face);
2101 *   saturation_fe_face_values.reinit(cell, face);
2102 *   assemble_saturation_rhs_boundary_term(saturation_fe_face_values,
2103 *   darcy_fe_face_values,
2104 *   local_dof_indices);
2105 *   }
2106 *   }
2107 *   }
2108 *  
2109 *  
2110 *  
2111 * @endcode
2112 *
2113 *
2114 * <a name="step_43-TwoPhaseFlowProblemdimassemble_saturation_rhs_cell_term"></a>
2115 * <h4>TwoPhaseFlowProblem<dim>::assemble_saturation_rhs_cell_term</h4>
2116 *
2117
2118 *
2119 * This function takes care of integrating the cell terms of the right hand
2120 * side of the saturation equation, and then assembling it into the global
2121 * right hand side vector. Given the discussion in the introduction, the
2122 * form of these contributions is clear. The only tricky part is getting the
2123 * artificial viscosity and all that is necessary to compute it. The first
2124 * half of the function is devoted to this task.
2125 *
2126
2127 *
2128 * The last part of the function is copying the local contributions into the
2129 * global vector with position specified in local_dof_indices.
2130 *
2131 * @code
2132 *   template <int dim>
2133 *   void TwoPhaseFlowProblem<dim>::assemble_saturation_rhs_cell_term(
2134 *   const FEValues<dim> &saturation_fe_values,
2135 *   const FEValues<dim> &darcy_fe_values,
2136 *   const double global_max_u_F_prime,
2137 *   const double global_S_variation,
2138 *   const std::vector<types::global_dof_index> &local_dof_indices)
2139 *   {
2140 *   const unsigned int dofs_per_cell = saturation_fe_values.dofs_per_cell;
2141 *   const unsigned int n_q_points = saturation_fe_values.n_quadrature_points;
2142 *  
2143 *   std::vector<double> old_saturation_solution_values(n_q_points);
2144 *   std::vector<double> old_old_saturation_solution_values(n_q_points);
2145 *   std::vector<Tensor<1, dim>> old_grad_saturation_solution_values(n_q_points);
2146 *   std::vector<Tensor<1, dim>> old_old_grad_saturation_solution_values(
2147 *   n_q_points);
2148 *   std::vector<Vector<double>> present_darcy_solution_values(
2149 *   n_q_points, Vector<double>(dim + 1));
2150 *  
2151 *   saturation_fe_values.get_function_values(old_saturation_solution,
2152 *   old_saturation_solution_values);
2153 *   saturation_fe_values.get_function_values(
2154 *   old_old_saturation_solution, old_old_saturation_solution_values);
2155 *   saturation_fe_values.get_function_gradients(
2156 *   old_saturation_solution, old_grad_saturation_solution_values);
2157 *   saturation_fe_values.get_function_gradients(
2158 *   old_old_saturation_solution, old_old_grad_saturation_solution_values);
2159 *   darcy_fe_values.get_function_values(darcy_solution,
2160 *   present_darcy_solution_values);
2161 *  
2162 *   const double nu =
2163 *   compute_viscosity(old_saturation_solution_values,
2164 *   old_old_saturation_solution_values,
2165 *   old_grad_saturation_solution_values,
2166 *   old_old_grad_saturation_solution_values,
2167 *   present_darcy_solution_values,
2168 *   global_max_u_F_prime,
2169 *   global_S_variation,
2170 *   saturation_fe_values.get_cell()->diameter());
2171 *  
2172 *   Vector<double> local_rhs(dofs_per_cell);
2173 *  
2174 *   for (unsigned int q = 0; q < n_q_points; ++q)
2175 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
2176 *   {
2177 *   const double old_s = old_saturation_solution_values[q];
2178 *   Tensor<1, dim> present_u;
2179 *   for (unsigned int d = 0; d < dim; ++d)
2180 *   present_u[d] = present_darcy_solution_values[q](d);
2181 *  
2182 *   const double phi_i_s = saturation_fe_values.shape_value(i, q);
2183 *   const Tensor<1, dim> grad_phi_i_s =
2184 *   saturation_fe_values.shape_grad(i, q);
2185 *  
2186 *   local_rhs(i) +=
2187 *   (time_step * fractional_flow(old_s, viscosity) * present_u *
2188 *   grad_phi_i_s -
2189 *   time_step * nu * old_grad_saturation_solution_values[q] *
2190 *   grad_phi_i_s +
2191 *   porosity * old_s * phi_i_s) *
2192 *   saturation_fe_values.JxW(q);
2193 *   }
2194 *  
2195 *   saturation_constraints.distribute_local_to_global(local_rhs,
2196 *   local_dof_indices,
2197 *   saturation_rhs);
2198 *   }
2199 *  
2200 *  
2201 * @endcode
2202 *
2203 *
2204 * <a name="step_43-TwoPhaseFlowProblemdimassemble_saturation_rhs_boundary_term"></a>
2205 * <h4>TwoPhaseFlowProblem<dim>::assemble_saturation_rhs_boundary_term</h4>
2206 *
2207
2208 *
2209 * The next function is responsible for the boundary integral terms in the
2210 * right hand side form of the saturation equation. For these, we have to
2211 * compute the upwinding flux on the global boundary faces, i.e. we impose
2212 * Dirichlet boundary conditions weakly only on inflow parts of the global
2213 * boundary. As before, this has been described in @ref step_21 "step-21" so we refrain
2214 * from giving more descriptions about that.
2215 *
2216 * @code
2217 *   template <int dim>
2218 *   void TwoPhaseFlowProblem<dim>::assemble_saturation_rhs_boundary_term(
2219 *   const FEFaceValues<dim> &saturation_fe_face_values,
2220 *   const FEFaceValues<dim> &darcy_fe_face_values,
2221 *   const std::vector<types::global_dof_index> &local_dof_indices)
2222 *   {
2223 *   const unsigned int dofs_per_cell = saturation_fe_face_values.dofs_per_cell;
2224 *   const unsigned int n_face_q_points =
2225 *   saturation_fe_face_values.n_quadrature_points;
2226 *  
2227 *   Vector<double> local_rhs(dofs_per_cell);
2228 *  
2229 *   std::vector<double> old_saturation_solution_values_face(n_face_q_points);
2230 *   std::vector<Vector<double>> present_darcy_solution_values_face(
2231 *   n_face_q_points, Vector<double>(dim + 1));
2232 *   std::vector<double> neighbor_saturation(n_face_q_points);
2233 *  
2234 *   saturation_fe_face_values.get_function_values(
2235 *   old_saturation_solution, old_saturation_solution_values_face);
2236 *   darcy_fe_face_values.get_function_values(
2237 *   darcy_solution, present_darcy_solution_values_face);
2238 *  
2239 *   SaturationBoundaryValues<dim> saturation_boundary_values;
2240 *   saturation_boundary_values.value_list(
2241 *   saturation_fe_face_values.get_quadrature_points(), neighbor_saturation);
2242 *  
2243 *   for (unsigned int q = 0; q < n_face_q_points; ++q)
2244 *   {
2245 *   Tensor<1, dim> present_u_face;
2246 *   for (unsigned int d = 0; d < dim; ++d)
2247 *   present_u_face[d] = present_darcy_solution_values_face[q](d);
2248 *  
2249 *   const double normal_flux =
2250 *   present_u_face * saturation_fe_face_values.normal_vector(q);
2251 *  
2252 *   const bool is_outflow_q_point = (normal_flux >= 0);
2253 *  
2254 *   for (unsigned int i = 0; i < dofs_per_cell; ++i)
2255 *   local_rhs(i) -=
2256 *   time_step * normal_flux *
2257 *   fractional_flow((is_outflow_q_point == true ?
2258 *   old_saturation_solution_values_face[q] :
2259 *   neighbor_saturation[q]),
2260 *   viscosity) *
2261 *   saturation_fe_face_values.shape_value(i, q) *
2262 *   saturation_fe_face_values.JxW(q);
2263 *   }
2264 *   saturation_constraints.distribute_local_to_global(local_rhs,
2265 *   local_dof_indices,
2266 *   saturation_rhs);
2267 *   }
2268 *  
2269 *  
2270 * @endcode
2271 *
2272 *
2273 * <a name="step_43-TwoPhaseFlowProblemdimsolve"></a>
2274 * <h3>TwoPhaseFlowProblem<dim>::solve</h3>
2275 *
2276
2277 *
2278 * This function implements the operator splitting algorithm, i.e. in each
2279 * time step it either re-computes the solution of the Darcy system or
2280 * extrapolates velocity/pressure from previous time steps, then determines
2281 * the size of the time step, and then updates the saturation variable. The
2282 * implementation largely follows similar code in @ref step_31 "step-31". It is, next to
2283 * the run() function, the central one in this program.
2284 *
2285
2286 *
2287 * At the beginning of the function, we ask whether to solve the
2288 * pressure-velocity part by evaluating the a posteriori criterion (see the
2289 * following function). If necessary, we will solve the pressure-velocity
2290 * part using the GMRES solver with the Schur complement block
2291 * preconditioner as is described in the introduction.
2292 *
2293 * @code
2294 *   template <int dim>
2295 *   void TwoPhaseFlowProblem<dim>::solve()
2296 *   {
2297 *   const bool solve_for_pressure_and_velocity =
2298 *   determine_whether_to_solve_for_pressure_and_velocity();
2299 *  
2300 *   if (solve_for_pressure_and_velocity == true)
2301 *   {
2302 *   std::cout << " Solving Darcy (pressure-velocity) system..."
2303 *   << std::endl;
2304 *  
2305 *   assemble_darcy_system();
2306 *   build_darcy_preconditioner();
2307 *  
2308 *   {
2309 *   const LinearSolvers::InverseMatrix<TrilinosWrappers::SparseMatrix,
2310 *   TrilinosWrappers::PreconditionIC>
2311 *   mp_inverse(darcy_preconditioner_matrix.block(1, 1),
2312 *   *bottom_right_preconditioner);
2313 *  
2314 *   const LinearSolvers::BlockSchurPreconditioner<
2315 *   TrilinosWrappers::PreconditionIC,
2316 *   TrilinosWrappers::PreconditionIC>
2317 *   preconditioner(darcy_matrix, mp_inverse, *top_left_preconditioner);
2318 *  
2319 *   SolverControl solver_control(darcy_matrix.m(),
2320 *   1e-16 * darcy_rhs.l2_norm());
2321 *  
2322 *   SolverGMRES<TrilinosWrappers::MPI::BlockVector> gmres(
2323 *   solver_control,
2324 *   SolverGMRES<TrilinosWrappers::MPI::BlockVector>::AdditionalData(
2325 *   100));
2326 *  
2327 *   for (unsigned int i = 0; i < darcy_solution.size(); ++i)
2328 *   if (darcy_constraints.is_constrained(i))
2329 *   darcy_solution(i) = 0;
2330 *  
2331 *   gmres.solve(darcy_matrix, darcy_solution, darcy_rhs, preconditioner);
2332 *  
2333 *   darcy_constraints.distribute(darcy_solution);
2334 *  
2335 *   std::cout << " ..." << solver_control.last_step()
2336 *   << " GMRES iterations." << std::endl;
2337 *   }
2338 *  
2339 *   {
2340 *   second_last_computed_darcy_solution = last_computed_darcy_solution;
2341 *   last_computed_darcy_solution = darcy_solution;
2342 *  
2343 *   saturation_matching_last_computed_darcy_solution =
2344 *   saturation_solution;
2345 *   }
2346 *   }
2347 * @endcode
2348 *
2349 * On the other hand, if we have decided that we don't want to compute the
2350 * solution of the Darcy system for the current time step, then we need to
2351 * simply extrapolate the previous two Darcy solutions to the same time as
2352 * we would have computed the velocity/pressure at. We do a simple linear
2353 * extrapolation, i.e. given the current length @f$dt@f$ of the macro time
2354 * step from the time when we last computed the Darcy solution to now
2355 * (given by <code>current_macro_time_step</code>), and @f$DT@f$ the length of
2356 * the last macro time step (given by <code>old_macro_time_step</code>),
2357 * then we get @f$u^\ast = u_p + dt \frac{u_p-u_{pp}}{DT} = (1+dt/DT)u_p -
2358 * dt/DT u_{pp}@f$, where @f$u_p@f$ and @f$u_{pp}@f$ are the last two computed Darcy
2359 * solutions. We can implement this formula using just two lines of code.
2360 *
2361
2362 *
2363 * Note that the algorithm here only works if we have at least two
2364 * previously computed Darcy solutions from which we can extrapolate to
2365 * the current time, and this is ensured by requiring re-computation of
2366 * the Darcy solution for the first 2 time steps.
2367 *
2368 * @code
2369 *   else
2370 *   {
2371 *   darcy_solution = last_computed_darcy_solution;
2372 *   darcy_solution.sadd(1 + current_macro_time_step / old_macro_time_step,
2373 *   -current_macro_time_step / old_macro_time_step,
2374 *   second_last_computed_darcy_solution);
2375 *   }
2376 *  
2377 *  
2378 * @endcode
2379 *
2380 * With the so computed velocity vector, compute the optimal time step
2381 * based on the CFL criterion discussed in the introduction...
2382 *
2383 * @code
2384 *   {
2385 *   old_time_step = time_step;
2386 *  
2387 *   const double max_u_F_prime = get_max_u_F_prime();
2388 *   if (max_u_F_prime > 0)
2389 *   time_step = porosity * GridTools::minimal_cell_diameter(triangulation) /
2390 *   saturation_degree / max_u_F_prime / 50;
2391 *   else
2392 *   time_step = end_time - time;
2393 *   }
2394 *  
2395 *  
2396 *  
2397 * @endcode
2398 *
2399 * ...and then also update the length of the macro time steps we use while
2400 * we're dealing with time step sizes. In particular, this involves: (i)
2401 * If we have just recomputed the Darcy solution, then the length of the
2402 * previous macro time step is now fixed and the length of the current
2403 * macro time step is, up to now, simply the length of the current (micro)
2404 * time step. (ii) If we have not recomputed the Darcy solution, then the
2405 * length of the current macro time step has just grown by
2406 * <code>time_step</code>.
2407 *
2408 * @code
2409 *   if (solve_for_pressure_and_velocity == true)
2410 *   {
2411 *   old_macro_time_step = current_macro_time_step;
2412 *   current_macro_time_step = time_step;
2413 *   }
2414 *   else
2415 *   current_macro_time_step += time_step;
2416 *  
2417 * @endcode
2418 *
2419 * The last step in this function is to recompute the saturation solution
2420 * based on the velocity field we've just obtained. This naturally happens
2421 * in every time step, and we don't skip any of these computations. At the
2422 * end of computing the saturation, we project back into the allowed
2423 * interval @f$[0,1]@f$ to make sure our solution remains physical.
2424 *
2425 * @code
2426 *   {
2427 *   std::cout << " Solving saturation transport equation..." << std::endl;
2428 *  
2429 *   assemble_saturation_system();
2430 *  
2431 *   SolverControl solver_control(saturation_matrix.m(),
2432 *   1e-16 * saturation_rhs.l2_norm());
2433 *   SolverCG<TrilinosWrappers::MPI::Vector> cg(solver_control);
2434 *  
2435 *   TrilinosWrappers::PreconditionIC preconditioner;
2436 *   preconditioner.initialize(saturation_matrix);
2437 *  
2438 *   cg.solve(saturation_matrix,
2439 *   saturation_solution,
2440 *   saturation_rhs,
2441 *   preconditioner);
2442 *  
2443 *   saturation_constraints.distribute(saturation_solution);
2444 *   project_back_saturation();
2445 *  
2446 *   std::cout << " ..." << solver_control.last_step()
2447 *   << " CG iterations." << std::endl;
2448 *   }
2449 *   }
2450 *  
2451 *  
2452 * @endcode
2453 *
2454 *
2455 * <a name="step_43-TwoPhaseFlowProblemdimrefine_mesh"></a>
2456 * <h3>TwoPhaseFlowProblem<dim>::refine_mesh</h3>
2457 *
2458
2459 *
2460 * The next function does the refinement and coarsening of the mesh. It does
2461 * its work in three blocks: (i) Compute refinement indicators by looking at
2462 * the gradient of a solution vector extrapolated linearly from the previous
2463 * two using the respective sizes of the time step (or taking the only
2464 * solution we have if this is the first time step). (ii) Flagging those
2465 * cells for refinement and coarsening where the gradient is larger or
2466 * smaller than a certain threshold, preserving minimal and maximal levels
2467 * of mesh refinement. (iii) Transferring the solution from the old to the
2468 * new mesh. None of this is particularly difficult.
2469 *
2470 * @code
2471 *   template <int dim>
2472 *   void TwoPhaseFlowProblem<dim>::refine_mesh(const unsigned int min_grid_level,
2473 *   const unsigned int max_grid_level)
2474 *   {
2475 *   Vector<double> refinement_indicators(triangulation.n_active_cells());
2476 *   {
2477 *   const QMidpoint<dim> quadrature_formula;
2478 *   FEValues<dim> fe_values(saturation_fe,
2479 *   quadrature_formula,
2480 *   update_gradients);
2481 *   std::vector<Tensor<1, dim>> grad_saturation(1);
2482 *  
2483 *   TrilinosWrappers::MPI::Vector extrapolated_saturation_solution(
2484 *   saturation_solution);
2485 *   if (timestep_number != 0)
2486 *   extrapolated_saturation_solution.sadd((1. + time_step / old_time_step),
2487 *   time_step / old_time_step,
2488 *   old_saturation_solution);
2489 *  
2490 *   for (const auto &cell : saturation_dof_handler.active_cell_iterators())
2491 *   {
2492 *   const unsigned int cell_no = cell->active_cell_index();
2493 *   fe_values.reinit(cell);
2494 *   fe_values.get_function_gradients(extrapolated_saturation_solution,
2495 *   grad_saturation);
2496 *  
2497 *   refinement_indicators(cell_no) = grad_saturation[0].norm();
2498 *   }
2499 *   }
2500 *  
2501 *   {
2502 *   for (const auto &cell : saturation_dof_handler.active_cell_iterators())
2503 *   {
2504 *   const unsigned int cell_no = cell->active_cell_index();
2505 *   cell->clear_coarsen_flag();
2506 *   cell->clear_refine_flag();
2507 *  
2508 *   if ((static_cast<unsigned int>(cell->level()) < max_grid_level) &&
2509 *   (std::fabs(refinement_indicators(cell_no)) >
2510 *   saturation_refinement_threshold))
2511 *   cell->set_refine_flag();
2512 *   else if ((static_cast<unsigned int>(cell->level()) >
2513 *   min_grid_level) &&
2514 *   (std::fabs(refinement_indicators(cell_no)) <
2515 *   0.5 * saturation_refinement_threshold))
2516 *   cell->set_coarsen_flag();
2517 *   }
2518 *   }
2519 *  
2520 *   triangulation.prepare_coarsening_and_refinement();
2521 *  
2522 *   {
2523 *   std::vector<TrilinosWrappers::MPI::Vector> x_saturation(3);
2524 *   x_saturation[0] = saturation_solution;
2525 *   x_saturation[1] = old_saturation_solution;
2526 *   x_saturation[2] = saturation_matching_last_computed_darcy_solution;
2527 *  
2528 *   std::vector<TrilinosWrappers::MPI::BlockVector> x_darcy(2);
2529 *   x_darcy[0] = last_computed_darcy_solution;
2530 *   x_darcy[1] = second_last_computed_darcy_solution;
2531 *  
2532 *   SolutionTransfer<dim, TrilinosWrappers::MPI::Vector> saturation_soltrans(
2533 *   saturation_dof_handler);
2534 *  
2535 *   SolutionTransfer<dim, TrilinosWrappers::MPI::BlockVector> darcy_soltrans(
2536 *   darcy_dof_handler);
2537 *  
2538 *  
2539 *   triangulation.prepare_coarsening_and_refinement();
2540 *   saturation_soltrans.prepare_for_coarsening_and_refinement(x_saturation);
2541 *  
2542 *   darcy_soltrans.prepare_for_coarsening_and_refinement(x_darcy);
2543 *  
2544 *   triangulation.execute_coarsening_and_refinement();
2545 *   setup_dofs();
2546 *  
2547 *   std::vector<TrilinosWrappers::MPI::Vector> tmp_saturation(3);
2548 *   tmp_saturation[0].reinit(saturation_solution);
2549 *   tmp_saturation[1].reinit(saturation_solution);
2550 *   tmp_saturation[2].reinit(saturation_solution);
2551 *   saturation_soltrans.interpolate(tmp_saturation);
2552 *  
2553 *   saturation_solution = tmp_saturation[0];
2554 *   old_saturation_solution = tmp_saturation[1];
2555 *   saturation_matching_last_computed_darcy_solution = tmp_saturation[2];
2556 *  
2557 *   saturation_constraints.distribute(saturation_solution);
2558 *   saturation_constraints.distribute(old_saturation_solution);
2559 *   saturation_constraints.distribute(
2560 *   saturation_matching_last_computed_darcy_solution);
2561 *  
2562 *   std::vector<TrilinosWrappers::MPI::BlockVector> tmp_darcy(2);
2563 *   tmp_darcy[0].reinit(darcy_solution);
2564 *   tmp_darcy[1].reinit(darcy_solution);
2565 *   darcy_soltrans.interpolate(tmp_darcy);
2566 *  
2567 *   last_computed_darcy_solution = tmp_darcy[0];
2568 *   second_last_computed_darcy_solution = tmp_darcy[1];
2569 *  
2570 *   darcy_constraints.distribute(last_computed_darcy_solution);
2571 *   darcy_constraints.distribute(second_last_computed_darcy_solution);
2572 *  
2573 *   rebuild_saturation_matrix = true;
2574 *   }
2575 *   }
2576 *  
2577 *  
2578 *  
2579 * @endcode
2580 *
2581 *
2582 * <a name="step_43-TwoPhaseFlowProblemdimoutput_results"></a>
2583 * <h3>TwoPhaseFlowProblem<dim>::output_results</h3>
2584 *
2585
2586 *
2587 * This function generates graphical output. It is in essence a copy of the
2588 * implementation in @ref step_31 "step-31".
2589 *
2590 * @code
2591 *   template <int dim>
2592 *   void TwoPhaseFlowProblem<dim>::output_results() const
2593 *   {
2594 *   const FESystem<dim> joint_fe(darcy_fe, 1, saturation_fe, 1);
2595 *   DoFHandler<dim> joint_dof_handler(triangulation);
2596 *   joint_dof_handler.distribute_dofs(joint_fe);
2597 *   Assert(joint_dof_handler.n_dofs() ==
2598 *   darcy_dof_handler.n_dofs() + saturation_dof_handler.n_dofs(),
2599 *   ExcInternalError());
2600 *  
2601 *   Vector<double> joint_solution(joint_dof_handler.n_dofs());
2602 *  
2603 *   {
2604 *   std::vector<types::global_dof_index> local_joint_dof_indices(
2605 *   joint_fe.n_dofs_per_cell());
2606 *   std::vector<types::global_dof_index> local_darcy_dof_indices(
2607 *   darcy_fe.n_dofs_per_cell());
2608 *   std::vector<types::global_dof_index> local_saturation_dof_indices(
2609 *   saturation_fe.n_dofs_per_cell());
2610 *  
2611 *   auto joint_cell = joint_dof_handler.begin_active();
2612 *   const auto joint_endc = joint_dof_handler.end();
2613 *   auto darcy_cell = darcy_dof_handler.begin_active();
2614 *   auto saturation_cell = saturation_dof_handler.begin_active();
2615 *  
2616 *   for (; joint_cell != joint_endc;
2617 *   ++joint_cell, ++darcy_cell, ++saturation_cell)
2618 *   {
2619 *   joint_cell->get_dof_indices(local_joint_dof_indices);
2620 *   darcy_cell->get_dof_indices(local_darcy_dof_indices);
2621 *   saturation_cell->get_dof_indices(local_saturation_dof_indices);
2622 *  
2623 *   for (unsigned int i = 0; i < joint_fe.n_dofs_per_cell(); ++i)
2624 *   if (joint_fe.system_to_base_index(i).first.first == 0)
2625 *   {
2626 *   Assert(joint_fe.system_to_base_index(i).second <
2627 *   local_darcy_dof_indices.size(),
2628 *   ExcInternalError());
2629 *   joint_solution(local_joint_dof_indices[i]) = darcy_solution(
2630 *   local_darcy_dof_indices[joint_fe.system_to_base_index(i)
2631 *   .second]);
2632 *   }
2633 *   else
2634 *   {
2635 *   Assert(joint_fe.system_to_base_index(i).first.first == 1,
2636 *   ExcInternalError());
2637 *   Assert(joint_fe.system_to_base_index(i).second <
2638 *   local_darcy_dof_indices.size(),
2639 *   ExcInternalError());
2640 *   joint_solution(local_joint_dof_indices[i]) =
2641 *   saturation_solution(
2642 *   local_saturation_dof_indices
2643 *   [joint_fe.system_to_base_index(i).second]);
2644 *   }
2645 *   }
2646 *   }
2647 *   std::vector<std::string> joint_solution_names(dim, "velocity");
2648 *   joint_solution_names.emplace_back("pressure");
2649 *   joint_solution_names.emplace_back("saturation");
2650 *  
2651 *   std::vector<DataComponentInterpretation::DataComponentInterpretation>
2652 *   data_component_interpretation(
2653 *   dim, DataComponentInterpretation::component_is_part_of_vector);
2654 *   data_component_interpretation.push_back(
2655 *   DataComponentInterpretation::component_is_scalar);
2656 *   data_component_interpretation.push_back(
2657 *   DataComponentInterpretation::component_is_scalar);
2658 *  
2659 *   DataOut<dim> data_out;
2660 *  
2661 *   data_out.attach_dof_handler(joint_dof_handler);
2662 *   data_out.add_data_vector(joint_solution,
2663 *   joint_solution_names,
2664 *   DataOut<dim>::type_dof_data,
2665 *   data_component_interpretation);
2666 *  
2667 *   data_out.build_patches();
2668 *  
2669 *   std::string filename =
2670 *   "solution-" + Utilities::int_to_string(timestep_number, 5) + ".vtu";
2671 *   std::ofstream output(filename);
2672 *   data_out.write_vtu(output);
2673 *   }
2674 *  
2675 *  
2676 *  
2677 * @endcode
2678 *
2679 *
2680 * <a name="step_43-Toolfunctions"></a>
2681 * <h3>Tool functions</h3>
2682 *
2683
2684 *
2685 *
2686 * <a name="step_43-TwoPhaseFlowProblemdimdetermine_whether_to_solve_for_pressure_and_velocity"></a>
2687 * <h4>TwoPhaseFlowProblem<dim>::determine_whether_to_solve_for_pressure_and_velocity</h4>
2688 *
2689
2690 *
2691 * This function implements the a posteriori criterion for adaptive operator
2692 * splitting. The function is relatively straightforward given the way we
2693 * have implemented other functions above and given the formula for the
2694 * criterion derived in the paper.
2695 *
2696
2697 *
2698 * If one decides that one wants the original IMPES method in which the
2699 * Darcy equation is solved in every time step, then this can be achieved by
2700 * setting the threshold value <code>AOS_threshold</code> (with a default of
2701 * @f$5.0@f$) to zero, thereby forcing the function to always return true.
2702 *
2703
2704 *
2705 * Finally, note that the function returns true unconditionally for the
2706 * first two time steps to ensure that we have always solved the Darcy
2707 * system at least twice when skipping its solution, thereby allowing us to
2708 * extrapolate the velocity from the last two solutions in
2709 * <code>solve()</code>.
2710 *
2711 * @code
2712 *   template <int dim>
2713 *   bool TwoPhaseFlowProblem<
2714 *   dim>::determine_whether_to_solve_for_pressure_and_velocity() const
2715 *   {
2716 *   if (timestep_number <= 2)
2717 *   return true;
2718 *  
2719 *   const QGauss<dim> quadrature_formula(saturation_degree + 2);
2720 *   const unsigned int n_q_points = quadrature_formula.size();
2721 *  
2722 *   FEValues<dim> fe_values(saturation_fe,
2723 *   quadrature_formula,
2724 *   update_values | update_quadrature_points);
2725 *  
2726 *   std::vector<double> old_saturation_after_solving_pressure(n_q_points);
2727 *   std::vector<double> present_saturation(n_q_points);
2728 *  
2729 *   std::vector<Tensor<2, dim>> k_inverse_values(n_q_points);
2730 *  
2731 *   double max_global_aop_indicator = 0.0;
2732 *  
2733 *   for (const auto &cell : saturation_dof_handler.active_cell_iterators())
2734 *   {
2735 *   double max_local_mobility_reciprocal_difference = 0.0;
2736 *   double max_local_permeability_inverse_l1_norm = 0.0;
2737 *  
2738 *   fe_values.reinit(cell);
2739 *   fe_values.get_function_values(
2740 *   saturation_matching_last_computed_darcy_solution,
2741 *   old_saturation_after_solving_pressure);
2742 *   fe_values.get_function_values(saturation_solution, present_saturation);
2743 *  
2744 *   k_inverse.value_list(fe_values.get_quadrature_points(),
2745 *   k_inverse_values);
2746 *  
2747 *   for (unsigned int q = 0; q < n_q_points; ++q)
2748 *   {
2749 *   const double mobility_reciprocal_difference = std::fabs(
2750 *   mobility_inverse(present_saturation[q], viscosity) -
2751 *   mobility_inverse(old_saturation_after_solving_pressure[q],
2752 *   viscosity));
2753 *  
2754 *   max_local_mobility_reciprocal_difference =
2755 *   std::max(max_local_mobility_reciprocal_difference,
2756 *   mobility_reciprocal_difference);
2757 *  
2758 *   max_local_permeability_inverse_l1_norm =
2759 *   std::max(max_local_permeability_inverse_l1_norm,
2760 *   l1_norm(k_inverse_values[q]));
2761 *   }
2762 *  
2763 *   max_global_aop_indicator =
2764 *   std::max(max_global_aop_indicator,
2765 *   (max_local_mobility_reciprocal_difference *
2766 *   max_local_permeability_inverse_l1_norm));
2767 *   }
2768 *  
2769 *   return (max_global_aop_indicator > AOS_threshold);
2770 *   }
2771 *  
2772 *  
2773 *  
2774 * @endcode
2775 *
2776 *
2777 * <a name="step_43-TwoPhaseFlowProblemdimproject_back_saturation"></a>
2778 * <h4>TwoPhaseFlowProblem<dim>::project_back_saturation</h4>
2779 *
2780
2781 *
2782 * The next function simply makes sure that the saturation values always
2783 * remain within the physically reasonable range of @f$[0,1]@f$. While the
2784 * continuous equations guarantee that this is so, the discrete equations
2785 * don't. However, if we allow the discrete solution to escape this range we
2786 * get into trouble because terms like @f$F(S)@f$ and @f$F'(S)@f$ will produce
2787 * unreasonable results (e.g. @f$F'(S)<0@f$ for @f$S<0@f$, which would imply that
2788 * the wetting fluid phase flows <i>against</i> the direction of the bulk
2789 * fluid velocity)). Consequently, at the end of each time step, we simply
2790 * project the saturation field back into the physically reasonable region.
2791 *
2792 * @code
2793 *   template <int dim>
2794 *   void TwoPhaseFlowProblem<dim>::project_back_saturation()
2795 *   {
2796 *   for (unsigned int i = 0; i < saturation_solution.size(); ++i)
2797 *   if (saturation_solution(i) < 0.2)
2798 *   saturation_solution(i) = 0.2;
2799 *   else if (saturation_solution(i) > 1)
2800 *   saturation_solution(i) = 1;
2801 *   }
2802 *  
2803 *  
2804 *  
2805 * @endcode
2806 *
2807 *
2808 * <a name="step_43-TwoPhaseFlowProblemdimget_max_u_F_prime"></a>
2809 * <h4>TwoPhaseFlowProblem<dim>::get_max_u_F_prime</h4>
2810 *
2811
2812 *
2813 * Another simpler helper function: Compute the maximum of the total
2814 * velocity times the derivative of the fraction flow function, i.e.,
2815 * compute @f$\|\mathbf{u} F'(S)\|_{L_\infty(\Omega)}@f$. This term is used in
2816 * both the computation of the time step as well as in normalizing the
2817 * entropy-residual term in the artificial viscosity.
2818 *
2819 * @code
2820 *   template <int dim>
2821 *   double TwoPhaseFlowProblem<dim>::get_max_u_F_prime() const
2822 *   {
2823 *   const QGauss<dim> quadrature_formula(darcy_degree + 2);
2824 *   const unsigned int n_q_points = quadrature_formula.size();
2825 *  
2826 *   FEValues<dim> darcy_fe_values(darcy_fe, quadrature_formula, update_values);
2827 *   FEValues<dim> saturation_fe_values(saturation_fe,
2828 *   quadrature_formula,
2829 *   update_values);
2830 *  
2831 *   std::vector<Vector<double>> darcy_solution_values(n_q_points,
2832 *   Vector<double>(dim + 1));
2833 *   std::vector<double> saturation_values(n_q_points);
2834 *  
2835 *   double max_velocity_times_dF_dS = 0;
2836 *  
2837 *   auto cell = darcy_dof_handler.begin_active();
2838 *   const auto endc = darcy_dof_handler.end();
2839 *   auto saturation_cell = saturation_dof_handler.begin_active();
2840 *   for (; cell != endc; ++cell, ++saturation_cell)
2841 *   {
2842 *   darcy_fe_values.reinit(cell);
2843 *   saturation_fe_values.reinit(saturation_cell);
2844 *  
2845 *   darcy_fe_values.get_function_values(darcy_solution,
2846 *   darcy_solution_values);
2847 *   saturation_fe_values.get_function_values(old_saturation_solution,
2848 *   saturation_values);
2849 *  
2850 *   for (unsigned int q = 0; q < n_q_points; ++q)
2851 *   {
2852 *   Tensor<1, dim> velocity;
2853 *   for (unsigned int i = 0; i < dim; ++i)
2854 *   velocity[i] = darcy_solution_values[q](i);
2855 *  
2856 *   const double dF_dS =
2857 *   fractional_flow_derivative(saturation_values[q], viscosity);
2858 *  
2859 *   max_velocity_times_dF_dS =
2860 *   std::max(max_velocity_times_dF_dS, velocity.norm() * dF_dS);
2861 *   }
2862 *   }
2863 *  
2864 *   return max_velocity_times_dF_dS;
2865 *   }
2866 *  
2867 *  
2868 * @endcode
2869 *
2870 *
2871 * <a name="step_43-TwoPhaseFlowProblemdimget_extrapolated_saturation_range"></a>
2872 * <h4>TwoPhaseFlowProblem<dim>::get_extrapolated_saturation_range</h4>
2873 *
2874
2875 *
2876 * For computing the stabilization term, we need to know the range of the
2877 * saturation variable. Unlike in @ref step_31 "step-31", this range is trivially bounded
2878 * by the interval @f$[0,1]@f$ but we can do a bit better by looping over a
2879 * collection of quadrature points and seeing what the values are there. If
2880 * we can, i.e., if there are at least two timesteps around, we can even
2881 * take the values extrapolated to the next time step.
2882 *
2883
2884 *
2885 * As before, the function is taken with minimal modifications from @ref step_31 "step-31".
2886 *
2887 * @code
2888 *   template <int dim>
2889 *   std::pair<double, double>
2890 *   TwoPhaseFlowProblem<dim>::get_extrapolated_saturation_range() const
2891 *   {
2892 *   const QGauss<dim> quadrature_formula(saturation_degree + 2);
2893 *   const unsigned int n_q_points = quadrature_formula.size();
2894 *  
2895 *   FEValues<dim> fe_values(saturation_fe, quadrature_formula, update_values);
2896 *   std::vector<double> old_saturation_values(n_q_points);
2897 *   std::vector<double> old_old_saturation_values(n_q_points);
2898 *  
2899 *   if (timestep_number != 0)
2900 *   {
2901 *   double min_saturation = std::numeric_limits<double>::max(),
2902 *   max_saturation = -std::numeric_limits<double>::max();
2903 *  
2904 *   for (const auto &cell : saturation_dof_handler.active_cell_iterators())
2905 *   {
2906 *   fe_values.reinit(cell);
2907 *   fe_values.get_function_values(old_saturation_solution,
2908 *   old_saturation_values);
2909 *   fe_values.get_function_values(old_old_saturation_solution,
2910 *   old_old_saturation_values);
2911 *  
2912 *   for (unsigned int q = 0; q < n_q_points; ++q)
2913 *   {
2914 *   const double saturation =
2915 *   (1. + time_step / old_time_step) * old_saturation_values[q] -
2916 *   time_step / old_time_step * old_old_saturation_values[q];
2917 *  
2918 *   min_saturation = std::min(min_saturation, saturation);
2919 *   max_saturation = std::max(max_saturation, saturation);
2920 *   }
2921 *   }
2922 *  
2923 *   return std::make_pair(min_saturation, max_saturation);
2924 *   }
2925 *   else
2926 *   {
2927 *   double min_saturation = std::numeric_limits<double>::max(),
2928 *   max_saturation = -std::numeric_limits<double>::max();
2929 *  
2930 *   for (const auto &cell : saturation_dof_handler.active_cell_iterators())
2931 *   {
2932 *   fe_values.reinit(cell);
2933 *   fe_values.get_function_values(old_saturation_solution,
2934 *   old_saturation_values);
2935 *  
2936 *   for (unsigned int q = 0; q < n_q_points; ++q)
2937 *   {
2938 *   const double saturation = old_saturation_values[q];
2939 *  
2940 *   min_saturation = std::min(min_saturation, saturation);
2941 *   max_saturation = std::max(max_saturation, saturation);
2942 *   }
2943 *   }
2944 *  
2945 *   return std::make_pair(min_saturation, max_saturation);
2946 *   }
2947 *   }
2948 *  
2949 *  
2950 *  
2951 * @endcode
2952 *
2953 *
2954 * <a name="step_43-TwoPhaseFlowProblemdimcompute_viscosity"></a>
2955 * <h4>TwoPhaseFlowProblem<dim>::compute_viscosity</h4>
2956 *
2957
2958 *
2959 * The final tool function is used to compute the artificial viscosity on a
2960 * given cell. This isn't particularly complicated if you have the formula
2961 * for it in front of you, and looking at the implementation in @ref step_31 "step-31". The
2962 * major difference to that tutorial program is that the velocity here is
2963 * not simply @f$\mathbf u@f$ but @f$\mathbf u F'(S)@f$ and some of the formulas
2964 * need to be adjusted accordingly.
2965 *
2966 * @code
2967 *   template <int dim>
2968 *   double TwoPhaseFlowProblem<dim>::compute_viscosity(
2969 *   const std::vector<double> &old_saturation,
2970 *   const std::vector<double> &old_old_saturation,
2971 *   const std::vector<Tensor<1, dim>> &old_saturation_grads,
2972 *   const std::vector<Tensor<1, dim>> &old_old_saturation_grads,
2973 *   const std::vector<Vector<double>> &present_darcy_values,
2974 *   const double global_max_u_F_prime,
2975 *   const double global_S_variation,
2976 *   const double cell_diameter) const
2977 *   {
2978 *   const double beta = .4 * dim;
2979 *   const double alpha = 1;
2980 *  
2981 *   if (global_max_u_F_prime == 0)
2982 *   return 5e-3 * cell_diameter;
2983 *  
2984 *   const unsigned int n_q_points = old_saturation.size();
2985 *  
2986 *   double max_residual = 0;
2987 *   double max_velocity_times_dF_dS = 0;
2988 *  
2989 *   const bool use_dF_dS = true;
2990 *  
2991 *   for (unsigned int q = 0; q < n_q_points; ++q)
2992 *   {
2993 *   Tensor<1, dim> u;
2994 *   for (unsigned int d = 0; d < dim; ++d)
2995 *   u[d] = present_darcy_values[q](d);
2996 *  
2997 *   const double dS_dt = porosity *
2998 *   (old_saturation[q] - old_old_saturation[q]) /
2999 *   old_time_step;
3000 *  
3001 *   const double dF_dS = fractional_flow_derivative(
3002 *   (old_saturation[q] + old_old_saturation[q]) / 2.0, viscosity);
3003 *  
3004 *   const double u_grad_S =
3005 *   u * dF_dS * (old_saturation_grads[q] + old_old_saturation_grads[q]) /
3006 *   2.0;
3007 *  
3008 *   const double residual =
3009 *   std::abs((dS_dt + u_grad_S) *
3010 *   std::pow((old_saturation[q] + old_old_saturation[q]) / 2,
3011 *   alpha - 1.));
3012 *  
3013 *   max_residual = std::max(residual, max_residual);
3014 *   max_velocity_times_dF_dS =
3015 *   std::max(std::sqrt(u * u) * (use_dF_dS ? std::max(dF_dS, 1.) : 1),
3016 *   max_velocity_times_dF_dS);
3017 *   }
3018 *  
3019 *   const double c_R = 1.0;
3020 *   const double global_scaling = c_R * porosity *
3021 *   (global_max_u_F_prime)*global_S_variation /
3022 *   std::pow(global_Omega_diameter, alpha - 2.);
3023 *  
3024 *   return (beta *
3025 *   (max_velocity_times_dF_dS)*std::min(cell_diameter,
3026 *   std::pow(cell_diameter, alpha) *
3027 *   max_residual /
3028 *   global_scaling));
3029 *   }
3030 *  
3031 *  
3032 * @endcode
3033 *
3034 *
3035 * <a name="step_43-TwoPhaseFlowProblemdimrun"></a>
3036 * <h3>TwoPhaseFlowProblem<dim>::run</h3>
3037 *
3038
3039 *
3040 * This function is, besides <code>solve()</code>, the primary function of
3041 * this program as it controls the time iteration as well as when the
3042 * solution is written into output files and when to do mesh refinement.
3043 *
3044
3045 *
3046 * With the exception of the startup code that loops back to the beginning
3047 * of the function through the <code>goto start_time_iteration</code> label,
3048 * everything should be relatively straightforward. In any case, it mimics
3049 * the corresponding function in @ref step_31 "step-31".
3050 *
3051 * @code
3052 *   template <int dim>
3053 *   void TwoPhaseFlowProblem<dim>::run()
3054 *   {
3055 *   const unsigned int initial_refinement = (dim == 2 ? 5 : 2);
3056 *   const unsigned int n_pre_refinement_steps = (dim == 2 ? 3 : 2);
3057 *  
3058 *  
3059 *   GridGenerator::hyper_cube(triangulation, 0, 1);
3060 *   triangulation.refine_global(initial_refinement);
3061 *   global_Omega_diameter = GridTools::diameter(triangulation);
3062 *  
3063 *   setup_dofs();
3064 *  
3065 *   unsigned int pre_refinement_step = 0;
3066 *  
3067 *   start_time_iteration:
3068 *  
3069 *   VectorTools::project(saturation_dof_handler,
3070 *   saturation_constraints,
3071 *   QGauss<dim>(saturation_degree + 2),
3072 *   SaturationInitialValues<dim>(),
3073 *   old_saturation_solution);
3074 *  
3075 *   time_step = old_time_step = 0;
3076 *   current_macro_time_step = old_macro_time_step = 0;
3077 *  
3078 *   time = 0;
3079 *  
3080 *   do
3081 *   {
3082 *   std::cout << "Timestep " << timestep_number << ": t=" << time
3083 *   << ", dt=" << time_step << std::endl;
3084 *  
3085 *   solve();
3086 *  
3087 *   std::cout << std::endl;
3088 *  
3089 *   if (timestep_number % 200 == 0)
3090 *   output_results();
3091 *  
3092 *   if (timestep_number % 25 == 0)
3093 *   refine_mesh(initial_refinement,
3094 *   initial_refinement + n_pre_refinement_steps);
3095 *  
3096 *   if ((timestep_number == 0) &&
3097 *   (pre_refinement_step < n_pre_refinement_steps))
3098 *   {
3099 *   ++pre_refinement_step;
3100 *   goto start_time_iteration;
3101 *   }
3102 *  
3103 *   time += time_step;
3104 *   ++timestep_number;
3105 *  
3106 *   old_old_saturation_solution = old_saturation_solution;
3107 *   old_saturation_solution = saturation_solution;
3108 *   }
3109 *   while (time <= end_time);
3110 *   }
3111 *   } // namespace Step43
3112 *  
3113 *  
3114 *  
3115 * @endcode
3116 *
3117 *
3118 * <a name="step_43-Thecodemaincodefunction"></a>
3119 * <h3>The <code>main()</code> function</h3>
3120 *
3121
3122 *
3123 * The main function looks almost the same as in all other programs. The need
3124 * to initialize the MPI subsystem for a program that uses Trilinos -- even
3125 * for programs that do not actually run in parallel -- is explained in
3126 * @ref step_31 "step-31".
3127 *
3128 * @code
3129 *   int main(int argc, char *argv[])
3130 *   {
3131 *   try
3132 *   {
3133 *   using namespace dealii;
3134 *   using namespace Step43;
3135 *  
3136 *   Utilities::MPI::MPI_InitFinalize mpi_initialization(
3137 *   argc, argv, numbers::invalid_unsigned_int);
3138 *  
3139 * @endcode
3140 *
3141 * This program can only be run in serial. Otherwise, throw an exception.
3142 *
3143 * @code
3144 *   AssertThrow(Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD) == 1,
3145 *   ExcMessage(
3146 *   "This program can only be run in serial, use ./step-43"));
3147 *  
3148 *   TwoPhaseFlowProblem<2> two_phase_flow_problem(1);
3149 *   two_phase_flow_problem.run();
3150 *   }
3151 *   catch (std::exception &exc)
3152 *   {
3153 *   std::cerr << std::endl
3154 *   << std::endl
3155 *   << "----------------------------------------------------"
3156 *   << std::endl;
3157 *   std::cerr << "Exception on processing: " << std::endl
3158 *   << exc.what() << std::endl
3159 *   << "Aborting!" << std::endl
3160 *   << "----------------------------------------------------"
3161 *   << std::endl;
3162 *  
3163 *   return 1;
3164 *   }
3165 *   catch (...)
3166 *   {
3167 *   std::cerr << std::endl
3168 *   << std::endl
3169 *   << "----------------------------------------------------"
3170 *   << std::endl;
3171 *   std::cerr << "Unknown exception!" << std::endl
3172 *   << "Aborting!" << std::endl
3173 *   << "----------------------------------------------------"
3174 *   << std::endl;
3175 *   return 1;
3176 *   }
3177 *  
3178 *   return 0;
3179 *   }
3180 * @endcode
3181<a name="step_43-Results"></a><h1>Results</h1>
3182
3183
3184
3185The output of this program is not really much different from that of
3186@ref step_21 "step-21": it solves the same problem, after all. Of more importance are
3187quantitative metrics such as the accuracy of the solution as well as
3188the time needed to compute it. These are documented in detail in the
3189two publications listed at the top of this page and we won't repeat
3190them here.
3191
3192That said, no tutorial program is complete without a couple of good
3193pictures, so here is some output of a run in 3d:
3194
3195<table align="center" class="tutorial" cellspacing="3" cellpadding="3">
3196 <tr>
3197 <td align="center">
3198 <img src="https://www.dealii.org/images/steps/developer/step-43.3d.velocity.png" alt="">
3199 <p align="center">
3200 Velocity vectors of flow through the porous medium with random
3201 permeability model. Streaming paths of high permeability and resulting
3202 high velocity are clearly visible.
3203 </p>
3204 </td>
3205 <td align="center">
3206 <img src="https://www.dealii.org/images/steps/developer/step-43.3d.streamlines.png" alt="">
3207 <p align="center">
3208 Streamlines colored by the saturation along the streamline path. Blue
3209 streamlines indicate low saturations, i.e., the flow along these
3210 streamlines must be slow or else more fluid would have been
3211 transported along them. On the other hand, green paths indicate high
3212 velocities since the fluid front has already reached further into the
3213 domain.
3214 </p>
3215 </td>
3216 </tr>
3217 <tr>
3218 <td align="center">
3219 <img src="https://www.dealii.org/images/steps/developer/step-43.3d.saturation.png" alt="">
3220 <p align="center">
3221 Streamlines with a volume rendering of the saturation, showing how far
3222 the fluid front has advanced at this time.
3223 </p>
3224 </td>
3225 <td align="center">
3226 <img src="https://www.dealii.org/images/steps/developer/step-43.3d.mesh.png" alt="">
3227 <p align="center">
3228 Surface of the mesh showing the adaptive refinement along the front.
3229 </p>
3230 </td>
3231 </tr>
3232</table>
3233
3234
3235<a name="step-43-extensions"></a>
3236<a name="step_43-Possibilitiesforextensions"></a><h3>Possibilities for extensions</h3>
3237
3238
3239The primary objection one may have to this program is that it is still too
3240slow: 3d computations on reasonably fine meshes are simply too expensive to be
3241done routinely and with reasonably quick turn-around. This is similar to the
3242situation we were in when we wrote @ref step_31 "step-31", from which this program has taken
3243much inspiration. The solution is similar as it was there as well: We need to
3244parallelize the program in a way similar to how we derived @ref step_32 "step-32" out of
3245@ref step_31 "step-31". In fact, all of the techniques used in @ref step_32 "step-32" would be transferable
3246to this program as well, making the program run on dozens or hundreds of
3247processors immediately.
3248
3249A different direction is to make the program more relevant to many other
3250porous media applications. Specifically, one avenue is to go to the primary
3251user of porous media flow simulators, namely the oil industry. There,
3252applications in this area are dominated by multiphase flow (i.e., more than
3253the two phases we have here), and the reactions they may have with each other
3254(or any other way phases may exchange mass, such as through dissolution in and
3255bubbling out of gas from the oil phase). Furthermore, the presence of gas
3256often leads to compressibility effects of the fluid. Jointly, these effects
3257are typically formulated in the widely-used "black oil model". True reactions
3258between multiple phases also play a role in oil reservoir modeling when
3259considering controlled burns of oil in the reservoir to raise pressure and
3260temperature. These are much more complex problems, though, and left for future
3261projects.
3262
3263Finally, from a mathematical perspective, we have derived the
3264criterion for re-computing the velocity/pressure solution at a given
3265time step under the assumption that we want to compare the solution we
3266would get at the current time step with that computed the last time we
3267actually solved this system. However, in the program, whenever we did
3268not re-compute the solution, we didn't just use the previously
3269computed solution but instead extrapolated from the previous two times
3270we solved the system. Consequently, the criterion was pessimistically
3271stated: what we should really compare is the solution we would get at
3272the current time step with the extrapolated one. Re-stating the
3273theorem in this regard is left as an exercise.
3274
3275There are also other ways to extend the mathematical foundation of
3276this program; for example, one may say that it isn't the velocity we
3277care about, but in fact the saturation. Thus, one may ask whether the
3278criterion we use here to decide whether @f$\mathbf u@f$ needs to be
3279recomputed is appropriate; one may, for example, suggest that it is
3280also important to decide whether (and by how much) a wrong velocity
3281field in fact affects the solution of the saturation equation. This
3282would then naturally lead to a sensitivity analysis.
3283
3284From an algorithmic viewpoint, we have here used a criterion for refinement
3285that is often used in engineering, namely by looking at the gradient of
3286the solution. However, if you inspect the solution, you will find that
3287it quickly leads to refinement almost everywhere, even in regions where it
3288is clearly not necessary: frequently used therefore does not need to imply
3289that it is a useful criterion to begin with. On the other hand, replacing
3290this criterion by a different and better one should not be very difficult.
3291For example, the KellyErrorEstimator class used in many other programs
3292should certainly be applicable to the current problem as well.
3293 *
3294 *
3295<a name="step_43-PlainProg"></a>
3296<h1> The plain program</h1>
3297@include "step-43.cc"
3298*/
void condense(SparsityPattern &sparsity) const
virtual RangeNumberType value(const Point< dim > &p, const unsigned int component=0) const
virtual void vector_value(const Point< dim > &p, Vector< RangeNumberType > &values) const
Definition point.h:111
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< value_type > &values) const
Point< 2 > second
Definition grid_out.cc:4624
Point< 2 > first
Definition grid_out.cc:4623
#define AssertDimension(dim1, dim2)
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
void make_hanging_node_constraints(const DoFHandler< dim, spacedim > &dof_handler, AffineConstraints< number > &constraints)
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity_pattern, const AffineConstraints< number > &constraints={}, const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
void make_zero_boundary_constraints(const DoFHandler< dim, spacedim > &dof, const types::boundary_id boundary_id, AffineConstraints< number > &zero_boundary_constraints, const ComponentMask &component_mask={})
IndexSet complete_index_set(const IndexSet::size_type N)
Definition index_set.h:1204
std::vector< index_type > data
Definition mpi.cc:735
CGAL::Exact_predicates_exact_constructions_kernel_with_sqrt K
void component_wise(DoFHandler< dim, spacedim > &dof_handler, const std::vector< unsigned int > &target_component=std::vector< unsigned int >())
void Cuthill_McKee(DoFHandler< dim, spacedim > &dof_handler, const bool reversed_numbering=false, const bool use_constraints=false, const std::vector< types::global_dof_index > &starting_indices=std::vector< types::global_dof_index >())
void random(DoFHandler< dim, spacedim > &dof_handler)
std::vector< types::global_dof_index > count_dofs_per_fe_block(const DoFHandler< dim, spacedim > &dof, const std::vector< unsigned int > &target_block=std::vector< unsigned int >())
void extrapolate(const DoFHandler< dim, spacedim > &dof1, const InVector &z1, const DoFHandler< dim, spacedim > &dof2, OutVector &z2)
double minimal_cell_diameter(const Triangulation< dim, spacedim > &triangulation, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()))
double volume(const Triangulation< dim, spacedim > &tria)
@ matrix
Contents is actually a matrix.
void apply_boundary_values(const std::map< types::global_dof_index, number > &boundary_values, SparseMatrix< number > &matrix, Vector< number > &solution, Vector< number > &right_hand_side, const bool eliminate_columns=true)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition utilities.cc:191
std::string escape(const std::string &input, const PatternBase::OutputStyle style)
Definition patterns.cc:52
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
Tensor< 2, dim, Number > F(const Tensor< 2, dim, Number > &Grad_u)
VectorType::value_type * begin(VectorType &V)
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)
int(&) functions(const void *v1, const void *v2)
::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 > sin(const ::VectorizedArray< Number, width > &)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation