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