553 *
const double m = 0.5;
554 *
const double c1 = 0.;
555 *
const double c2 = 0.;
556 *
const double factor =
558 *
double result = 1.;
559 *
for (
unsigned int d = 0;
d < dim; ++
d)
570 * <a name=
"step_48-SineGordonProblemclass"></a>
571 * <h3>SineGordonProblem
class</h3>
575 * This is the main
class that builds on the class in @ref step_25
"step-25". However, we
577 * the geometry data. Also, we use a distributed
triangulation in
this
582 *
class SineGordonProblem
585 * SineGordonProblem();
591 *
void make_grid_and_dofs();
592 *
void output_results(
const unsigned int timestep_number);
594 * #ifdef DEAL_II_WITH_P4EST
612 *
const unsigned int n_global_refinements;
613 *
double time, time_step;
614 *
const double final_time;
615 *
const double cfl_number;
616 *
const unsigned int output_timestep_skip;
623 * <a name=
"step_48-SineGordonProblemSineGordonProblem"></a>
624 * <h4>SineGordonProblem::SineGordonProblem</h4>
628 * This is the constructor of the SineGordonProblem
class. The time interval
629 * and time step size are defined here. Moreover, we use the degree of the
630 * finite element that we defined at the top of the program to initialize a
631 *
FE_Q finite element based on Gauss-Lobatto support points. These points
632 * are convenient because in conjunction with a
QGaussLobatto quadrature
634 * compromising accuracy too much (note that the integration is inexact,
635 * though), see also the discussion in the introduction. Note that
FE_Q
636 * selects the Gauss-Lobatto nodal points by
default due to their improved
637 * conditioning versus equidistant points. To make things more
explicit, we
638 * state the selection of the nodal points nonetheless.
642 * SineGordonProblem<dim>::SineGordonProblem()
644 * #ifdef DEAL_II_WITH_P4EST
649 * , n_global_refinements(10 - 2 * dim)
653 * , cfl_number(.1 / fe_degree)
654 * , output_timestep_skip(200)
660 * <a name=
"step_48-SineGordonProblemmake_grid_and_dofs"></a>
661 * <h4>SineGordonProblem::make_grid_and_dofs</h4>
665 * As in @ref step_25
"step-25" this functions sets up a cube grid in <code>dim</code>
666 * dimensions of extent @f$[-15,15]@f$. We
refine the mesh more in the
center of
667 * the domain since the solution is concentrated there. We
first refine all
668 * cells whose
center is within a radius of 11, and then
refine once more
669 *
for a radius 6. This simple ad hoc refinement could be done better by
670 * adapting the mesh to the solution
using error estimators during the time
671 * stepping as done in other example programs, and
using
677 *
void SineGordonProblem<dim>::make_grid_and_dofs()
682 *
for (
const auto &cell :
triangulation.active_cell_iterators())
683 * if (cell->is_locally_owned())
685 * cell->set_refine_flag();
688 *
for (
const auto &cell :
triangulation.active_cell_iterators())
689 * if (cell->is_locally_owned())
691 * cell->set_refine_flag();
695 * pcout <<
" Number of global active cells: "
698 * dof_handler.distribute_dofs(fe);
700 * pcout <<
" Number of degrees of freedom: " << dof_handler.n_dofs()
706 * We generate hanging node constraints
for ensuring continuity of the
707 * solution. As in @ref step_40
"step-40", we need to equip the constraint
matrix with
708 * the
IndexSet of locally active and locally relevant degrees of freedom
709 * to avoid it consuming too much memory
for big problems. Next, the
710 *
MatrixFree object for the problem is
set up. Note that we specify a
711 * particular scheme
for shared-memory parallelization (hence one would
712 * use multithreading
for intra-node parallelism and not
MPI; we here
713 * choose the standard option —
if we wanted to disable shared
714 * memory parallelization even in
case where there is more than one TBB
715 * thread available in the program, we would choose
717 * instead of
using the
default QGauss quadrature argument, we supply a
719 * behavior. Finally, three solution vectors are initialized.
MatrixFree
720 * expects a particular layout of ghost indices (as it handles index
721 * access in
MPI-local
numbers that need to match between the vector and
722 *
MatrixFree), so we just ask it to initialize the vectors to be sure the
723 * ghost exchange is properly handled.
726 * locally_relevant_dofs =
728 * constraints.clear();
729 * constraints.reinit(dof_handler.locally_owned_dofs(), locally_relevant_dofs);
731 * constraints.close();
737 * matrix_free_data.
reinit(mapping,
743 * matrix_free_data.initialize_dof_vector(solution);
744 * old_solution.reinit(solution);
745 * old_old_solution.reinit(solution);
753 * <a name=
"step_48-SineGordonProblemoutput_results"></a>
754 * <h4>SineGordonProblem::output_results</h4>
758 * This function prints the
norm of the solution and writes the solution
759 * vector to a file. The
norm is standard (except
for the fact that we need
760 * to accumulate the norms over all processors
for the
parallel grid which
762 *
second is similar to what we did in @ref step_40
"step-40" or @ref step_37
"step-37". Note that we can
763 * use the same vector
for output as the one used during computations: The
765 * all locally owned cells (
this is what is needed in the local evaluations,
766 * too), including ghost vector entries on these cells. This is the only
768 * as well as in
DataOut. The only action to take at this point is to make
769 * sure that the vector updates its ghost values before we read from
770 * them, and to reset ghost values once done. This is a feature present only
772 *
PETSc and Trilinos, on the other hand, need to be copied to special
773 * vectors including ghost values (see the relevant section in @ref step_40 "step-40"). If
774 * we also wanted to access all degrees of freedom on ghost cells (e.g. when
775 * computing error estimators that use the jump of solution over cell
776 * boundaries), we would need more information and create a vector
777 * initialized with locally relevant dofs just as in @ref step_40 "step-40". Observe also
778 * that we need to distribute constraints for output - they are not filled
779 * during computations (rather, they are interpolated on the fly in the
785 * SineGordonProblem<dim>::output_results(const
unsigned int timestep_number)
787 * constraints.distribute(solution);
790 * solution.update_ghost_values();
798 *
const double solution_norm =
803 * pcout <<
" Time:" << std::setw(8) << std::setprecision(3) << time
804 * <<
", solution norm: " << std::setprecision(5) << std::setw(7)
805 * << solution_norm << std::endl;
810 * data_out.add_data_vector(solution,
"solution");
811 * data_out.build_patches(mapping);
813 * data_out.write_vtu_with_pvtu_record(
814 *
"./",
"solution", timestep_number, MPI_COMM_WORLD, 3);
816 * solution.zero_out_ghost_values();
823 * <a name=
"step_48-SineGordonProblemrun"></a>
824 * <h4>SineGordonProblem::run</h4>
828 * This function is called by the main function and steps into the
829 * subroutines of the
class.
833 * After printing some information about the
parallel setup, the
first
834 * action is to
set up the grid and the cell
operator. Then, the time step
835 * is computed from the CFL number given in the constructor and the finest
836 * mesh size. The finest mesh size is computed as the
diameter of the last
838 * the mesh. This is only possible
for meshes where all elements on a
level
839 * have the same size, otherwise, one needs to
loop over all cells. Note
840 * that we need to query all the processors
for their finest cell since
841 * not all processors might hold a region where the mesh is at the finest
842 *
level. Then, we readjust the time step a little to hit the
final time
847 *
void SineGordonProblem<dim>::run()
850 * pcout <<
"Number of MPI ranks: "
852 * pcout <<
"Number of threads on each rank: "
855 *
const unsigned int n_vect_bits = 8 *
sizeof(double) * n_vect_doubles;
856 * pcout <<
"Vectorization over " << n_vect_doubles
857 * <<
" doubles = " << n_vect_bits <<
" bits ("
862 * make_grid_and_dofs();
864 *
const double local_min_cell_diameter =
866 *
const double global_min_cell_diameter =
868 * time_step = cfl_number * global_min_cell_diameter;
869 * time_step = (final_time - time) / (
int((final_time - time) / time_step));
870 * pcout <<
" Time step size: " << time_step
871 * <<
", finest cell: " << global_min_cell_diameter << std::endl
876 * Next the
initial value is
set. Since we have a two-step time stepping
877 * method, we also need a
value of the solution at time-time_step. For
878 * accurate results, one would need to compute
this from the time
879 * derivative of the solution at
initial time, but here we ignore
this
885 * We then go on by writing the
initial state to file and collecting
886 * the two starting solutions in a <tt>std::vector</tt> of pointers that
887 * get later consumed by the SineGordonOperation::apply() function. Next,
888 * an instance of the <code> SineGordonOperation class </code> based on
889 * the finite element degree specified at the top of this file is set up.
894 * InitialCondition<dim>(1, time),
898 * InitialCondition<dim>(1, time - time_step),
903 * previous_solutions({&old_solution, &old_old_solution});
905 * SineGordonOperation<dim, fe_degree> sine_gordon_op(matrix_free_data,
910 * Now
loop over the time steps. In each iteration, we
shift the solution
911 * vectors by one and call the `apply` function of the
912 * `SineGordonOperator`
class. Then, we write the solution to a file. We
913 * clock the wall times
for the computational time needed as wall as the
914 * time needed to create the output and report the
numbers when the time
915 * stepping is finished.
919 * Note how
this shift is implemented: We simply call the
swap method on
920 * the two vectors which swaps only some pointers without the need to
copy
921 * data around, a relatively expensive operation within an
explicit time
922 * stepping method. Let us see what happens in more detail: First, we
923 * exchange <code>old_solution</code> with <code>old_old_solution</code>,
924 * which means that <code>old_old_solution</code> gets
925 * <code>old_solution</code>, which is what we expect. Similarly,
926 * <code>old_solution</code> gets the content from <code>solution</code>
927 * in the next step. After
this, <code>solution</code> holds
928 * <code>old_old_solution</code>, but that will be overwritten during
this
932 *
unsigned int timestep_number = 1;
936 *
double output_time = 0;
937 *
for (time += time_step; time <= final_time;
938 * time += time_step, ++timestep_number)
941 * old_old_solution.swap(old_solution);
942 * old_solution.swap(solution);
943 * sine_gordon_op.apply(solution, previous_solutions);
944 * wtime += timer.wall_time();
947 *
if (timestep_number % output_timestep_skip == 0)
948 * output_results(timestep_number / output_timestep_skip);
950 * output_time += timer.wall_time();
953 * output_results(timestep_number / output_timestep_skip + 1);
954 * output_time += timer.wall_time();
957 * <<
" Performed " << timestep_number <<
" time steps." << std::endl;
959 * pcout <<
" Average wallclock time per time step: "
960 * << wtime / timestep_number <<
's' << std::endl;
962 * pcout <<
" Spent " << output_time <<
"s on output and " << wtime
963 * <<
"s on computations." << std::endl;
972 * <a name=
"step_48-Thecodemaincodefunction"></a>
973 * <h3>The <code>main</code> function</h3>
977 * As in @ref step_40
"step-40", we initialize
MPI at the start of the program. Since we will
978 * in
general mix
MPI parallelization with threads, we also
set the third
979 * argument in MPI_InitFinalize that controls the number of threads to an
980 *
invalid number, which means that the TBB library chooses the number of
981 * threads automatically, typically to the number of available cores in the
982 * system. As an alternative, you can also
set this number manually
if you
983 * want to
set a specific number of threads (
e.g. when
MPI-only is required).
986 *
int main(
int argc,
char **argv)
988 *
using namespace Step48;
996 * SineGordonProblem<dimension> sg_problem;
999 *
catch (std::exception &exc)
1001 * std::cerr << std::endl
1003 * <<
"----------------------------------------------------"
1005 * std::cerr <<
"Exception on processing: " << std::endl
1006 * << exc.what() << std::endl
1007 * <<
"Aborting!" << std::endl
1008 * <<
"----------------------------------------------------"
1015 * std::cerr << std::endl
1017 * <<
"----------------------------------------------------"
1019 * std::cerr <<
"Unknown exception!" << std::endl
1020 * <<
"Aborting!" << std::endl
1021 * <<
"----------------------------------------------------"
1029<a name=
"step_48-Results"></a><h1>Results</h1>
1032<a name=
"step_48-Comparisonwithasparsematrix"></a><h3>Comparison with a sparse
matrix</h3>
1035In order to demonstrate the gain in
using the
MatrixFree class instead of
1036the standard <code>deal.II</code> assembly routines for evaluating the
1037information from old time steps, we study a simple
serial run of the code on a
1038nonadaptive mesh. Since much time is spent on evaluating the sine function, we
1039do not only show the
numbers of the full sine-Gordon equation but also for the
1040wave equation (the sine-term skipped from the sine-Gordon equation). We use
1041both
second and fourth order elements. The results are summarized in the
1044<table align=
"center" class=
"doxtable">
1047 <th colspan=
"3">wave equation</th>
1048 <th colspan=
"2">sine-Gordon</th>
1059 <td>2D, @f$\mathcal{Q}_2@f$</td>
1060 <td align=
"right"> 0.0106</td>
1061 <td align=
"right"> 0.00971</td>
1062 <td align=
"right"> 0.109</td>
1063 <td align=
"right"> 0.0243</td>
1064 <td align=
"right"> 0.124</td>
1067 <td>2D, @f$\mathcal{Q}_4@f$</td>
1068 <td align=
"right"> 0.0328</td>
1069 <td align=
"right"> 0.0706</td>
1070 <td align=
"right"> 0.528</td>
1071 <td align=
"right"> 0.0714</td>
1072 <td align=
"right"> 0.502</td>
1075 <td>3D, @f$\mathcal{Q}_2@f$</td>
1076 <td align=
"right"> 0.0151</td>
1077 <td align=
"right"> 0.0320</td>
1078 <td align=
"right"> 0.331</td>
1079 <td align=
"right"> 0.0376</td>
1080 <td align=
"right"> 0.364</td>
1083 <td>3D, @f$\mathcal{Q}_4@f$</td>
1084 <td align=
"right"> 0.0918</td>
1085 <td align=
"right"> 0.844</td>
1086 <td align=
"right"> 6.83</td>
1087 <td align=
"right"> 0.194</td>
1088 <td align=
"right"> 6.95</td>
1092It is apparent that the
matrix-
free code outperforms the standard assembly
1093routines in deal.II by far. In 3D and
for fourth order elements, one
operator
1094evaluation is also almost ten times as fast as a sparse
matrix-vector
1097<a name=
"step_48-Parallelrunin2Dand3D"></a><h3>Parallel
run in 2D and 3D</h3>
1100We start with the program output obtained on a workstation with 12 cores / 24
1101threads (one Intel Xeon E5-2687W v4 CPU running at 3.2 GHz, hyperthreading
1102enabled), running the program in release mode:
1105Number of
MPI ranks: 1
1106Number of threads on each rank: 24
1107Vectorization over 4 doubles = 256 bits (AVX)
1109 Number of global active cells: 15412
1110 Number of degrees of freedom: 249065
1111 Time step size: 0.00292997, finest cell: 0.117188
1113 Time: -10, solution
norm: 9.5599
1114 Time: -9.41, solution
norm: 17.678
1115 Time: -8.83, solution
norm: 23.504
1116 Time: -8.24, solution
norm: 27.5
1117 Time: -7.66, solution
norm: 29.513
1118 Time: -7.07, solution
norm: 29.364
1119 Time: -6.48, solution
norm: 27.23
1120 Time: -5.9, solution
norm: 23.527
1121 Time: -5.31, solution
norm: 18.439
1122 Time: -4.73, solution
norm: 11.935
1123 Time: -4.14, solution
norm: 5.5284
1124 Time: -3.55, solution
norm: 8.0354
1125 Time: -2.97, solution
norm: 14.707
1126 Time: -2.38, solution
norm: 20
1127 Time: -1.8, solution
norm: 22.834
1128 Time: -1.21, solution
norm: 22.771
1129 Time: -0.624, solution
norm: 20.488
1130 Time: -0.0381, solution
norm: 16.697
1131 Time: 0.548, solution
norm: 11.221
1132 Time: 1.13, solution
norm: 5.3912
1133 Time: 1.72, solution
norm: 8.4528
1134 Time: 2.31, solution
norm: 14.335
1135 Time: 2.89, solution
norm: 18.555
1136 Time: 3.48, solution
norm: 20.894
1137 Time: 4.06, solution
norm: 21.305
1138 Time: 4.65, solution
norm: 19.903
1139 Time: 5.24, solution
norm: 16.864
1140 Time: 5.82, solution
norm: 12.223
1141 Time: 6.41, solution
norm: 6.758
1142 Time: 6.99, solution
norm: 7.2423
1143 Time: 7.58, solution
norm: 12.888
1144 Time: 8.17, solution
norm: 17.273
1145 Time: 8.75, solution
norm: 19.654
1146 Time: 9.34, solution
norm: 19.838
1147 Time: 9.92, solution
norm: 17.964
1148 Time: 10, solution
norm: 17.595
1150 Performed 6826 time steps.
1151 Average wallclock time per time step: 0.0013453s
1152 Spent 14.976s on output and 9.1831s on computations.
1155In 3D, the respective output looks like
1158Number of
MPI ranks: 1
1159Number of threads on each rank: 24
1160Vectorization over 4 doubles = 256 bits (AVX)
1162 Number of global active cells: 17592
1163 Number of degrees of freedom: 1193881
1164 Time step size: 0.0117233, finest cell: 0.46875
1166 Time: -10, solution
norm: 29.558
1167 Time: -7.66, solution
norm: 129.13
1168 Time: -5.31, solution
norm: 67.753
1169 Time: -2.97, solution
norm: 79.245
1170 Time: -0.621, solution
norm: 123.52
1171 Time: 1.72, solution
norm: 43.525
1172 Time: 4.07, solution
norm: 93.285
1173 Time: 6.41, solution
norm: 97.722
1174 Time: 8.76, solution
norm: 36.734
1175 Time: 10, solution
norm: 94.115
1177 Performed 1706 time steps.
1178 Average wallclock time per time step: 0.0084542s
1179 Spent 16.766s on output and 14.423s on computations.
1182It takes 0.008 seconds
for one time step with more than a million
1183degrees of freedom (note that we would need many processors to reach such
1184numbers when solving linear systems).
1186If we replace the thread-parallelization by a pure
MPI parallelization, the
1189\$ mpirun -n 24 ./step-48
1190Number of
MPI ranks: 24
1191Number of threads on each rank: 1
1192Vectorization over 4 doubles = 256 bits (AVX)
1194 Performed 1706 time steps.
1195 Average wallclock time per time step: 0.0051747s
1196 Spent 2.0535s on output and 8.828s on computations.
1199We observe a dramatic speedup
for the output (which makes sense, given that
1200most code of the output is not parallelized via threads, whereas it is
for
1201MPI), but less than the theoretical factor of 12 we would expect from the
1202parallelism. More interestingly, the computations also get faster when
1203switching from the threads-only variant to the
MPI-only variant. This is a
12052019). The main reason is that the decisions regarding work on conflicting
1206cell batches made to enable execution in
parallel are overly pessimistic:
1207While they ensure that no work on neighboring cells is done on different
1208threads at the same time,
this conservative setting implies that data from
1209neighboring cells is also evicted from caches by the time neighbors get
1210touched. Furthermore, the current scheme is not able to provide a
constant
1211load
for all 24 threads
for the given mesh with 17,592 cells.
1213The current program allows to also mix
MPI parallelization with thread
1214parallelization. This is most beneficial when running programs on clusters
1215with multiple nodes,
using MPI for the inter-node parallelization and threads
1216for the intra-node parallelization. On the workstation used above, we can
run
1217threads in the hyperthreading region (i.e.,
using 2 threads
for each of the 12
1218MPI ranks). An important setting
for mixing
MPI with threads is to ensure
1219proper binning of tasks to CPUs. On many clusters the placing is either
1220automatically via the `mpirun/mpiexec` environment, or there can be manual
1221settings. Here, we simply report the
run times the plain version of the
1222program (noting that things could be improved towards the timings of the
1223MPI-only program when proper pinning is done):
1225\$ mpirun -n 12 ./step-48
1226Number of
MPI ranks: 12
1227Number of threads on each rank: 2
1228Vectorization over 4 doubles = 256 bits (AVX)
1230 Performed 1706 time steps.
1231 Average wallclock time per time step: 0.0056651s
1232 Spent 2.5175s on output and 9.6646s on computations.
1237<a name=
"step_48-Possibilitiesforextensions"></a><h3>Possibilities for extensions</h3>
1240There are several things in this program that could be improved to make it
1241even more efficient (besides improved boundary conditions and physical
1242stuff as discussed in @ref step_25
"step-25"):
1244<ul> <li> <
b>Faster evaluation of sine terms:</
b> As becomes obvious
1245 from the comparison of the plain wave equation and the sine-Gordon
1246 equation above, the evaluation of the sine terms dominates the total
1247 time for the finite element operator application. There are a few
1248 reasons for this: Firstly, the deal.II sine computation of a
1250 the operator application). This could be cured by handing the sine
1251 computation to a library with
vectorized sine computations like
1252 Intel
's math kernel library (MKL). By using the function
1253 <code>vdSin</code> in MKL, the program uses half the computing time
1254 in 2D and 40 percent less time in 3D. On the other hand, the sine
1255 computation is structurally much more complicated than the simple
1256 arithmetic operations like additions and multiplications in the rest
1257 of the local operation.
1259 <li> <b>Higher order time stepping:</b> While the implementation allows for
1260 arbitrary order in the spatial part (by adjusting the degree of the finite
1261 element), the time stepping scheme is a standard second-order leap-frog
1262 scheme. Since solutions in wave propagation problems are usually very
1263 smooth, the error is likely dominated by the time stepping part. Of course,
1264 this could be cured by using smaller time steps (at a fixed spatial
1265 resolution), but it would be more efficient to use higher order time
1266 stepping as well. While it would be straight-forward to do so for a
1267 first-order system (use some Runge–Kutta scheme of higher order,
1268 probably combined with adaptive time step selection like the <a
1269 href="http://en.wikipedia.org/wiki/Dormand%E2%80%93Prince_method">Dormand–Prince
1270 method</a>), it is more challenging for the second-order formulation. At
1271 least in the finite difference community, people usually use the PDE to find
1272 spatial correction terms that improve the temporal error.
1277<a name="step_48-PlainProg"></a>
1278<h1> The plain program</h1>
1279@include "step-48.cc"
void attach_dof_handler(const DoFHandler< dim, spacedim > &)
void reinit(const MappingType &mapping, const DoFHandler< dim > &dof_handler, const AffineConstraints< number2 > &constraint, const QuadratureType &quad, const AdditionalData &additional_data=AdditionalData())
static unsigned int n_threads()
unsigned int n_active_cells() const
void refine_global(const unsigned int times=1)
cell_iterator last() const
static constexpr std::size_t size()
virtual types::global_cell_index n_global_active_cells() const override
virtual void execute_coarsening_and_refinement() override
__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_hanging_node_constraints(const DoFHandler< dim, spacedim > &dof_handler, AffineConstraints< number > &constraints)
void hyper_cube(Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false)
void refine(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold, const unsigned int max_to_mark=numbers::invalid_unsigned_int)
@ matrix
Contents is actually a matrix.
@ diagonal
Matrix is diagonal.
@ general
No special properties.
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
std::vector< unsigned int > serial(const std::vector< unsigned int > &targets, const std::function< RequestType(const unsigned int)> &create_request, const std::function< AnswerType(const unsigned int, const RequestType &)> &answer_request, const std::function< void(const unsigned int, const AnswerType &)> &process_answer, const MPI_Comm comm)
unsigned int n_mpi_processes(const MPI_Comm mpi_communicator)
T max(const T &t, const MPI_Comm mpi_communicator)
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
std::string get_current_vectorization_level()
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)
void copy(const T *begin, const T *end, U *dest)
int(& functions)(const void *v1, const void *v2)
static const unsigned int invalid_unsigned_int
::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
void swap(SmartPointer< T, P > &t1, SmartPointer< T, Q > &t2)
TasksParallelScheme tasks_parallel_scheme