Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
step-52.h
Go to the documentation of this file.
1 ,
710  * const double tau,
711  * const Vector<double> &y)
712  * {
713  * SparseDirectUMFPACK inverse_mass_minus_tau_Jacobian;
714  *
715  * mass_minus_tau_Jacobian.copy_from(mass_matrix);
716  * mass_minus_tau_Jacobian.add(-tau, system_matrix);
717  *
718  * inverse_mass_minus_tau_Jacobian.initialize(mass_minus_tau_Jacobian);
719  *
720  * Vector<double> tmp(dof_handler.n_dofs());
721  * mass_matrix.vmult(tmp, y);
722  *
723  * Vector<double> result(y);
724  * inverse_mass_minus_tau_Jacobian.vmult(result, tmp);
725  *
726  * return result;
727  * }
728  *
729  *
730  *
731  * @endcode
732  *
733  *
734  * <a name="codeDiffusionoutput_resultscode"></a>
735  * <h4><code>Diffusion::output_results</code></h4>
736  *
737 
738  *
739  * The following function then outputs the solution in vtu files indexed by
740  * the number of the time step and the name of the time stepping method. Of
741  * course, the (exact) result should really be the same for all time
742  * stepping method, but the output here at least allows us to compare them.
743  *
744  * @code
745  * void Diffusion::output_results(const double time,
746  * const unsigned int time_step,
747  * TimeStepping::runge_kutta_method method) const
748  * {
749  * std::string method_name;
750  *
751  * switch (method)
752  * {
754  * {
755  * method_name = "forward_euler";
756  * break;
757  * }
759  * {
760  * method_name = "rk3";
761  * break;
762  * }
764  * {
765  * method_name = "rk4";
766  * break;
767  * }
769  * {
770  * method_name = "backward_euler";
771  * break;
772  * }
774  * {
775  * method_name = "implicit_midpoint";
776  * break;
777  * }
779  * {
780  * method_name = "sdirk";
781  * break;
782  * }
784  * {
785  * method_name = "heun_euler";
786  * break;
787  * }
789  * {
790  * method_name = "bocacki_shampine";
791  * break;
792  * }
793  * case TimeStepping::DOPRI:
794  * {
795  * method_name = "dopri";
796  * break;
797  * }
798  * case TimeStepping::FEHLBERG:
799  * {
800  * method_name = "fehlberg";
801  * break;
802  * }
804  * {
805  * method_name = "cash_karp";
806  * break;
807  * }
808  * default:
809  * {
810  * break;
811  * }
812  * }
813  *
814  * DataOut<2> data_out;
815  *
816  * data_out.attach_dof_handler(dof_handler);
817  * data_out.add_data_vector(solution, "solution");
818  *
819  * data_out.build_patches();
820  *
821  * data_out.set_flags(DataOutBase::VtkFlags(time, time_step));
822  *
823  * const std::string filename = "solution_" + method_name + "-" +
824  * Utilities::int_to_string(time_step, 3) +
825  * ".vtu";
826  * std::ofstream output(filename);
827  * data_out.write_vtu(output);
828  *
829  * static std::vector<std::pair<double, std::string>> times_and_names;
830  *
831  * static std::string method_name_prev = "";
832  * static std::string pvd_filename;
833  * if (method_name_prev != method_name)
834  * {
835  * times_and_names.clear();
836  * method_name_prev = method_name;
837  * pvd_filename = "solution_" + method_name + ".pvd";
838  * }
839  * times_and_names.emplace_back(time, filename);
840  * std::ofstream pvd_output(pvd_filename);
841  * DataOutBase::write_pvd_record(pvd_output, times_and_names);
842  * }
843  *
844  *
845  * @endcode
846  *
847  *
848  * <a name="codeDiffusionexplicit_methodcode"></a>
849  * <h4><code>Diffusion::explicit_method</code></h4>
850  *
851 
852  *
853  * This function is the driver for all the explicit methods. At the
854  * top it initializes the time stepping and the solution (by setting
855  * it to zero and then ensuring that boundary value and hanging node
856  * constraints are respected; of course, with the mesh we use here,
857  * hanging node constraints are not in fact an issue). It then calls
858  * <code>evolve_one_time_step</code> which performs one time step.
859  *
860 
861  *
862  * For explicit methods, <code>evolve_one_time_step</code> needs to
863  * evaluate @f$M^{-1}(f(t,y))@f$, i.e, it needs
864  * <code>evaluate_diffusion</code>. Because
865  * <code>evaluate_diffusion</code> is a member function, it needs to
866  * be bound to <code>this</code>. After each evolution step, we
867  * again apply the correct boundary values and hanging node
868  * constraints.
869  *
870 
871  *
872  * Finally, the solution is output
873  * every 10 time steps.
874  *
875  * @code
876  * void Diffusion::explicit_method(const TimeStepping::runge_kutta_method method,
877  * const unsigned int n_time_steps,
878  * const double initial_time,
879  * const double final_time)
880  * {
881  * const double time_step =
882  * (final_time - initial_time) / static_cast<double>(n_time_steps);
883  * double time = initial_time;
884  *
885  * solution = 0.;
886  * constraint_matrix.distribute(solution);
887  *
888  * TimeStepping::ExplicitRungeKutta<Vector<double>> explicit_runge_kutta(
889  * method);
890  * output_results(time, 0, method);
891  * for (unsigned int i = 0; i < n_time_steps; ++i)
892  * {
893  * time = explicit_runge_kutta.evolve_one_time_step(
894  * [this](const double time, const Vector<double> &y) {
895  * return this->evaluate_diffusion(time, y);
896  * },
897  * time,
898  * time_step,
899  * solution);
900  *
901  * constraint_matrix.distribute(solution);
902  *
903  * if ((i + 1) % 10 == 0)
904  * output_results(time, i + 1, method);
905  * }
906  * }
907  *
908  *
909  *
910  * @endcode
911  *
912  *
913  * <a name="codeDiffusionimplicit_methodcode"></a>
914  * <h4><code>Diffusion::implicit_method</code></h4>
915  * This function is equivalent to <code>explicit_method</code> but for
916  * implicit methods. When using implicit methods, we need to evaluate
917  * @f$M^{-1}(f(t,y))@f$ and @f$\left(I-\tau M^{-1} \frac{\partial f(t,y)}{\partial
918  * y}\right)^{-1}@f$ for which we use the two member functions previously
919  * introduced.
920  *
921  * @code
922  * void Diffusion::implicit_method(const TimeStepping::runge_kutta_method method,
923  * const unsigned int n_time_steps,
924  * const double initial_time,
925  * const double final_time)
926  * {
927  * const double time_step =
928  * (final_time - initial_time) / static_cast<double>(n_time_steps);
929  * double time = initial_time;
930  *
931  * solution = 0.;
932  * constraint_matrix.distribute(solution);
933  *
934  * TimeStepping::ImplicitRungeKutta<Vector<double>> implicit_runge_kutta(
935  * method);
936  * output_results(time, 0, method);
937  * for (unsigned int i = 0; i < n_time_steps; ++i)
938  * {
939  * time = implicit_runge_kutta.evolve_one_time_step(
940  * [this](const double time, const Vector<double> &y) {
941  * return this->evaluate_diffusion(time, y);
942  * },
943  * [this](const double time, const double tau, const Vector<double> &y) {
944  * return this->id_minus_tau_J_inverse(time, tau, y);
945  * },
946  * time,
947  * time_step,
948  * solution);
949  *
950  * constraint_matrix.distribute(solution);
951  *
952  * if ((i + 1) % 10 == 0)
953  * output_results(time, i + 1, method);
954  * }
955  * }
956  *
957  *
958  *
959  * @endcode
960  *
961  *
962  * <a name="codeDiffusionembedded_explicit_methodcode"></a>
963  * <h4><code>Diffusion::embedded_explicit_method</code></h4>
964  * This function is the driver for the embedded explicit methods. It requires
965  * more parameters:
966  * - coarsen_param: factor multiplying the current time step when the error
967  * is below the threshold.
968  * - refine_param: factor multiplying the current time step when the error
969  * is above the threshold.
970  * - min_delta: smallest time step acceptable.
971  * - max_delta: largest time step acceptable.
972  * - refine_tol: threshold above which the time step is refined.
973  * - coarsen_tol: threshold below which the time step is coarsen.
974  * Embedded methods use a guessed time step. If the error using this time step
975  * is too large, the time step will be reduced. If the error is below the
976  * threshold, a larger time step will be tried for the next time step.
977  * <code>delta_t_guess</code> is the guessed time step produced by the
978  * embedded method.
979  *
980  * @code
981  * unsigned int Diffusion::embedded_explicit_method(
982  * const TimeStepping::runge_kutta_method method,
983  * const unsigned int n_time_steps,
984  * const double initial_time,
985  * const double final_time)
986  * {
987  * double time_step =
988  * (final_time - initial_time) / static_cast<double>(n_time_steps);
989  * double time = initial_time;
990  * const double coarsen_param = 1.2;
991  * const double refine_param = 0.8;
992  * const double min_delta = 1e-8;
993  * const double max_delta = 10 * time_step;
994  * const double refine_tol = 1e-1;
995  * const double coarsen_tol = 1e-5;
996  *
997  * solution = 0.;
998  * constraint_matrix.distribute(solution);
999  *
1001  * embedded_explicit_runge_kutta(method,
1002  * coarsen_param,
1003  * refine_param,
1004  * min_delta,
1005  * max_delta,
1006  * refine_tol,
1007  * coarsen_tol);
1008  * output_results(time, 0, method);
1009  *
1010  * @endcode
1011  *
1012  * Now for the time loop. The last time step is chosen such that the final
1013  * time is exactly reached.
1014  *
1015  * @code
1016  * unsigned int n_steps = 0;
1017  * while (time < final_time)
1018  * {
1019  * if (time + time_step > final_time)
1020  * time_step = final_time - time;
1021  *
1022  * time = embedded_explicit_runge_kutta.evolve_one_time_step(
1023  * [this](const double time, const Vector<double> &y) {
1024  * return this->evaluate_diffusion(time, y);
1025  * },
1026  * time,
1027  * time_step,
1028  * solution);
1029  *
1030  * constraint_matrix.distribute(solution);
1031  *
1032  * if ((n_steps + 1) % 10 == 0)
1033  * output_results(time, n_steps + 1, method);
1034  *
1035  * time_step = embedded_explicit_runge_kutta.get_status().delta_t_guess;
1036  * ++n_steps;
1037  * }
1038  *
1039  * return n_steps;
1040  * }
1041  *
1042  *
1043  *
1044  * @endcode
1045  *
1046  *
1047  * <a name="codeDiffusionruncode"></a>
1048  * <h4><code>Diffusion::run</code></h4>
1049  *
1050 
1051  *
1052  * The following is the main function of the program. At the top, we create
1053  * the grid (a [0,5]x[0,5] square) and refine it four times to get a mesh
1054  * that has 16 by 16 cells, for a total of 256. We then set the boundary
1055  * indicator to 1 for those parts of the boundary where @f$x=0@f$ and @f$x=5@f$.
1056  *
1057  * @code
1058  * void Diffusion::run()
1059  * {
1061  * triangulation.refine_global(4);
1062  *
1063  * for (const auto &cell : triangulation.active_cell_iterators())
1064  * for (const auto &face : cell->face_iterators())
1065  * if (face->at_boundary())
1066  * {
1067  * if ((face->center()[0] == 0.) || (face->center()[0] == 5.))
1068  * face->set_boundary_id(1);
1069  * else
1070  * face->set_boundary_id(0);
1071  * }
1072  *
1073  * @endcode
1074  *
1075  * Next, we set up the linear systems and fill them with content so that
1076  * they can be used throughout the time stepping process:
1077  *
1078  * @code
1079  * setup_system();
1080  *
1081  * assemble_system();
1082  *
1083  * @endcode
1084  *
1085  * Finally, we solve the diffusion problem using several of the
1086  * Runge-Kutta methods implemented in namespace TimeStepping, each time
1087  * outputting the error at the end time. (As explained in the
1088  * introduction, since the exact solution is zero at the final time, the
1089  * error equals the numerical solution and can be computed by just taking
1090  * the @f$l_2@f$ norm of the solution vector.)
1091  *
1092  * @code
1093  * unsigned int n_steps = 0;
1094  * const unsigned int n_time_steps = 200;
1095  * const double initial_time = 0.;
1096  * const double final_time = 10.;
1097  *
1098  * std::cout << "Explicit methods:" << std::endl;
1099  * explicit_method(TimeStepping::FORWARD_EULER,
1100  * n_time_steps,
1101  * initial_time,
1102  * final_time);
1103  * std::cout << " Forward Euler: error=" << solution.l2_norm()
1104  * << std::endl;
1105  *
1106  * explicit_method(TimeStepping::RK_THIRD_ORDER,
1107  * n_time_steps,
1108  * initial_time,
1109  * final_time);
1110  * std::cout << " Third order Runge-Kutta: error=" << solution.l2_norm()
1111  * << std::endl;
1112  *
1113  * explicit_method(TimeStepping::RK_CLASSIC_FOURTH_ORDER,
1114  * n_time_steps,
1115  * initial_time,
1116  * final_time);
1117  * std::cout << " Fourth order Runge-Kutta: error=" << solution.l2_norm()
1118  * << std::endl;
1119  * std::cout << std::endl;
1120  *
1121  *
1122  * std::cout << "Implicit methods:" << std::endl;
1123  * implicit_method(TimeStepping::BACKWARD_EULER,
1124  * n_time_steps,
1125  * initial_time,
1126  * final_time);
1127  * std::cout << " Backward Euler: error=" << solution.l2_norm()
1128  * << std::endl;
1129  *
1130  * implicit_method(TimeStepping::IMPLICIT_MIDPOINT,
1131  * n_time_steps,
1132  * initial_time,
1133  * final_time);
1134  * std::cout << " Implicit Midpoint: error=" << solution.l2_norm()
1135  * << std::endl;
1136  *
1137  * implicit_method(TimeStepping::CRANK_NICOLSON,
1138  * n_time_steps,
1139  * initial_time,
1140  * final_time);
1141  * std::cout << " Crank-Nicolson: error=" << solution.l2_norm()
1142  * << std::endl;
1143  *
1144  * implicit_method(TimeStepping::SDIRK_TWO_STAGES,
1145  * n_time_steps,
1146  * initial_time,
1147  * final_time);
1148  * std::cout << " SDIRK: error=" << solution.l2_norm()
1149  * << std::endl;
1150  * std::cout << std::endl;
1151  *
1152  *
1153  * std::cout << "Embedded explicit methods:" << std::endl;
1154  * n_steps = embedded_explicit_method(TimeStepping::HEUN_EULER,
1155  * n_time_steps,
1156  * initial_time,
1157  * final_time);
1158  * std::cout << " Heun-Euler: error=" << solution.l2_norm()
1159  * << std::endl;
1160  * std::cout << " steps performed=" << n_steps << std::endl;
1161  *
1162  * n_steps = embedded_explicit_method(TimeStepping::BOGACKI_SHAMPINE,
1163  * n_time_steps,
1164  * initial_time,
1165  * final_time);
1166  * std::cout << " Bogacki-Shampine: error=" << solution.l2_norm()
1167  * << std::endl;
1168  * std::cout << " steps performed=" << n_steps << std::endl;
1169  *
1170  * n_steps = embedded_explicit_method(TimeStepping::DOPRI,
1171  * n_time_steps,
1172  * initial_time,
1173  * final_time);
1174  * std::cout << " Dopri: error=" << solution.l2_norm()
1175  * << std::endl;
1176  * std::cout << " steps performed=" << n_steps << std::endl;
1177  *
1178  * n_steps = embedded_explicit_method(TimeStepping::FEHLBERG,
1179  * n_time_steps,
1180  * initial_time,
1181  * final_time);
1182  * std::cout << " Fehlberg: error=" << solution.l2_norm()
1183  * << std::endl;
1184  * std::cout << " steps performed=" << n_steps << std::endl;
1185  *
1186  * n_steps = embedded_explicit_method(TimeStepping::CASH_KARP,
1187  * n_time_steps,
1188  * initial_time,
1189  * final_time);
1190  * std::cout << " Cash-Karp: error=" << solution.l2_norm()
1191  * << std::endl;
1192  * std::cout << " steps performed=" << n_steps << std::endl;
1193  * }
1194  * } // namespace Step52
1195  *
1196  *
1197  *
1198  * @endcode
1199  *
1200  *
1201  * <a name="Thecodemaincodefunction"></a>
1202  * <h3>The <code>main()</code> function</h3>
1203  *
1204 
1205  *
1206  * The following <code>main</code> function is similar to previous examples
1207  * and need not be commented on.
1208  *
1209  * @code
1210  * int main()
1211  * {
1212  * try
1213  * {
1214  * Step52::Diffusion diffusion;
1215  * diffusion.run();
1216  * }
1217  * catch (std::exception &exc)
1218  * {
1219  * std::cerr << std::endl
1220  * << std::endl
1221  * << "----------------------------------------------------"
1222  * << std::endl;
1223  * std::cerr << "Exception on processing: " << std::endl
1224  * << exc.what() << std::endl
1225  * << "Aborting!" << std::endl
1226  * << "----------------------------------------------------"
1227  * << std::endl;
1228  * return 1;
1229  * }
1230  * catch (...)
1231  * {
1232  * std::cerr << std::endl
1233  * << std::endl
1234  * << "----------------------------------------------------"
1235  * << std::endl;
1236  * std::cerr << "Unknown exception!" << std::endl
1237  * << "Aborting!" << std::endl
1238  * << "----------------------------------------------------"
1239  * << std::endl;
1240  * return 1;
1241  * };
1242  *
1243  * return 0;
1244  * }
1245  * @endcode
1246 <a name="Results"></a><h1>Results</h1>
1247 
1248 
1249 The point of this program is less to show particular results, but instead to
1250 show how it is done. This we have already demonstrated simply by discussing
1251 the code above. Consequently, the output the program yields is relatively
1252 sparse and consists only of the console output and the solutions given in VTU
1253 format for visualization.
1254 
1255 The console output contains both errors and, for some of the methods, the
1256 number of steps they performed:
1257 @code
1258 Explicit methods:
1259 Forward Euler: error=1.00883
1260 Third order Runge-Kutta: error=0.000227982
1261 Fourth order Runge-Kutta: error=1.90541e-06
1262 
1263 Implicit methods:
1264 Backward Euler: error=1.03428
1265 Implicit Midpoint: error=0.00862702
1266 Crank-Nicolson: error=0.00862675
1267 SDIRK: error=0.0042349
1268 
1269 Embedded %explicit methods:
1270 Heun-Euler: error=0.0073012
1271  steps performed=284
1272 Bogacki-Shampine: error=0.000403281
1273  steps performed=181
1274 Dopri: error=0.0165485
1275  steps performed=119
1276 Fehlberg: error=0.00104926
1277  steps performed=106
1278 Cash-Karp: error=8.59366e-07
1279  steps performed=107
1280 
1281 @endcode
1282 
1283 As expected the higher order methods give (much) more accurate solutions. We
1284 also see that the (rather inaccurate) Heun-Euler method increased the number of
1285 time steps in order to satisfy the tolerance. On the other hand, the other
1286 embedded methods used a lot less time steps than what was prescribed.
1287  *
1288  *
1289 <a name="PlainProg"></a>
1290 <h1> The plain program</h1>
1291 @include "step-52.cc"
1292 */
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:443
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
void write_pvd_record(std::ostream &out, const std::vector< std::pair< double, std::string >> &times_and_names)
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim >>> &Du)
Definition: divergence.h:548
static const types::blas_int one
virtual void build_patches(const unsigned int n_subdivisions=0)
Definition: data_out.cc:1071
void vmult(Vector< double > &dst, const Vector< double > &src) const
void mass_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
Definition: l2.h:63
auto apply(F &&fn, Tuple &&t) -> decltype(apply_impl(std::forward< F >(fn), std::forward< Tuple >(t), std_cxx14::make_index_sequence< std::tuple_size< typename std::remove_reference< Tuple >::type >::value >()))
Definition: tuple.h:40
void add_data_vector(const VectorType &data, const std::vector< std::string > &names, const DataVectorType type=type_automatic, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation=std::vector< DataComponentInterpretation::DataComponentInterpretation >())
VectorType::value_type * end(VectorType &V)
void coarsen(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold)
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:474
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition: utilities.cc:188
void initialize(const SparsityPattern &sparsity_pattern)
void hyper_cube(Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false)
void run(const std::vector< std::vector< Iterator >> &colored_iterators, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length=2 *MultithreadInfo::n_threads(), const unsigned int chunk_size=8)
Definition: work_stream.h:1185
static const types::blas_int zero
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
static const bool value
void refine(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold, const unsigned int max_to_mark=numbers::invalid_unsigned_int)
int(&) functions(const void *v1, const void *v2)