xref: /petsc/src/vec/is/sf/utils/sfcoord.c (revision 8b5d2d36b1bd7331337e6600e2fff187f080efc8)
1*8b5d2d36SJed Brown #include <petsc/private/sfimpl.h> /*I  "petscsf.h"   I*/
2*8b5d2d36SJed Brown 
3*8b5d2d36SJed Brown static PetscErrorCode GetBoundingBox_Internal(PetscInt npoints, PetscInt dim, const PetscReal *coords, PetscReal *bbox)
4*8b5d2d36SJed Brown {
5*8b5d2d36SJed Brown   PetscFunctionBegin;
6*8b5d2d36SJed Brown   for (PetscInt d = 0; d < dim; d++) {
7*8b5d2d36SJed Brown     bbox[0 * dim + d] = PETSC_MAX_REAL;
8*8b5d2d36SJed Brown     bbox[1 * dim + d] = PETSC_MIN_REAL;
9*8b5d2d36SJed Brown   }
10*8b5d2d36SJed Brown   for (PetscInt i = 0; i < npoints; i++) {
11*8b5d2d36SJed Brown     for (PetscInt d = 0; d < dim; d++) {
12*8b5d2d36SJed Brown       bbox[0 * dim + d] = PetscMin(bbox[0 * dim + d], coords[i * dim + d]);
13*8b5d2d36SJed Brown       bbox[1 * dim + d] = PetscMax(bbox[1 * dim + d], coords[i * dim + d]);
14*8b5d2d36SJed Brown     }
15*8b5d2d36SJed Brown   }
16*8b5d2d36SJed Brown   PetscFunctionReturn(0);
17*8b5d2d36SJed Brown }
18*8b5d2d36SJed Brown 
19*8b5d2d36SJed Brown static PetscBool IntersectBoundingBox_Internal(PetscInt dim, const PetscReal *a, const PetscReal *b, PetscReal tol)
20*8b5d2d36SJed Brown {
21*8b5d2d36SJed Brown   for (PetscInt d = 0; d < dim; d++) {
22*8b5d2d36SJed Brown     if (a[1 * dim + d] + tol < b[0 * dim + d] || b[1 * dim + d] + tol < a[0 * dim + d]) return PETSC_FALSE;
23*8b5d2d36SJed Brown   }
24*8b5d2d36SJed Brown   return PETSC_TRUE;
25*8b5d2d36SJed Brown }
26*8b5d2d36SJed Brown 
27*8b5d2d36SJed Brown static PetscBool InBoundingBox_Internal(PetscInt dim, const PetscReal *x, const PetscReal *bbox, PetscReal tol)
28*8b5d2d36SJed Brown {
29*8b5d2d36SJed Brown   for (PetscInt d = 0; d < dim; d++) {
30*8b5d2d36SJed Brown     if (x[d] + tol < bbox[0 * dim + d] || bbox[1 * dim + d] + tol < x[d]) return PETSC_FALSE;
31*8b5d2d36SJed Brown   }
32*8b5d2d36SJed Brown   return PETSC_TRUE;
33*8b5d2d36SJed Brown }
34*8b5d2d36SJed Brown 
35*8b5d2d36SJed Brown /*@
36*8b5d2d36SJed Brown   PetscSFSetGraphFromCoordinates - Create SF by fuzzy matching leaf coordinates to root coordinates
37*8b5d2d36SJed Brown 
38*8b5d2d36SJed Brown   Collective
39*8b5d2d36SJed Brown 
40*8b5d2d36SJed Brown   Input Parameters:
41*8b5d2d36SJed Brown + sf - PetscSF to set graph on
42*8b5d2d36SJed Brown . nroots - number of root coordinates
43*8b5d2d36SJed Brown . nleaves - number of leaf coordinates
44*8b5d2d36SJed Brown . dim - spatial dimension of coordinates
45*8b5d2d36SJed Brown . tol - positive tolerance for matching
46*8b5d2d36SJed Brown . rootcoords - array of root coordinates in which root i component d is [i*dim+d]
47*8b5d2d36SJed Brown - leafcoords - array of root coordinates in which leaf i component d is [i*dim+d]
48*8b5d2d36SJed Brown 
49*8b5d2d36SJed Brown   Notes:
50*8b5d2d36SJed Brown   The tolerance typically represents the rounding error incurred by numerically computing coordinates via
51*8b5d2d36SJed Brown   possibly-different procedures. Passing anything from `PETSC_SMALL` to `100 * PETSC_MACHINE_EPSILON` is appropriate for
52*8b5d2d36SJed Brown   most use cases.
53*8b5d2d36SJed Brown 
54*8b5d2d36SJed Brown   Example:
55*8b5d2d36SJed Brown   As a motivating example, consider fluid flow in the x direction with y (distance from a wall). The spanwise direction,
56*8b5d2d36SJed Brown   z, has periodic boundary conditions and needs some spanwise length to allow turbulent structures to develop. The
57*8b5d2d36SJed Brown   distribution is stationary with respect to z, so you want to average turbulence variables (like Reynolds stress) over
58*8b5d2d36SJed Brown   the z direction. It is complicated in a 3D simulation with arbitrary partitioner to uniquely number the nodes or
59*8b5d2d36SJed Brown   quadrature point coordinates to average these quantities into a 2D plane where they will be visualized, but it's easy
60*8b5d2d36SJed Brown   to compute the projection of each 3D point into the 2D plane.
61*8b5d2d36SJed Brown 
62*8b5d2d36SJed Brown   Suppose a 2D target mesh and 3D source mesh (logically an extrusion of the 2D, though perhaps not created in that way)
63*8b5d2d36SJed Brown   are distributed independently on a communicator. Each rank passes its 2D target points as root coordinates and the 2D
64*8b5d2d36SJed Brown   projection of its 3D source points as leaf coordinates. Calling `PetscSFReduceBegin()`/`PetscSFReduceEnd()` on the
65*8b5d2d36SJed Brown   result will sum data from the 3D sources to the 2D targets.
66*8b5d2d36SJed Brown 
67*8b5d2d36SJed Brown   As a concrete example, consider three MPI ranks with targets (roots)
68*8b5d2d36SJed Brown .vb
69*8b5d2d36SJed Brown Rank 0: (0, 0), (0, 1)
70*8b5d2d36SJed Brown Rank 1: (0.1, 0), (0.1, 1)
71*8b5d2d36SJed Brown Rank 2: (0.2, 0), (0.2, 1)
72*8b5d2d36SJed Brown .ve
73*8b5d2d36SJed Brown   Note that targets must be uniquely owned. Suppose also that we identify the following leaf coordinates (perhaps via projection from a 3D space).
74*8b5d2d36SJed Brown .vb
75*8b5d2d36SJed Brown Rank 0: (0, 0), (0.1, 0), (0, 1), (0.1, 1)
76*8b5d2d36SJed Brown Rank 1: (0, 0), (0.1, 0), (0.2, 0), (0, 1), (0.1, 1)
77*8b5d2d36SJed Brown Rank 2: (0.1, 0), (0.2, 0), (0.1, 1), (0.2, 1)
78*8b5d2d36SJed Brown .ve
79*8b5d2d36SJed Brown   Leaf coordinates may be repeated, both on a rank and between ranks. This example yields the following `PetscSF` capable of reducing from sources to targets.
80*8b5d2d36SJed Brown .vb
81*8b5d2d36SJed Brown Roots by rank
82*8b5d2d36SJed Brown [0]  0:   0.0000e+00   0.0000e+00   0.0000e+00   1.0000e+00
83*8b5d2d36SJed Brown [1]  0:   1.0000e-01   0.0000e+00   1.0000e-01   1.0000e+00
84*8b5d2d36SJed Brown [2]  0:   2.0000e-01   0.0000e+00   2.0000e-01   1.0000e+00
85*8b5d2d36SJed Brown Leaves by rank
86*8b5d2d36SJed Brown [0]  0:   0.0000e+00   0.0000e+00   1.0000e-01   0.0000e+00   0.0000e+00
87*8b5d2d36SJed Brown [0]  5:   1.0000e+00   1.0000e-01   1.0000e+00
88*8b5d2d36SJed Brown [1]  0:   0.0000e+00   0.0000e+00   1.0000e-01   0.0000e+00   2.0000e-01
89*8b5d2d36SJed Brown [1]  5:   0.0000e+00   0.0000e+00   1.0000e+00   1.0000e-01   1.0000e+00
90*8b5d2d36SJed Brown [1] 10:   2.0000e-01   1.0000e+00
91*8b5d2d36SJed Brown [2]  0:   1.0000e-01   0.0000e+00   2.0000e-01   0.0000e+00   1.0000e-01
92*8b5d2d36SJed Brown [2]  5:   1.0000e+00   2.0000e-01   1.0000e+00
93*8b5d2d36SJed Brown PetscSF Object: 3 MPI processes
94*8b5d2d36SJed Brown   type: basic
95*8b5d2d36SJed Brown   [0] Number of roots=2, leaves=4, remote ranks=2
96*8b5d2d36SJed Brown   [0] 0 <- (0,0)
97*8b5d2d36SJed Brown   [0] 1 <- (1,0)
98*8b5d2d36SJed Brown   [0] 2 <- (0,1)
99*8b5d2d36SJed Brown   [0] 3 <- (1,1)
100*8b5d2d36SJed Brown   [1] Number of roots=2, leaves=6, remote ranks=3
101*8b5d2d36SJed Brown   [1] 0 <- (0,0)
102*8b5d2d36SJed Brown   [1] 1 <- (1,0)
103*8b5d2d36SJed Brown   [1] 2 <- (2,0)
104*8b5d2d36SJed Brown   [1] 3 <- (0,1)
105*8b5d2d36SJed Brown   [1] 4 <- (1,1)
106*8b5d2d36SJed Brown   [1] 5 <- (2,1)
107*8b5d2d36SJed Brown   [2] Number of roots=2, leaves=4, remote ranks=2
108*8b5d2d36SJed Brown   [2] 0 <- (1,0)
109*8b5d2d36SJed Brown   [2] 1 <- (2,0)
110*8b5d2d36SJed Brown   [2] 2 <- (1,1)
111*8b5d2d36SJed Brown   [2] 3 <- (2,1)
112*8b5d2d36SJed Brown .ve
113*8b5d2d36SJed Brown 
114*8b5d2d36SJed Brown   Level: advanced
115*8b5d2d36SJed Brown 
116*8b5d2d36SJed Brown .seealso: `PetscSFCreate()`, `PetscSFSetGraph()`, `PetscSFCreateByMatchingIndices()`
117*8b5d2d36SJed Brown @*/
118*8b5d2d36SJed Brown PetscErrorCode PetscSFSetGraphFromCoordinates(PetscSF sf, PetscInt nroots, PetscInt nleaves, PetscInt dim, PetscReal tol, const PetscReal *rootcoords, const PetscReal *leafcoords)
119*8b5d2d36SJed Brown {
120*8b5d2d36SJed Brown   PetscReal    bbox[6], *bboxes, *target_coords;
121*8b5d2d36SJed Brown   PetscMPIInt  size, *ranks_needed, num_ranks;
122*8b5d2d36SJed Brown   PetscInt    *root_sizes, *root_starts;
123*8b5d2d36SJed Brown   PetscSFNode *premote, *lremote;
124*8b5d2d36SJed Brown   PetscSF      psf;
125*8b5d2d36SJed Brown   MPI_Datatype unit;
126*8b5d2d36SJed Brown   MPI_Comm     comm;
127*8b5d2d36SJed Brown 
128*8b5d2d36SJed Brown   PetscFunctionBegin;
129*8b5d2d36SJed Brown   PetscCall(PetscObjectGetComm((PetscObject)sf, &comm));
130*8b5d2d36SJed Brown   PetscCall(GetBoundingBox_Internal(nroots, dim, rootcoords, bbox));
131*8b5d2d36SJed Brown   PetscCallMPI(MPI_Comm_size(comm, &size));
132*8b5d2d36SJed Brown   PetscCall(PetscMalloc1(size * 2 * dim, &bboxes));
133*8b5d2d36SJed Brown   PetscCallMPI(MPI_Allgather(bbox, 2 * dim, MPIU_REAL, bboxes, 2 * dim, MPIU_REAL, comm));
134*8b5d2d36SJed Brown   PetscCall(GetBoundingBox_Internal(nleaves, dim, leafcoords, bbox));
135*8b5d2d36SJed Brown   PetscCall(PetscMalloc1(size, &root_sizes));
136*8b5d2d36SJed Brown   PetscCallMPI(MPI_Allgather(&nroots, 1, MPIU_INT, root_sizes, 1, MPIU_INT, comm));
137*8b5d2d36SJed Brown 
138*8b5d2d36SJed Brown   PetscCall(PetscMalloc2(size, &ranks_needed, size + 1, &root_starts));
139*8b5d2d36SJed Brown   root_starts[0] = 0;
140*8b5d2d36SJed Brown   num_ranks      = 0;
141*8b5d2d36SJed Brown   for (PetscMPIInt r = 0; r < size; r++) {
142*8b5d2d36SJed Brown     if (IntersectBoundingBox_Internal(dim, bbox, &bboxes[2 * dim * r], tol)) {
143*8b5d2d36SJed Brown       ranks_needed[num_ranks++] = r;
144*8b5d2d36SJed Brown       root_starts[num_ranks]    = root_starts[num_ranks - 1] + root_sizes[r];
145*8b5d2d36SJed Brown     }
146*8b5d2d36SJed Brown   }
147*8b5d2d36SJed Brown   PetscCall(PetscFree(root_sizes));
148*8b5d2d36SJed Brown   PetscCall(PetscMalloc1(root_starts[num_ranks], &premote));
149*8b5d2d36SJed Brown   for (PetscInt i = 0; i < num_ranks; i++) {
150*8b5d2d36SJed Brown     for (PetscInt j = root_starts[i]; j < root_starts[i + 1]; j++) {
151*8b5d2d36SJed Brown       premote[j].rank  = ranks_needed[i];
152*8b5d2d36SJed Brown       premote[j].index = j - root_starts[i];
153*8b5d2d36SJed Brown     }
154*8b5d2d36SJed Brown   }
155*8b5d2d36SJed Brown   PetscCall(PetscSFCreate(comm, &psf));
156*8b5d2d36SJed Brown   PetscCall(PetscSFSetGraph(psf, nroots, root_starts[num_ranks], NULL, PETSC_USE_POINTER, premote, PETSC_USE_POINTER));
157*8b5d2d36SJed Brown   PetscCall(PetscMalloc1(root_starts[num_ranks] * dim, &target_coords));
158*8b5d2d36SJed Brown   PetscCallMPI(MPI_Type_contiguous(dim, MPIU_REAL, &unit));
159*8b5d2d36SJed Brown   PetscCallMPI(MPI_Type_commit(&unit));
160*8b5d2d36SJed Brown   PetscCall(PetscSFBcastBegin(psf, unit, rootcoords, target_coords, MPI_REPLACE));
161*8b5d2d36SJed Brown   PetscCall(PetscSFBcastEnd(psf, unit, rootcoords, target_coords, MPI_REPLACE));
162*8b5d2d36SJed Brown   PetscCallMPI(MPI_Type_free(&unit));
163*8b5d2d36SJed Brown   PetscCall(PetscSFDestroy(&psf));
164*8b5d2d36SJed Brown 
165*8b5d2d36SJed Brown   // Condense targets to only those that lie within our bounding box
166*8b5d2d36SJed Brown   PetscInt num_targets = 0;
167*8b5d2d36SJed Brown   for (PetscInt i = 0; i < root_starts[num_ranks]; i++) {
168*8b5d2d36SJed Brown     if (InBoundingBox_Internal(dim, &target_coords[i * dim], bbox, tol)) {
169*8b5d2d36SJed Brown       premote[num_targets] = premote[i];
170*8b5d2d36SJed Brown       for (PetscInt d = 0; d < dim; d++) target_coords[num_targets * dim + d] = target_coords[i * dim + d];
171*8b5d2d36SJed Brown       num_targets++;
172*8b5d2d36SJed Brown     }
173*8b5d2d36SJed Brown   }
174*8b5d2d36SJed Brown   PetscCall(PetscFree(bboxes));
175*8b5d2d36SJed Brown   PetscCall(PetscFree2(ranks_needed, root_starts));
176*8b5d2d36SJed Brown 
177*8b5d2d36SJed Brown   PetscCall(PetscMalloc1(nleaves, &lremote));
178*8b5d2d36SJed Brown   for (PetscInt i = 0; i < nleaves; i++) {
179*8b5d2d36SJed Brown     for (PetscInt j = 0; j < num_targets; j++) {
180*8b5d2d36SJed Brown       PetscReal sum = 0;
181*8b5d2d36SJed Brown       for (PetscInt d = 0; d < dim; d++) sum += PetscSqr(leafcoords[i * dim + d] - target_coords[j * dim + d]);
182*8b5d2d36SJed Brown       if (sum < tol * tol) {
183*8b5d2d36SJed Brown         lremote[i] = premote[j];
184*8b5d2d36SJed Brown         goto matched;
185*8b5d2d36SJed Brown       }
186*8b5d2d36SJed Brown     }
187*8b5d2d36SJed Brown     switch (dim) {
188*8b5d2d36SJed Brown     case 1:
189*8b5d2d36SJed Brown       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "No target found for leaf coordinate %g", (double)leafcoords[i * dim + 0]);
190*8b5d2d36SJed Brown     case 2:
191*8b5d2d36SJed Brown       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "No target found for leaf coordinate (%g, %g)", (double)leafcoords[i * dim + 0], (double)leafcoords[i * dim + 1]);
192*8b5d2d36SJed Brown     case 3:
193*8b5d2d36SJed Brown       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "No target found for leaf coordinate (%g, %g, %g)", (double)leafcoords[i * dim + 0], (double)leafcoords[i * dim + 1], (double)leafcoords[i * dim + 2]);
194*8b5d2d36SJed Brown     }
195*8b5d2d36SJed Brown   matched:;
196*8b5d2d36SJed Brown   }
197*8b5d2d36SJed Brown   PetscCall(PetscFree(premote));
198*8b5d2d36SJed Brown   PetscCall(PetscFree(target_coords));
199*8b5d2d36SJed Brown   PetscCall(PetscSFSetGraph(sf, nroots, nleaves, NULL, PETSC_USE_POINTER, lremote, PETSC_OWN_POINTER));
200*8b5d2d36SJed Brown   PetscFunctionReturn(0);
201*8b5d2d36SJed Brown }
202