deal.II version 9.7.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-87.h
Go to the documentation of this file.
1);
1005 *  
1006 * @endcode
1007 *
1008 * For the iterative solution procedure of the closest-point projection,
1009 * the maximum number of iterations and the tolerance for the maximum
1010 * absolute acceptable change in the distance in one iteration are set.
1011 *
1012 * @code
1013 *   pcout << " - determine closest point iteratively" << std::endl;
1014 *   constexpr int max_iter = 30;
1015 *   constexpr double tol_distance = 1e-6;
1016 *  
1017 * @endcode
1018 *
1019 * Now, we are ready to perform the algorithm by setting an initial guess
1020 * for the projection points simply corresponding to the collected support
1021 * points. We collect the global indices of the support points and the
1022 * total number of points that need to be processed and do not
1023 * fulfill the required tolerance. Those will be gradually reduced
1024 * upon the iterative process.
1025 *
1026 * @code
1027 *   std::vector<Point<dim>> closest_points = support_points; // initial guess
1028 *  
1029 *   std::vector<unsigned int> unmatched_points_idx(closest_points.size());
1030 *   std::iota(unmatched_points_idx.begin(), unmatched_points_idx.end(), 0);
1031 *  
1032 *   int n_unmatched_points =
1033 *   Utilities::MPI::sum(unmatched_points_idx.size(), MPI_COMM_WORLD);
1034 *  
1035 * @endcode
1036 *
1037 * Now, we create a Utilities::MPI::RemotePointEvaluation cache object and
1038 * start the loop for the fix-point iteration. We update the list of points
1039 * that still need to be processed and subsequently pass this information
1040 * to the Utilities::MPI::RemotePointEvaluation object. For the sake of
1041 * illustration, we export the coordinates of the points to be updated for
1042 * each iteration to a CSV file. Next, we can evaluate the signed distance
1043 * function and the gradient at those points to update the current solution
1044 * for the closest points. We perform the update only if the signed
1045 * distance of the closest point is not already within the tolerance
1046 * and register those points that still need to be processed.
1047 *
1048 * @code
1050 *  
1051 *   for (int it = 0; it < max_iter && n_unmatched_points > 0; ++it)
1052 *   {
1053 *   pcout << " - iteration " << it << ": " << n_unmatched_points;
1054 *  
1055 *   std::vector<Point<dim>> unmatched_points(unmatched_points_idx.size());
1056 *   for (unsigned int i = 0; i < unmatched_points_idx.size(); ++i)
1057 *   unmatched_points[i] = closest_points[unmatched_points_idx[i]];
1058 *  
1059 *   const auto all_unmatched_points =
1061 *   unmatched_points, MPI_COMM_WORLD, [](const auto &a, const auto &b) {
1062 *   auto result = a;
1063 *   result.insert(result.end(), b.begin(), b.end());
1064 *   return result;
1065 *   });
1066 *  
1067 *   if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
1068 *   {
1069 *   std::ofstream file("example_2_" + std::to_string(it) + ".csv");
1070 *   for (const auto &p : all_unmatched_points)
1071 *   file << p << std::endl;
1072 *   file.close();
1073 *   }
1074 *  
1075 *   rpe.reinit(unmatched_points, tria, mapping);
1076 *  
1077 *   AssertThrow(rpe.all_points_found(),
1078 *   ExcMessage("Processed point is outside domain."));
1079 *  
1080 *   const auto eval_values =
1081 *   VectorTools::point_values<1>(rpe, dof_handler, signed_distance);
1082 *  
1083 *   const auto eval_gradient =
1084 *   VectorTools::point_gradients<1>(rpe, dof_handler, signed_distance);
1085 *  
1086 *   std::vector<unsigned int> unmatched_points_idx_next;
1087 *  
1088 *   for (unsigned int i = 0; i < unmatched_points_idx.size(); ++i)
1089 *   if (std::abs(eval_values[i]) > tol_distance)
1090 *   {
1091 *   closest_points[unmatched_points_idx[i]] -=
1092 *   eval_values[i] * eval_gradient[i];
1093 *  
1094 *   unmatched_points_idx_next.emplace_back(unmatched_points_idx[i]);
1095 *   }
1096 *  
1097 *   unmatched_points_idx.swap(unmatched_points_idx_next);
1098 *  
1099 *   n_unmatched_points =
1100 *   Utilities::MPI::sum(unmatched_points_idx.size(), MPI_COMM_WORLD);
1101 *  
1102 *   pcout << " -> " << n_unmatched_points << std::endl;
1103 *   }
1104 *  
1105 * @endcode
1106 *
1107 * We print a warning message if we exceed the maximum number of allowed
1108 * iterations and if there are still projection points with a distance
1109 * value exceeding the tolerance.
1110 *
1111 * @code
1112 *   if (n_unmatched_points > 0)
1113 *   pcout << "WARNING: The tolerance of " << n_unmatched_points
1114 *   << " points is not yet attained." << std::endl;
1115 *  
1116 * @endcode
1117 *
1118 * As a result, we obtain a list of support points and corresponding
1119 * closest points at the zero-isosurface level set. This information
1120 * can be used to update the signed distance function, i.e., the
1121 * reinitialization the values of the level-set function to maintain
1122 * the signed distance property @cite henri2022geometrical.
1123 *
1124 * @code
1125 *   pcout << " - determine distance in narrow band" << std::endl;
1126 *   LinearAlgebra::distributed::Vector<double> solution_distance;
1127 *   solution_distance.reinit(solution);
1128 *  
1129 *   for (unsigned int i = 0; i < closest_points.size(); ++i)
1130 *   solution_distance[support_points_idx[i]] =
1131 *   support_points[i].distance(closest_points[i]);
1132 *  
1133 * @endcode
1134 *
1135 * In addition, we use the information of the closest point to
1136 * extrapolate values from the interface, i.e., the zero-level
1137 * set isosurface, to the support points within the narrow band.
1138 * This might be helpful to improve accuracy, e.g., for
1139 * diffuse interface fluxes where certain quantities are only
1140 * accurately determined at the interface (e.g. curvature
1141 * for surface tension @cite coquerelle2016fourth).
1142 *
1143 * @code
1144 *   pcout << " - perform extrapolation in narrow band" << std::endl;
1145 *   rpe.reinit(closest_points, tria, mapping);
1146 *   solution.update_ghost_values();
1147 *   const auto vals = VectorTools::point_values<1>(rpe, dof_handler, solution);
1148 *  
1149 *   LinearAlgebra::distributed::Vector<double> solution_extrapolated;
1150 *   solution_extrapolated.reinit(solution);
1151 *  
1152 *   for (unsigned int i = 0; i < closest_points.size(); ++i)
1153 *   solution_extrapolated[support_points_idx[i]] = vals[i];
1154 *  
1155 * @endcode
1156 *
1157 * Finally, we output the results to a VTU file.
1158 *
1159 * @code
1160 *   pcout << " - output results" << std::endl;
1161 *   DataOut<dim> data_out;
1162 *   data_out.add_data_vector(dof_handler, signed_distance, "signed_distance");
1163 *   data_out.add_data_vector(dof_handler, solution, "solution");
1164 *   data_out.add_data_vector(dof_handler,
1165 *   solution_distance,
1166 *   "solution_distance");
1167 *   data_out.add_data_vector(dof_handler,
1168 *   solution_extrapolated,
1169 *   "solution_extrapolated");
1170 *   data_out.build_patches(mapping);
1171 *   data_out.write_vtu_in_parallel("example_2.vtu", MPI_COMM_WORLD);
1172 *  
1173 *   pcout << std::endl;
1174 *   }
1175 *  
1176 * @endcode
1177 *
1178 *
1179 * <a name="step_87-Miniexample3Sharpinterfacemethodontheexampleofsurfacetensionforfronttracking"></a>
1180 * <h3>Mini example 3: Sharp interface method on the example of surface tension for front tracking</h3>
1181 *
1182
1183 *
1184 * The final mini example presents a basic implementation of
1185 * front tracking @cite peskin1977numerical, @cite unverdi1992front
1186 * of a surface mesh @f$\mathbb{T}_\Gamma@f$ immersed
1187 * in a Eulerian background fluid mesh @f$\mathbb{T}_\Omega@f$.
1188 *
1189
1190 *
1191 * We assume that the immersed surface is transported according to a
1192 * prescribed velocity field from the background mesh. Subsequently,
1193 * we perform a sharp computation of the surface-tension force:
1194 * @f[
1195 * (\boldsymbol v_i (\boldsymbol{x}), \boldsymbol F_S
1196 * (\boldsymbol{x}))_{\Omega}
1197 * =
1198 * \left( \boldsymbol{v}_i (\boldsymbol{x}), \sigma (\boldsymbol{x}) \kappa
1199 * (\boldsymbol{x}) \boldsymbol{n} (\boldsymbol{x}) \right)_\Gamma \approx
1200 * \sum_{q\in\mathbb{T}_\Gamma} \boldsymbol{v}_i^T (\boldsymbol{x}_q)
1201 * \sigma (\boldsymbol{x}_q) \kappa (\boldsymbol{x}_q) \boldsymbol{n}
1202 * (\boldsymbol{x}_q) |J(\boldsymbol{x}_q)| w(\boldsymbol{x}_q) \quad \forall
1203 * i\in\mathbb{T}_\Omega
1204 * .
1205 * @f]
1206 * We decompose this operation into two steps. In the first step, we evaluate
1207 * the force contributions @f$\sigma (\boldsymbol{x}_q) \kappa
1208 * (\boldsymbol{x}_q) \boldsymbol{n}
1209 * (\boldsymbol{x}_q)@f$ at the quadrature points defined on the immersed mesh
1210 * and multiply them with the mapped quadrature weight @f$|J(\boldsymbol{x}_q)|
1211 * w_q@f$:
1212 * @f[
1213 * \boldsymbol{F}_S (\boldsymbol{x}_q) \gets \sigma (\boldsymbol{x}_q) \kappa
1214 * (\boldsymbol{x}_q) \boldsymbol{n} (\boldsymbol{x}_q) |J(\boldsymbol{x}_q)|
1215 * w_q \quad \forall q\in\mathbb{T}_\Gamma.
1216 * @f]
1217 * In the second step, we compute the discretized weak form by multiplying
1218 * with test functions on the background mesh:
1219 * @f[
1220 * (\boldsymbol v_i (\boldsymbol{x}), \boldsymbol F_S
1221 * (\boldsymbol{x}))_{\Omega} \gets \sum_q \boldsymbol{v}_i^T
1222 * (\boldsymbol{x}_q) \boldsymbol{F}_S
1223 * (\boldsymbol{x}_q)
1224 * \quad \forall i\in\mathbb{T}_\Omega
1225 * .
1226 * @f]
1227 * Obviously, we need to communicate between the two steps. The second step
1228 * can be handled completely by Utilities::MPI::RemotePointEvaluation, which
1229 * provides the function
1230 * Utilities::MPI::RemotePointEvaluation::process_and_evaluate() for this
1231 * purpose.
1232 *
1233
1234 *
1235 * We start with setting the parameters consisting of the polynomial degree of
1236 * the shape functions, the dimension of the background mesh, the time-step
1237 * size to be considered for transporting the surface mesh and the number of
1238 * time steps.
1239 *
1240
1241 *
1242 *
1243 * @code
1244 *   void example_3()
1245 *   {
1246 *   constexpr unsigned int degree = 3;
1247 *   constexpr unsigned int dim = 2;
1248 *   const double dt = 0.01;
1249 *   const unsigned int n_time_steps = 200;
1250 *  
1251 * @endcode
1252 *
1253 * This program is intended to be executed in 2D or 3D.
1254 *
1255 * @code
1256 *   static_assert(dim == 2 || dim == 3, "Only implemented for 2D or 3D.");
1257 *  
1258 *   ConditionalOStream pcout(std::cout,
1259 *   Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) ==
1260 *   0);
1261 *  
1262 *   pcout << "Running: example 3" << std::endl;
1263 *  
1264 * @endcode
1265 *
1266 * Next, we create the standard objects necessary for the finite element
1267 * representation of the background mesh
1268 *
1269 * @code
1270 *   pcout << " - creating background mesh" << std::endl;
1271 *   DistributedTriangulation<dim> tria_background(MPI_COMM_WORLD);
1272 *   GridGenerator::hyper_cube(tria_background);
1273 *   tria_background.refine_global(5);
1274 *  
1275 *   const MappingQ1<dim> mapping_background;
1276 *   const FESystem<dim> fe_background(FE_Q<dim>(degree), dim);
1277 *   DoFHandler<dim> dof_handler_background(tria_background);
1278 *   dof_handler_background.distribute_dofs(fe_background);
1279 *  
1280 * @endcode
1281 *
1282 * and, similarly, for the immersed surface mesh.
1283 * We use a sphere with radius @f$r=0.75@f$ which is
1284 * placed in the center of the top half of the cubic background domain.
1285 *
1286 * @code
1287 *   pcout << " - creating immersed mesh" << std::endl;
1288 *   const Point<dim> center((dim == 2) ? Point<dim>(0.5, 0.75) :
1289 *   Point<dim>(0.5, 0.75, 0.5));
1290 *   const double radius = 0.15;
1291 *  
1292 *   DistributedTriangulation<dim - 1, dim> tria_immersed(MPI_COMM_WORLD);
1293 *   GridGenerator::hyper_sphere(tria_immersed, center, radius);
1294 *   tria_immersed.refine_global(4);
1295 *  
1296 * @endcode
1297 *
1298 * Two different mappings are considered for the immersed
1299 * surface mesh: one valid for the initial configuration and one
1300 * that is updated in every time step according to the nodal
1301 * displacements. Two types of finite elements are used to
1302 * represent scalar and vector-valued DoF values.
1303 *
1304 * @code
1305 *   const MappingQ<dim - 1, dim> mapping_immersed_base(3);
1306 *   MappingQCache<dim - 1, dim> mapping_immersed(3);
1307 *   mapping_immersed.initialize(mapping_immersed_base, tria_immersed);
1308 *   const QGauss<dim - 1> quadrature_immersed(degree + 1);
1309 *  
1310 *   const FE_Q<dim - 1, dim> fe_scalar_immersed(degree);
1311 *   const FESystem<dim - 1, dim> fe_immersed(fe_scalar_immersed, dim);
1312 *   DoFHandler<dim - 1, dim> dof_handler_immersed(tria_immersed);
1313 *   dof_handler_immersed.distribute_dofs(fe_immersed);
1314 *  
1315 * @endcode
1316 *
1317 * We renumber the DoFs related to the vector-valued problem to
1318 * simplify access to the individual components.
1319 *
1320 * @code
1321 *   DoFRenumbering::support_point_wise(dof_handler_immersed);
1322 *  
1323 * @endcode
1324 *
1325 * We fill a DoF vector on the background mesh with an analytical
1326 * velocity field considering the Rayleigh-Kothe vortex flow and
1327 * initialize a DoF vector for the weak form of the surface-tension force.
1328 *
1329 * @code
1331 *   velocity.reinit(dof_handler_background.locally_owned_dofs(),
1333 *   dof_handler_background),
1334 *   MPI_COMM_WORLD);
1336 *  
1338 *   dof_handler_background.locally_owned_dofs(),
1339 *   DoFTools::extract_locally_active_dofs(dof_handler_background),
1340 *   MPI_COMM_WORLD);
1341 *  
1342 * @endcode
1343 *
1344 * Next, we collect the real positions @f$\boldsymbol{x}_q@f$ of the quadrature
1345 * points of the surface mesh in a vector.
1346 *
1347 * @code
1348 *   LinearAlgebra::distributed::Vector<double> immersed_support_points;
1349 *   collect_support_points(mapping_immersed,
1350 *   dof_handler_immersed,
1351 *   immersed_support_points);
1352 *  
1353 * @endcode
1354 *
1355 * We initialize a Utilities::MPI::RemotePointEvaluation object and start
1356 * the time loop. For any other step than the initial one, we first move the
1357 * support points of the surface mesh according to the fluid velocity of the
1358 * background mesh. Thereto, we first update the time of the velocity
1359 * function. Then, we update the internal data structures of the
1360 * Utilities::MPI::RemotePointEvaluation object with the collected support
1361 * points of the immersed mesh. We throw an exception if one of the points
1362 * cannot be found within the domain of the background mesh. Next, we
1363 * evaluate the velocity at the surface-mesh support points and compute the
1364 * resulting update of the coordinates. Finally, we update the mapping of
1365 * the immersed surface mesh to the current position.
1366 *
1367 * @code
1369 *   double time = 0.0;
1370 *   for (unsigned int it = 0; it <= n_time_steps; ++it, time += dt)
1371 *   {
1372 *   pcout << "time: " << time << std::endl;
1373 *   if (it > 0)
1374 *   {
1375 *   pcout << " - move support points (immersed mesh)" << std::endl;
1376 *   vortex.set_time(time);
1377 *   VectorTools::interpolate(mapping_background,
1378 *   dof_handler_background,
1379 *   vortex,
1380 *   velocity);
1381 *   rpe.reinit(convert<dim>(immersed_support_points),
1382 *   tria_background,
1383 *   mapping_background);
1384 *  
1385 *   AssertThrow(rpe.all_points_found(),
1386 *   ExcMessage(
1387 *   "Immersed domain leaves background domain!"));
1388 *  
1389 *   velocity.update_ghost_values();
1390 *   const auto immersed_velocity =
1392 *   dof_handler_background,
1393 *   velocity);
1394 *   velocity.zero_out_ghost_values();
1395 *  
1396 *   for (unsigned int i = 0, c = 0;
1397 *   i < immersed_support_points.locally_owned_size() / dim;
1398 *   ++i)
1399 *   for (unsigned int d = 0; d < dim; ++d, ++c)
1400 *   immersed_support_points.local_element(c) +=
1401 *   immersed_velocity[i][d] * dt;
1402 *  
1403 *   mapping_immersed.initialize(mapping_immersed_base,
1404 *   dof_handler_immersed,
1405 *   immersed_support_points,
1406 *   false);
1407 *   }
1408 *  
1409 * @endcode
1410 *
1411 * Next, we loop over all locally owned cells of the immersed mesh and
1412 * its quadrature points to compute the value for the local surface
1413 * tension force contribution @f$\boldsymbol{F}_S(\boldsymbol{x}_q)@f$. We
1414 * store the real coordinates of the quadrature points and the
1415 * corresponding force contributions in two individual vectors. For
1416 * computation of the latter, the normal vector
1417 * @f$\boldsymbol{n}(\boldsymbol{x}_q)@f$ can be directly extracted from the
1418 * surface mesh via FEValues and, for the curvature, we use the
1419 * following approximation:
1420 * @f[
1421 * \kappa(\boldsymbol{x}_q)
1422 * =
1423 * \nabla \cdot \boldsymbol{n}(\boldsymbol{x}_q)
1424 * =
1425 * \text{tr}\left({\nabla \boldsymbol{n}(\boldsymbol{x}_q)}\right)
1426 * \approx
1427 * \text{tr}\left({\nabla \sum_i \boldsymbol{N}_i (\boldsymbol{x}_q)
1428 * \boldsymbol n_i}\right)
1429 * =
1430 * \sum_i\text{tr}\left({\nabla \boldsymbol{N}_i (\boldsymbol{x}_q)
1431 * \boldsymbol n_i}\right)
1432 * \;\text{with}\; i\in[0,n_{\text{dofs\_per\_cell}}),
1433 * @f]
1434 * which we can apply since the immersed mesh is consistently
1435 * orientated. The surface tension coefficient is set to 1 for the
1436 * sake of demonstration.
1437 *
1438 * @code
1439 *   pcout << " - compute to be tested values (immersed mesh)" << std::endl;
1440 *   using value_type = Tensor<1, dim, double>;
1441 *  
1442 *   std::vector<Point<dim>> integration_points;
1443 *   std::vector<value_type> integration_values;
1444 *  
1445 *   FEValues<dim - 1, dim> fe_values(mapping_immersed,
1446 *   fe_immersed,
1447 *   quadrature_immersed,
1451 *  
1452 *   FEValues<dim - 1, dim> fe_values_co(
1453 *   mapping_immersed,
1454 *   fe_scalar_immersed,
1455 *   fe_scalar_immersed.get_unit_support_points(),
1457 *  
1458 *   std::vector<unsigned int> component_to_system_index(
1459 *   fe_immersed.n_dofs_per_cell());
1460 *  
1461 *   for (unsigned int i = 0, c = 0;
1462 *   i < fe_scalar_immersed.n_dofs_per_cell();
1463 *   ++i)
1464 *   for (unsigned int d = 0; d < dim; ++d, ++c)
1465 *   component_to_system_index[c] =
1466 *   fe_immersed.component_to_system_index(d, i);
1467 *  
1468 *   for (const auto &cell : tria_immersed.active_cell_iterators() |
1470 *   {
1471 *   fe_values.reinit(cell);
1472 *   fe_values_co.reinit(cell);
1473 *  
1474 *   for (const auto &q : fe_values.quadrature_point_indices())
1475 *   {
1476 *   integration_points.emplace_back(fe_values.quadrature_point(q));
1477 *  
1478 *   const auto sigma = 1.0; // surface tension coefficient
1479 *  
1480 *   const auto normal = fe_values.normal_vector(q);
1481 *   double curvature = 0;
1482 *   for (unsigned int i = 0, c = 0;
1483 *   i < fe_scalar_immersed.n_dofs_per_cell();
1484 *   ++i)
1485 *   for (unsigned int d = 0; d < dim; ++d, ++c)
1486 *   curvature += fe_values.shape_grad_component(
1487 *   component_to_system_index[c], q, d)[d] *
1488 *   fe_values_co.normal_vector(i)[d];
1489 *  
1490 *   const auto FxJxW =
1491 *   sigma * curvature * normal * fe_values.JxW(q);
1492 *  
1493 *   integration_values.emplace_back(FxJxW);
1494 *   }
1495 *   }
1496 *  
1497 * @endcode
1498 *
1499 * Before we evaluate the weak form of the surface-tension force, the
1500 * communication pattern of Utilities::MPI::RemotePointEvaluation is
1501 * set up from the quadrature points of the immersed mesh, determining
1502 * the surrounding cells on the background mesh.
1503 *
1504 * @code
1505 *   pcout << " - test values (background mesh)" << std::endl;
1506 *  
1507 *   rpe.reinit(integration_points, tria_background, mapping_background);
1508 *  
1509 * @endcode
1510 *
1511 * In preparation for utilizing
1513 * performs the
1514 * multiplication with the test function, we set up a callback function
1515 * that contains the operation on the intersected cells of the
1516 * background mesh. Within this function, we initialize a
1517 * FEPointEvaluation object that allows us to integrate values at
1518 * arbitrary points within a cell. We loop over the cells that surround
1519 * quadrature points of the immersed mesh -- provided by the callback
1520 * function. From the provided CellData object, we retrieve the unit
1521 * points, i.e., the quadrature points of the immersed mesh that lie
1522 * within the current cell and a pointer to the stored values on the
1523 * current cell (local surface-tension force) for convenience. We
1524 * reinitialize the data structure of FEPointEvaluation on every cell
1525 * according to the unit points. Next, we loop over the quadrature
1526 * points attributed to the cell and submit the local surface-tension
1527 * force to the FEPointEvaluation object. Via
1528 * FEPointEvaluation::test_and_sum(), the submitted values are
1529 * multiplied by the values of the test function and a summation over
1530 * all given points is performed. Subsequently, the contributions are
1531 * assembled into the global vector containing the weak form of the
1532 * surface-tension force.
1533 *
1534 * @code
1535 *   const auto integration_function = [&](const auto &values,
1536 *   const auto &cell_data) {
1537 *   FEPointEvaluation<dim, dim> phi_force(mapping_background,
1538 *   fe_background,
1539 *   update_values);
1540 *  
1541 *   std::vector<double> local_values;
1542 *   std::vector<types::global_dof_index> local_dof_indices;
1543 *  
1544 *   for (const auto cell : cell_data.cell_indices())
1545 *   {
1546 *   const auto cell_dofs =
1547 *   cell_data.get_active_cell_iterator(cell)
1548 *   ->as_dof_handler_iterator(dof_handler_background);
1549 *  
1550 *   const auto unit_points = cell_data.get_unit_points(cell);
1551 *   const auto FxJxW = cell_data.get_data_view(cell, values);
1552 *  
1553 *   phi_force.reinit(cell_dofs, unit_points);
1554 *  
1555 *   for (const auto q : phi_force.quadrature_point_indices())
1556 *   phi_force.submit_value(FxJxW[q], q);
1557 *  
1558 *   local_values.resize(cell_dofs->get_fe().n_dofs_per_cell());
1559 *   phi_force.test_and_sum(local_values, EvaluationFlags::values);
1560 *  
1561 *   local_dof_indices.resize(cell_dofs->get_fe().n_dofs_per_cell());
1562 *   cell_dofs->get_dof_indices(local_dof_indices);
1564 *   local_values, local_dof_indices, force_vector);
1565 *   }
1566 *   };
1567 *  
1568 * @endcode
1569 *
1570 * The callback function is passed together with the vector holding the
1571 * surface-tension force contribution at each quadrature point of the
1572 * immersed mesh to
1574 * missing step is to compress the distributed force vector.
1575 *
1576 * @code
1577 *   rpe.process_and_evaluate<value_type>(integration_values,
1578 *   integration_function);
1579 *   force_vector.compress(VectorOperation::add);
1580 *  
1581 * @endcode
1582 *
1583 * After every tenth step or at the beginning/end of the time loop, we
1584 * output the force vector and the velocity of the background mesh to
1585 * a VTU file. In addition, we also export the geometry of the
1586 * (deformed) immersed surface mesh to a separate VTU file.
1587 *
1588 * @code
1589 *   if (it % 10 == 0 || it == n_time_steps)
1590 *   {
1591 *   std::vector<
1593 *   vector_component_interpretation(
1595 *   pcout << " - write data (background mesh)" << std::endl;
1596 *   DataOut<dim> data_out_background;
1597 *   DataOutBase::VtkFlags flags_background;
1598 *   flags_background.write_higher_order_cells = true;
1599 *   data_out_background.set_flags(flags_background);
1600 *   data_out_background.add_data_vector(
1601 *   dof_handler_background,
1602 *   force_vector,
1603 *   std::vector<std::string>(dim, "force"),
1604 *   vector_component_interpretation);
1605 *   data_out_background.add_data_vector(
1606 *   dof_handler_background,
1607 *   velocity,
1608 *   std::vector<std::string>(dim, "velocity"),
1609 *   vector_component_interpretation);
1610 *   data_out_background.build_patches(mapping_background, 3);
1611 *   data_out_background.write_vtu_in_parallel("example_3_background_" +
1612 *   std::to_string(it) +
1613 *   ".vtu",
1614 *   MPI_COMM_WORLD);
1615 *  
1616 *   pcout << " - write mesh (immersed mesh)" << std::endl;
1617 *   DataOut<dim - 1, dim> data_out_immersed;
1618 *   data_out_immersed.attach_triangulation(tria_immersed);
1619 *   data_out_immersed.build_patches(mapping_immersed, 3);
1620 *   data_out_immersed.write_vtu_in_parallel("example_3_immersed_" +
1621 *   std::to_string(it) +
1622 *   ".vtu",
1623 *   MPI_COMM_WORLD);
1624 *   }
1625 *   pcout << std::endl;
1626 *   }
1627 *   }
1628 *  
1629 * @endcode
1630 *
1631 *
1632 * <a name="step_87-Utilityfunctionsdefinition"></a>
1633 * <h3>Utility functions (definition)</h3>
1634 *
1635 * @code
1636 *   template <int spacedim>
1637 *   std::vector<Point<spacedim>>
1638 *   create_points_along_line(const Point<spacedim> &p0,
1639 *   const Point<spacedim> &p1,
1640 *   const unsigned int n_subdivisions)
1641 *   {
1642 *   Assert(n_subdivisions >= 1, ExcInternalError());
1643 *  
1644 *   std::vector<Point<spacedim>> points;
1645 *   points.reserve(n_subdivisions + 1);
1646 *  
1647 *   points.emplace_back(p0);
1648 *   for (unsigned int i = 1; i < n_subdivisions; ++i)
1649 *   points.emplace_back(p0 + (p1 - p0) * static_cast<double>(i) /
1650 *   static_cast<double>(n_subdivisions));
1651 *   points.emplace_back(p1);
1652 *  
1653 *   return points;
1654 *   }
1655 *  
1656 *   template <int spacedim, typename T0, typename T1>
1657 *   void print_along_line(const std::string &file_name,
1658 *   const std::vector<Point<spacedim>> &points,
1659 *   const std::vector<T0> &values_0,
1660 *   const std::vector<T1> &values_1)
1661 *   {
1662 *   AssertThrow(points.size() == values_0.size() &&
1663 *   (values_1.size() == points.size() || values_1.empty()),
1664 *   ExcMessage("The provided vectors must have the same length."));
1665 *  
1666 *   std::ofstream file(file_name);
1667 *  
1668 *   for (unsigned int i = 0; i < points.size(); ++i)
1669 *   {
1670 *   file << std::fixed << std::right << std::setw(5) << std::setprecision(3)
1671 *   << points[0].distance(points[i]);
1672 *  
1673 *   for (unsigned int d = 0; d < spacedim; ++d)
1674 *   file << std::fixed << std::right << std::setw(10)
1675 *   << std::setprecision(3) << points[i][d];
1676 *  
1677 *   file << std::fixed << std::right << std::setw(10)
1678 *   << std::setprecision(3) << values_0[i];
1679 *  
1680 *   if (!values_1.empty())
1681 *   file << std::fixed << std::right << std::setw(10)
1682 *   << std::setprecision(3) << values_1[i];
1683 *   file << std::endl;
1684 *   }
1685 *   }
1686 *  
1687 *   template <int dim, int spacedim>
1688 *   void collect_support_points(
1689 *   const Mapping<dim, spacedim> &mapping,
1690 *   const DoFHandler<dim, spacedim> &dof_handler,
1691 *   LinearAlgebra::distributed::Vector<double> &support_points)
1692 *   {
1693 *   support_points.reinit(dof_handler.locally_owned_dofs(),
1694 *   DoFTools::extract_locally_active_dofs(dof_handler),
1695 *   dof_handler.get_mpi_communicator());
1696 *  
1697 *   const auto &fe = dof_handler.get_fe();
1698 *  
1699 *   FEValues<dim, spacedim> fe_values(
1700 *   mapping,
1701 *   fe,
1702 *   fe.base_element(0).get_unit_support_points(),
1704 *  
1705 *   std::vector<types::global_dof_index> local_dof_indices(
1706 *   fe.n_dofs_per_cell());
1707 *  
1708 *   std::vector<unsigned int> component_to_system_index(
1709 *   fe_values.n_quadrature_points * spacedim);
1710 *  
1711 *   for (unsigned int q = 0, c = 0; q < fe_values.n_quadrature_points; ++q)
1712 *   for (unsigned int d = 0; d < spacedim; ++d, ++c)
1713 *   component_to_system_index[c] = fe.component_to_system_index(d, q);
1714 *  
1715 *   for (const auto &cell : dof_handler.active_cell_iterators() |
1717 *   {
1718 *   fe_values.reinit(cell);
1719 *   cell->get_dof_indices(local_dof_indices);
1720 *  
1721 *   for (unsigned int q = 0, c = 0; q < fe_values.n_quadrature_points; ++q)
1722 *   for (unsigned int d = 0; d < spacedim; ++d, ++c)
1723 *   support_points[local_dof_indices[component_to_system_index[c]]] =
1724 *   fe_values.quadrature_point(q)[d];
1725 *   }
1726 *   }
1727 *  
1728 *   template <int dim, int spacedim, typename T>
1729 *   std::tuple<std::vector<Point<spacedim>>, std::vector<types::global_dof_index>>
1730 *   collect_support_points_with_narrow_band(
1731 *   const Mapping<dim, spacedim> &mapping,
1732 *   const DoFHandler<dim, spacedim> &dof_handler_signed_distance,
1733 *   const LinearAlgebra::distributed::Vector<T> &signed_distance,
1734 *   const DoFHandler<dim, spacedim> &dof_handler_support_points,
1735 *   const double narrow_band_threshold)
1736 *   {
1737 *   AssertThrow(narrow_band_threshold >= 0,
1738 *   ExcMessage("The narrow band threshold"
1739 *   " must be larger than or equal to 0."));
1740 *   const auto &tria = dof_handler_signed_distance.get_triangulation();
1741 *   const Quadrature<dim> quad(dof_handler_support_points.get_fe()
1742 *   .base_element(0)
1743 *   .get_unit_support_points());
1744 *  
1745 *   FEValues<dim> distance_values(mapping,
1746 *   dof_handler_signed_distance.get_fe(),
1747 *   quad,
1748 *   update_values);
1749 *  
1750 *   FEValues<dim> req_values(mapping,
1751 *   dof_handler_support_points.get_fe(),
1752 *   quad,
1754 *  
1755 *   std::vector<T> temp_distance(quad.size());
1756 *   std::vector<types::global_dof_index> local_dof_indices(
1757 *   dof_handler_support_points.get_fe().n_dofs_per_cell());
1758 *  
1759 *   std::vector<Point<dim>> support_points;
1760 *   std::vector<types::global_dof_index> support_points_idx;
1761 *  
1762 *   const bool has_ghost_elements = signed_distance.has_ghost_elements();
1763 *  
1764 *   const auto &locally_owned_dofs_req =
1765 *   dof_handler_support_points.locally_owned_dofs();
1766 *   std::vector<bool> flags(locally_owned_dofs_req.n_elements(), false);
1767 *  
1768 *   if (has_ghost_elements == false)
1769 *   signed_distance.update_ghost_values();
1770 *  
1771 *   for (const auto &cell :
1772 *   tria.active_cell_iterators() | IteratorFilters::LocallyOwnedCell())
1773 *   {
1774 *   const auto cell_distance =
1775 *   cell->as_dof_handler_iterator(dof_handler_signed_distance);
1776 *   distance_values.reinit(cell_distance);
1777 *   distance_values.get_function_values(signed_distance, temp_distance);
1778 *  
1779 *   const auto cell_req =
1780 *   cell->as_dof_handler_iterator(dof_handler_support_points);
1781 *   req_values.reinit(cell_req);
1782 *   cell_req->get_dof_indices(local_dof_indices);
1783 *  
1784 *   for (const auto q : req_values.quadrature_point_indices())
1785 *   if (std::abs(temp_distance[q]) < narrow_band_threshold)
1786 *   {
1787 *   const auto idx = local_dof_indices[q];
1788 *  
1789 *   if (locally_owned_dofs_req.is_element(idx) == false ||
1790 *   flags[locally_owned_dofs_req.index_within_set(idx)])
1791 *   continue;
1792 *  
1793 *   flags[locally_owned_dofs_req.index_within_set(idx)] = true;
1794 *  
1795 *   support_points_idx.emplace_back(idx);
1796 *   support_points.emplace_back(req_values.quadrature_point(q));
1797 *   }
1798 *   }
1799 *  
1800 *   if (has_ghost_elements == false)
1801 *   signed_distance.zero_out_ghost_values();
1802 *  
1803 *   return {support_points, support_points_idx};
1804 *   }
1805 *  
1806 *   template <int spacedim>
1807 *   std::vector<Point<spacedim>> convert(
1808 *   const LinearAlgebra::distributed::Vector<double> &support_points_unrolled)
1809 *   {
1810 *   const unsigned int n_points =
1811 *   support_points_unrolled.locally_owned_size() / spacedim;
1812 *  
1813 *   std::vector<Point<spacedim>> points(n_points);
1814 *  
1815 *   for (unsigned int i = 0, c = 0; i < n_points; ++i)
1816 *   for (unsigned int d = 0; d < spacedim; ++d, ++c)
1817 *   points[i][d] = support_points_unrolled.local_element(c);
1818 *  
1819 *   return points;
1820 *   }
1821 *  
1822 *   } // namespace Step87
1823 *  
1824 * @endcode
1825 *
1826 *
1827 * <a name="step_87-Driver"></a>
1828 * <h3>Driver</h3>
1829 *
1830
1831 *
1832 * Finally, the driver of the program executes the four mini examples.
1833 *
1834 * @code
1835 *   int main(int argc, char **argv)
1836 *   {
1837 *   using namespace dealii;
1838 *   Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
1839 *   std::cout.precision(5);
1840 *  
1841 *   if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
1842 *   Step87::example_0(); // only run on root process
1843 *  
1844 *   Step87::example_1();
1845 *   Step87::example_2();
1846 *   Step87::example_3();
1847 *   }
1848 * @endcode
1849<a name="step_87-Results"></a><h1>Results</h1>
1850
1851
1852<a name="step_87-Miniexample0"></a><h3>Mini example 0</h3>
1853
1854
1855We present a part of the terminal output. It shows, for each point, the
1856determined cell and reference position. Also, one can see that
1857the values evaluated with FEValues, FEPointEvaluation, and
1858VectorTools::point_values() are identical, as expected.
1859
1860@verbatim
1861Running: example 0
1862 - Found point with real coordinates: 0 0.5
1863 - in cell with vertices: (0 0.4) (0.2 0.4) (0 0.6) (0.2 0.6)
1864 - with coordinates on the unit cell: (0 0.5)
1865 - Values at point:
1866 - 0.25002 (w. FEValues)
1867 - 0.25002 (w. FEPointEvaluation)
1868 - 0.25002 (w. VectorTools::point_value())
1869
1870 - Found point with real coordinates: 0.05 0.5
1871 - in cell with vertices: (0 0.4) (0.2 0.4) (0 0.6) (0.2 0.6)
1872 - with coordinates on the unit cell: (0.25 0.5)
1873 - Values at point:
1874 - 0.20003 (w. FEValues)
1875 - 0.20003 (w. FEPointEvaluation)
1876 - 0.20003 (w. VectorTools::point_value())
1877
1878...
1879
1880 - Found point with real coordinates: 1 0.5
1881 - in cell with vertices: (0.8 0.4) (1 0.4) (0.8 0.6) (1 0.6)
1882 - with coordinates on the unit cell: (1 0.5)
1883 - Values at point:
1884 - 0.25002 (w. FEValues)
1885 - 0.25002 (w. FEPointEvaluation)
1886 - 0.25002 (w. VectorTools::point_value())
1887
1888 - writing csv file
1889@endverbatim
1890
1891The CSV output is as follows and contains, in the
1892first column, the distances with respect to the first point,
1893the second and the third column represent the coordinates
1894of the points and the fourth column the evaluated solution
1895values at those points.
1896
1897@verbatim
18980.000 0.000 0.500 0.250
18990.050 0.050 0.500 0.200
19000.100 0.100 0.500 0.150
19010.150 0.150 0.500 0.100
19020.200 0.200 0.500 0.050
19030.250 0.250 0.500 0.000
19040.300 0.300 0.500 -0.050
19050.350 0.350 0.500 -0.099
19060.400 0.400 0.500 -0.148
19070.450 0.450 0.500 -0.195
19080.500 0.500 0.500 -0.211
19090.550 0.550 0.500 -0.195
19100.600 0.600 0.500 -0.148
19110.650 0.650 0.500 -0.099
19120.700 0.700 0.500 -0.050
19130.750 0.750 0.500 0.000
19140.800 0.800 0.500 0.050
19150.850 0.850 0.500 0.100
19160.900 0.900 0.500 0.150
19170.950 0.950 0.500 0.200
19181.000 1.000 0.500 0.250
1919@endverbatim
1920
1921<a name="step_87-Miniexample1"></a><h3>Mini example 1</h3>
1922
1923
1924We show the terminal output.
1925
1926@verbatim
1927Running: example 1
1928 - writing csv file
1929@endverbatim
1930
1931The CSV output is as follows and identical to the results
1932of the serial case presented in mini example 0.
1933The fifth column shows the
1934user quantity evaluated additionally in this mini example.
1935
1936@verbatim
19370.000 0.000 0.500 0.250 0.000
19380.050 0.050 0.500 0.200 0.050
19390.100 0.100 0.500 0.150 0.100
19400.150 0.150 0.500 0.100 0.150
19410.200 0.200 0.500 0.050 0.200
19420.250 0.250 0.500 0.000 0.250
19430.300 0.300 0.500 -0.050 0.300
19440.350 0.350 0.500 -0.100 0.350
19450.400 0.400 0.500 -0.149 0.400
19460.450 0.450 0.500 -0.200 0.450
19470.500 0.500 0.500 -0.222 0.500
19480.550 0.550 0.500 -0.200 0.550
19490.600 0.600 0.500 -0.149 0.600
19500.650 0.650 0.500 -0.100 0.650
19510.700 0.700 0.500 -0.050 0.700
19520.750 0.750 0.500 0.000 0.750
19530.800 0.800 0.500 0.050 0.800
19540.850 0.850 0.500 0.100 0.850
19550.900 0.900 0.500 0.150 0.900
19560.950 0.950 0.500 0.200 0.950
19571.000 1.000 0.500 0.250 1.000
1958@endverbatim
1959
1960
1961<a name="step_87-Miniexample2"></a><h3>Mini example 2</h3>
1962
1963
1964We show the terminal output.
1965@verbatim
1966Running: example 2
1967 - create system
1968 - determine narrow band
1969 - determine closest point iteratively
1970 - iteration 0: 7076 -> 7076
1971 - iteration 1: 7076 -> 104
1972 - iteration 2: 104 -> 0
1973 - determine distance in narrow band
1974 - perform extrapolation in narrow band
1975 - output results
1976@endverbatim
1977
1978The following three plots, representing the performed iterations of the
1979closest-point projection, show the current position of the closest
1980points exceeding the required tolerance of the discrete interface
1981of the circle and still need to
1982be corrected.
1983It can be seen that the numbers of points to be processed decrease
1984from iteration to iteration.
1985<table align="center" class="doxtable">
1986 <tr>
1987 <td>
1988 @image html https://www.dealii.org/images/steps/developer/step_87_ex_2_p_0.png
1989 </td>
1990 <td>
1991 @image html https://www.dealii.org/images/steps/developer/step_87_ex_2_p_1.png
1992 </td>
1993 <td>
1994 @image html https://www.dealii.org/images/steps/developer/step_87_ex_2_p_2.png
1995 </td>
1996 </tr>
1997</table>
1998
1999The output visualized in Paraview looks like the following: On the
2000left, the original distance function is shown as the light gray surface.
2001In addition, the contour values refer to the distance values determined
2002from calculation of the distance to the closest points at the interface
2003in the narrow band. It can be seen that the two functions coincide.
2004Similarly, on the right, the original solution and the extrapolated
2005solution from the interface is shown.
2006
2007<table align="center" class="doxtable">
2008 <tr>
2009 <td>
2010 @image html https://www.dealii.org/images/steps/developer/step_87_ex_2_res_0.png
2011 </td>
2012 <td>
2013 @image html https://www.dealii.org/images/steps/developer/step_87_ex_2_res_1.png
2014 </td>
2015 </tr>
2016</table>
2017
2018<a name="step_87-Miniexample3"></a><h3>Mini example 3</h3>
2019
2020
2021We show a shortened version of the terminal output.
2022
2023@verbatim
2024Running: example 3
2025 - creating background mesh
2026 - creating immersed mesh
2027time: 0
2028 - compute to be tested values (immersed mesh)
2029 - test values (background mesh)
2030 - write data (background mesh)
2031 - write mesh (immersed mesh)
2032
2033time: 0.01
2034 - move support points (immersed mesh)
2035 - compute to be tested values (immersed mesh)
2036 - test values (background mesh)
2037
2038time: 0.02
2039 - move support points (immersed mesh)
2040 - compute to be tested values (immersed mesh)
2041 - test values (background mesh)
2042
2043...
2044
2045time: 2
2046 - move support points (immersed mesh)
2047 - compute to be tested values (immersed mesh)
2048 - test values (background mesh)
2049 - write data (background mesh)
2050 - write mesh (immersed mesh)
2051@endverbatim
2052
2053The output visualized in Paraview looks like the following: The deformation of
2054the immersed mesh by the reversible vortex flow can be seen. Due to
2055discretization errors, the shape is not exactly circular at the end, illustrated
2056in the right figure. The sharp nature of the surface-tension force vector, shown
2057as vector plots, can be seen by its restriction to cells that are intersected by
2058the immersed mesh.
2059
2060<table align="center" class="doxtable">
2061 <tr>
2062 <td>
2063 @image html https://www.dealii.org/images/steps/developer/step_87_ex_3_force.0000.png
2064 </td>
2065 <td>
2066 @image html https://www.dealii.org/images/steps/developer/step_87_ex_3_force.0010.png
2067 </td>
2068 <td>
2069 @image html https://www.dealii.org/images/steps/developer/step_87_ex_3_force.0020.png
2070 </td>
2071 </tr>
2072</table>
2073
2074<a name="step_87-Possibilitiesforextension"></a><h3>Possibilities for extension</h3>
2075
2076
2077This program highlights some of the main capabilities
2078of the distributed evaluation routines in deal.II. However, there are many
2079related topics worth mentioning:
2080- Performing a distributed search is an expensive step. That is why we suggest
2081to provide hints to Utilities::MPI::RemotePointEvaluation and to reuse
2082Utilities::MPI::RemotePointEvaluation
2083instances in the case that the communication pattern has not changed.
2084Furthermore, there are instances where no search is needed and the points are
2085already sorted into the right cells. This is the case if the points are
2086generated on the cell level (see @ref step_85 "step-85"; CutFEM) or the points are
2087automatically sorted into the correct (neighboring) cell (see @ref step_68 "step-68"; PIC with
2088Particles::ParticleHandler). Having said that, the
2089Particles::ParticleHandler::insert_global_particles() function uses
2090the described infrastructure to perform the initial sorting of particles into
2091cells.
2092- We concentrated on parallelization aspects in this tutorial. However, we would
2093like to point out the need for fast evaluation on cell level.
2094The task for this in deal.II is FEPointEvaluation. It exploits the structure of
2095@f[
2096\hat{u}(\hat{\boldsymbol{x}}) = \sum_i \hat{N}_i(\hat{\boldsymbol{x}}) \hat{u}_i
2097@f]
2098to derive fast implementations, e.g., for tensor-product elements
2099@f[
2100\hat{u}(\hat{x}_0, \hat{x}_1, \hat{x}_2) =
2101\sum_k \hat{N}^{\text{1D}}_k(\hat{x}_2)
2102\sum_j \hat{N}^{\text{1D}}_j(\hat{x}_1)
2103\sum_i \hat{N}^{\text{1D}}_i(\hat{x}_0)
2104\hat{u}_{ijk}.
2105@f]
2106Since only 1D shape functions are queried and are re-used in re-occurring terms,
2107this formulation is faster than without exploitation of the structure.
2108- Utilities::MPI::RemotePointEvaluation is used in multiple places in deal.II.
2109The class DataOutResample allows to output results on a different mesh than
2110the computational mesh. This is useful if one wants to output the results
2111on a coarser mesh or one does not want to output 3D results but instead 2D
2112slices. In addition, MGTwoLevelTransferNonNested allows to prolongate solutions
2113and restrict residuals between two independent meshes. By passing a sequence
2114of such two-level transfer operators to MGTransferMF and, finally, to Multigrid,
2115non-nested multigrid can be computed.
2116- Utilities::MPI::RemotePointEvaluation can be used to couple non-matching
2117grids via surfaces (example: fluid-structure interaction, independently created
2118grids). The evaluation points can come from any side (pointwise interpolation)
2119or are defined on intersected meshes (Nitsche-type mortaring
2120@cite heinz2022high). Concerning the creation of such intersected meshes and the
2121corresponding evaluation points, see
2122GridTools::internal::distributed\_compute_intersection_locations().
2123- Alternatively to the coupling via Utilities::MPI::RemotePointEvaluation,
2124preCICE @cite bungartz2016precice @cite chourdakis2021precice can be used. The
2125code-gallery program "Laplace equation coupled to an external simulation
2126program" shows how to use this library with deal.II.
2127 *
2128 *
2129<a name="step_87-PlainProg"></a>
2130<h1> The plain program</h1>
2131@include "step-87.cc"
2132*/
void distribute_local_to_global(const InVector &local_vector, const std::vector< size_type > &local_dof_indices, OutVector &global_vector) const
void attach_triangulation(const Triangulation< dim, spacedim > &)
void add_data_vector(const VectorType &data, const std::vector< std::string > &names, const DataVectorType type=type_automatic, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation={})
void test_and_sum(const StridedArrayView< ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags, const bool sum_into_values=false)
Definition fe_q.h:554
void reinit(const size_type size, const bool omit_zeroing_entries=false)
Abstract base class for mapping classes.
Definition mapping.h:320
Definition point.h:113
Communicate values between a mesh and arbitrary points.
void process_and_evaluate(const std::vector< DataType > &input, std::vector< DataType > &buffer, const std::function< void(const ArrayView< const DataType > &, const CellData &)> &evaluation_function, const bool sort_data=true) const
void reinit(const std::vector< Point< spacedim > > &points, const Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping)
const unsigned int DoFAccessor< structdim, dim, spacedim, level_dof_access >::dimension
#define Assert(cond, exc)
#define AssertThrow(cond, exc)
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())
Definition loop.h:564
const bool IsBlockVector< VectorType >::value
@ update_values
Shape function values.
@ update_normal_vectors
Normal vectors.
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
@ update_quadrature_points
Transformed quadrature points.
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
Definition mapping_q1.h:104
const Event initial
Definition event.cc:68
void support_point_wise(DoFHandler< dim, spacedim > &dof_handler)
IndexSet extract_locally_active_dofs(const DoFHandler< dim, spacedim > &dof_handler)
void extrapolate(const DoFHandler< dim, spacedim > &dof1, const InVector &z1, const DoFHandler< dim, spacedim > &dof2, OutVector &z2)
void hyper_sphere(Triangulation< spacedim - 1, spacedim > &tria, const Point< spacedim > &center=Point< spacedim >(), const double radius=1.)
void hyper_cube(Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false)
@ valid
Iterator points to a valid object.
constexpr char N
constexpr char T
constexpr types::blas_int zero
constexpr types::blas_int one
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition utilities.cc:193
Tensor< 2, dim, Number > w(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
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)
Tensor< 2, dim, Number > F(const Tensor< 2, dim, Number > &Grad_u)
VectorType::value_type * end(VectorType &V)
T sum(const T &t, const MPI_Comm mpi_communicator)
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
Definition mpi.cc:120
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)
std::string compress(const std::string &input)
Definition utilities.cc:383
void interpolate(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Function< spacedim, typename VectorType::value_type > &function, VectorType &vec, const ComponentMask &component_mask={}, const unsigned int level=numbers::invalid_unsigned_int)
std::vector< typename FEPointEvaluation< n_components, dim, spacedim, typename VectorType::value_type >::value_type > point_values(const Mapping< dim > &mapping, const MeshType< dim, spacedim > &mesh, const VectorType &vector, const std::vector< Point< spacedim > > &evaluation_points, Utilities::MPI::RemotePointEvaluation< dim, spacedim > &cache, const EvaluationFlags::EvaluationFlags flags=EvaluationFlags::avg, const unsigned int first_selected_component=0)
std::vector< typename FEPointEvaluation< n_components, dim, spacedim, typename VectorType::value_type >::gradient_type > point_gradients(const Mapping< dim > &mapping, const MeshType< dim, spacedim > &mesh, const VectorType &vector, const std::vector< Point< spacedim > > &evaluation_points, Utilities::MPI::RemotePointEvaluation< dim, spacedim > &cache, const EvaluationFlags::EvaluationFlags flags=EvaluationFlags::avg, const unsigned int first_selected_component=0)
int(& functions)(const void *v1, const void *v2)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
Definition types.h:32