Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/d6/d6c/generation_8h_source.html
DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
generation.h
1// Copyright (C) 2005-2023 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 <array>
13#include <cfloat>
14#include <concepts>
15#include <cstddef>
16#include <limits>
17#include <mpi.h>
18#include <vector>
19
20namespace dolfinx::mesh
21{
23enum class DiagonalType
24{
25 left,
26 right,
27 crossed,
28 shared_facet,
29 left_right,
30 right_left
31};
32
33namespace impl
34{
35template <std::floating_point T>
36Mesh<T> build_tri(MPI_Comm comm, std::array<std::array<double, 2>, 2> p,
37 std::array<std::size_t, 2> n,
38 const CellPartitionFunction& partitioner,
39 DiagonalType diagonal);
40
41template <std::floating_point T>
42Mesh<T> build_quad(MPI_Comm comm, const std::array<std::array<double, 2>, 2> p,
43 std::array<std::size_t, 2> n,
44 const CellPartitionFunction& partitioner);
45
46template <std::floating_point T>
47std::vector<T> create_geom(MPI_Comm comm,
48 std::array<std::array<double, 3>, 2> p,
49 std::array<std::size_t, 3> n);
50
51template <std::floating_point T>
52Mesh<T> build_tet(MPI_Comm comm, MPI_Comm subcomm,
53 std::array<std::array<double, 3>, 2> p,
54 std::array<std::size_t, 3> n,
55 const CellPartitionFunction& partitioner);
56
57template <std::floating_point T>
58Mesh<T> build_hex(MPI_Comm comm, MPI_Comm subcomm,
59 std::array<std::array<double, 3>, 2> p,
60 std::array<std::size_t, 3> n,
61 const CellPartitionFunction& partitioner);
62
63template <std::floating_point T>
64Mesh<T> build_prism(MPI_Comm comm, MPI_Comm subcomm,
65 std::array<std::array<double, 3>, 2> p,
66 std::array<std::size_t, 3> n,
67 const CellPartitionFunction& partitioner);
68} // namespace impl
69
90template <std::floating_point T = double>
91Mesh<T> create_box(MPI_Comm comm, MPI_Comm subcomm,
92 std::array<std::array<double, 3>, 2> p,
93 std::array<std::size_t, 3> n, CellType celltype,
94 CellPartitionFunction partitioner = nullptr)
95{
96 if (!partitioner and dolfinx::MPI::size(comm) > 1)
97 partitioner = create_cell_partitioner();
98
99 switch (celltype)
100 {
101 case CellType::tetrahedron:
102 return impl::build_tet<T>(comm, subcomm, p, n, partitioner);
103 case CellType::hexahedron:
104 return impl::build_hex<T>(comm, subcomm, p, n, partitioner);
105 case CellType::prism:
106 return impl::build_prism<T>(comm, subcomm, p, n, partitioner);
107 default:
108 throw std::runtime_error("Generate box mesh. Wrong cell type");
109 }
110}
111
128template <std::floating_point T = double>
129Mesh<T> create_box(MPI_Comm comm, std::array<std::array<double, 3>, 2> p,
130 std::array<std::size_t, 3> n, CellType celltype,
131 const CellPartitionFunction& partitioner = nullptr)
132{
133 return create_box<T>(comm, comm, p, n, celltype, partitioner);
134}
135
152template <std::floating_point T = double>
153Mesh<T> create_rectangle(MPI_Comm comm, std::array<std::array<double, 2>, 2> p,
154 std::array<std::size_t, 2> n, CellType celltype,
155 CellPartitionFunction partitioner,
156 DiagonalType diagonal = DiagonalType::right)
157{
158 if (!partitioner and dolfinx::MPI::size(comm) > 1)
159 partitioner = create_cell_partitioner();
160
161 switch (celltype)
162 {
163 case CellType::triangle:
164 return impl::build_tri<T>(comm, p, n, partitioner, diagonal);
165 case CellType::quadrilateral:
166 return impl::build_quad<T>(comm, p, n, partitioner);
167 default:
168 throw std::runtime_error("Generate rectangle mesh. Wrong cell type");
169 }
170}
171
186template <std::floating_point T = double>
187Mesh<T> create_rectangle(MPI_Comm comm, std::array<std::array<double, 2>, 2> p,
188 std::array<std::size_t, 2> n, CellType celltype,
189 DiagonalType diagonal = DiagonalType::right)
190{
191 return create_rectangle<T>(comm, p, n, celltype, nullptr, diagonal);
192}
193
206template <std::floating_point T = double>
207Mesh<T> create_interval(MPI_Comm comm, std::size_t nx, std::array<double, 2> p,
208 CellPartitionFunction partitioner = nullptr)
209{
210 if (!partitioner and dolfinx::MPI::size(comm) > 1)
211 partitioner = create_cell_partitioner();
212
213 fem::CoordinateElement<T> element(CellType::interval, 1);
214 std::vector<T> x;
215 std::vector<std::int64_t> cells;
216 if (dolfinx::MPI::rank(comm) == 0)
217 {
218 const T a = p[0];
219 const T b = p[1];
220 const T ab = (b - a) / static_cast<T>(nx);
221
222 if (std::abs(a - b) < std::numeric_limits<double>::epsilon())
223 {
224 throw std::runtime_error(
225 "Length of interval is zero. Check your dimensions.");
226 }
227
228 if (b < a)
229 {
230 throw std::runtime_error(
231 "Interval length is negative. Check order of arguments.");
232 }
233
234 if (nx < 1)
235 throw std::runtime_error(
236 "Number of points on interval must be at least 1");
237
238 // Create vertices
239 x.resize(nx + 1);
240 for (std::size_t ix = 0; ix <= nx; ix++)
241 x[ix] = a + ab * static_cast<T>(ix);
242
243 // Create intervals
244 cells.resize(nx * 2);
245 for (std::size_t ix = 0; ix < nx; ++ix)
246 for (std::size_t j = 0; j < 2; ++j)
247 cells[2 * ix + j] = ix + j;
248
249 return create_mesh(comm, MPI_COMM_SELF, cells, element, MPI_COMM_SELF, x,
250 {x.size(), 1}, partitioner);
251 }
252 else
253 {
254 return create_mesh(comm, MPI_COMM_NULL, {}, element, MPI_COMM_NULL, x,
255 {x.size(), 1}, partitioner);
256 }
257}
258
259namespace impl
260{
261template <std::floating_point T>
262std::vector<T> create_geom(MPI_Comm comm,
263 std::array<std::array<double, 3>, 2> p,
264 std::array<std::size_t, 3> n)
265{
266 // Extract data
267 const std::array<double, 3> p0 = p[0];
268 const std::array<double, 3> p1 = p[1];
269 std::int64_t nx = n[0];
270 std::int64_t ny = n[1];
271 std::int64_t nz = n[2];
272
273 const std::int64_t n_points = (nx + 1) * (ny + 1) * (nz + 1);
274 std::array range_p = dolfinx::MPI::local_range(
275 dolfinx::MPI::rank(comm), n_points, dolfinx::MPI::size(comm));
276
277 // Extract minimum and maximum coordinates
278 const double x0 = std::min(p0[0], p1[0]);
279 const double x1 = std::max(p0[0], p1[0]);
280 const double y0 = std::min(p0[1], p1[1]);
281 const double y1 = std::max(p0[1], p1[1]);
282 const double z0 = std::min(p0[2], p1[2]);
283 const double z1 = std::max(p0[2], p1[2]);
284
285 const T a = x0;
286 const T b = x1;
287 const T ab = (b - a) / static_cast<T>(nx);
288 const T c = y0;
289 const T d = y1;
290 const T cd = (d - c) / static_cast<T>(ny);
291 const T e = z0;
292 const T f = z1;
293 const T ef = (f - e) / static_cast<T>(nz);
294
295 if (std::abs(x0 - x1) < 2.0 * std::numeric_limits<double>::epsilon()
296 or std::abs(y0 - y1) < 2.0 * std::numeric_limits<double>::epsilon()
297 or std::abs(z0 - z1) < 2.0 * std::numeric_limits<double>::epsilon())
298 {
299 throw std::runtime_error(
300 "Box seems to have zero width, height or depth. Check dimensions");
301 }
302
303 if (nx < 1 or ny < 1 or nz < 1)
304 {
305 throw std::runtime_error(
306 "BoxMesh has non-positive number of vertices in some dimension");
307 }
308
309 std::vector<T> geom;
310 geom.reserve((range_p[1] - range_p[0]) * 3);
311 const std::int64_t sqxy = (nx + 1) * (ny + 1);
312 for (std::int64_t v = range_p[0]; v < range_p[1]; ++v)
313 {
314 const std::int64_t iz = v / sqxy;
315 const std::int64_t p = v % sqxy;
316 const std::int64_t iy = p / (nx + 1);
317 const std::int64_t ix = p % (nx + 1);
318 const T z = e + ef * static_cast<T>(iz);
319 const T y = c + cd * static_cast<T>(iy);
320 const T x = a + ab * static_cast<T>(ix);
321 geom.insert(geom.end(), {x, y, z});
322 }
323
324 return geom;
325}
326
327template <std::floating_point T>
328Mesh<T> build_tet(MPI_Comm comm, MPI_Comm subcomm,
329 std::array<std::array<double, 3>, 2> p,
330 std::array<std::size_t, 3> n,
331 const CellPartitionFunction& partitioner)
332{
333 common::Timer timer("Build BoxMesh (tetrahedra)");
334 std::vector<T> x;
335 std::vector<std::int64_t> cells;
336 if (subcomm != MPI_COMM_NULL)
337 {
338 x = create_geom<T>(subcomm, p, n);
339
340 const std::int64_t nx = n[0];
341 const std::int64_t ny = n[1];
342 const std::int64_t nz = n[2];
343 const std::int64_t n_cells = nx * ny * nz;
344
345 std::array range_c = dolfinx::MPI::local_range(
346 dolfinx::MPI::rank(subcomm), n_cells, dolfinx::MPI::size(subcomm));
347 cells.reserve(6 * (range_c[1] - range_c[0]) * 4);
348
349 // Create tetrahedra
350 for (std::int64_t i = range_c[0]; i < range_c[1]; ++i)
351 {
352 const std::size_t iz = i / (nx * ny);
353 const std::size_t j = i % (nx * ny);
354 const std::size_t iy = j / nx;
355 const std::size_t ix = j % nx;
356 const std::int64_t v0 = iz * (nx + 1) * (ny + 1) + iy * (nx + 1) + ix;
357 const std::int64_t v1 = v0 + 1;
358 const std::int64_t v2 = v0 + (nx + 1);
359 const std::int64_t v3 = v1 + (nx + 1);
360 const std::int64_t v4 = v0 + (nx + 1) * (ny + 1);
361 const std::int64_t v5 = v1 + (nx + 1) * (ny + 1);
362 const std::int64_t v6 = v2 + (nx + 1) * (ny + 1);
363 const std::int64_t v7 = v3 + (nx + 1) * (ny + 1);
364
365 // Note that v0 < v1 < v2 < v3 < vmid
366 cells.insert(cells.end(),
367 {v0, v1, v3, v7, v0, v1, v7, v5, v0, v5, v7, v4,
368 v0, v3, v2, v7, v0, v6, v4, v7, v0, v2, v6, v7});
369 }
370 }
371
372 fem::CoordinateElement<T> element(CellType::tetrahedron, 1);
373 return create_mesh(comm, subcomm, cells, element, subcomm, x,
374 {x.size() / 3, 3}, partitioner);
375}
376
377template <std::floating_point T>
378mesh::Mesh<T> build_hex(MPI_Comm comm, MPI_Comm subcomm,
379 std::array<std::array<double, 3>, 2> p,
380 std::array<std::size_t, 3> n,
381 const CellPartitionFunction& partitioner)
382{
383 common::Timer timer("Build BoxMesh (hexahedra)");
384 std::vector<T> x;
385 std::vector<std::int64_t> cells;
386 if (subcomm != MPI_COMM_NULL)
387 {
388 x = create_geom<T>(subcomm, p, n);
389
390 // Create cuboids
391 const std::int64_t nx = n[0];
392 const std::int64_t ny = n[1];
393 const std::int64_t nz = n[2];
394 const std::int64_t n_cells = nx * ny * nz;
395 std::array range_c = dolfinx::MPI::local_range(
396 dolfinx::MPI::rank(subcomm), n_cells, dolfinx::MPI::size(subcomm));
397 cells.reserve((range_c[1] - range_c[0]) * 8);
398 for (std::int64_t i = range_c[0]; i < range_c[1]; ++i)
399 {
400 const std::int64_t iz = i / (nx * ny);
401 const std::int64_t j = i % (nx * ny);
402 const std::int64_t iy = j / nx;
403 const std::int64_t ix = j % nx;
404
405 const std::int64_t v0 = (iz * (ny + 1) + iy) * (nx + 1) + ix;
406 const std::int64_t v1 = v0 + 1;
407 const std::int64_t v2 = v0 + (nx + 1);
408 const std::int64_t v3 = v1 + (nx + 1);
409 const std::int64_t v4 = v0 + (nx + 1) * (ny + 1);
410 const std::int64_t v5 = v1 + (nx + 1) * (ny + 1);
411 const std::int64_t v6 = v2 + (nx + 1) * (ny + 1);
412 const std::int64_t v7 = v3 + (nx + 1) * (ny + 1);
413 cells.insert(cells.end(), {v0, v1, v2, v3, v4, v5, v6, v7});
414 }
415 }
416
417 fem::CoordinateElement<T> element(CellType::hexahedron, 1);
418 return create_mesh(comm, subcomm, cells, element, subcomm, x,
419 {x.size() / 3, 3}, partitioner);
420}
421
422template <std::floating_point T>
423Mesh<T> build_prism(MPI_Comm comm, MPI_Comm subcomm,
424 std::array<std::array<double, 3>, 2> p,
425 std::array<std::size_t, 3> n,
426 const CellPartitionFunction& partitioner)
427{
428 std::vector<T> x;
429 std::vector<std::int64_t> cells;
430 if (subcomm != MPI_COMM_NULL)
431 {
432 x = create_geom<T>(subcomm, p, n);
433
434 const std::int64_t nx = n[0];
435 const std::int64_t ny = n[1];
436 const std::int64_t nz = n[2];
437 const std::int64_t n_cells = nx * ny * nz;
438 std::array range_c = dolfinx::MPI::local_range(
439 dolfinx::MPI::rank(comm), n_cells, dolfinx::MPI::size(comm));
440 const std::size_t cell_range = range_c[1] - range_c[0];
441
442 // Create cuboids
443
444 cells.reserve(2 * cell_range * 6);
445 for (std::int64_t i = range_c[0]; i < range_c[1]; ++i)
446 {
447 const std::int64_t iz = i / (nx * ny);
448 const std::int64_t j = i % (nx * ny);
449 const std::int64_t iy = j / nx;
450 const std::int64_t ix = j % nx;
451
452 const std::int64_t v0 = (iz * (ny + 1) + iy) * (nx + 1) + ix;
453 const std::int64_t v1 = v0 + 1;
454 const std::int64_t v2 = v0 + (nx + 1);
455 const std::int64_t v3 = v1 + (nx + 1);
456 const std::int64_t v4 = v0 + (nx + 1) * (ny + 1);
457 const std::int64_t v5 = v1 + (nx + 1) * (ny + 1);
458 const std::int64_t v6 = v2 + (nx + 1) * (ny + 1);
459 const std::int64_t v7 = v3 + (nx + 1) * (ny + 1);
460 cells.insert(cells.end(), {v0, v1, v2, v4, v5, v6});
461 cells.insert(cells.end(), {v1, v2, v3, v5, v6, v7});
462 }
463 }
464
465 fem::CoordinateElement<T> element(CellType::prism, 1);
466 return create_mesh(comm, subcomm, cells, element, subcomm, x,
467 {x.size() / 3, 3}, partitioner);
468}
469
470template <std::floating_point T>
471Mesh<T> build_tri(MPI_Comm comm, std::array<std::array<double, 2>, 2> p,
472 std::array<std::size_t, 2> n,
473 const CellPartitionFunction& partitioner,
474 DiagonalType diagonal)
475{
476 fem::CoordinateElement<T> element(CellType::triangle, 1);
477 std::vector<T> x;
478 std::vector<std::int64_t> cells;
479 if (dolfinx::MPI::rank(comm) == 0)
480 {
481 const std::array<double, 2> p0 = p[0];
482 const std::array<double, 2> p1 = p[1];
483
484 const std::size_t nx = n[0];
485 const std::size_t ny = n[1];
486
487 // Extract minimum and maximum coordinates
488 const T x0 = std::min(p0[0], p1[0]);
489 const T x1 = std::max(p0[0], p1[0]);
490 const T y0 = std::min(p0[1], p1[1]);
491 const T y1 = std::max(p0[1], p1[1]);
492
493 const T a = x0;
494 const T b = x1;
495 const T ab = (b - a) / static_cast<T>(nx);
496 const T c = y0;
497 const T d = y1;
498 const T cd = (d - c) / static_cast<T>(ny);
499
500 if (std::abs(x0 - x1) < std::numeric_limits<double>::epsilon()
501 or std::abs(y0 - y1) < std::numeric_limits<double>::epsilon())
502 {
503 throw std::runtime_error("Rectangle seems to have zero width, height or "
504 "depth. Check dimensions");
505 }
506
507 if (nx < 1 or ny < 1)
508 {
509 throw std::runtime_error(
510 "Rectangle has non-positive number of vertices in some dimension: "
511 "number of vertices must be at least 1 in each dimension");
512 }
513
514 // Create vertices and cells
515 std::size_t nv, nc;
516 switch (diagonal)
517 {
518 case DiagonalType::crossed:
519 nv = (nx + 1) * (ny + 1) + nx * ny;
520 nc = 4 * nx * ny;
521 break;
522 default:
523 nv = (nx + 1) * (ny + 1);
524 nc = 2 * nx * ny;
525 }
526
527 x.reserve(nv * 2);
528 cells.reserve(nc * 3);
529
530 // Create main vertices
531 std::size_t vertex = 0;
532 for (std::size_t iy = 0; iy <= ny; iy++)
533 {
534 const T x1 = c + cd * static_cast<T>(iy);
535 for (std::size_t ix = 0; ix <= nx; ix++)
536 x.insert(x.end(), {a + ab * static_cast<T>(ix), x1});
537 }
538
539 // Create midpoint vertices if the mesh type is crossed
540 switch (diagonal)
541 {
542 case DiagonalType::crossed:
543 for (std::size_t iy = 0; iy < ny; iy++)
544 {
545 const T x1 = c + cd * (static_cast<T>(iy) + 0.5);
546 for (std::size_t ix = 0; ix < nx; ix++)
547 x.insert(x.end(), {a + ab * (static_cast<T>(ix) + 0.5), x1});
548 }
549 break;
550 default:
551 break;
552 }
553
554 // Create triangles
555 switch (diagonal)
556 {
557 case DiagonalType::crossed:
558 {
559 for (std::size_t iy = 0; iy < ny; iy++)
560 {
561 for (std::size_t ix = 0; ix < nx; ix++)
562 {
563 const std::size_t v0 = iy * (nx + 1) + ix;
564 const std::size_t v1 = v0 + 1;
565 const std::size_t v2 = v0 + (nx + 1);
566 const std::size_t v3 = v1 + (nx + 1);
567 const std::size_t vmid = (nx + 1) * (ny + 1) + iy * nx + ix;
568
569 // Note that v0 < v1 < v2 < v3 < vmid
570 cells.insert(cells.end(), {v0, v1, vmid, v0, v2, vmid, v1, v3, vmid,
571 v2, v3, vmid});
572 }
573 }
574 break;
575 }
576 default:
577 {
578 DiagonalType local_diagonal = diagonal;
579 for (std::size_t iy = 0; iy < ny; iy++)
580 {
581 // Set up alternating diagonal
582 switch (diagonal)
583 {
584 case DiagonalType::right_left:
585 if (iy % 2)
586 local_diagonal = DiagonalType::right;
587 else
588 local_diagonal = DiagonalType::left;
589 break;
590 case DiagonalType::left_right:
591 if (iy % 2)
592 local_diagonal = DiagonalType::left;
593 else
594 local_diagonal = DiagonalType::right;
595 break;
596 default:
597 break;
598 }
599 for (std::size_t ix = 0; ix < nx; ix++)
600 {
601 const std::size_t v0 = iy * (nx + 1) + ix;
602 const std::size_t v1 = v0 + 1;
603 const std::size_t v2 = v0 + (nx + 1);
604 const std::size_t v3 = v1 + (nx + 1);
605
606 std::size_t offset = iy * nx + ix;
607 switch (local_diagonal)
608 {
609 case DiagonalType::left:
610 {
611 cells.insert(cells.end(), {v0, v1, v2, v1, v2, v3});
612 if (diagonal == DiagonalType::right_left
613 or diagonal == DiagonalType::left_right)
614 {
615 local_diagonal = DiagonalType::right;
616 }
617 break;
618 }
619 default:
620 {
621 cells.insert(cells.end(), {v0, v1, v3, v0, v2, v3});
622 if (diagonal == DiagonalType::right_left
623 or diagonal == DiagonalType::left_right)
624 {
625 local_diagonal = DiagonalType::left;
626 }
627 }
628 }
629 }
630 }
631 }
632 }
633
634 return create_mesh(comm, MPI_COMM_SELF, cells, element, MPI_COMM_SELF, x,
635 {x.size() / 2, 2}, partitioner);
636 }
637 else
638 {
639 return create_mesh(comm, MPI_COMM_NULL, cells, element, MPI_COMM_NULL, x,
640 {x.size() / 2, 2}, partitioner);
641 }
642}
643
644template <std::floating_point T>
645Mesh<T> build_quad(MPI_Comm comm, const std::array<std::array<double, 2>, 2> p,
646 std::array<std::size_t, 2> n,
647 const CellPartitionFunction& partitioner)
648{
649 fem::CoordinateElement<T> element(CellType::quadrilateral, 1);
650 std::vector<std::int64_t> cells;
651 std::vector<T> x;
652 if (dolfinx::MPI::rank(comm) == 0)
653 {
654 const std::size_t nx = n[0];
655 const std::size_t ny = n[1];
656 const T a = p[0][0];
657 const T b = p[1][0];
658 const T ab = (b - a) / static_cast<T>(nx);
659 const T c = p[0][1];
660 const T d = p[1][1];
661 const T cd = (d - c) / static_cast<T>(ny);
662
663 // Create vertices
664 x.reserve((nx + 1) * (ny + 1) * 2);
665 std::size_t vertex = 0;
666 for (std::size_t ix = 0; ix <= nx; ix++)
667 {
668 T x0 = a + ab * static_cast<T>(ix);
669 for (std::size_t iy = 0; iy <= ny; iy++)
670 x.insert(x.end(), {x0, c + cd * static_cast<T>(iy)});
671 }
672
673 // Create rectangles
674 cells.reserve(nx * ny * 4);
675 for (std::size_t ix = 0; ix < nx; ix++)
676 {
677 for (std::size_t iy = 0; iy < ny; iy++)
678 {
679 std::size_t i0 = ix * (ny + 1);
680 cells.insert(cells.end(), {i0 + iy, i0 + iy + 1, i0 + iy + ny + 1,
681 i0 + iy + ny + 2});
682 }
683 }
684
685 return create_mesh(comm, MPI_COMM_SELF, cells, element, MPI_COMM_SELF, x,
686 {x.size() / 2, 2}, partitioner);
687 }
688 else
689 {
690 return create_mesh(comm, MPI_COMM_NULL, cells, element, MPI_COMM_NULL, x,
691 {x.size() / 2, 2}, partitioner);
692 }
693}
694} // namespace impl
695} // namespace dolfinx::mesh
Definition CoordinateElement.h:38
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:15
Mesh data structures and algorithms on meshes.
Definition DofMap.h:32
Mesh< T > create_rectangle(MPI_Comm comm, std::array< std::array< double, 2 >, 2 > p, std::array< std::size_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:153
Mesh< T > create_box(MPI_Comm comm, MPI_Comm subcomm, std::array< std::array< double, 3 >, 2 > p, std::array< std::size_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:91
DiagonalType
Enum for different diagonal types.
Definition generation.h:24
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:785
std::function< graph::AdjacencyList< std::int32_t >( MPI_Comm comm, int nparts, CellType cell_type, const graph::AdjacencyList< std::int64_t > &cells)> CellPartitionFunction
Signature for the cell partitioning function. The function should compute the destination rank for ce...
Definition utils.h:178
Mesh< T > create_interval(MPI_Comm comm, std::size_t nx, std::array< double, 2 > p, CellPartitionFunction partitioner=nullptr)
Interval mesh of the 1D line [a, b].
Definition generation.h:207
CellType
Cell type identifier.
Definition cell_types.h:22
CellPartitionFunction create_cell_partitioner(mesh::GhostMode ghost_mode=mesh::GhostMode::none, const graph::partition_fn &partfn=&graph::partition_graph)
Definition utils.cpp:87