DOLFINx 0.10.0.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
generation.h
1// Copyright (C) 2005-2024 Anders Logg and Garth N. Wells
2//
3// This file is part of DOLFINx (https://www.fenicsproject.org)
4//
5// SPDX-License-Identifier: LGPL-3.0-or-later
6
7#pragma once
8
9#include "Mesh.h"
10#include "cell_types.h"
11#include "utils.h"
12#include <algorithm>
13#include <array>
14#include <cmath>
15#include <concepts>
16#include <cstddef>
17#include <cstdint>
18#include <limits>
19#include <mpi.h>
20#include <stdexcept>
21#include <utility>
22#include <vector>
23
24namespace dolfinx::mesh
25{
27enum class DiagonalType
28{
29 left,
30 right,
31 crossed,
32 shared_facet,
33 left_right,
34 right_left
35};
36
37namespace impl
38{
39template <std::floating_point T>
40std::tuple<std::vector<T>, std::vector<std::int64_t>>
41create_interval_cells(std::array<T, 2> p, std::int64_t n);
42
43template <std::floating_point T>
44Mesh<T> build_tri(MPI_Comm comm, std::array<std::array<T, 2>, 2> p,
45 std::array<std::int64_t, 2> n,
46 const CellPartitionFunction& partitioner,
47 DiagonalType diagonal);
48
49template <std::floating_point T>
50Mesh<T> build_quad(MPI_Comm comm, const std::array<std::array<T, 2>, 2> p,
51 std::array<std::int64_t, 2> n,
52 const CellPartitionFunction& partitioner);
53
54template <std::floating_point T>
55std::vector<T> create_geom(MPI_Comm comm, std::array<std::array<T, 3>, 2> p,
56 std::array<std::int64_t, 3> n);
57
58template <std::floating_point T>
59Mesh<T> build_tet(MPI_Comm comm, MPI_Comm subcomm,
60 std::array<std::array<T, 3>, 2> p,
61 std::array<std::int64_t, 3> n,
62 const CellPartitionFunction& partitioner);
63
64template <std::floating_point T>
65Mesh<T> build_hex(MPI_Comm comm, MPI_Comm subcomm,
66 std::array<std::array<T, 3>, 2> p,
67 std::array<std::int64_t, 3> n,
68 const CellPartitionFunction& partitioner);
69
70template <std::floating_point T>
71Mesh<T> build_prism(MPI_Comm comm, MPI_Comm subcomm,
72 std::array<std::array<T, 3>, 2> p,
73 std::array<std::int64_t, 3> n,
74 const CellPartitionFunction& partitioner);
75} // namespace impl
76
97template <std::floating_point T = double>
98Mesh<T> create_box(MPI_Comm comm, MPI_Comm subcomm,
99 std::array<std::array<T, 3>, 2> p,
100 std::array<std::int64_t, 3> n, CellType celltype,
101 CellPartitionFunction partitioner = nullptr)
102{
103 if (std::ranges::any_of(n, [](auto e) { return e < 1; }))
104 throw std::runtime_error("At least one cell per dimension is required");
105
106 for (int32_t i = 0; i < 3; i++)
107 {
108 if (p[0][i] >= p[1][i])
109 throw std::runtime_error("It must hold p[0] < p[1].");
110 }
111
112 if (!partitioner and dolfinx::MPI::size(comm) > 1)
113 partitioner = create_cell_partitioner();
114
115 switch (celltype)
116 {
117 case CellType::tetrahedron:
118 return impl::build_tet<T>(comm, subcomm, p, n, partitioner);
119 case CellType::hexahedron:
120 return impl::build_hex<T>(comm, subcomm, p, n, partitioner);
121 case CellType::prism:
122 return impl::build_prism<T>(comm, subcomm, p, n, partitioner);
123 default:
124 throw std::runtime_error("Generate box mesh. Wrong cell type");
125 }
126}
127
144template <std::floating_point T = double>
145Mesh<T> create_box(MPI_Comm comm, std::array<std::array<T, 3>, 2> p,
146 std::array<std::int64_t, 3> n, CellType celltype,
147 const CellPartitionFunction& partitioner = nullptr)
148{
149 return create_box<T>(comm, comm, p, n, celltype, partitioner);
150}
151
168template <std::floating_point T = double>
169Mesh<T> create_rectangle(MPI_Comm comm, std::array<std::array<T, 2>, 2> p,
170 std::array<std::int64_t, 2> n, CellType celltype,
171 CellPartitionFunction partitioner,
172 DiagonalType diagonal = DiagonalType::right)
173{
174 if (std::ranges::any_of(n, [](auto e) { return e < 1; }))
175 throw std::runtime_error("At least one cell per dimension is required");
176
177 for (int32_t i = 0; i < 2; i++)
178 {
179 if (p[0][i] >= p[1][i])
180 throw std::runtime_error("It must hold p[0] < p[1].");
181 }
182
183 if (!partitioner and dolfinx::MPI::size(comm) > 1)
184 partitioner = create_cell_partitioner();
185
186 switch (celltype)
187 {
188 case CellType::triangle:
189 return impl::build_tri<T>(comm, p, n, partitioner, diagonal);
190 case CellType::quadrilateral:
191 return impl::build_quad<T>(comm, p, n, partitioner);
192 default:
193 throw std::runtime_error("Generate rectangle mesh. Wrong cell type");
194 }
195}
196
211template <std::floating_point T = double>
212Mesh<T> create_rectangle(MPI_Comm comm, std::array<std::array<T, 2>, 2> p,
213 std::array<std::int64_t, 2> n, CellType celltype,
214 DiagonalType diagonal = DiagonalType::right)
215{
216 return create_rectangle<T>(comm, p, n, celltype, nullptr, diagonal);
217}
218
232template <std::floating_point T = double>
233Mesh<T> create_interval(MPI_Comm comm, std::int64_t n, std::array<T, 2> p,
234 mesh::GhostMode ghost_mode = mesh::GhostMode::none,
235 CellPartitionFunction partitioner = nullptr)
236{
237 if (n < 1)
238 throw std::runtime_error("At least one cell is required.");
239
240 const auto [a, b] = p;
241 if (a >= b)
242 throw std::runtime_error("It must hold p[0] < p[1].");
243 if (std::abs(a - b) < std::numeric_limits<T>::epsilon())
244 {
245 throw std::runtime_error(
246 "Length of interval is zero. Check your dimensions.");
247 }
248
249 if (!partitioner and dolfinx::MPI::size(comm) > 1)
250 partitioner = create_cell_partitioner(ghost_mode);
251
252 fem::CoordinateElement<T> element(CellType::interval, 1);
253 if (dolfinx::MPI::rank(comm) == 0)
254 {
255 auto [x, cells] = impl::create_interval_cells<T>(p, n);
256 return create_mesh(comm, MPI_COMM_SELF, cells, element, MPI_COMM_SELF, x,
257 {x.size(), 1}, partitioner);
258 }
259 else
260 {
261 return create_mesh(comm, MPI_COMM_NULL, {}, element, MPI_COMM_NULL,
262 std::vector<T>{}, {0, 1}, partitioner);
263 }
264}
265
266namespace impl
267{
268
269template <std::floating_point T>
270std::tuple<std::vector<T>, std::vector<std::int64_t>>
271create_interval_cells(std::array<T, 2> p, std::int64_t n)
272{
273 const auto [a, b] = p;
274
275 const T h = (b - a) / static_cast<T>(n);
276
277 // Create vertices
278 std::vector<T> x(n + 1);
279 std::ranges::generate(x, [i = std::int64_t(0), a, h]() mutable
280 { return a + h * static_cast<T>(i++); });
281
282 // Create intervals -> cells=[0, 1, 1, ..., n-1, n-1, n]
283 std::vector<std::int64_t> cells(2 * n);
284 for (std::size_t ix = 0; ix < cells.size() / 2; ++ix)
285 {
286 cells[2 * ix] = ix;
287 cells[2 * ix + 1] = ix + 1;
288 }
289
290 return {std::move(x), std::move(cells)};
291}
292
293template <std::floating_point T>
294std::vector<T> create_geom(MPI_Comm comm, std::array<std::array<T, 3>, 2> p,
295 std::array<std::int64_t, 3> n)
296{
297 // Extract data
298 auto [p0, p1] = p;
299 const auto [nx, ny, nz] = n;
300
301 assert(std::ranges::all_of(n, [](auto e) { return e >= 1; }));
302 for (std::int64_t i = 0; i < 3; i++)
303 assert(p0[i] < p1[i]);
304
305 // Structured grid cuboid extents
306 const std::array<T, 3> extents = {
307 (p1[0] - p0[0]) / static_cast<T>(nx),
308 (p1[1] - p0[1]) / static_cast<T>(ny),
309 (p1[2] - p0[2]) / static_cast<T>(nz),
310 };
311
312 if (std::ranges::any_of(
313 extents, [](auto e)
314 { return std::abs(e) < 2.0 * std::numeric_limits<T>::epsilon(); }))
315 {
316 throw std::runtime_error(
317 "Box seems to have zero width, height or depth. Check dimensions");
318 }
319
320 const std::int64_t n_points = (nx + 1) * (ny + 1) * (nz + 1);
321 const auto [range_begin, range_end] = dolfinx::MPI::local_range(
322 dolfinx::MPI::rank(comm), n_points, dolfinx::MPI::size(comm));
323
324 std::vector<T> geom;
325 geom.reserve((range_end - range_begin) * 3);
326 const std::int64_t sqxy = (nx + 1) * (ny + 1);
327 for (std::int64_t v = range_begin; v < range_end; ++v)
328 {
329 // lexiographic index to spatial index
330 const std::int64_t p = v % sqxy;
331 std::array<std::int64_t, 3> idx = {p % (nx + 1), p / (nx + 1), v / sqxy};
332
333 // vertex = p0 + idx * extents (elementwise)
334 for (std::size_t i = 0; i < idx.size(); i++)
335 geom.push_back(p0[i] + static_cast<T>(idx[i]) * extents[i]);
336 }
337
338 return geom;
339}
340
341template <std::floating_point T>
342Mesh<T> build_tet(MPI_Comm comm, MPI_Comm subcomm,
343 std::array<std::array<T, 3>, 2> p,
344 std::array<std::int64_t, 3> n,
345 const CellPartitionFunction& partitioner)
346{
347 common::Timer timer("Build BoxMesh (tetrahedra)");
348 std::vector<T> x;
349 std::vector<std::int64_t> cells;
350 fem::CoordinateElement<T> element(CellType::tetrahedron, 1);
351 if (subcomm != MPI_COMM_NULL)
352 {
353 x = create_geom<T>(subcomm, p, n);
354
355 const auto [nx, ny, nz] = n;
356 const std::int64_t n_cells = nx * ny * nz;
357
358 std::array range_c = dolfinx::MPI::local_range(
359 dolfinx::MPI::rank(subcomm), n_cells, dolfinx::MPI::size(subcomm));
360 cells.reserve(6 * (range_c[1] - range_c[0]) * 4);
361
362 // Create tetrahedra
363 for (std::int64_t i = range_c[0]; i < range_c[1]; ++i)
364 {
365 const std::int64_t iz = i / (nx * ny);
366 const std::int64_t j = i % (nx * ny);
367 const std::int64_t iy = j / nx;
368 const std::int64_t ix = j % nx;
369 const std::int64_t v0 = iz * (nx + 1) * (ny + 1) + iy * (nx + 1) + ix;
370 const std::int64_t v1 = v0 + 1;
371 const std::int64_t v2 = v0 + (nx + 1);
372 const std::int64_t v3 = v1 + (nx + 1);
373 const std::int64_t v4 = v0 + (nx + 1) * (ny + 1);
374 const std::int64_t v5 = v1 + (nx + 1) * (ny + 1);
375 const std::int64_t v6 = v2 + (nx + 1) * (ny + 1);
376 const std::int64_t v7 = v3 + (nx + 1) * (ny + 1);
377
378 // Note that v0 < v1 < v2 < v3 < vmid
379 cells.insert(cells.end(),
380 {v0, v1, v3, v7, v0, v1, v7, v5, v0, v5, v7, v4,
381 v0, v3, v2, v7, v0, v6, v4, v7, v0, v2, v6, v7});
382 }
383 }
384
385 return create_mesh(comm, subcomm, cells, element, subcomm, x,
386 {x.size() / 3, 3}, partitioner);
387}
388
389template <std::floating_point T>
390mesh::Mesh<T> build_hex(MPI_Comm comm, MPI_Comm subcomm,
391 std::array<std::array<T, 3>, 2> p,
392 std::array<std::int64_t, 3> n,
393 const CellPartitionFunction& partitioner)
394{
395 common::Timer timer("Build BoxMesh (hexahedra)");
396 std::vector<T> x;
397 std::vector<std::int64_t> cells;
398 fem::CoordinateElement<T> element(CellType::hexahedron, 1);
399 if (subcomm != MPI_COMM_NULL)
400 {
401 x = create_geom<T>(subcomm, p, n);
402
403 // Create cuboids
404 const auto [nx, ny, nz] = n;
405 const std::int64_t n_cells = nx * ny * nz;
406 std::array range_c = dolfinx::MPI::local_range(
407 dolfinx::MPI::rank(subcomm), n_cells, dolfinx::MPI::size(subcomm));
408 cells.reserve((range_c[1] - range_c[0]) * 8);
409 for (std::int64_t i = range_c[0]; i < range_c[1]; ++i)
410 {
411 const std::int64_t iz = i / (nx * ny);
412 const std::int64_t j = i % (nx * ny);
413 const std::int64_t iy = j / nx;
414 const std::int64_t ix = j % nx;
415
416 const std::int64_t v0 = (iz * (ny + 1) + iy) * (nx + 1) + ix;
417 const std::int64_t v1 = v0 + 1;
418 const std::int64_t v2 = v0 + (nx + 1);
419 const std::int64_t v3 = v1 + (nx + 1);
420 const std::int64_t v4 = v0 + (nx + 1) * (ny + 1);
421 const std::int64_t v5 = v1 + (nx + 1) * (ny + 1);
422 const std::int64_t v6 = v2 + (nx + 1) * (ny + 1);
423 const std::int64_t v7 = v3 + (nx + 1) * (ny + 1);
424 cells.insert(cells.end(), {v0, v1, v2, v3, v4, v5, v6, v7});
425 }
426 }
427
428 return create_mesh(comm, subcomm, cells, element, subcomm, x,
429 {x.size() / 3, 3}, partitioner);
430}
431
432template <std::floating_point T>
433Mesh<T> build_prism(MPI_Comm comm, MPI_Comm subcomm,
434 std::array<std::array<T, 3>, 2> p,
435 std::array<std::int64_t, 3> n,
436 const CellPartitionFunction& partitioner)
437{
438 std::vector<T> x;
439 std::vector<std::int64_t> cells;
440 fem::CoordinateElement<T> element(CellType::prism, 1);
441 if (subcomm != MPI_COMM_NULL)
442 {
443 x = create_geom<T>(subcomm, p, n);
444
445 const std::int64_t nx = n[0];
446 const std::int64_t ny = n[1];
447 const std::int64_t nz = n[2];
448 const std::int64_t n_cells = nx * ny * nz;
449 std::array range_c = dolfinx::MPI::local_range(
450 dolfinx::MPI::rank(comm), n_cells, dolfinx::MPI::size(comm));
451 const std::int64_t cell_range = range_c[1] - range_c[0];
452
453 // Create cuboids
454 cells.reserve(2 * cell_range * 6);
455 for (std::int64_t i = range_c[0]; i < range_c[1]; ++i)
456 {
457 const std::int64_t iz = i / (nx * ny);
458 const std::int64_t j = i % (nx * ny);
459 const std::int64_t iy = j / nx;
460 const std::int64_t ix = j % nx;
461
462 const std::int64_t v0 = (iz * (ny + 1) + iy) * (nx + 1) + ix;
463 const std::int64_t v1 = v0 + 1;
464 const std::int64_t v2 = v0 + (nx + 1);
465 const std::int64_t v3 = v1 + (nx + 1);
466 const std::int64_t v4 = v0 + (nx + 1) * (ny + 1);
467 const std::int64_t v5 = v1 + (nx + 1) * (ny + 1);
468 const std::int64_t v6 = v2 + (nx + 1) * (ny + 1);
469 const std::int64_t v7 = v3 + (nx + 1) * (ny + 1);
470 cells.insert(cells.end(), {v0, v1, v2, v4, v5, v6});
471 cells.insert(cells.end(), {v1, v2, v3, v5, v6, v7});
472 }
473 }
474
475 return create_mesh(comm, subcomm, cells, element, subcomm, x,
476 {x.size() / 3, 3}, partitioner);
477}
478
479template <std::floating_point T>
480Mesh<T> build_tri(MPI_Comm comm, std::array<std::array<T, 2>, 2> p,
481 std::array<std::int64_t, 2> n,
482 const CellPartitionFunction& partitioner,
483 DiagonalType diagonal)
484{
485 fem::CoordinateElement<T> element(CellType::triangle, 1);
486 if (dolfinx::MPI::rank(comm) == 0)
487 {
488 const auto [p0, p1] = p;
489 const auto [nx, ny] = n;
490
491 const auto [a, c] = p0;
492 const auto [b, d] = p1;
493
494 const T ab = (b - a) / static_cast<T>(nx);
495 const T cd = (d - c) / static_cast<T>(ny);
496 if (std::abs(b - a) < std::numeric_limits<T>::epsilon()
497 or std::abs(d - c) < std::numeric_limits<T>::epsilon())
498 {
499 throw std::runtime_error("Rectangle seems to have zero width, height or "
500 "depth. Check dimensions");
501 }
502
503 // Create vertices and cells
504 std::int64_t nv, nc;
505 switch (diagonal)
506 {
507 case DiagonalType::crossed:
508 nv = (nx + 1) * (ny + 1) + nx * ny;
509 nc = 4 * nx * ny;
510 break;
511 default:
512 nv = (nx + 1) * (ny + 1);
513 nc = 2 * nx * ny;
514 }
515
516 std::vector<T> x;
517 x.reserve(nv * 2);
518 std::vector<std::int64_t> cells;
519 cells.reserve(nc * 3);
520
521 // Create main vertices
522 std::int64_t vertex = 0;
523 for (std::int64_t iy = 0; iy <= ny; iy++)
524 {
525 T x1 = c + cd * static_cast<T>(iy);
526 for (std::int64_t ix = 0; ix <= nx; ix++)
527 x.insert(x.end(), {a + ab * static_cast<T>(ix), x1});
528 }
529
530 // Create midpoint vertices if the mesh type is crossed
531 switch (diagonal)
532 {
533 case DiagonalType::crossed:
534 for (std::int64_t iy = 0; iy < ny; iy++)
535 {
536 T x1 = c + cd * (static_cast<T>(iy) + 0.5);
537 for (std::int64_t ix = 0; ix < nx; ix++)
538 {
539 T x0 = a + ab * (static_cast<T>(ix) + 0.5);
540 x.insert(x.end(), {x0, x1});
541 }
542 }
543 break;
544 default:
545 break;
546 }
547
548 // Create triangles
549 switch (diagonal)
550 {
551 case DiagonalType::crossed:
552 {
553 for (std::int64_t iy = 0; iy < ny; iy++)
554 {
555 for (std::int64_t ix = 0; ix < nx; ix++)
556 {
557 std::int64_t v0 = iy * (nx + 1) + ix;
558 std::int64_t v1 = v0 + 1;
559 std::int64_t v2 = v0 + (nx + 1);
560 std::int64_t v3 = v1 + (nx + 1);
561 std::int64_t vmid = (nx + 1) * (ny + 1) + iy * nx + ix;
562
563 // Note that v0 < v1 < v2 < v3 < vmid
564 cells.insert(cells.end(), {v0, v1, vmid, v0, v2, vmid, v1, v3, vmid,
565 v2, v3, vmid});
566 }
567 }
568 break;
569 }
570 default:
571 {
572 DiagonalType local_diagonal = diagonal;
573 for (std::int64_t iy = 0; iy < ny; iy++)
574 {
575 // Set up alternating diagonal
576 switch (diagonal)
577 {
578 case DiagonalType::right_left:
579 if (iy % 2)
580 local_diagonal = DiagonalType::right;
581 else
582 local_diagonal = DiagonalType::left;
583 break;
584 case DiagonalType::left_right:
585 if (iy % 2)
586 local_diagonal = DiagonalType::left;
587 else
588 local_diagonal = DiagonalType::right;
589 break;
590 default:
591 break;
592 }
593 for (std::int64_t ix = 0; ix < nx; ix++)
594 {
595 std::int64_t v0 = iy * (nx + 1) + ix;
596 std::int64_t v1 = v0 + 1;
597 std::int64_t v2 = v0 + (nx + 1);
598 std::int64_t v3 = v1 + (nx + 1);
599 switch (local_diagonal)
600 {
601 case DiagonalType::left:
602 {
603 cells.insert(cells.end(), {v0, v1, v2, v1, v2, v3});
604 if (diagonal == DiagonalType::right_left
605 or diagonal == DiagonalType::left_right)
606 {
607 local_diagonal = DiagonalType::right;
608 }
609 break;
610 }
611 default:
612 {
613 cells.insert(cells.end(), {v0, v1, v3, v0, v2, v3});
614 if (diagonal == DiagonalType::right_left
615 or diagonal == DiagonalType::left_right)
616 {
617 local_diagonal = DiagonalType::left;
618 }
619 }
620 }
621 }
622 }
623 }
624 }
625
626 return create_mesh(comm, MPI_COMM_SELF, cells, element, MPI_COMM_SELF, x,
627 {x.size() / 2, 2}, partitioner);
628 }
629 else
630 {
631 return create_mesh(comm, MPI_COMM_NULL, {}, element, MPI_COMM_NULL,
632 std::vector<T>{}, {0, 2}, partitioner);
633 }
634}
635
636template <std::floating_point T>
637Mesh<T> build_quad(MPI_Comm comm, const std::array<std::array<T, 2>, 2> p,
638 std::array<std::int64_t, 2> n,
639 const CellPartitionFunction& partitioner)
640{
641 fem::CoordinateElement<T> element(CellType::quadrilateral, 1);
642 if (dolfinx::MPI::rank(comm) == 0)
643 {
644 const auto [nx, ny] = n;
645 const auto [a, c] = p[0];
646 const auto [b, d] = p[1];
647
648 const T ab = (b - a) / static_cast<T>(nx);
649 const T cd = (d - c) / static_cast<T>(ny);
650
651 // Create vertices
652 std::vector<T> x;
653 x.reserve((nx + 1) * (ny + 1) * 2);
654 std::int64_t vertex = 0;
655 for (std::int64_t ix = 0; ix <= nx; ix++)
656 {
657 T x0 = a + ab * static_cast<T>(ix);
658 for (std::int64_t iy = 0; iy <= ny; iy++)
659 x.insert(x.end(), {x0, c + cd * static_cast<T>(iy)});
660 }
661
662 // Create rectangles
663 std::vector<std::int64_t> cells;
664 cells.reserve(nx * ny * 4);
665 for (std::int64_t ix = 0; ix < nx; ix++)
666 {
667 for (std::int64_t iy = 0; iy < ny; iy++)
668 {
669 std::int64_t i0 = ix * (ny + 1);
670 cells.insert(cells.end(), {i0 + iy, i0 + iy + 1, i0 + iy + ny + 1,
671 i0 + iy + ny + 2});
672 }
673 }
674
675 return create_mesh(comm, MPI_COMM_SELF, cells, element, MPI_COMM_SELF, x,
676 {x.size() / 2, 2}, partitioner);
677 }
678 else
679 {
680 return create_mesh(comm, MPI_COMM_NULL, {}, element, MPI_COMM_NULL,
681 std::vector<T>{}, {0, 2}, partitioner);
682 }
683}
684} // namespace impl
685} // namespace dolfinx::mesh
Definition XDMFFile.h:27
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition Mesh.h:23
Functions supporting mesh operations.
int size(MPI_Comm comm)
Definition MPI.cpp:72
int rank(MPI_Comm comm)
Return process rank for the communicator.
Definition MPI.cpp:64
constexpr std::array< std::int64_t, 2 > local_range(int rank, std::int64_t N, int size)
Return local range for the calling process, partitioning the global [0, N - 1] range across all ranks...
Definition MPI.h:91
void cells(la::SparsityPattern &pattern, std::array< std::span< const std::int32_t >, 2 > cells, std::array< std::reference_wrapper< const DofMap >, 2 > dofmaps)
Iterate over cells and insert entries into sparsity pattern.
Definition sparsitybuild.cpp:16
Mesh data structures and algorithms on meshes.
Definition DofMap.h:32
Mesh< T > create_box(MPI_Comm comm, MPI_Comm subcomm, std::array< std::array< T, 3 >, 2 > p, std::array< std::int64_t, 3 > n, CellType celltype, CellPartitionFunction partitioner=nullptr)
Create a uniform mesh::Mesh over rectangular prism spanned by the two points p.
Definition generation.h:98
GhostMode
Enum for different partitioning ghost modes.
Definition utils.h:36
DiagonalType
Enum for different diagonal types.
Definition generation.h:28
Mesh< typename std::remove_reference_t< typename U::value_type > > create_mesh(MPI_Comm comm, MPI_Comm commt, std::span< const std::int64_t > cells, const fem::CoordinateElement< typename std::remove_reference_t< typename U::value_type > > &element, MPI_Comm commg, const U &x, std::array< std::size_t, 2 > xshape, const CellPartitionFunction &partitioner)
Create a distributed mesh from mesh data using a provided graph partitioning function for determining...
Definition utils.h:783
std::vector< T > h(const Mesh< T > &mesh, std::span< const std::int32_t > entities, int dim)
Compute greatest distance between any two vertices of the mesh entities (h).
Definition utils.h:211
Mesh< T > create_rectangle(MPI_Comm comm, std::array< std::array< T, 2 >, 2 > p, std::array< std::int64_t, 2 > n, CellType celltype, CellPartitionFunction partitioner, DiagonalType diagonal=DiagonalType::right)
Create a uniform mesh::Mesh over the rectangle spanned by the two points p.
Definition generation.h:169
Mesh< T > create_interval(MPI_Comm comm, std::int64_t n, std::array< T, 2 > p, mesh::GhostMode ghost_mode=mesh::GhostMode::none, CellPartitionFunction partitioner=nullptr)
Interval mesh of the 1D line [a, b].
Definition generation.h:233
CellType
Cell type identifier.
Definition cell_types.h:22
std::function< graph::AdjacencyList< std::int32_t >( MPI_Comm comm, int nparts, const std::vector< CellType > &cell_types, const std::vector< std::span< const std::int64_t > > &cells)> CellPartitionFunction
Signature for the cell partitioning function. The function should compute the destination rank for ce...
Definition utils.h:185
CellPartitionFunction create_cell_partitioner(mesh::GhostMode ghost_mode=mesh::GhostMode::none, const graph::partition_fn &partfn=&graph::partition_graph)
Definition utils.cpp:85