660 *
const unsigned int component = 0)
const override
663 *
Assert(component == 0, ExcIndexRange(component, 0, 1));
671 *
class InitialValuesV :
public Function<dim>
675 *
const unsigned int component = 0) const override
687 * Secondly, we have the right hand side forcing term. Boring as we are, we
688 * choose zero here as well:
692 *
class RightHandSide :
public Function<dim>
696 *
const unsigned int component = 0) const override
708 * Finally, we have boundary values
for @f$u@f$ and @f$v@f$. They are as described
709 * in the introduction, one being the time derivative of the other:
713 *
class BoundaryValuesU :
public Function<dim>
717 *
const unsigned int component = 0) const override
722 *
if ((this->
get_time() <= 0.5) && (p[0] < 0) && (p[1] < 1. / 3) &&
733 *
class BoundaryValuesV :
public Function<dim>
737 *
const unsigned int component = 0) const override
742 *
if ((this->
get_time() <= 0.5) && (p[0] < 0) && (p[1] < 1. / 3) &&
755 * <a name=
"step_23-ImplementationofthecodeWaveEquationcodeclass"></a>
756 * <h3>Implementation of the <code>WaveEquation</code>
class</h3>
760 * The implementation of the actual logic is actually fairly short, since we
761 * relegate things like assembling the matrices and right hand side vectors
762 * to the library. The rest boils down to not much more than 130 lines of
763 * actual code, a significant fraction of which is boilerplate code that can
764 * be taken from previous example programs (
e.g. the functions that solve
765 * linear systems, or that generate output).
769 * Let
's start with the constructor (for an explanation of the choice of
770 * time step, see the section on Courant, Friedrichs, and Lewy in the
775 * WaveEquation<dim>::WaveEquation()
777 * , dof_handler(triangulation)
778 * , time_step(1. / 64)
780 * , timestep_number(1)
788 * <a name="step_23-WaveEquationsetup_system"></a>
789 * <h4>WaveEquation::setup_system</h4>
793 * The next function is the one that sets up the mesh, DoFHandler, and
794 * matrices and vectors at the beginning of the program, i.e. before the
795 * first time step. The first few lines are pretty much standard if you've
796 * read through the tutorial programs at least up to @ref step_6
"step-6":
800 *
void WaveEquation<dim>::setup_system()
808 * dof_handler.distribute_dofs(fe);
810 * std::cout <<
"Number of degrees of freedom: " << dof_handler.n_dofs()
816 * sparsity_pattern.copy_from(dsp);
820 * Then comes a block where we have to initialize the 3 matrices we need
821 * in the course of the program: the mass
matrix, the Laplace
matrix, and
822 * the
matrix @f$M+k^2\theta^2A@f$ used when solving
for @f$U^n@f$ in each time
827 * When setting up these matrices, note that they all make use of the same
828 * sparsity pattern
object. Finally, the reason why matrices and sparsity
829 * patterns are separate objects in deal.II (unlike in many other finite
830 * element or linear algebra classes) becomes clear: in a significant
831 * fraction of applications, one has to hold several matrices that happen
832 * to have the same sparsity pattern, and there is no reason
for them not
833 * to share
this information, rather than re-building and wasting memory
834 * on it several times.
838 * After initializing all of these matrices, we call library
functions
839 * that build the Laplace and mass matrices. All they need is a
DoFHandler
840 *
object and a quadrature formula
object that is to be used
for numerical
841 * integration. Note that in many respects these
functions are better than
842 * what we would usually
do in application programs,
for example because
843 * they automatically parallelize building the matrices
if multiple
844 * processors are available in a machine:
for more information see the
846 * @ref threads
"Parallel computing with multiple processors"
847 * topic. The matrices
for solving linear systems will be filled in the
848 *
run() method because we need to re-apply boundary conditions every time
852 * mass_matrix.reinit(sparsity_pattern);
853 * laplace_matrix.reinit(sparsity_pattern);
854 * matrix_u.reinit(sparsity_pattern);
855 * matrix_v.reinit(sparsity_pattern);
858 *
QGauss<dim>(fe.degree + 1),
861 *
QGauss<dim>(fe.degree + 1),
866 * The rest of the function is spent on setting vector sizes to the
867 * correct value. The final line closes the hanging node constraints
868 *
object. Since we work on a uniformly refined mesh, no constraints exist
869 * or have been computed (i.e. there was no need to call
870 *
DoFTools::make_hanging_node_constraints as in other programs), but we
871 * need a constraints
object in one place further down below anyway.
874 * solution_u.reinit(dof_handler.n_dofs());
875 * solution_v.reinit(dof_handler.n_dofs());
876 * old_solution_u.reinit(dof_handler.n_dofs());
877 * old_solution_v.reinit(dof_handler.n_dofs());
878 * system_rhs.reinit(dof_handler.n_dofs());
880 * constraints.close();
888 * <a name="step_23-WaveEquationsolve_uandWaveEquationsolve_v"></a>
889 * <h4>WaveEquation::solve_u and WaveEquation::solve_v</h4>
893 * The next two functions deal with solving the linear systems associated
894 * with the equations for @f$U^n@f$ and @f$V^n@f$. Both are not particularly
895 * interesting as they pretty much follow the scheme used in all the
896 * previous tutorial programs.
900 * One can make little experiments with preconditioners for the two matrices
901 * we have to
invert. As it turns out, however, for the matrices at hand
902 * here, using Jacobi or SSOR preconditioners reduces the number of
903 * iterations necessary to solve the linear system slightly, but due to the
904 * cost of applying the preconditioner it is no win in terms of run-time. It
905 * is not much of a loss either, but let's keep it simple and just do
910 *
void WaveEquation<dim>::solve_u()
912 *
SolverControl solver_control(1000, 1e-8 * system_rhs.l2_norm());
917 * std::cout <<
" u-equation: " << solver_control.last_step()
918 * <<
" CG iterations." << std::endl;
924 *
void WaveEquation<dim>::solve_v()
926 *
SolverControl solver_control(1000, 1e-8 * system_rhs.l2_norm());
931 * std::cout <<
" v-equation: " << solver_control.last_step()
932 * <<
" CG iterations." << std::endl;
940 * <a name=
"step_23-WaveEquationoutput_results"></a>
941 * <h4>WaveEquation::output_results</h4>
945 * Likewise, the following function is pretty much what we
've done
946 * before. The only thing worth mentioning is how here we generate a string
947 * representation of the time step number padded with leading zeros to 3
948 * character length using the Utilities::int_to_string function's
second
953 *
void WaveEquation<dim>::output_results() const
958 * data_out.add_data_vector(solution_u,
"U");
959 * data_out.add_data_vector(solution_v,
"V");
961 * data_out.build_patches();
963 *
const std::string filename =
967 * Like @ref step_15
"step-15", since we write output at every time step (and the system
968 * we have to solve is relatively easy), we instruct
DataOut to use the
969 * zlib compression algorithm that is optimized
for speed instead of disk
970 * usage since otherwise plotting the output becomes a bottleneck:
975 * data_out.set_flags(vtk_flags);
976 * std::ofstream output(filename);
977 * data_out.write_vtu(output);
985 * <a name=
"step_23-WaveEquationrun"></a>
986 * <h4>WaveEquation::run</h4>
990 * The following is really the only interesting function of the program. It
991 * contains the
loop over all time steps, but before we get to that we have
992 * to
set up the grid,
DoFHandler, and matrices. In addition, we have to
993 * somehow get started with
initial values. To
this end, we use the
995 * continuous function and computes the @f$L^2@f$ projection of
this function
996 * onto the finite element space described by the
DoFHandler object. Can
't
997 * be any simpler than that:
1000 * template <int dim>
1001 * void WaveEquation<dim>::run()
1005 * VectorTools::project(dof_handler,
1007 * QGauss<dim>(fe.degree + 1),
1008 * InitialValuesU<dim>(),
1010 * VectorTools::project(dof_handler,
1012 * QGauss<dim>(fe.degree + 1),
1013 * InitialValuesV<dim>(),
1018 * The next thing is to loop over all the time steps until we reach the
1019 * end time (@f$T=5@f$ in this case). In each time step, we first have to
1020 * solve for @f$U^n@f$, using the equation @f$(M^n + k^2\theta^2 A^n)U^n =@f$
1021 * @f$(M^{n,n-1} - k^2\theta(1-\theta) A^{n,n-1})U^{n-1} + kM^{n,n-1}V^{n-1}
1022 * +@f$ @f$k\theta \left[k \theta F^n + k(1-\theta) F^{n-1} \right]@f$. Note
1023 * that we use the same mesh for all time steps, so that @f$M^n=M^{n,n-1}=M@f$
1024 * and @f$A^n=A^{n,n-1}=A@f$. What we therefore have to do first is to add up
1025 * @f$MU^{n-1} - k^2\theta(1-\theta) AU^{n-1} + kMV^{n-1}@f$ and the forcing
1026 * terms, and put the result into the <code>system_rhs</code> vector. (For
1027 * these additions, we need a temporary vector that we declare before the
1028 * loop to avoid repeated memory allocations in each time step.)
1032 * The one thing to realize here is how we communicate the time variable
1033 * to the object describing the right hand side: each object derived from
1034 * the Function class has a time field that can be set using the
1035 * Function::set_time and read by Function::get_time. In essence, using
1036 * this mechanism, all functions of space and time are therefore
1037 * considered functions of space evaluated at a particular time. This
1038 * matches well what we typically need in finite element programs, where
1039 * we almost always work on a single time step at a time, and where it
1040 * never happens that, for example, one would like to evaluate a
1041 * space-time function for all times at any given spatial location.
1044 * Vector<double> tmp(solution_u.size());
1045 * Vector<double> forcing_terms(solution_u.size());
1047 * for (; time <= 5; time += time_step, ++timestep_number)
1049 * std::cout << "Time step " << timestep_number << " at t=" << time
1052 * mass_matrix.vmult(system_rhs, old_solution_u);
1054 * mass_matrix.vmult(tmp, old_solution_v);
1055 * system_rhs.add(time_step, tmp);
1057 * laplace_matrix.vmult(tmp, old_solution_u);
1058 * system_rhs.add(-theta * (1 - theta) * time_step * time_step, tmp);
1060 * RightHandSide<dim> rhs_function;
1061 * rhs_function.set_time(time);
1062 * VectorTools::create_right_hand_side(dof_handler,
1063 * QGauss<dim>(fe.degree + 1),
1066 * forcing_terms = tmp;
1067 * forcing_terms *= theta * time_step;
1069 * rhs_function.set_time(time - time_step);
1070 * VectorTools::create_right_hand_side(dof_handler,
1071 * QGauss<dim>(fe.degree + 1),
1075 * forcing_terms.add((1 - theta) * time_step, tmp);
1077 * system_rhs.add(theta * time_step, forcing_terms);
1081 * After so constructing the right hand side vector of the first
1082 * equation, all we have to do is apply the correct boundary
1083 * values. As for the right hand side, this is a space-time function
1084 * evaluated at a particular time, which we interpolate at boundary
1085 * nodes and then use the result to apply boundary values as we
1086 * usually do. The result is then handed off to the solve_u()
1091 * BoundaryValuesU<dim> boundary_values_u_function;
1092 * boundary_values_u_function.set_time(time);
1094 * std::map<types::global_dof_index, double> boundary_values;
1095 * VectorTools::interpolate_boundary_values(dof_handler,
1097 * boundary_values_u_function,
1102 * The matrix for solve_u() is the same in every time steps, so one
1103 * could think that it is enough to do this only once at the
1104 * beginning of the simulation. However, since we need to apply
1105 * boundary values to the linear system (which eliminate some matrix
1106 * rows and columns and give contributions to the right hand side),
1107 * we have to refill the matrix in every time steps before we
1108 * actually apply boundary data. The actual content is very simple:
1109 * it is the sum of the mass matrix and a weighted Laplace matrix:
1112 * matrix_u.copy_from(mass_matrix);
1113 * matrix_u.add(theta * theta * time_step * time_step, laplace_matrix);
1114 * MatrixTools::apply_boundary_values(boundary_values,
1124 * The second step, i.e. solving for @f$V^n@f$, works similarly, except
1125 * that this time the matrix on the left is the mass matrix (which we
1126 * copy again in order to be able to apply boundary conditions, and
1127 * the right hand side is @f$MV^{n-1} - k\left[ \theta A U^n +
1128 * (1-\theta) AU^{n-1}\right]@f$ plus forcing terms. Boundary values
1129 * are applied in the same way as before, except that now we have to
1130 * use the BoundaryValuesV class:
1133 * laplace_matrix.vmult(system_rhs, solution_u);
1134 * system_rhs *= -theta * time_step;
1136 * mass_matrix.vmult(tmp, old_solution_v);
1137 * system_rhs += tmp;
1139 * laplace_matrix.vmult(tmp, old_solution_u);
1140 * system_rhs.add(-time_step * (1 - theta), tmp);
1142 * system_rhs += forcing_terms;
1145 * BoundaryValuesV<dim> boundary_values_v_function;
1146 * boundary_values_v_function.set_time(time);
1148 * std::map<types::global_dof_index, double> boundary_values;
1149 * VectorTools::interpolate_boundary_values(dof_handler,
1151 * boundary_values_v_function,
1153 * matrix_v.copy_from(mass_matrix);
1154 * MatrixTools::apply_boundary_values(boundary_values,
1163 * Finally, after both solution components have been computed, we
1164 * output the result, compute the energy in the solution, and go on to
1165 * the next time step after shifting the present solution into the
1166 * vectors that hold the solution at the previous time step. Note the
1167 * function SparseMatrix::matrix_norm_square that can compute
1168 * @f$\left<V^n,MV^n\right>@f$ and @f$\left<U^n,AU^n\right>@f$ in one step,
1169 * saving us the expense of a temporary vector and several lines of
1175 * std::cout << " Total energy: "
1176 * << (mass_matrix.matrix_norm_square(solution_v) +
1177 * laplace_matrix.matrix_norm_square(solution_u)) /
1181 * old_solution_u = solution_u;
1182 * old_solution_v = solution_v;
1185 * } // namespace Step23
1191 * <a name="step_23-Thecodemaincodefunction"></a>
1192 * <h3>The <code>main</code> function</h3>
1196 * What remains is the main function of the program. There is nothing here
1197 * that hasn't been shown in several of the previous programs:
1204 *
using namespace Step23;
1206 * WaveEquation<2> wave_equation_solver;
1207 * wave_equation_solver.run();
1209 *
catch (std::exception &exc)
1211 * std::cerr << std::endl
1213 * <<
"----------------------------------------------------"
1215 * std::cerr <<
"Exception on processing: " << std::endl
1216 * << exc.what() << std::endl
1217 * <<
"Aborting!" << std::endl
1218 * <<
"----------------------------------------------------"
1225 * std::cerr << std::endl
1227 * <<
"----------------------------------------------------"
1229 * std::cerr <<
"Unknown exception!" << std::endl
1230 * <<
"Aborting!" << std::endl
1231 * <<
"----------------------------------------------------"
1239<a name=
"step_23-Results"></a><h1>Results</h1>
1242When the program is
run, it produces the following output:
1244Number of active cells: 16384
1245Number of degrees of freedom: 16641
1247Time step 1 at t=0.015625
1248 u-equation: 8 CG iterations.
1249 v-equation: 22 CG iterations.
1250 Total energy: 1.17887
1251Time step 2 at t=0.03125
1252 u-equation: 8 CG iterations.
1253 v-equation: 20 CG iterations.
1254 Total energy: 2.9655
1255Time step 3 at t=0.046875
1256 u-equation: 8 CG iterations.
1257 v-equation: 21 CG iterations.
1258 Total energy: 4.33761
1259Time step 4 at t=0.0625
1260 u-equation: 7 CG iterations.
1261 v-equation: 21 CG iterations.
1262 Total energy: 5.35499
1263Time step 5 at t=0.078125
1264 u-equation: 7 CG iterations.
1265 v-equation: 21 CG iterations.
1266 Total energy: 6.18652
1267Time step 6 at t=0.09375
1268 u-equation: 7 CG iterations.
1269 v-equation: 20 CG iterations.
1270 Total energy: 6.6799
1274Time step 31 at t=0.484375
1275 u-equation: 7 CG iterations.
1276 v-equation: 20 CG iterations.
1277 Total energy: 21.9068
1278Time step 32 at t=0.5
1279 u-equation: 7 CG iterations.
1280 v-equation: 20 CG iterations.
1281 Total energy: 23.3394
1282Time step 33 at t=0.515625
1283 u-equation: 7 CG iterations.
1284 v-equation: 20 CG iterations.
1285 Total energy: 23.1019
1289Time step 319 at t=4.98438
1290 u-equation: 7 CG iterations.
1291 v-equation: 20 CG iterations.
1292 Total energy: 23.1019
1294 u-equation: 7 CG iterations.
1295 v-equation: 20 CG iterations.
1296 Total energy: 23.1019
1299What we see immediately is that the energy is a
constant at least after
1300@f$t=\frac 12@f$ (until which the boundary source term @f$g@f$ is nonzero, injecting
1301energy into the system).
1303In addition to the screen output, the program writes the solution of each time
1304step to an output file. If we process them adequately and paste them into a
1305movie, we get the following:
1307<img src=
"https://www.dealii.org/images/steps/developer/step-23.movie.gif" alt=
"Animation of the solution of step 23.">
1309The movie shows the generated wave nice traveling through the domain and back,
1310being reflected at the clamped boundary. Some numerical noise is trailing the
1311wave, an artifact of a too-large mesh size that can be reduced by reducing the
1312mesh width and the time step.
1315<a name=
"step-23-extensions"></a>
1316<a name=
"step_23-Possibilitiesforextensions"></a><h3>Possibilities
for extensions</h3>
1319If you want to explore a bit,
try out some of the following things:
1321 <li>Varying @f$\theta@f$. This gives different time stepping schemes, some of
1322 which are stable
while others are not. Take a look at how the energy
1325 <li>Different
initial and boundary conditions, right hand sides.
1327 <li>More complicated domains or more refined meshes. Remember that the time
1328 step needs to be bounded by the mesh width, so changing the mesh should
1329 always involve also changing the time step. We will come back to
this issue
1330 in @ref step_24
"step-24".
1332 <li>Variable coefficients: In real media, the wave speed is often
1333 variable. In particular, the
"real" wave equation in realistic media would
1336 \rho(x) \frac{\partial^2 u}{\partial t^2}
1341 where @f$\rho(x)@f$ is the density of the material, and @f$a(x)@f$ is related to the
1342 stiffness coefficient. The wave speed is then @f$c=\
sqrt{a/\rho}@f$.
1344 To make such a change, we would have to compute the mass and Laplace
1345 matrices with a variable coefficient. Fortunately,
this isn
't too hard: the
1346 functions MatrixCreator::create_laplace_matrix and
1347 MatrixCreator::create_mass_matrix have additional default parameters that can
1348 be used to pass non-constant coefficient functions to them. The required
1349 changes are therefore relatively small. On the other hand, care must be
1350 taken again to make sure the time step is within the allowed range.
1352 <li>In the in-code comments, we discussed the fact that the matrices for
1353 solving for @f$U^n@f$ and @f$V^n@f$ need to be reset in every time because of
1354 boundary conditions, even though the actual content does not change. It is
1355 possible to avoid copying by not eliminating columns in the linear systems,
1356 which is implemented by appending a @p false argument to the call:
1358 MatrixTools::apply_boundary_values(boundary_values,
1365 <li>deal.II being a library that supports adaptive meshes it would of course be
1366 nice if this program supported change the mesh every few time steps. Given the
1367 structure of the solution — a wave that travels through the domain —
1368 it would seem appropriate if we only refined the mesh where the wave currently is,
1369 and not simply everywhere. It is intuitively clear that we should be able to
1370 save a significant amount of cells this way. (Though upon further thought one
1371 realizes that this is really only the case in the initial stages of the simulation.
1372 After some time, for wave phenomena, the domain is filled with reflections of
1373 the initial wave going in every direction and filling every corner of the domain.
1374 At this point, there is in general little one can gain using local mesh
1377 To make adaptively changing meshes possible, there are basically two routes.
1378 The "correct" way would be to go back to the weak form we get using Rothe's
1379 method. For example, the
first of the two equations to be solved in each time
1380 step looked like
this:
1382 (u^n,\varphi) + k^2\theta^2(\nabla u^n,\nabla \varphi) &=&
1383 (u^{n-1},\varphi) - k^2\theta(1-\theta)(\nabla u^{n-1},\nabla \varphi)
1388 \theta (f^n,\varphi) + (1-\theta) (f^{n-1},\varphi)
1391 Now, note that we solve
for @f$u^n@f$ on mesh @f${\mathbb T}^n@f$, and
1392 consequently the test
functions @f$\varphi@f$ have to be from the space
1393 @f$V_h^n@f$ as well. As discussed in the introduction, terms like
1394 @f$(u^{n-1},\varphi)@f$ then require us to integrate the solution of the
1395 previous step (which may have been computed on a different mesh
1396 @f${\mathbb T}^{n-1}@f$) against the test functions of the current mesh,
1397 leading to a matrix @f$M^{n,n-1}@f$. This process of integrating shape
1398 functions from different meshes is, at best, awkward. It can be done
1399 but because it is difficult to ensure that @f${\mathbb T}^{n-1}@f$ and
1400 @f${\mathbb T}^{n}@f$ differ by at most one
level of refinement, one
1401 has to recursively match cells from both meshes. It is feasible to
1402 do this, but it leads to lengthy and not entirely obvious code.
1404 The
second approach is the following: whenever we change the mesh,
1405 we simply
interpolate the solution from the last time step on the old
1407 instead of the equation above, we would solve
1409 (u^n,\varphi) + k^2\theta^2(\nabla u^n,\nabla \varphi) &=&
1410 (I^n u^{n-1},\varphi) - k^2\theta(1-\theta)(\nabla I^n u^{n-1},\nabla \varphi)
1412 k(I^n v^{n-1},\varphi)
1415 \theta (f^n,\varphi) + (1-\theta) (f^{n-1},\varphi)
1418 where @f$I^n@f$ interpolates a given function onto mesh @f${\mathbb T}^n@f$.
1419 This is a much simpler approach because, in each time step, we no
1420 longer have to worry whether @f$u^{n-1},v^{n-1}@f$ were computed on the
1421 same mesh as we are
using now or on a different mesh. Consequently,
1422 the only changes to the code necessary are the addition of a function
1423 that computes the error, marks cells
for refinement, sets up a
1425 rebuilds matrices and right hand side vectors on the
new mesh. Neither
1426 the
functions building the matrices and right hand sides, nor the
1427 solvers need to be changed.
1429 While
this second approach is, strictly speaking,
1430 not quite correct in the Rothe framework (it introduces an addition source
1431 of error, namely the interpolation), it is nevertheless what
1432 almost everyone solving time dependent equations does. We will use
this
1433 method in @ref step_31
"step-31",
for example.
1437<a name=
"step_23-PlainProg"></a>
1438<h1> The plain program</h1>
1439@include
"step-23.cc"
void attach_dof_handler(const DoFHandler< dim, spacedim > &)
numbers::NumberTraits< double >::real_type get_time() const
virtual RangeNumberType value(const Point< dim > &p, const unsigned int component=0) const
unsigned int n_active_cells() const
void refine_global(const unsigned int times=1)
__global__ void set(Number *val, const Number s, const size_type N)
#define Assert(cond, exc)
static ::ExceptionBase & ExcIndexRange(std::size_t arg1, std::size_t arg2, std::size_t arg3)
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)
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.
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
VectorType::value_type * end(VectorType &V)
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
void run(const Iterator &begin, const std_cxx20::type_identity_t< Iterator > &end, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length, const unsigned int chunk_size)
int(& functions)(const void *v1, const void *v2)
static constexpr double PI
::VectorizedArray< Number, width > cos(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sin(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
DataOutBase::CompressionLevel compression_level
DEAL_II_HOST constexpr SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)