Reference documentation for deal.II version 9.5.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\}}\)
Loading...
Searching...
No Matches
step-8.h
Go to the documentation of this file.
1 false);
746 *   sparsity_pattern.copy_from(dsp);
747 *  
748 *   system_matrix.reinit(sparsity_pattern);
749 *   }
750 *  
751 *  
752 * @endcode
753 *
754 *
755 * <a name="ElasticProblemassemble_system"></a>
756 * <h4>ElasticProblem::assemble_system</h4>
757 *
758
759 *
760 * The big changes in this program are in the creation of matrix and right
761 * hand side, since they are problem-dependent. We will go through that
762 * process step-by-step, since it is a bit more complicated than in previous
763 * examples.
764 *
765
766 *
767 * The first parts of this function are the same as before, however: setting
768 * up a suitable quadrature formula, initializing an FEValues object for the
769 * (vector-valued) finite element we use as well as the quadrature object,
770 * and declaring a number of auxiliary arrays. In addition, we declare the
771 * ever same two abbreviations: <code>n_q_points</code> and
772 * <code>dofs_per_cell</code>. The number of degrees of freedom per cell we
773 * now obviously ask from the composed finite element rather than from the
774 * underlying scalar Q1 element. Here, it is <code>dim</code> times the
775 * number of degrees of freedom per cell of the Q1 element, though this is
776 * not explicit knowledge we need to care about:
777 *
778 * @code
779 *   template <int dim>
780 *   void ElasticProblem<dim>::assemble_system()
781 *   {
782 *   QGauss<dim> quadrature_formula(fe.degree + 1);
783 *  
784 *   FEValues<dim> fe_values(fe,
785 *   quadrature_formula,
788 *  
789 *   const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
790 *   const unsigned int n_q_points = quadrature_formula.size();
791 *  
792 *   FullMatrix<double> cell_matrix(dofs_per_cell, dofs_per_cell);
793 *   Vector<double> cell_rhs(dofs_per_cell);
794 *  
795 *   std::vector<types::global_dof_index> local_dof_indices(dofs_per_cell);
796 *  
797 * @endcode
798 *
799 * As was shown in previous examples as well, we need a place where to
800 * store the values of the coefficients at all the quadrature points on a
801 * cell. In the present situation, we have two coefficients, lambda and
802 * mu.
803 *
804 * @code
805 *   std::vector<double> lambda_values(n_q_points);
806 *   std::vector<double> mu_values(n_q_points);
807 *  
808 * @endcode
809 *
810 * Well, we could as well have omitted the above two arrays since we will
811 * use constant coefficients for both lambda and mu, which can be declared
812 * like this. They both represent functions always returning the constant
813 * value 1.0. Although we could omit the respective factors in the
814 * assemblage of the matrix, we use them here for purpose of
815 * demonstration.
816 *
817 * @code
819 *  
820 * @endcode
821 *
822 * Like the two constant functions above, we will call the function
823 * right_hand_side just once per cell to make things simpler.
824 *
825 * @code
826 *   std::vector<Tensor<1, dim>> rhs_values(n_q_points);
827 *  
828 * @endcode
829 *
830 * Now we can begin with the loop over all cells:
831 *
832 * @code
833 *   for (const auto &cell : dof_handler.active_cell_iterators())
834 *   {
835 *   cell_matrix = 0;
836 *   cell_rhs = 0;
837 *  
838 *   fe_values.reinit(cell);
839 *  
840 * @endcode
841 *
842 * Next we get the values of the coefficients at the quadrature
843 * points. Likewise for the right hand side:
844 *
845 * @code
846 *   lambda.value_list(fe_values.get_quadrature_points(), lambda_values);
847 *   mu.value_list(fe_values.get_quadrature_points(), mu_values);
848 *   right_hand_side(fe_values.get_quadrature_points(), rhs_values);
849 *  
850 * @endcode
851 *
852 * Then assemble the entries of the local @ref GlossStiffnessMatrix "stiffness matrix" and right
853 * hand side vector. This follows almost one-to-one the pattern
854 * described in the introduction of this example. One of the few
855 * comments in place is that we can compute the number
856 * <code>comp(i)</code>, i.e. the index of the only nonzero vector
857 * component of shape function <code>i</code> using the
858 * <code>fe.system_to_component_index(i).first</code> function call
859 * below.
860 *
861
862 *
863 * (By accessing the <code>first</code> variable of the return value
864 * of the <code>system_to_component_index</code> function, you might
865 * already have guessed that there is more in it. In fact, the
866 * function returns a <code>std::pair@<unsigned int, unsigned
867 * int@></code>, of which the first element is <code>comp(i)</code>
868 * and the second is the value <code>base(i)</code> also noted in the
869 * introduction, i.e. the index of this shape function within all the
870 * shape functions that are nonzero in this component,
871 * i.e. <code>base(i)</code> in the diction of the introduction. This
872 * is not a number that we are usually interested in, however.)
873 *
874
875 *
876 * With this knowledge, we can assemble the local matrix
877 * contributions:
878 *
879 * @code
880 *   for (const unsigned int i : fe_values.dof_indices())
881 *   {
882 *   const unsigned int component_i =
883 *   fe.system_to_component_index(i).first;
884 *  
885 *   for (const unsigned int j : fe_values.dof_indices())
886 *   {
887 *   const unsigned int component_j =
888 *   fe.system_to_component_index(j).first;
889 *  
890 *   for (const unsigned int q_point :
891 *   fe_values.quadrature_point_indices())
892 *   {
893 *   cell_matrix(i, j) +=
894 * @endcode
895 *
896 * The first term is @f$(\lambda \partial_i u_i, \partial_j
897 * v_j) + (\mu \partial_i u_j, \partial_j v_i)@f$. Note
898 * that <code>shape_grad(i,q_point)</code> returns the
899 * gradient of the only nonzero component of the i-th
900 * shape function at quadrature point q_point. The
901 * component <code>comp(i)</code> of the gradient, which
902 * is the derivative of this only nonzero vector
903 * component of the i-th shape function with respect to
904 * the comp(i)th coordinate is accessed by the appended
905 * brackets.
906 *
907 * @code
908 *   (
909 *   (fe_values.shape_grad(i, q_point)[component_i] *
910 *   fe_values.shape_grad(j, q_point)[component_j] *
911 *   lambda_values[q_point])
912 *   +
913 *   (fe_values.shape_grad(i, q_point)[component_j] *
914 *   fe_values.shape_grad(j, q_point)[component_i] *
915 *   mu_values[q_point])
916 *   +
917 * @endcode
918 *
919 * The second term is @f$(\mu \nabla u_i, \nabla
920 * v_j)@f$. We need not access a specific component of
921 * the gradient, since we only have to compute the
922 * scalar product of the two gradients, of which an
923 * overloaded version of <tt>operator*</tt> takes
924 * care, as in previous examples.
925 *
926
927 *
928 * Note that by using the <tt>?:</tt> operator, we only
929 * do this if <tt>component_i</tt> equals
930 * <tt>component_j</tt>, otherwise a zero is added
931 * (which will be optimized away by the compiler).
932 *
933 * @code
934 *   ((component_i == component_j) ?
935 *   (fe_values.shape_grad(i, q_point) *
936 *   fe_values.shape_grad(j, q_point) *
937 *   mu_values[q_point]) :
938 *   0)
939 *   ) *
940 *   fe_values.JxW(q_point);
941 *   }
942 *   }
943 *   }
944 *  
945 * @endcode
946 *
947 * Assembling the right hand side is also just as discussed in the
948 * introduction:
949 *
950 * @code
951 *   for (const unsigned int i : fe_values.dof_indices())
952 *   {
953 *   const unsigned int component_i =
954 *   fe.system_to_component_index(i).first;
955 *  
956 *   for (const unsigned int q_point :
957 *   fe_values.quadrature_point_indices())
958 *   cell_rhs(i) += fe_values.shape_value(i, q_point) *
959 *   rhs_values[q_point][component_i] *
960 *   fe_values.JxW(q_point);
961 *   }
962 *  
963 * @endcode
964 *
965 * The transfer from local degrees of freedom into the global matrix
966 * and right hand side vector does not depend on the equation under
967 * consideration, and is thus the same as in all previous
968 * examples.
969 *
970 * @code
971 *   cell->get_dof_indices(local_dof_indices);
972 *   constraints.distribute_local_to_global(
973 *   cell_matrix, cell_rhs, local_dof_indices, system_matrix, system_rhs);
974 *   }
975 *   }
976 *  
977 *  
978 *  
979 * @endcode
980 *
981 *
982 * <a name="ElasticProblemsolve"></a>
983 * <h4>ElasticProblem::solve</h4>
984 *
985
986 *
987 * The solver does not care about where the system of equations comes, as
988 * long as it stays positive definite and symmetric (which are the
989 * requirements for the use of the CG solver), which the system indeed
990 * is. Therefore, we need not change anything.
991 *
992 * @code
993 *   template <int dim>
994 *   void ElasticProblem<dim>::solve()
995 *   {
996 *   SolverControl solver_control(1000, 1e-12);
997 *   SolverCG<Vector<double>> cg(solver_control);
998 *  
999 *   PreconditionSSOR<SparseMatrix<double>> preconditioner;
1000 *   preconditioner.initialize(system_matrix, 1.2);
1001 *  
1002 *   cg.solve(system_matrix, solution, system_rhs, preconditioner);
1003 *  
1004 *   constraints.distribute(solution);
1005 *   }
1006 *  
1007 *  
1008 * @endcode
1009 *
1010 *
1011 * <a name="ElasticProblemrefine_grid"></a>
1012 * <h4>ElasticProblem::refine_grid</h4>
1013 *
1014
1015 *
1016 * The function that does the refinement of the grid is the same as in the
1017 * @ref step_6 "step-6" example. The quadrature formula is adapted to the linear elements
1018 * again. Note that the error estimator by default adds up the estimated
1019 * obtained from all components of the finite element solution, i.e., it
1020 * uses the displacement in all directions with the same weight. If we would
1021 * like the grid to be adapted to the x-displacement only, we could pass the
1022 * function an additional parameter which tells it to do so and do not
1023 * consider the displacements in all other directions for the error
1024 * indicators. However, for the current problem, it seems appropriate to
1025 * consider all displacement components with equal weight.
1026 *
1027 * @code
1028 *   template <int dim>
1029 *   void ElasticProblem<dim>::refine_grid()
1030 *   {
1031 *   Vector<float> estimated_error_per_cell(triangulation.n_active_cells());
1032 *  
1033 *   KellyErrorEstimator<dim>::estimate(dof_handler,
1034 *   QGauss<dim - 1>(fe.degree + 1),
1035 *   {},
1036 *   solution,
1037 *   estimated_error_per_cell);
1038 *  
1040 *   estimated_error_per_cell,
1041 *   0.3,
1042 *   0.03);
1043 *  
1044 *   triangulation.execute_coarsening_and_refinement();
1045 *   }
1046 *  
1047 *  
1048 * @endcode
1049 *
1050 *
1051 * <a name="ElasticProblemoutput_results"></a>
1052 * <h4>ElasticProblem::output_results</h4>
1053 *
1054
1055 *
1056 * The output happens mostly as has been shown in previous examples
1057 * already. The only difference is that the solution function is vector
1058 * valued. The DataOut class takes care of this automatically, but we have
1059 * to give each component of the solution vector a different name.
1060 *
1061
1062 *
1063 * To do this, the DataOut::add_vector() function wants a vector of
1064 * strings. Since the number of components is the same as the number
1065 * of dimensions we are working in, we use the <code>switch</code>
1066 * statement below.
1067 *
1068
1069 *
1070 * We note that some graphics programs have restriction on what
1071 * characters are allowed in the names of variables. deal.II therefore
1072 * supports only the minimal subset of these characters that is supported
1073 * by all programs. Basically, these are letters, numbers, underscores,
1074 * and some other characters, but in particular no whitespace and
1075 * minus/hyphen. The library will throw an exception otherwise, at least
1076 * if in debug mode.
1077 *
1078
1079 *
1080 * After listing the 1d, 2d, and 3d case, it is good style to let the
1081 * program die if we run upon a case which we did not consider. Remember
1082 * that the Assert macro generates an exception if the condition in the
1083 * first parameter is not satisfied. Of course, the condition
1084 * <code>false</code> can never be satisfied, so the program will always
1085 * abort whenever it gets to the default statement:
1086 *
1087 * @code
1088 *   template <int dim>
1089 *   void ElasticProblem<dim>::output_results(const unsigned int cycle) const
1090 *   {
1091 *   DataOut<dim> data_out;
1092 *   data_out.attach_dof_handler(dof_handler);
1093 *  
1094 *   std::vector<std::string> solution_names;
1095 *   switch (dim)
1096 *   {
1097 *   case 1:
1098 *   solution_names.emplace_back("displacement");
1099 *   break;
1100 *   case 2:
1101 *   solution_names.emplace_back("x_displacement");
1102 *   solution_names.emplace_back("y_displacement");
1103 *   break;
1104 *   case 3:
1105 *   solution_names.emplace_back("x_displacement");
1106 *   solution_names.emplace_back("y_displacement");
1107 *   solution_names.emplace_back("z_displacement");
1108 *   break;
1109 *   default:
1110 *   Assert(false, ExcNotImplemented());
1111 *   }
1112 *  
1113 * @endcode
1114 *
1115 * After setting up the names for the different components of the
1116 * solution vector, we can add the solution vector to the list of
1117 * data vectors scheduled for output. Note that the following
1118 * function takes a vector of strings as second argument, whereas
1119 * the one which we have used in all previous examples accepted a
1120 * string there. (In fact, the function we had used before would
1121 * convert the single string into a vector with only one element
1122 * and forwards that to the other function.)
1123 *
1124 * @code
1125 *   data_out.add_data_vector(solution, solution_names);
1126 *   data_out.build_patches();
1127 *  
1128 *   std::ofstream output("solution-" + std::to_string(cycle) + ".vtk");
1129 *   data_out.write_vtk(output);
1130 *   }
1131 *  
1132 *  
1133 *  
1134 * @endcode
1135 *
1136 *
1137 * <a name="ElasticProblemrun"></a>
1138 * <h4>ElasticProblem::run</h4>
1139 *
1140
1141 *
1142 * The <code>run</code> function does the same things as in @ref step_6 "step-6", for
1143 * example. This time, we use the square [-1,1]^d as domain, and we refine
1144 * it globally four times before starting the first iteration.
1145 *
1146
1147 *
1148 * The reason for refining is a bit accidental: we use the QGauss
1149 * quadrature formula with two points in each direction for integration of the
1150 * right hand side; that means that there are four quadrature points on each
1151 * cell (in 2d). If we only refine the initial grid once globally, then there
1152 * will be only four quadrature points in each direction on the
1153 * domain. However, the right hand side function was chosen to be rather
1154 * localized and in that case, by pure chance, it happens that all quadrature
1155 * points lie at points where the right hand side function is zero (in
1156 * mathematical terms, the quadrature points happen to be at points outside
1157 * the <i>support</i> of the right hand side function). The right hand side
1158 * vector computed with quadrature will then contain only zeroes (even though
1159 * it would of course be nonzero if we had computed the right hand side vector
1160 * exactly using the integral) and the solution of the system of
1161 * equations is the zero vector, i.e., a finite element function that is zero
1162 * everywhere. In a sense, we
1163 * should not be surprised that this is happening since we have chosen
1164 * an initial grid that is totally unsuitable for the problem at hand.
1165 *
1166
1167 *
1168 * The unfortunate thing is that if the discrete solution is constant, then
1169 * the error indicators computed by the KellyErrorEstimator class are zero
1170 * for each cell as well, and the call to
1171 * Triangulation::refine_and_coarsen_fixed_number() will not flag any cells
1172 * for refinement (why should it if the indicated error is zero for each
1173 * cell?). The grid in the next iteration will therefore consist of four
1174 * cells only as well, and the same problem occurs again.
1175 *
1176
1177 *
1178 * The conclusion needs to be: while of course we will not choose the
1179 * initial grid to be well-suited for the accurate solution of the problem,
1180 * we must at least choose it such that it has the chance to capture the
1181 * important features of the solution. In this case, it needs to be able to
1182 * see the right hand side. Thus, we refine globally four times. (Any larger
1183 * number of global refinement steps would of course also work.)
1184 *
1185 * @code
1186 *   template <int dim>
1187 *   void ElasticProblem<dim>::run()
1188 *   {
1189 *   for (unsigned int cycle = 0; cycle < 8; ++cycle)
1190 *   {
1191 *   std::cout << "Cycle " << cycle << ':' << std::endl;
1192 *  
1193 *   if (cycle == 0)
1194 *   {
1196 *   triangulation.refine_global(4);
1197 *   }
1198 *   else
1199 *   refine_grid();
1200 *  
1201 *   std::cout << " Number of active cells: "
1202 *   << triangulation.n_active_cells() << std::endl;
1203 *  
1204 *   setup_system();
1205 *  
1206 *   std::cout << " Number of degrees of freedom: " << dof_handler.n_dofs()
1207 *   << std::endl;
1208 *  
1209 *   assemble_system();
1210 *   solve();
1211 *   output_results(cycle);
1212 *   }
1213 *   }
1214 *   } // namespace Step8
1215 *  
1216 * @endcode
1217 *
1218 *
1219 * <a name="Thecodemaincodefunction"></a>
1220 * <h3>The <code>main</code> function</h3>
1221 *
1222
1223 *
1224 * After closing the <code>Step8</code> namespace in the last line above, the
1225 * following is the main function of the program and is again exactly like in
1226 * @ref step_6 "step-6" (apart from the changed class names, of course).
1227 *
1228 * @code
1229 *   int main()
1230 *   {
1231 *   try
1232 *   {
1233 *   Step8::ElasticProblem<2> elastic_problem_2d;
1234 *   elastic_problem_2d.run();
1235 *   }
1236 *   catch (std::exception &exc)
1237 *   {
1238 *   std::cerr << std::endl
1239 *   << std::endl
1240 *   << "----------------------------------------------------"
1241 *   << std::endl;
1242 *   std::cerr << "Exception on processing: " << std::endl
1243 *   << exc.what() << std::endl
1244 *   << "Aborting!" << std::endl
1245 *   << "----------------------------------------------------"
1246 *   << std::endl;
1247 *  
1248 *   return 1;
1249 *   }
1250 *   catch (...)
1251 *   {
1252 *   std::cerr << std::endl
1253 *   << std::endl
1254 *   << "----------------------------------------------------"
1255 *   << std::endl;
1256 *   std::cerr << "Unknown exception!" << std::endl
1257 *   << "Aborting!" << std::endl
1258 *   << "----------------------------------------------------"
1259 *   << std::endl;
1260 *   return 1;
1261 *   }
1262 *  
1263 *   return 0;
1264 *   }
1265 * @endcode
1266<a name="Results"></a><h1>Results</h1>
1267
1268
1269
1270There is not much to be said about the results of this program, other than
1271that they look nice. All images were made using VisIt from the
1272output files that the program wrote to disk. The first two pictures show
1273the @f$x@f$- and @f$y@f$-displacements as scalar components:
1274
1275<table width="100%">
1276<tr>
1277<td>
1278<img src="https://www.dealii.org/images/steps/developer/step-8.x.png" alt="">
1279</td>
1280<td>
1281<img src="https://www.dealii.org/images/steps/developer/step-8.y.png" alt="">
1282</td>
1283</tr>
1284</table>
1285
1286
1287You can clearly see the sources of @f$x@f$-displacement around @f$x=0.5@f$ and
1288@f$x=-0.5@f$, and of @f$y@f$-displacement at the origin.
1289
1290What one frequently would like to do is to show the displacement as a vector
1291field, i.e., vectors that for each point illustrate the direction and magnitude
1292of displacement. Unfortunately, that's a bit more involved. To understand why
1293this is so, remember that we have just defined our finite element as a
1294collection of two components (in <code>dim=2</code> dimensions). Nowhere have
1295we said that this is not just a pressure and a concentration (two scalar
1296quantities) but that the two components actually are the parts of a
1297vector-valued quantity, namely the displacement. Absent this knowledge, the
1298DataOut class assumes that all individual variables we print are separate
1299scalars, and VisIt and Paraview then faithfully assume that this is indeed what it is. In
1300other words, once we have written the data as scalars, there is nothing in
1301these programs that allows us to paste these two scalar fields back together as a
1302vector field. Where we would have to attack this problem is at the root,
1303namely in <code>ElasticProblem::output_results()</code>. We won't do so here but
1304instead refer the reader to the @ref step_22 "step-22" program where we show how to do this
1305for a more general situation. That said, we couldn't help generating the data
1306anyway that would show how this would look if implemented as discussed in
1307@ref step_22 "step-22". The vector field then looks like this (VisIt and Paraview
1308randomly select a few
1309hundred vertices from which to draw the vectors; drawing them from each
1310individual vertex would make the picture unreadable):
1311
1312<img src="https://www.dealii.org/images/steps/developer/step-8.vectors.png" alt="">
1313
1314
1315We note that one may have intuitively expected the
1316solution to be symmetric about the @f$x@f$- and @f$y@f$-axes since the @f$x@f$- and
1317@f$y@f$-forces are symmetric with respect to these axes. However, the force
1318considered as a vector is not symmetric and consequently neither is
1319the solution.
1320 *
1321 *
1322<a name="PlainProg"></a>
1323<h1> The plain program</h1>
1324@include "step-8.cc"
1325*/
void attach_dof_handler(const DoFHandler< dim, spacedim > &)
static void estimate(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, typename InputVector::value_type > * > &neumann_bc, const InputVector &solution, Vector< float > &error, const ComponentMask &component_mask=ComponentMask(), const Function< spacedim > *coefficients=nullptr, const unsigned int n_threads=numbers::invalid_unsigned_int, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id, const types::material_id material_id=numbers::invalid_material_id, const Strategy strategy=cell_diameter_over_24)
void initialize(const MatrixType &A, const AdditionalData &parameters=AdditionalData())
Point< 2 > second
Definition grid_out.cc:4616
Point< 2 > first
Definition grid_out.cc:4615
#define Assert(cond, exc)
void loop(ITERATOR begin, std_cxx20::type_identity_t< ITERATOR > end, DOFINFO &dinfo, INFOBOX &info, const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &cell_worker, const std::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &boundary_worker, const std::function< void(DOFINFO &, DOFINFO &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, ASSEMBLER &assembler, const LoopControl &lctrl=LoopControl())
Definition loop.h:439
@ update_values
Shape function values.
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
const Event initial
Definition event.cc:65
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)
void refine_and_coarsen_fixed_number(Triangulation< dim, spacedim > &triangulation, const Vector< Number > &criteria, const double top_fraction_of_cells, const double bottom_fraction_of_cells, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max())
@ matrix
Contents is actually a matrix.
@ symmetric
Matrix is symmetric.
@ general
No special properties.
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const ArrayView< const std::vector< double > > &velocity, const double factor=1.)
Definition advection.h:75
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition utilities.cc:189
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void call(const std::function< RT()> &function, internal::return_value< RT > &ret_val)
VectorType::value_type * begin(VectorType &V)
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)
void assemble(const MeshWorker::DoFInfoBox< dim, DOFINFO > &dinfo, A *assembler)
Definition loop.h:71
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation