513 *
const double t = this->
get_time();
519 *
const double m = 0.5;
520 *
const double c1 = 0.;
521 *
const double c2 = 0.;
530 *
const double lambda = 1.;
531 *
const double a0 = 1.;
532 *
const double s = 1.;
533 *
const double arg = p[0] *
std::cos(theta) +
543 *
const double tau = 1.;
544 *
const double c0 = 1.;
545 *
const double s = 1.;
546 *
const double arg = p[0] *
std::cos(theta) +
562 * In the
second part of
this section, we provide the
initial conditions. We
563 * are lazy (and cautious) and don
't want to implement the same functions as
564 * above a second time. Rather, if we are queried for initial conditions, we
565 * create an object <code>ExactSolution</code>, set it to the correct time,
566 * and let it compute whatever values the exact solution has at that time:
570 * class InitialValues : public Function<dim>
573 * InitialValues(const unsigned int n_components = 1, const double time = 0.)
574 * : Function<dim>(n_components, time)
577 * virtual double value(const Point<dim> &p,
578 * const unsigned int component = 0) const override
580 * return ExactSolution<dim>(1, this->get_time()).value(p, component);
588 * <a name="step_25-ImplementationofthecodeSineGordonProblemcodeclass"></a>
589 * <h3>Implementation of the <code>SineGordonProblem</code> class</h3>
593 * Let's move on to the implementation of the main
class, as it implements
594 * the algorithm outlined in the introduction.
599 * <a name=
"step_25-SineGordonProblemSineGordonProblem"></a>
600 * <h4>SineGordonProblem::SineGordonProblem</h4>
604 * This is the constructor of the <code>SineGordonProblem</code>
class. It
605 * specifies the desired polynomial degree of the finite elements,
607 * object (just as in the example programs @ref step_3
"step-3" and @ref step_4
"step-4"), initializes
608 * the current or
initial time, the
final time, the time step size, and the
609 *
value of @f$\theta@f$
for the time stepping scheme. Since the solutions we
610 * compute here are time-periodic, the actual
value of the start-time
611 * doesn
't matter, and we choose it so that we start at an interesting time.
615 * Note that if we were to chose the explicit Euler time stepping scheme
616 * (@f$\theta = 0@f$), then we must pick a time step @f$k \le h@f$, otherwise the
617 * scheme is not stable and oscillations might arise in the solution. The
618 * Crank-Nicolson scheme (@f$\theta = \frac{1}{2}@f$) and the implicit Euler
619 * scheme (@f$\theta=1@f$) do not suffer from this deficiency, since they are
620 * unconditionally stable. However, even then the time step should be chosen
621 * to be on the order of @f$h@f$ in order to obtain a good solution. Since we
622 * know that our mesh results from the uniform subdivision of a rectangle,
623 * we can compute that time step easily; if we had a different domain, the
624 * technique in @ref step_24 "step-24" using GridTools::minimal_cell_diameter would work as
629 * SineGordonProblem<dim>::SineGordonProblem()
631 * , dof_handler(triangulation)
632 * , n_global_refinements(6)
634 * , final_time(2.7207)
635 * , time_step(10 * 1. / std::pow(2., 1. * n_global_refinements))
637 * , output_timestep_skip(1)
643 * <a name="step_25-SineGordonProblemmake_grid_and_dofs"></a>
644 * <h4>SineGordonProblem::make_grid_and_dofs</h4>
648 * This function creates a rectangular grid in <code>dim</code> dimensions
649 * and refines it several times. Also, all matrix and vector members of the
650 * <code>SineGordonProblem</code> class are initialized to their appropriate
651 * sizes once the degrees of freedom have been assembled. Like @ref step_24 "step-24", we
652 * use <code>MatrixCreator</code> functions to generate a mass matrix @f$M@f$
653 * and a Laplace matrix @f$A@f$ and store them in the appropriate variables for
654 * the remainder of the program's life.
658 *
void SineGordonProblem<dim>::make_grid_and_dofs()
668 * dof_handler.distribute_dofs(fe);
670 * std::cout <<
" Number of degrees of freedom: " << dof_handler.n_dofs()
675 * sparsity_pattern.copy_from(dsp);
677 * system_matrix.reinit(sparsity_pattern);
679 * laplace_matrix.reinit(sparsity_pattern);
688 * solution.reinit(dof_handler.n_dofs());
689 * solution_update.reinit(dof_handler.n_dofs());
690 * old_solution.reinit(dof_handler.n_dofs());
691 * M_x_velocity.reinit(dof_handler.n_dofs());
692 * system_rhs.reinit(dof_handler.n_dofs());
698 * <a name=
"step_25-SineGordonProblemassemble_system"></a>
699 * <h4>SineGordonProblem::assemble_system</h4>
703 * This function assembles the system
matrix and right-hand side vector
for
704 * each iteration of Newton
's method. The reader should refer to the
705 * Introduction for the explicit formulas for the system matrix and
710 * Note that during each time step, we have to add up the various
711 * contributions to the matrix and right hand sides. In contrast to @ref step_23 "step-23"
712 * and @ref step_24 "step-24", this requires assembling a few more terms, since they depend
713 * on the solution of the previous time step or previous nonlinear step. We
714 * use the functions <code>compute_nl_matrix</code> and
715 * <code>compute_nl_term</code> to do this, while the present function
716 * provides the top-level logic.
720 * void SineGordonProblem<dim>::assemble_system()
724 * First we assemble the Jacobian matrix @f$F'_h(U^{n,
l})@f$, where @f$U^{n,
l}@f$
725 * is stored in the vector <code>solution</code>
for convenience.
728 * system_matrix.copy_from(mass_matrix);
733 * compute_nl_matrix(old_solution, solution, tmp_matrix);
738 * Next we compute the right-hand side vector. This is just the
739 * combination of
matrix-vector products implied by the description of
740 * @f$-F_h(U^{n,
l})@f$ in the introduction.
748 * laplace_matrix.vmult(tmp_vector, solution);
752 * system_rhs.add(-1.0, tmp_vector);
753 * laplace_matrix.vmult(tmp_vector, old_solution);
757 * system_rhs.add(-time_step, M_x_velocity);
759 * compute_nl_term(old_solution, solution, tmp_vector);
768 * <a name=
"step_25-SineGordonProblemcompute_nl_term"></a>
769 * <h4>SineGordonProblem::compute_nl_term</h4>
773 * This function computes the vector @f$S(\cdot,\cdot)@f$, which appears in the
774 * nonlinear term in both equations of the
split formulation. This
775 * function not only simplifies the repeated computation of
this term, but
776 * it is also a fundamental part of the nonlinear iterative solver that we
777 * use when the time stepping is implicit (i.e. @f$\theta\ne 0@f$). Moreover, we
778 * must allow the function to receive as input an
"old" and a
"new"
779 * solution. These may not be the actual solutions of the problem stored in
780 * <code>old_solution</code> and <code>solution</code>, but are simply the
781 * two
functions we linearize about. For the purposes of
this function, let
782 * us call the
first two arguments @f$w_{\mathrm{old}}@f$ and @f$w_{\mathrm{
new}}@f$
783 * in the documentation of
this class below, respectively.
787 * As a side-note, it is perhaps worth investigating what order quadrature
788 * formula is best suited
for this type of integration. Since @f$\sin(\cdot)@f$
789 * is not a polynomial, there are probably no quadrature formulas that can
790 * integrate these terms exactly. It is usually sufficient to just make sure
791 * that the right hand side is integrated up to the same order of accuracy
792 * as the discretization scheme is, but it may be possible to improve on the
793 *
constant in the asymptotic statement of convergence by choosing a more
794 * accurate quadrature formula.
798 *
void SineGordonProblem<dim>::compute_nl_term(
const Vector<double> &old_data,
803 *
const QGauss<dim> quadrature_formula(fe.degree + 1);
805 * quadrature_formula,
809 *
const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
810 *
const unsigned int n_q_points = quadrature_formula.size();
813 * std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
814 * std::vector<double> old_data_values(n_q_points);
815 * std::vector<double> new_data_values(n_q_points);
817 *
for (
const auto &cell : dof_handler.active_cell_iterators())
822 * Once we re-initialize our <code>
FEValues</code> instantiation to
823 * the current cell, we make use of the
824 * <code>get_function_values</code> routine to get the values of the
825 *
"old" data (presumably at @f$t=t_{n-1}@f$) and the
"new" data
826 * (presumably at @f$t=t_n@f$) at the nodes of the chosen quadrature
831 * fe_values.get_function_values(old_data, old_data_values);
832 * fe_values.get_function_values(new_data, new_data_values);
836 * Now, we can evaluate @f$\int_K \sin\left[\theta w_{\mathrm{
new}} +
837 * (1-\theta) w_{\mathrm{old}}\right] \,\varphi_j\,\mathrm{
d}x@f$
using
838 * the desired quadrature formula.
841 *
for (
unsigned int q_point = 0; q_point < n_q_points; ++q_point)
842 *
for (
unsigned int i = 0; i < dofs_per_cell; ++i)
843 * local_nl_term(i) +=
844 * (
std::sin(theta * new_data_values[q_point] +
845 * (1 - theta) * old_data_values[q_point]) *
846 * fe_values.shape_value(i, q_point) * fe_values.JxW(q_point));
850 * We conclude by adding up the contributions of the integrals over
851 * the cells to the global integral.
854 * cell->get_dof_indices(local_dof_indices);
856 *
for (
unsigned int i = 0; i < dofs_per_cell; ++i)
857 * nl_term(local_dof_indices[i]) += local_nl_term(i);
864 * <a name=
"step_25-SineGordonProblemcompute_nl_matrix"></a>
865 * <h4>SineGordonProblem::compute_nl_matrix</h4>
869 * This is the
second function dealing with the nonlinear scheme. It
870 * computes the
matrix @f$N(\cdot,\cdot)@f$, which appears in the nonlinear
871 * term in the Jacobian of @f$F(\cdot)@f$. Just as <code>compute_nl_term</code>,
872 * we must allow
this function to receive as input an
"old" and a
"new"
873 * solution, which we again call @f$w_{\mathrm{old}}@f$ and @f$w_{\mathrm{
new}}@f$
874 * below, respectively.
878 *
void SineGordonProblem<dim>::compute_nl_matrix(
883 *
const QGauss<dim> quadrature_formula(fe.degree + 1);
885 * quadrature_formula,
889 *
const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
890 *
const unsigned int n_q_points = quadrature_formula.size();
893 * std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
894 * std::vector<double> old_data_values(n_q_points);
895 * std::vector<double> new_data_values(n_q_points);
897 *
for (
const auto &cell : dof_handler.active_cell_iterators())
899 * local_nl_matrix = 0;
903 * instantiation to the current cell.
907 * fe_values.get_function_values(old_data, old_data_values);
908 * fe_values.get_function_values(new_data, new_data_values);
912 * Then, we evaluate @f$\int_K \cos\left[\theta w_{\mathrm{
new}} +
913 * (1-\theta) w_{\mathrm{old}}\right]\, \varphi_i\,
914 * \varphi_j\,\mathrm{
d}x@f$
using the desired quadrature formula.
917 *
for (
unsigned int q_point = 0; q_point < n_q_points; ++q_point)
918 *
for (
unsigned int i = 0; i < dofs_per_cell; ++i)
919 *
for (
unsigned int j = 0; j < dofs_per_cell; ++j)
920 * local_nl_matrix(i, j) +=
921 * (
std::cos(theta * new_data_values[q_point] +
922 * (1 - theta) * old_data_values[q_point]) *
923 * fe_values.shape_value(i, q_point) *
924 * fe_values.shape_value(j, q_point) * fe_values.JxW(q_point));
928 * Finally, we add up the contributions of the integrals over the
929 * cells to the global integral.
932 * cell->get_dof_indices(local_dof_indices);
934 *
for (
unsigned int i = 0; i < dofs_per_cell; ++i)
935 *
for (
unsigned int j = 0; j < dofs_per_cell; ++j)
936 * nl_matrix.add(local_dof_indices[i],
937 * local_dof_indices[j],
938 * local_nl_matrix(i, j));
947 * <a name=
"step_25-SineGordonProblemsolve"></a>
948 * <h4>SineGordonProblem::solve</h4>
952 * As discussed in the Introduction,
this function uses the CG iterative
953 * solver on the linear system of equations resulting from the finite
954 * element spatial discretization of each iteration of Newton
's method for
955 * the (nonlinear) first equation of the split formulation. The solution to
956 * the system is, in fact, @f$\delta U^{n,l}@f$ so it is stored in
957 * <code>solution_update</code> and used to update <code>solution</code> in
958 * the <code>run</code> function.
962 * Note that we re-set the solution update to zero before solving for
963 * it. This is not necessary: iterative solvers can start from any point and
964 * converge to the correct solution. If one has a good estimate about the
965 * solution of a linear system, it may be worthwhile to start from that
966 * vector, but as a general observation it is a fact that the starting point
967 * doesn't matter very much: it has to be a very, very good guess to
reduce
968 * the number of iterations by more than a few. It turns out that
for this
969 * problem,
using the previous nonlinear update as a starting
point actually
970 * hurts convergence and increases the number of iterations needed, so we
971 * simply
set it to zero.
975 * The function returns the number of iterations it took to converge to a
976 * solution. This number will later be used to generate output on the screen
977 * showing how many iterations were needed in each nonlinear iteration.
981 *
unsigned int SineGordonProblem<dim>::solve()
983 *
SolverControl solver_control(1000, 1e-12 * system_rhs.l2_norm());
987 * preconditioner.
initialize(system_matrix, 1.2);
989 * cg.solve(system_matrix, solution_update, system_rhs, preconditioner);
991 *
return solver_control.last_step();
997 * <a name=
"step_25-SineGordonProblemoutput_results"></a>
998 * <h4>SineGordonProblem::output_results</h4>
1002 * This function outputs the results to a file. It is pretty much identical
1003 * to the respective
functions in @ref step_23
"step-23" and @ref step_24
"step-24":
1006 *
template <
int dim>
1007 *
void SineGordonProblem<dim>::output_results(
1008 *
const unsigned int timestep_number)
const
1013 * data_out.add_data_vector(solution,
"u");
1014 * data_out.build_patches();
1016 *
const std::string filename =
1020 * data_out.set_flags(vtk_flags);
1021 * std::ofstream output(filename);
1022 * data_out.write_vtu(output);
1028 * <a name=
"step_25-SineGordonProblemrun"></a>
1029 * <h4>SineGordonProblem::run</h4>
1033 * This function has the top-
level control over everything: it runs the
1034 * (outer) time-stepping loop, the (inner) nonlinear-solver
loop, and
1035 * outputs the solution after each time step.
1038 *
template <
int dim>
1039 *
void SineGordonProblem<dim>::run()
1041 * make_grid_and_dofs();
1045 * To acknowledge the
initial condition, we must use the function @f$u_0(x)@f$
1046 * to compute @f$U^0@f$. To
this end, below we will create an
object of type
1047 * <code>InitialValues</code>; note that when we create
this object (which
1049 * time variable to @f$t_0@f$, to indicate that the
initial condition is a
1050 * function of space and time evaluated at @f$t=t_0@f$.
1054 * Then we produce @f$U^0@f$ by projecting @f$u_0(x)@f$ onto the grid
using
1057 * function
requires a hanging node constraints object, but to be used we
1058 *
first need to close it:
1063 * constraints.
close();
1067 * InitialValues<dim>(1, time),
1073 * For completeness, we output the zeroth time step to a file just like
1074 * any other time step.
1077 * output_results(0);
1081 * Now we perform the time stepping: at every time step we solve the
1082 *
matrix equation(s) corresponding to the finite element discretization
1083 * of the problem, and then
advance our solution according to the time
1084 * stepping formulas we discussed in the Introduction.
1087 *
unsigned int timestep_number = 1;
1088 *
for (time += time_step; time <= final_time;
1089 * time += time_step, ++timestep_number)
1091 * old_solution = solution;
1093 * std::cout << std::endl
1094 * <<
"Time step #" << timestep_number <<
"; "
1095 * <<
"advancing to t = " << time <<
'.' << std::endl;
1099 * At the beginning of each time step we must solve the nonlinear
1100 * equation in the
split formulation via Newton
's method ---
1101 * i.e. solve for @f$\delta U^{n,l}@f$ then compute @f$U^{n,l+1}@f$ and so
1102 * on. The stopping criterion for this nonlinear iteration is that
1103 * @f$\|F_h(U^{n,l})\|_2 \le 10^{-6} \|F_h(U^{n,0})\|_2@f$. Consequently,
1104 * we need to record the norm of the residual in the first iteration.
1108 * At the end of each iteration, we output to the console how many
1109 * linear solver iterations it took us. When the loop below is done,
1110 * we have (an approximation of) @f$U^n@f$.
1113 * double initial_rhs_norm = 0.;
1114 * bool first_iteration = true;
1117 * assemble_system();
1119 * if (first_iteration == true)
1120 * initial_rhs_norm = system_rhs.l2_norm();
1122 * const unsigned int n_iterations = solve();
1124 * solution += solution_update;
1126 * if (first_iteration == true)
1127 * std::cout << " " << n_iterations;
1129 * std::cout << '+
' << n_iterations;
1130 * first_iteration = false;
1132 * while (system_rhs.l2_norm() > 1e-6 * initial_rhs_norm);
1134 * std::cout << " CG iterations per nonlinear step." << std::endl;
1138 * Upon obtaining the solution to the first equation of the problem at
1139 * @f$t=t_n@f$, we must update the auxiliary velocity variable
1140 * @f$V^n@f$. However, we do not compute and store @f$V^n@f$ since it is not a
1141 * quantity we use directly in the problem. Hence, for simplicity, we
1142 * update @f$MV^n@f$ directly:
1145 * Vector<double> tmp_vector(solution.size());
1146 * laplace_matrix.vmult(tmp_vector, solution);
1147 * M_x_velocity.add(-time_step * theta, tmp_vector);
1149 * laplace_matrix.vmult(tmp_vector, old_solution);
1150 * M_x_velocity.add(-time_step * (1 - theta), tmp_vector);
1152 * compute_nl_term(old_solution, solution, tmp_vector);
1153 * M_x_velocity.add(-time_step, tmp_vector);
1157 * Oftentimes, in particular for fine meshes, we must pick the time
1158 * step to be quite small in order for the scheme to be
1159 * stable. Therefore, there are a lot of time steps during which
1160 * "nothing interesting happens" in the solution. To improve overall
1161 * efficiency -- in particular, speed up the program and save disk
1162 * space -- we only output the solution every
1163 * <code>output_timestep_skip</code> time steps:
1166 * if (timestep_number % output_timestep_skip == 0)
1167 * output_results(timestep_number);
1170 * } // namespace Step25
1175 * <a name="step_25-Thecodemaincodefunction"></a>
1176 * <h3>The <code>main</code> function</h3>
1180 * This is the main function of the program. It creates an object of top-level
1181 * class and calls its principal function. If exceptions are thrown during the
1182 * execution of the run method of the <code>SineGordonProblem</code> class, we
1183 * catch and report them here. For more information about exceptions the
1184 * reader should consult @ref step_6 "step-6".
1191 * using namespace Step25;
1193 * SineGordonProblem<1> sg_problem;
1196 * catch (std::exception &exc)
1198 * std::cerr << std::endl
1200 * << "----------------------------------------------------"
1202 * std::cerr << "Exception on processing: " << std::endl
1203 * << exc.what() << std::endl
1204 * << "Aborting!" << std::endl
1205 * << "----------------------------------------------------"
1212 * std::cerr << std::endl
1214 * << "----------------------------------------------------"
1216 * std::cerr << "Unknown exception!" << std::endl
1217 * << "Aborting!" << std::endl
1218 * << "----------------------------------------------------"
1226<a name="step_25-Results"></a><h1>Results</h1>
1228The explicit Euler time stepping scheme (@f$\theta=0@f$) performs adequately for the problems we wish to solve. Unfortunately, a rather small time step has to be chosen due to stability issues --- @f$k\sim h/10@f$ appears to work for most the simulations we performed. On the other hand, the Crank-Nicolson scheme (@f$\theta=\frac{1}{2}@f$) is unconditionally stable, and (at least for the case of the 1D breather) we can pick the time step to be as large as @f$25h@f$ without any ill effects on the solution. The implicit Euler scheme (@f$\theta=1@f$) is "exponentially damped," so it is not a good choice for solving the sine-Gordon equation, which is conservative. However, some of the damped schemes in the continuum that is offered by the @f$\theta@f$-method were useful for eliminating spurious oscillations due to boundary effects.
1230In the simulations below, we solve the sine-Gordon equation on the interval @f$\Omega =
1231[-10,10]@f$ in 1D and on the square @f$\Omega = [-10,10]\times [-10,10]@f$ in 2D. In
1232each case, the respective grid is refined uniformly 6 times, i.e. @f$h\sim
1235<a name="step_25-An11dSolution"></a><h3>An (1+1)-d Solution</h3>
1237The first example we discuss is the so-called 1D (stationary) breather
1238solution of the sine-Gordon equation. The breather has the following
1239closed-form expression, as mentioned in the Introduction:
1241u_{\mathrm{breather}}(x,t) = -4\arctan \left(\frac{m}{\sqrt{1-m^2}} \frac{\sin\left(\sqrt{1-m^2}t +c_2\right)}{\cosh(mx+c_1)} \right),
1243where @f$c_1@f$, @f$c_2@f$ and @f$m<1@f$ are constants. In the simulation below, we have chosen @f$c_1=0@f$, @f$c_2=0@f$, @f$m=0.5@f$. Moreover, it is know that the period of oscillation of the breather is @f$2\pi\sqrt{1-m^2}@f$, hence we have chosen @f$t_0=-5.4414@f$ and @f$t_f=2.7207@f$ so that we can observe three oscillations of the solution. Then, taking @f$u_0(x) = u_{\mathrm{breather}}(x,t_0)@f$, @f$\theta=0@f$ and @f$k=h/10@f$, the program computed the following solution.
1245<img src="https://www.dealii.org/images/steps/developer/step-25.1d-breather.gif" alt="Animation of the 1D stationary breather.">
1247Though not shown how to do this in the program, another way to visualize the
1248(1+1)-d solution is to use output generated by the DataOutStack class; it
1249allows to "stack" the solutions of individual time steps, so that we get
12502D space-time graphs from 1D time-dependent
1251solutions. This produces the space-time plot below instead of the animation
1254<img src="https://www.dealii.org/images/steps/developer/step-25.1d-breather_stp.png" alt="A space-time plot of the 1D stationary breather.">
1256Furthermore, since the breather is an analytical solution of the sine-Gordon
1257equation, we can use it to validate our code, although we have to assume that
1258the error introduced by our choice of Neumann boundary conditions is small
1259compared to the numerical error. Under this assumption, one could use the
1260VectorTools::integrate_difference function to compute the difference between
1261the numerical solution and the function described by the
1262<code>ExactSolution</code> class of this program. For the
1263simulation shown in the two images above, the @f$L^2@f$ norm of the error in the
1264finite element solution at each time step remained on the order of
1265@f$10^{-2}@f$. Hence, we can conclude that the numerical method has been
1266implemented correctly in the program.
1269<a name="step_25-Afew21DSolutions"></a><h3>A few (2+1)D Solutions</h3>
1272The only analytical solution to the sine-Gordon equation in (2+1)D that can be found in the literature is the so-called kink solitary wave. It has the following closed-form expression:
1274 u(x,y,t) = 4 \arctan \left[a_0 e^{s\xi}\right]
1278 \xi = x \cos\vartheta + \sin(\vartheta) (y\cosh\lambda + t\sinh \lambda)
1280where @f$a_0@f$, @f$\vartheta@f$ and @f$\lambda@f$ are constants. In the simulation below
1281we have chosen @f$a_0=\lambda=1@f$. Notice that if @f$\vartheta=\pi@f$ the kink is
1282stationary, hence it would make a good solution against which we can
1283validate the program in 2D because no reflections off the boundary of the
1286The simulation shown below was performed with @f$u_0(x) = u_{\mathrm{kink}}(x,t_0)@f$, @f$\theta=\frac{1}{2}@f$, @f$k=20h@f$, @f$t_0=1@f$ and @f$t_f=500@f$. The @f$L^2@f$ norm of the error of the finite element solution at each time step remained on the order of @f$10^{-2}@f$, showing that the program is working correctly in 2D, as well as 1D. Unfortunately, the solution is not very interesting, nonetheless we have included a snapshot of it below for completeness.
1288<img src="https://www.dealii.org/images/steps/developer/step-25.2d-kink.png" alt="Stationary 2D kink.">
1290Now that we have validated the code in 1D and 2D, we move to a problem where the analytical solution is unknown.
1292To this end, we rotate the kink solution discussed above about the @f$z@f$
1293axis: we let @f$\vartheta=\frac{\pi}{4}@f$. The latter results in a
1294solitary wave that is not aligned with the grid, so reflections occur
1295at the boundaries of the domain immediately. For the simulation shown
1296below, we have taken @f$u_0(x)=u_{\mathrm{kink}}(x,t_0)@f$,
1297@f$\theta=\frac{2}{3}@f$, @f$k=20h@f$, @f$t_0=0@f$ and @f$t_f=20@f$. Moreover, we had
1298to pick @f$\theta=\frac{2}{3}@f$ because for any @f$\theta\le\frac{1}{2}@f$
1299oscillations arose at the boundary, which are likely due to the scheme
1300and not the equation, thus picking a value of @f$\theta@f$ a good bit into
1301the "exponentially damped" spectrum of the time stepping schemes
1302assures these oscillations are not created.
1304<img src="https://www.dealii.org/images/steps/developer/step-25.2d-angled_kink.gif" alt="Animation of a moving 2D kink, at 45 degrees to the axes of the grid, showing boundary effects.">
1306Another interesting solution to the sine-Gordon equation (which cannot be
1307obtained analytically) can be produced by using two 1D breathers to construct
1308the following separable 2D initial condition:
1311 u_{\mathrm{pseudobreather}}(x,t_0) =
1313 \frac{m}{\sqrt{1-m^2}}
1314 \frac{\sin\left(\sqrt{1-m^2}t_0\right)}{\cosh(mx_1)} \right)
1316 \frac{m}{\sqrt{1-m^2}}
1317 \frac{\sin\left(\sqrt{1-m^2}t_0\right)}{\cosh(mx_2)} \right),
1319where @f$x=(x_1,x_2)\in{R}^2@f$, @f$m=0.5<1@f$ as in the 1D case we discussed
1320above. For the simulation shown below, we have chosen @f$\theta=\frac{1}{2}@f$,
1321@f$k=10h@f$, @f$t_0=-5.4414@f$ and @f$t_f=2.7207@f$. The solution is pretty interesting
1322--- it acts like a breather (as far as the pictures are concerned); however,
1323it appears to break up and reassemble, rather than just oscillate.
1325<img src="https://www.dealii.org/images/steps/developer/step-25.2d-pseudobreather.gif" alt="Animation of a 2D pseudobreather.">
1328<a name="step-25-extensions"></a>
1329<a name="step_25-Possibilitiesforextensions"></a><h3>Possibilities for extensions</h3>
1332It is instructive to change the initial conditions. Most choices will not lead
1333to solutions that stay localized (in the soliton community, such
1334solutions are called "stationary", though the solution does change
1335with time), but lead to solutions where the wave-like
1336character of the equation dominates and a wave travels away from the location
1337of a localized initial condition. For example, it is worth playing around with
1338the <code>InitialValues</code> class, by replacing the call to the
1339<code>ExactSolution</code> class by something like this function:
1341 u_0(x,y) = \cos\left(\frac x2\right)\cos\left(\frac y2\right)
1343if @f$|x|,|y|\le \frac\pi 2@f$, and @f$u_0(x,y)=0@f$ outside this region.
1345A second area would be to investigate whether the scheme is
1346energy-preserving. For the pure wave equation, discussed in @ref
1347step_23 "step-23", this is the case if we choose the time stepping
1348parameter such that we get the Crank-Nicolson scheme. One could do a
1349similar thing here, noting that the energy in the sine-Gordon solution
1352 E(t) = \frac 12 \int_\Omega \left(\frac{\partial u}{\partial
1354 + \left(\nabla u\right)^2 + 2 (1-\cos u) \; dx.
1356(We use @f$1-\cos u@f$ instead of @f$-\cos u@f$ in the formula to ensure that all
1357contributions to the energy are positive, and so that decaying solutions have
1358finite energy on unbounded domains.)
1360Beyond this, there are two obvious areas:
1362- Clearly, adaptivity (i.e. time-adaptive grids) would be of interest
1363 to problems like these. Their complexity leads us to leave this out
1364 of this program again, though the general comments in the
1365 introduction of @ref step_23 "step-23" remain true.
1367- Faster schemes to solve this problem. While computers today are
1368 plenty fast enough to solve 2d and, frequently, even 3d stationary
1369 problems within not too much time, time dependent problems present
1370 an entirely different class of problems. We address this topic in
1371 @ref step_48 "step-48" where we show how to solve this problem in parallel and
1372 without assembling or inverting any matrix at all.
1375<a name="step_25-PlainProg"></a>
1376<h1> The plain program</h1>
1377@include "step-25.cc"
void attach_dof_handler(const DoFHandler< dim, spacedim > &)
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
void initialize(const MatrixType &A, const AdditionalData ¶meters=AdditionalData())
unsigned int n_active_cells() const
void refine_global(const unsigned int times=1)
unsigned int n_cells() const
__global__ void set(Number *val, const Number s, const size_type N)
void loop(IteratorType begin, std_cxx20::type_identity_t< IteratorType > end, DOFINFO &dinfo, INFOBOX &info, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &)> &cell_worker, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &)> &boundary_worker, const std::function< void(std_cxx20::type_identity_t< DOFINFO > &, std_cxx20::type_identity_t< DOFINFO > &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, AssemblerType &assembler, const LoopControl &lctrl=LoopControl())
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity_pattern, const AffineConstraints< number > &constraints={}, const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
@ update_values
Shape function values.
@ update_JxW_values
Transformed quadrature weights.
@ update_quadrature_points
Transformed quadrature points.
#define DEAL_II_NOT_IMPLEMENTED()
std::vector< value_type > split(const typename ::Triangulation< dim, spacedim >::cell_iterator &parent, const value_type parent_value)
void hyper_cube(Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false)
@ matrix
Contents is actually a matrix.
void mass_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const double factor=1.)
void create_mass_matrix(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &q, MatrixType &matrix, const Function< spacedim, typename MatrixType::value_type > *const a=nullptr, const AffineConstraints< typename MatrixType::value_type > &constraints=AffineConstraints< typename MatrixType::value_type >())
void create_laplace_matrix(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &q, MatrixType &matrix, const Function< spacedim, typename MatrixType::value_type > *const a=nullptr, const AffineConstraints< typename MatrixType::value_type > &constraints=AffineConstraints< typename MatrixType::value_type >())
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
VectorType::value_type * end(VectorType &V)
T reduce(const T &local_value, const MPI_Comm comm, const std::function< T(const T &, const T &)> &combiner, const unsigned int root_process=0)
constexpr T fixed_power(const T t)
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
int(& functions)(const void *v1, const void *v2)
static constexpr double PI
::VectorizedArray< Number, width > exp(const ::VectorizedArray< Number, width > &)
inline ::VectorizedArray< Number, width > sinh(const ::VectorizedArray< Number, width > &x)
::VectorizedArray< Number, width > cos(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sin(const ::VectorizedArray< Number, width > &)
inline ::VectorizedArray< Number, width > cosh(const ::VectorizedArray< Number, width > &x)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
inline ::VectorizedArray< Number, width > atan(const ::VectorizedArray< Number, width > &x)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
DataOutBase::CompressionLevel compression_level
void advance(std::tuple< I1, I2 > &t, const unsigned int n)