PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
Macros | Functions
walkingsearch.c File Reference
#include "walkingsearch.h"
#include "logging.h"
#include <petsc.h>
#include <stdbool.h>
#include <math.h>
Include dependency graph for walkingsearch.c:

Go to the source code of this file.

Macros

#define MAX_TRAVERSAL   1000
 
#define DISTANCE_THRESHOLD   1e-11
 
#define REPEAT_COUNT_THRESHOLD   5
 
#define __FUNCT__   "GetCellCharacteristicSize"
 
#define __FUNCT__   "ComputeSignedDistanceToPlane"
 
#define __FUNCT__   "CalculateDistancesToCellFaces"
 
#define __FUNCT__   "DeterminePointPosition"
 
#define __FUNCT__   "GetCellVerticesFromGrid"
 
#define __FUNCT__   "InitializeTraversalParameters"
 
#define __FUNCT__   "CheckCellWithinLocalGrid"
 
#define __FUNCT__   "RetrieveCurrentCell"
 
#define __FUNCT__   "EvaluateParticlePosition"
 
#define __FUNCT__   "UpdateCellIndicesBasedOnDistances"
 
#define __FUNCT__   "FinalizeTraversal"
 
#define __FUNCT__   "FindOwnerOfCell"
 
#define __FUNCT__   "LocateParticleOrFindMigrationTarget"
 
#define __FUNCT__   "ReportSearchOutcome"
 

Functions

PetscErrorCode GetCellCharacteristicSize (const Cell *cell, PetscReal *cellSize)
 Implementation of GetCellCharacteristicSize().
 
PetscErrorCode ComputeSignedDistanceToPlane (const Cmpnts v1, const Cmpnts v2, const Cmpnts v3, const Cmpnts v4, const Cmpnts cell_centroid, const Cmpnts p_target, PetscReal *d_signed, const PetscReal threshold)
 Internal helper implementation: ComputeSignedDistanceToPlane().
 
PetscErrorCode CalculateDistancesToCellFaces (const Cmpnts p, const Cell *cell, PetscReal *d, const PetscReal threshold)
 Implementation of CalculateDistancesToCellFaces().
 
PetscErrorCode DeterminePointPosition (PetscReal *d, PetscInt *result)
 Implementation of DeterminePointPosition().
 
PetscErrorCode GetCellVerticesFromGrid (Cmpnts ***coor, PetscInt idx, PetscInt idy, PetscInt idz, Cell *cell)
 Implementation of GetCellVerticesFromGrid().
 
PetscErrorCode InitializeTraversalParameters (UserCtx *user, Particle *particle, PetscInt *idx, PetscInt *idy, PetscInt *idz, PetscInt *traversal_steps)
 Implementation of InitializeTraversalParameters().
 
PetscErrorCode CheckCellWithinLocalGrid (UserCtx *user, PetscInt idx, PetscInt idy, PetscInt idz, PetscBool *is_within)
 Implementation of CheckCellWithinLocalGrid().
 
PetscErrorCode RetrieveCurrentCell (UserCtx *user, PetscInt idx, PetscInt idy, PetscInt idz, Cell *cell)
 Implementation of RetrieveCurrentCell().
 
PetscErrorCode EvaluateParticlePosition (const Cell *cell, PetscReal *d, const Cmpnts p, PetscInt *position, const PetscReal threshold)
 Implementation of EvaluateParticlePosition().
 
PetscErrorCode UpdateCellIndicesBasedOnDistances (PetscReal d[NUM_FACES], PetscInt *idx, PetscInt *idy, PetscInt *idz)
 Internal helper implementation: UpdateCellIndicesBasedOnDistances().
 
PetscErrorCode FinalizeTraversal (UserCtx *user, Particle *particle, PetscInt traversal_steps, PetscBool cell_found, PetscInt idx, PetscInt idy, PetscInt idz)
 Implementation of FinalizeTraversal().
 
PetscErrorCode FindOwnerOfCell (UserCtx *user, PetscInt i, PetscInt j, PetscInt k, PetscMPIInt *owner_rank)
 Internal helper implementation: FindOwnerOfCell().
 
PetscErrorCode LocateParticleOrFindMigrationTarget (UserCtx *user, Particle *particle, ParticleLocationStatus *status_out)
 Implementation of LocateParticleOrFindMigrationTarget().
 
PetscErrorCode ReportSearchOutcome (const Particle *particle, ParticleLocationStatus status, PetscInt traversal_steps)
 Implementation of ReportSearchOutcome().
 

Macro Definition Documentation

◆ MAX_TRAVERSAL

#define MAX_TRAVERSAL   1000

Definition at line 10 of file walkingsearch.c.

◆ DISTANCE_THRESHOLD

#define DISTANCE_THRESHOLD   1e-11

Definition at line 11 of file walkingsearch.c.

◆ REPEAT_COUNT_THRESHOLD

#define REPEAT_COUNT_THRESHOLD   5

Definition at line 12 of file walkingsearch.c.

◆ __FUNCT__ [1/14]

#define __FUNCT__   "GetCellCharacteristicSize"

Definition at line 15 of file walkingsearch.c.

◆ __FUNCT__ [2/14]

#define __FUNCT__   "ComputeSignedDistanceToPlane"

Definition at line 15 of file walkingsearch.c.

◆ __FUNCT__ [3/14]

#define __FUNCT__   "CalculateDistancesToCellFaces"

Definition at line 15 of file walkingsearch.c.

◆ __FUNCT__ [4/14]

#define __FUNCT__   "DeterminePointPosition"

Definition at line 15 of file walkingsearch.c.

◆ __FUNCT__ [5/14]

#define __FUNCT__   "GetCellVerticesFromGrid"

Definition at line 15 of file walkingsearch.c.

◆ __FUNCT__ [6/14]

#define __FUNCT__   "InitializeTraversalParameters"

Definition at line 15 of file walkingsearch.c.

◆ __FUNCT__ [7/14]

#define __FUNCT__   "CheckCellWithinLocalGrid"

Definition at line 15 of file walkingsearch.c.

◆ __FUNCT__ [8/14]

#define __FUNCT__   "RetrieveCurrentCell"

Definition at line 15 of file walkingsearch.c.

◆ __FUNCT__ [9/14]

#define __FUNCT__   "EvaluateParticlePosition"

Definition at line 15 of file walkingsearch.c.

◆ __FUNCT__ [10/14]

#define __FUNCT__   "UpdateCellIndicesBasedOnDistances"

Definition at line 15 of file walkingsearch.c.

◆ __FUNCT__ [11/14]

#define __FUNCT__   "FinalizeTraversal"

Definition at line 15 of file walkingsearch.c.

◆ __FUNCT__ [12/14]

#define __FUNCT__   "FindOwnerOfCell"

Definition at line 15 of file walkingsearch.c.

◆ __FUNCT__ [13/14]

#define __FUNCT__   "LocateParticleOrFindMigrationTarget"

Definition at line 15 of file walkingsearch.c.

◆ __FUNCT__ [14/14]

#define __FUNCT__   "ReportSearchOutcome"

Definition at line 15 of file walkingsearch.c.

Function Documentation

◆ GetCellCharacteristicSize()

PetscErrorCode GetCellCharacteristicSize ( const Cell cell,
PetscReal *  cellSize 
)

Implementation of GetCellCharacteristicSize().

Estimates a characteristic length of the cell for threshold scaling.

Full API contract (arguments, ownership, side effects) is documented with the header declaration in include/walkingsearch.h.

See also
GetCellCharacteristicSize()

Definition at line 22 of file walkingsearch.c.

23{
24 PetscErrorCode ierr = 0;
25 PetscFunctionBeginUser;
27 if (!cell || !cellSize) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL,
28 "GetCellCharacteristicSize: Null pointer(s).");
29
30 // A simple approach: compute the average of the distances between
31 // all pairs of adjacent vertices that share an edge. That gives
32 // a typical measure of cell dimension. For a uniform grid, you
33 // could do something simpler.
34
35 PetscInt edges[12][2] = {
36 {0,1},{1,2},{2,3},{3,0}, // bottom face edges
37 {4,5},{5,6},{6,7},{7,4}, // top face edges
38 {0,7},{1,6},{2,5},{3,4} // vertical edges
39 };
40
41 PetscReal totalEdgeLen = 0.0;
42 for (PetscInt i=0; i<12; i++) {
43 PetscInt vA = edges[i][0];
44 PetscInt vB = edges[i][1];
45 Cmpnts A = cell->vertices[vA];
46 Cmpnts B = cell->vertices[vB];
47 PetscReal dx = B.x - A.x;
48 PetscReal dy = B.y - A.y;
49 PetscReal dz = B.z - A.z;
50 PetscReal edgeLen = sqrt(dx*dx + dy*dy + dz*dz);
51 totalEdgeLen += edgeLen;
52 }
53
54 // Average edge length
55 *cellSize = totalEdgeLen / 12.0;
56
58 PetscFunctionReturn(ierr);
59}
#define PROFILE_FUNCTION_END
Marks the end of a profiled code block.
Definition logging.h:779
#define PROFILE_FUNCTION_BEGIN
Marks the beginning of a profiled code block (typically a function).
Definition logging.h:770
PetscScalar x
Definition variables.h:101
PetscScalar z
Definition variables.h:101
PetscScalar y
Definition variables.h:101
Cmpnts vertices[8]
Coordinates of the eight vertices of the cell.
Definition variables.h:161
A 3D point or vector with PetscScalar components.
Definition variables.h:100
Here is the caller graph for this function:

◆ ComputeSignedDistanceToPlane()

PetscErrorCode ComputeSignedDistanceToPlane ( const Cmpnts  v1,
const Cmpnts  v2,
const Cmpnts  v3,
const Cmpnts  v4,
const Cmpnts  cell_centroid,
const Cmpnts  p_target,
PetscReal *  d_signed,
const PetscReal  threshold 
)

Internal helper implementation: ComputeSignedDistanceToPlane().

Computes the signed distance from a point to the plane approximating a quadrilateral face.

Local to this translation unit.

Definition at line 68 of file walkingsearch.c.

71{
72 PetscErrorCode ierr = 0;
73 PetscMPIInt rank;
74
75 PetscFunctionBeginUser;
77
78 if (d_signed == NULL) {
79 ierr = MPI_Comm_rank(MPI_COMM_WORLD, &rank); CHKERRQ(ierr);
80 LOG_ALLOW(LOCAL, LOG_ERROR, "Output pointer 'd_signed' is NULL on rank %d.\n", rank);
81 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Output pointer 'd_signed' must not be NULL.");
82 }
83
84 LOG_ALLOW(LOCAL, LOG_VERBOSE, "Target point: (%.6e, %.6e, %.6e)\n", p_target.x, p_target.y, p_target.z);
85 LOG_ALLOW(LOCAL, LOG_VERBOSE, " Cell Centroid: (%.6e, %.6e, %.6e)\n", cell_centroid.x, cell_centroid.y, cell_centroid.z);
86 LOG_ALLOW(LOCAL, LOG_VERBOSE, " Face Vertices:\n");
87 LOG_ALLOW(LOCAL, LOG_VERBOSE, " v1: (%.6e, %.6e, %.6e)\n", v1.x, v1.y, v1.z);
88 LOG_ALLOW(LOCAL, LOG_VERBOSE, " v2: (%.6e, %.6e, %.6e)\n", v2.x, v2.y, v2.z);
89 LOG_ALLOW(LOCAL, LOG_VERBOSE, " v3: (%.6e, %.6e, %.6e)\n", v3.x, v3.y, v3.z);
90 LOG_ALLOW(LOCAL, LOG_VERBOSE, " v4: (%.6e, %.6e, %.6e)\n", v4.x, v4.y, v4.z);
91
92 // --- Calculate Edge Vectors for Initial Normal Computation ---
93 PetscReal edge1_x = v2.x - v1.x;
94 PetscReal edge1_y = v2.y - v1.y;
95 PetscReal edge1_z = v2.z - v1.z;
96
97 PetscReal edge2_x = v4.x - v1.x;
98 PetscReal edge2_y = v4.y - v1.y;
99 PetscReal edge2_z = v4.z - v1.z;
100 LOG_ALLOW(LOCAL, LOG_VERBOSE, " Edge1 (v2-v1): (%.6e, %.6e, %.6e)\n", edge1_x, edge1_y, edge1_z);
101 LOG_ALLOW(LOCAL, LOG_VERBOSE, " Edge2 (v4-v1): (%.6e, %.6e, %.6e)\n", edge2_x, edge2_y, edge2_z);
102
103 // --- Compute Initial Normal Vector (Cross Product: edge1 x edge2) ---
104 PetscReal normal_x_initial = edge1_y * edge2_z - edge1_z * edge2_y;
105 PetscReal normal_y_initial = edge1_z * edge2_x - edge1_x * edge2_z;
106 PetscReal normal_z_initial = edge1_x * edge2_y - edge1_y * edge2_x;
107 LOG_ALLOW(LOCAL, LOG_VERBOSE, " Initial Raw Normal (edge1 x edge2): (%.6e, %.6e, %.6e)\n", normal_x_initial, normal_y_initial, normal_z_initial);
108
109 PetscReal normal_magnitude = sqrt(normal_x_initial * normal_x_initial +
110 normal_y_initial * normal_y_initial +
111 normal_z_initial * normal_z_initial);
112 LOG_ALLOW(LOCAL, LOG_VERBOSE, " Initial Normal Magnitude: %.6e\n", normal_magnitude);
113
114 if (normal_magnitude < 1.0e-12) {
115 ierr = MPI_Comm_rank(MPI_COMM_WORLD, &rank); CHKERRQ(ierr);
117 "Degenerate plane detected on rank %d. Normal magnitude (%.3e) is too small.\n",
118 rank, normal_magnitude);
119 LOG_ALLOW(LOCAL, LOG_ERROR, " Offending vertices for normal: v1(%.3e,%.3e,%.3e), v2(%.3e,%.3e,%.3e), v4(%.3e,%.3e,%.3e)\n",
120 v1.x,v1.y,v1.z, v2.x,v2.y,v2.z, v4.x,v4.y,v4.z);
121 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Degenerate plane detected (normal vector is near zero).");
122 }
123
124 // --- Calculate the Centroid of the Four Face Vertices (v1, v2, v3, v4) ---
125 PetscReal face_centroid_x = 0.25 * (v1.x + v2.x + v3.x + v4.x);
126 PetscReal face_centroid_y = 0.25 * (v1.y + v2.y + v3.y + v4.y);
127 PetscReal face_centroid_z = 0.25 * (v1.z + v2.z + v3.z + v4.z);
128 LOG_ALLOW(LOCAL, LOG_VERBOSE, " Face Centroid: (%.6e, %.6e, %.6e)\n", face_centroid_x, face_centroid_y, face_centroid_z);
129
130 // --- Orient the Normal to Point Outward from the Cell ---
131 // Vector from face centroid to cell centroid
132 PetscReal vec_fc_to_cc_x = cell_centroid.x - face_centroid_x;
133 PetscReal vec_fc_to_cc_y = cell_centroid.y - face_centroid_y;
134 PetscReal vec_fc_to_cc_z = cell_centroid.z - face_centroid_z;
135 LOG_ALLOW(LOCAL, LOG_VERBOSE, " Vec (FaceCentroid -> CellCentroid): (%.6e, %.6e, %.6e)\n", vec_fc_to_cc_x, vec_fc_to_cc_y, vec_fc_to_cc_z);
136
137 // Dot product of initial normal with vector from face centroid to cell centroid
138 PetscReal dot_prod_orientation = normal_x_initial * vec_fc_to_cc_x +
139 normal_y_initial * vec_fc_to_cc_y +
140 normal_z_initial * vec_fc_to_cc_z;
141 LOG_ALLOW(LOCAL, LOG_VERBOSE, " Dot Product for Orientation (N_initial . Vec_FC_to_CC): %.6e\n", dot_prod_orientation);
142
143 PetscReal normal_x = normal_x_initial;
144 PetscReal normal_y = normal_y_initial;
145 PetscReal normal_z = normal_z_initial;
146
147 // If dot_prod_orientation > 0, initial normal points towards cell interior, so flip it.
148 if (dot_prod_orientation > 1.0e-9) { // Use a small epsilon to avoid issues if dot product is extremely close to zero
149 normal_x = -normal_x_initial;
150 normal_y = -normal_y_initial;
151 normal_z = -normal_z_initial;
152 LOG_ALLOW(LOCAL, LOG_VERBOSE, " Initial normal was inward (dot_prod > 0). Flipped normal.\n");
153 } else if (dot_prod_orientation == 0.0 && normal_magnitude > 1e-12) {
154 // This case is ambiguous or face plane contains cell centroid.
155 // This might happen for highly symmetric cells or if face_centroid IS cell_centroid (e.g. 2D cell).
156 // For now, we keep the original normal direction based on v1,v2,v4 ordering.
157 // A more robust solution for this edge case might be needed if it occurs often.
158 ierr = MPI_Comm_rank(MPI_COMM_WORLD, &rank); CHKERRQ(ierr);
159 LOG_ALLOW(LOCAL, LOG_WARNING, "Rank %d: Dot product for normal orientation is zero. Normal direction might be ambiguous. Keeping initial normal direction from (v2-v1)x(v4-v1).\n", rank);
160 }
161 LOG_ALLOW(LOCAL, LOG_VERBOSE, " Oriented Raw Normal: (%.6e, %.6e, %.6e)\n", normal_x, normal_y, normal_z);
162
163
164 // --- Normalize the (Now Outward-Pointing) Normal Vector ---
165 // Note: normal_magnitude was calculated from initial normal.
166 // If we flipped, magnitude is the same.
167 normal_x /= normal_magnitude;
168 normal_y /= normal_magnitude;
169 normal_z /= normal_magnitude;
170 LOG_ALLOW(LOCAL, LOG_VERBOSE, " Normalized Outward Normal: (%.6f, %.6f, %.6f)\n", normal_x, normal_y, normal_z);
171
172
173 // --- Compute Vector from Target Point to Face Centroid ---
174 PetscReal vec_p_to_fc_x = face_centroid_x - p_target.x;
175 PetscReal vec_p_to_fc_y = face_centroid_y - p_target.y;
176 PetscReal vec_p_to_fc_z = face_centroid_z - p_target.z;
177
178 // --- Compute Signed Distance ---
179 *d_signed = vec_p_to_fc_x * normal_x +
180 vec_p_to_fc_y * normal_y +
181 vec_p_to_fc_z * normal_z;
182
183 LOG_ALLOW(LOCAL, LOG_VERBOSE, " Raw Signed Distance (using outward normal): %.15e\n", *d_signed);
184
185 // --- Apply Threshold ---
186 if (fabs(*d_signed) < threshold) {
187 LOG_ALLOW(LOCAL, LOG_VERBOSE, " Distance %.15e is less than threshold %.1e. Setting to 0.0.\n", *d_signed, threshold);
188 *d_signed = 0.0;
189 }
190
191 LOG_ALLOW(LOCAL, LOG_VERBOSE, "Final Signed Distance: %.15e\n", *d_signed);
192
194 PetscFunctionReturn(0);
195}
#define LOCAL
Logging scope definitions for controlling message output.
Definition logging.h:44
#define LOG_ALLOW(scope, level, fmt,...)
Logging macro that checks both the log level and whether the calling function is in the allowed-funct...
Definition logging.h:199
@ LOG_ERROR
Critical errors that may halt the program.
Definition logging.h:28
@ LOG_WARNING
Non-critical issues that warrant attention.
Definition logging.h:29
@ LOG_VERBOSE
Extremely detailed logs, typically for development use only.
Definition logging.h:33
Here is the caller graph for this function:

◆ CalculateDistancesToCellFaces()

PetscErrorCode CalculateDistancesToCellFaces ( const Cmpnts  p,
const Cell cell,
PetscReal *  d,
const PetscReal  threshold 
)

Implementation of CalculateDistancesToCellFaces().

Computes the signed distances from a point to each face of a cubic cell.

Full API contract (arguments, ownership, side effects) is documented with the header declaration in include/walkingsearch.h.

See also
CalculateDistancesToCellFaces()

Definition at line 206 of file walkingsearch.c.

207{
208 PetscErrorCode ierr;
209 PetscFunctionBeginUser;
211
212 // Validate that the 'cell' pointer is not NULL to prevent dereferencing a null pointer.
213 if (cell == NULL) {
214 LOG_ALLOW(LOCAL,LOG_ERROR, "'cell' is NULL.\n");
215 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "'cell' is NULL.");
216 }
217
218 // Validate that the 'd' pointer is not NULL to ensure there is memory allocated for distance storage.
219 if (d == NULL) {
220 LOG_ALLOW(LOCAL,LOG_ERROR, "'d' is NULL.\n");
221 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "'d' is NULL.");
222 }
223
224 // Compute the centroid of the entire cell
225 Cmpnts cell_centroid = {0.0, 0.0, 0.0};
226 for (int i = 0; i < 8; ++i) {
227 cell_centroid.x += cell->vertices[i].x;
228 cell_centroid.y += cell->vertices[i].y;
229 cell_centroid.z += cell->vertices[i].z;
230 }
231 cell_centroid.x /= 8.0;
232 cell_centroid.y /= 8.0;
233 cell_centroid.z /= 8.0;
234
235 LOG_ALLOW(LOCAL,LOG_DEBUG, "Cell Centroid: (%.6e, %.6e, %.6e)\n",
236 cell_centroid.x, cell_centroid.y, cell_centroid.z);
237
238
239 // Compute the signed distance from point 'p' to the BACK face of the cell.
240 // The BACK face is defined by vertices 0, 3, 2, and 1, with its normal vector pointing in the -z direction.
242 cell->vertices[0], // Vertex 0
243 cell->vertices[3], // Vertex 3
244 cell->vertices[2], // Vertex 2
245 cell->vertices[1], // Vertex 1
246 cell_centroid, // Cell centroid
247 p, // Target point
248 &d[BACK], // Storage location for the BACK face distance
249 threshold // Threshold for zero distance
250 ); CHKERRQ(ierr);
251
252 // Compute the signed distance from point 'p' to the FRONT face of the cell.
253 // The FRONT face is defined by vertices 4, 7, 6, and 5, with its normal vector pointing in the +z direction.
255 cell->vertices[4], // Vertex 4
256 cell->vertices[7], // Vertex 7
257 cell->vertices[6], // Vertex 6
258 cell->vertices[5], // Vertex 5
259 cell_centroid, // Cell centroid
260 p, // Target point
261 &d[FRONT], // Storage location for the FRONT face distance
262 threshold // Threshold for zero distance
263 ); CHKERRQ(ierr);
264
265 // Compute the signed distance from point 'p' to the BOTTOM face of the cell.
266 // The BOTTOM face is defined by vertices 0, 1, 6, and 7, with its normal vector pointing in the -y direction.
268 cell->vertices[0], // Vertex 0
269 cell->vertices[1], // Vertex 1
270 cell->vertices[6], // Vertex 6
271 cell->vertices[7], // Vertex 7
272 cell_centroid, // Cell centroid
273 p, // Target point
274 &d[BOTTOM], // Storage location for the BOTTOM face distance
275 threshold // Threshold for zero distance
276 ); CHKERRQ(ierr);
277
278 // Compute the signed distance from point 'p' to the TOP face of the cell.
279 // The TOP face is defined by vertices 3, 4, 5, and 2, with its normal vector pointing in the +y direction.
281 cell->vertices[3], // Vertex 3
282 cell->vertices[4], // Vertex 4
283 cell->vertices[5], // Vertex 5
284 cell->vertices[2], // Vertex 2
285 cell_centroid, // Cell centroid
286 p, // Target point
287 &d[TOP], // Storage location for the TOP face distance
288 threshold // Threshold for zero distance
289 ); CHKERRQ(ierr);
290
291 // Compute the signed distance from point 'p' to the LEFT face of the cell.
292 // The LEFT face is defined by vertices 0, 7, 4, and 3, with its normal vector pointing in the -x direction.
294 cell->vertices[0], // Vertex 0
295 cell->vertices[7], // Vertex 7
296 cell->vertices[4], // Vertex 4
297 cell->vertices[3], // Vertex 3
298 cell_centroid, // Cell centroid
299 p, // Target point
300 &d[LEFT], // Storage location for the LEFT face distance
301 threshold // Threshold for zero distance
302 ); CHKERRQ(ierr);
303
304 // Compute the signed distance from point 'p' to the RIGHT face of the cell.
305 // The RIGHT face is defined by vertices 1, 2, 5, and 6, with its normal vector pointing in the +x direction.
307 cell->vertices[1], // Vertex 1
308 cell->vertices[2], // Vertex 2
309 cell->vertices[5], // Vertex 5
310 cell->vertices[6], // Vertex 6
311 cell_centroid, // Cell centroid
312 p, // Target point
313 &d[RIGHT], // Storage location for the RIGHT face distance
314 threshold // Threshold for zero distance
315 ); CHKERRQ(ierr);
316
317 LOG_ALLOW(LOCAL, LOG_DEBUG, "Populated d: "
318 "d[LEFT=%d]=%.3e, d[RIGHT=%d]=%.3e, d[BOTTOM=%d]=%.3e, "
319 "d[TOP=%d]=%.3e, d[FRONT=%d]=%.3e, d[BACK=%d]=%.3e\n",
320 LEFT, d[LEFT], RIGHT, d[RIGHT], BOTTOM, d[BOTTOM],
321 TOP, d[TOP], FRONT, d[FRONT], BACK, d[BACK]);
322 LOG_ALLOW(LOCAL, LOG_DEBUG, "Raw d: "
323 "[%.3e, %.3e, %.3e, %.3e, %.3e, %.3e]\n",
324 d[0], d[1], d[2], d[3], d[4], d[5]);
325
327 PetscFunctionReturn(0); // Indicate successful execution of the function.
328}
@ LOG_DEBUG
Detailed debugging information.
Definition logging.h:31
@ TOP
Definition variables.h:145
@ FRONT
Definition variables.h:145
@ BOTTOM
Definition variables.h:145
@ BACK
Definition variables.h:145
@ LEFT
Definition variables.h:145
@ RIGHT
Definition variables.h:145
PetscErrorCode ComputeSignedDistanceToPlane(const Cmpnts v1, const Cmpnts v2, const Cmpnts v3, const Cmpnts v4, const Cmpnts cell_centroid, const Cmpnts p_target, PetscReal *d_signed, const PetscReal threshold)
Internal helper implementation: ComputeSignedDistanceToPlane().
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DeterminePointPosition()

PetscErrorCode DeterminePointPosition ( PetscReal *  d,
PetscInt *  result 
)

Implementation of DeterminePointPosition().

Classifies a point based on precomputed face distances.

Full API contract (arguments, ownership, side effects) is documented with the header declaration in include/walkingsearch.h.

See also
DeterminePointPosition()

Definition at line 338 of file walkingsearch.c.

339{
340
341 PetscFunctionBeginUser;
343 // Validate input pointers
344 if (d == NULL) {
345 LOG_ALLOW(LOCAL,LOG_ERROR, "'d' is NULL.\n");
346 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Input parameter 'd' is NULL.");
347 }
348 if (result == NULL) {
349 LOG_ALLOW(LOCAL,LOG_ERROR, "'result' is NULL.\n");
350 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Output parameter 'result' is NULL.");
351 }
352
353 // Initialize flags
354 PetscBool isInside = PETSC_TRUE;
355 PetscBool isOnBoundary = PETSC_FALSE;
356 PetscInt IntersectionCount = 0; // Counts the number of intersections of the point with various planes of the cell.
357
358 // Analyze distances to determine position
359 for(int i = 0; i < NUM_FACES; i++) {
360 if(d[i] < 0.0) {
361 isInside = PETSC_FALSE; // Point is outside in at least one direction
362 }
363 if(d[i] == 0.0) {
364 isOnBoundary = PETSC_TRUE; // Point is exactly at least one face
365 IntersectionCount++;
366 }
367 }
368
369 // Set the result based on flags
370 if(isInside && isOnBoundary) {
371 if(IntersectionCount == 1){
372 *result = 1; // on face
373 LOG_ALLOW(LOCAL,LOG_VERBOSE, "Particle is on a face of the cell.\n");
374 }
375 else if(IntersectionCount == 2){
376 *result = 2; // on edge
377 LOG_ALLOW(LOCAL,LOG_VERBOSE, "Particle is on an edge of the cell.\n");
378 }
379 else if(IntersectionCount >= 3){
380 *result = 3; // on corner
381 LOG_ALLOW(LOCAL,LOG_VERBOSE, "Particle is on a corner of the cell.\n");
382 }
383 }
384 else if(isInside) {
385 *result = 0; // Inside the cell
386 LOG_ALLOW(LOCAL,LOG_VERBOSE, "Particle is inside the cell.\n");
387 }
388 else {
389 *result = -1; // Outside the cell
390 LOG_ALLOW(LOCAL,LOG_VERBOSE, "Particle is outside the cell.\n");
391 }
392
394 PetscFunctionReturn(0); // Indicate successful execution
395}
@ NUM_FACES
Definition variables.h:145
Here is the caller graph for this function:

◆ GetCellVerticesFromGrid()

PetscErrorCode GetCellVerticesFromGrid ( Cmpnts ***  coor,
PetscInt  idx,
PetscInt  idy,
PetscInt  idz,
Cell cell 
)

Implementation of GetCellVerticesFromGrid().

Retrieves the coordinates of the eight vertices of a cell based on grid indices.

Full API contract (arguments, ownership, side effects) is documented with the header declaration in include/walkingsearch.h.

See also
GetCellVerticesFromGrid()

Definition at line 405 of file walkingsearch.c.

407{
408
409 PetscFunctionBeginUser;
411 // Validate input pointers
412 if (coor == NULL) {
413 LOG_ALLOW(LOCAL,LOG_ERROR, "'coor' is NULL.\n");
414 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Input array 'coor' is NULL.");
415 }
416 if (cell == NULL) {
417 LOG_ALLOW(LOCAL,LOG_ERROR, "'cell' is NULL.\n");
418 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Output parameter 'cell' is NULL.");
419 }
420
421 // Assign vertices based on grid indices matching the figure
422 // Vertex numbering follows the convention depicted in the figure
423 cell->vertices[0] = coor[idz][idy][idx]; // Vertex 0: (i, j, k)
424 cell->vertices[1] = coor[idz][idy][idx+1]; // Vertex 1: (i+1, j, k)
425 cell->vertices[2] = coor[idz][idy+1][idx+1]; // Vertex 2: (i+1, j+1, k)
426 cell->vertices[3] = coor[idz][idy+1][idx]; // Vertex 3: (i, j+1, k)
427 cell->vertices[4] = coor[idz+1][idy+1][idx]; // Vertex 4: (i, j+1, k+1)
428 cell->vertices[5] = coor[idz+1][idy+1][idx+1]; // Vertex 5: (i+1, j+1, k+1)
429 cell->vertices[6] = coor[idz+1][idy][idx+1]; // Vertex 6: (i+1, j, k+1)
430 cell->vertices[7] = coor[idz+1][idy][idx]; // Vertex 7: (i, j, k+1)
431
432 LOG_ALLOW(LOCAL,LOG_VERBOSE, "Retrieved vertices for cell (%d, %d, %d).\n", idx, idy, idz);
433
435 PetscFunctionReturn(0); // Indicate successful execution
436}
Here is the caller graph for this function:

◆ InitializeTraversalParameters()

PetscErrorCode InitializeTraversalParameters ( UserCtx user,
Particle particle,
PetscInt *  idx,
PetscInt *  idy,
PetscInt *  idz,
PetscInt *  traversal_steps 
)

Implementation of InitializeTraversalParameters().

Initializes traversal parameters for locating a particle.

Full API contract (arguments, ownership, side effects) is documented with the header declaration in include/walkingsearch.h.

See also
InitializeTraversalParameters()

Definition at line 446 of file walkingsearch.c.

447{
448 PetscErrorCode ierr;
449 DMDALocalInfo info;
450
451 PetscFunctionBeginUser;
453
454 // Validate input pointers
455 if (user == NULL || particle == NULL || idx == NULL || idy == NULL || idz == NULL || traversal_steps == NULL) {
456 LOG_ALLOW(LOCAL,LOG_ERROR, "One or more input pointers are NULL.\n");
457 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "One or more input pointers are NULL.");
458 }
459
460 // Get grid information (needed for fallback start and potentially validation)
461 ierr = DMDAGetLocalInfo(user->da, &info); CHKERRQ(ierr);
462
463 // --- Check if the particle has a valid previous cell ID ---
464 // Assuming particle->cell stores the *global* indices
465 if (particle->cell[0] >= 0 && particle->cell[1] >= 0 && particle->cell[2] >= 0) {
466
467 // It sees a valid cell ID exists.
468 // Before using it, it explicitly checks if this cell is accessible.
469 PetscBool is_handoff_cell_valid;
470 ierr = CheckCellWithinLocalGrid(user, particle->cell[0], particle->cell[1], particle->cell[2], &is_handoff_cell_valid); CHKERRQ(ierr);
471
472 if (is_handoff_cell_valid) {
473 // FAST PATH: The check passed. The cell is safe. Use it.
474 *idx = particle->cell[0];
475 *idy = particle->cell[1];
476 *idz = particle->cell[2];
477
478 LOG_ALLOW(LOCAL,LOG_DEBUG, "Particle %lld has previous cell (%d, %d, %d). Starting search there.\n",
479 (long long)particle->PID, *idx, *idy, *idz); // Cast id to long long for printing PetscInt64
480
481 } else {
482 // SAFE FALLBACK: The check failed! The handoff cell is NOT safe.
483 // Discard the unsafe handoff cell and revert to starting at the
484 // guaranteed-safe local corner.
485 *idx = info.xs;
486 *idy = info.ys;
487 *idz = info.zs;
488 }
489 } else {
490 // No valid previous cell ID (e.g., -1,-1,-1 or first time).
491 // IMPROVED: Use distance-based multi-seed approach to find the closest strategic cell.
492 // Sample the 6 face centers + volume center (7 total candidates) and start from the closest.
493
494 Cmpnts ***cent;
495
496 // Get local cell center array (lCent)
497 // For cell (i,j,k), the center is at cent[k+1][j+1][i+1]
498 ierr = DMDAVecGetArrayRead(user->fda, user->lCent, &cent); CHKERRQ(ierr);
499
500 // Define 7 strategic sampling points: 6 face centers + 1 volume center
501 PetscInt candidates[7][3];
502 const char* candidate_names[7] = {"X-min face", "X-max face", "Y-min face",
503 "Y-max face", "Z-min face", "Z-max face", "Volume center"};
504
505 // Calculate middle indices for each dimension
506 PetscInt mid_i = info.xs + (info.xm - 1) / 2;
507 PetscInt mid_j = info.ys + (info.ym - 1) / 2;
508 PetscInt mid_k = info.zs + (info.zm - 1) / 2;
509 PetscInt max_i = info.xs + info.xm - 2;
510 PetscInt max_j = info.ys + info.ym - 2;
511 PetscInt max_k = info.zs + info.zm - 2;
512
513 // Clamp max indices to ensure valid cells
514 if (max_i < info.xs) max_i = info.xs;
515 if (max_j < info.ys) max_j = info.ys;
516 if (max_k < info.zs) max_k = info.zs;
517
518 // X-min face center (i=min, j=mid, k=mid)
519 candidates[0][0] = info.xs; candidates[0][1] = mid_j; candidates[0][2] = mid_k;
520 // X-max face center (i=max, j=mid, k=mid)
521 candidates[1][0] = max_i; candidates[1][1] = mid_j; candidates[1][2] = mid_k;
522 // Y-min face center (i=mid, j=min, k=mid)
523 candidates[2][0] = mid_i; candidates[2][1] = info.ys; candidates[2][2] = mid_k;
524 // Y-max face center (i=mid, j=max, k=mid)
525 candidates[3][0] = mid_i; candidates[3][1] = max_j; candidates[3][2] = mid_k;
526 // Z-min face center (i=mid, j=mid, k=min)
527 candidates[4][0] = mid_i; candidates[4][1] = mid_j; candidates[4][2] = info.zs;
528 // Z-max face center (i=mid, j=mid, k=max)
529 candidates[5][0] = mid_i; candidates[5][1] = mid_j; candidates[5][2] = max_k;
530 // Volume center (i=mid, j=mid, k=mid)
531 candidates[6][0] = mid_i; candidates[6][1] = mid_j; candidates[6][2] = mid_k;
532
533 // Find the closest candidate cell
534 PetscReal min_distance = PETSC_MAX_REAL;
535 PetscInt best_candidate = 6; // Default to volume center
536
537 for (PetscInt i = 0; i < 7; i++) {
538 PetscInt ci = candidates[i][0];
539 PetscInt cj = candidates[i][1];
540 PetscInt ck = candidates[i][2];
541
542 // Get cell center coordinates: cent[k+1][j+1][i+1] for cell (i,j,k)
543 Cmpnts cell_center = cent[ck+1][cj+1][ci+1];
544
545 // Calculate Euclidean distance from particle to cell center
546 PetscReal dx = particle->loc.x - cell_center.x;
547 PetscReal dy = particle->loc.y - cell_center.y;
548 PetscReal dz = particle->loc.z - cell_center.z;
549 PetscReal distance = sqrt(dx*dx + dy*dy + dz*dz);
550
551 LOG_ALLOW(LOCAL, LOG_DEBUG, " Candidate %d (%s) at cell (%d,%d,%d): center=(%.3e,%.3e,%.3e), distance=%.3e\n",
552 i, candidate_names[i], ci, cj, ck, cell_center.x, cell_center.y, cell_center.z, distance);
553
554 if (distance < min_distance) {
555 min_distance = distance;
556 best_candidate = i;
557 *idx = ci;
558 *idy = cj;
559 *idz = ck;
560 }
561 }
562
563 // Restore array
564 ierr = DMDAVecRestoreArrayRead(user->fda, user->lCent, &cent); CHKERRQ(ierr);
565
566 LOG_ALLOW(LOCAL, LOG_INFO, "Particle %lld has no valid previous cell. Multi-seed search selected %s at cell (%d,%d,%d) with distance %.3e.\n",
567 (long long)particle->PID, candidate_names[best_candidate], *idx, *idy, *idz, min_distance);
568 }
569
570 // Initialize traversal step counter
571 *traversal_steps = 0;
572
573 // Log the chosen starting point
574 LOG_ALLOW(LOCAL,LOG_INFO, "Traversal for particle %lld initialized to start at cell (%d, %d, %d).\n",
575 (long long)particle->PID, *idx, *idy, *idz);
576
577 PetscFunctionReturn(0);
578}
@ LOG_INFO
Informational messages about program execution.
Definition logging.h:30
Vec lCent
Definition variables.h:858
PetscInt cell[3]
Definition variables.h:167
Cmpnts loc
Definition variables.h:168
PetscInt64 PID
Definition variables.h:166
PetscErrorCode CheckCellWithinLocalGrid(UserCtx *user, PetscInt idx, PetscInt idy, PetscInt idz, PetscBool *is_within)
Implementation of CheckCellWithinLocalGrid().
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CheckCellWithinLocalGrid()

PetscErrorCode CheckCellWithinLocalGrid ( UserCtx user,
PetscInt  idx,
PetscInt  idy,
PetscInt  idz,
PetscBool *  is_within 
)

Implementation of CheckCellWithinLocalGrid().

Checks if the current cell indices are within the local grid boundaries.

Full API contract (arguments, ownership, side effects) is documented with the header declaration in include/walkingsearch.h.

See also
CheckCellWithinLocalGrid()

Definition at line 589 of file walkingsearch.c.

590{
591 PetscErrorCode ierr;
592 DMDALocalInfo info_nodes; // Node information from the DMDA that defines ghost regions (user->fda)
593
594 PetscFunctionBeginUser; // Assuming this is part of your PETSc style
596
597 // Validate inputs
598 if (user == NULL || is_within == NULL) {
599 LOG_ALLOW(LOCAL, LOG_ERROR, "Input pointer is NULL.\n");
600 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Input pointer is NULL in CheckCellWithinLocalGrid.");
601 }
602 if (user->fda == NULL) {
603 LOG_ALLOW(LOCAL, LOG_ERROR, "user->fda is NULL.Cannot get ghost info.\n");
604 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "user->fda is NULL. Cannot get ghost info.");
605 }
606
607 // Get node info from user->fda (this DMDA has the ghost layer information for nodes)
608 ierr = DMDAGetLocalInfo(user->fda, &info_nodes); CHKERRQ(ierr);
609
610 // Determine the range of GLOBAL CELL INDICES that are covered by this rank's ghosted NODAL region.
611 // A cell C(i,j,k) has origin node N(i,j,k).
612 // The ghosted nodal region starts at global node index info_nodes.gxs and has info_nodes.gxm nodes.
613
614 // Global starting index of the first cell whose origin node is within the ghosted nodal region.
615 PetscInt gxs_cell_global_start = info_nodes.gxs;
616 PetscInt gys_cell_global_start = info_nodes.gys;
617 PetscInt gzs_cell_global_start = info_nodes.gzs;
618
619 // Number of cells that can be formed starting from nodes within the ghosted nodal region.
620 // If there are N ghosted nodes (info_nodes.gxm), they can be origins for N-1 cells.
621 PetscInt gxm_cell_local_count = (info_nodes.gxm > 0) ? info_nodes.gxm - 1 : 0;
622 PetscInt gym_cell_local_count = (info_nodes.gym > 0) ? info_nodes.gym - 1 : 0;
623 PetscInt gzm_cell_local_count = (info_nodes.gzm > 0) ? info_nodes.gzm - 1 : 0;
624
625 // Global exclusive end index for cells whose origins are in the ghosted nodal region.
626 PetscInt gxe_cell_global_end_exclusive = gxs_cell_global_start + gxm_cell_local_count;
627 PetscInt gye_cell_global_end_exclusive = gys_cell_global_start + gym_cell_local_count;
628 PetscInt gze_cell_global_end_exclusive = gzs_cell_global_start + gzm_cell_local_count;
629
630 // Check if the given global cell index (idx, idy, idz) falls within this range.
631 // This means the origin node of cell (idx,idy,idz) is within the rank's accessible ghosted node region,
632 // and that node can indeed serve as a cell origin (i.e., it's not the very last node in the ghosted region).
633 if (idx >= gxs_cell_global_start && idx < gxe_cell_global_end_exclusive &&
634 idy >= gys_cell_global_start && idy < gye_cell_global_end_exclusive &&
635 idz >= gzs_cell_global_start && idz < gze_cell_global_end_exclusive) {
636 *is_within = PETSC_TRUE;
637 } else {
638 *is_within = PETSC_FALSE;
639 }
640
641 LOG_ALLOW(LOCAL, LOG_DEBUG, "Cell (origin node glob idx) (%d, %d, %d) is %s the ghosted local grid (covered cell origins x:[%d..%d), y:[%d..%d), z:[%d..%d)).\n",
642 idx, idy, idz, (*is_within) ? "within" : "outside",
643 gxs_cell_global_start, gxe_cell_global_end_exclusive,
644 gys_cell_global_start, gye_cell_global_end_exclusive,
645 gzs_cell_global_start, gze_cell_global_end_exclusive);
646
648 PetscFunctionReturn(0);
649}
Here is the caller graph for this function:

◆ RetrieveCurrentCell()

PetscErrorCode RetrieveCurrentCell ( UserCtx user,
PetscInt  idx,
PetscInt  idy,
PetscInt  idz,
Cell cell 
)

Implementation of RetrieveCurrentCell().

Retrieves the coordinates of the eight vertices of the current cell.

Full API contract (arguments, ownership, side effects) is documented with the header declaration in include/walkingsearch.h.

See also
RetrieveCurrentCell()

Definition at line 659 of file walkingsearch.c.

660{
661 PetscErrorCode ierr;
662 Vec Coor;
663 Cmpnts ***coor_array;
664 PetscMPIInt rank;
665 DMDALocalInfo info_nodes;
666
667 PetscFunctionBeginUser;
669
670 // Validate input pointers
671 if (user == NULL || cell == NULL) {
672 LOG_ALLOW(LOCAL,LOG_ERROR, "One or more input pointers are NULL.\n");
673 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "One or more input pointers are NULL.");
674 }
675
676 // Get local coordinates
677 ierr = DMGetCoordinatesLocal(user->da, &Coor); CHKERRQ(ierr);
678 ierr = DMDAVecGetArrayRead(user->fda, Coor, &coor_array); CHKERRQ(ierr);
679
680 // Get the local grid information FOR THE GHOSTED ARRAY from user->fda
681 // This info contains the mapping between global and local indices.
682 ierr = DMDAGetLocalInfo(user->fda, &info_nodes); CHKERRQ(ierr);
683
684 PetscInt idx_local = idx; // - info_nodes.gxs;
685 PetscInt idy_local = idy; // - info_nodes.gys;
686 PetscInt idz_local = idz; // - info_nodes.gzs;
687
688 LOG_ALLOW(LOCAL,LOG_VERBOSE," [Rank %d] Getting vertex coordinates for cell: %d,%d,%d whose local coordinates are %d,%d,%d.\n",rank,idx,idy,idz,idx_local,idy_local,idz_local);
689
690 // Get the current cell's vertices
691 ierr = GetCellVerticesFromGrid(coor_array, idx_local, idy_local, idz_local, cell); CHKERRQ(ierr);
692
693 // Restore array
694 ierr = DMDAVecRestoreArrayRead(user->fda, Coor, &coor_array); CHKERRQ(ierr);
695
696 // Get MPI rank for debugging
697 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
698
699 // Debug: Print cell vertices
700 LOG_ALLOW(LOCAL,LOG_VERBOSE, "Cell (%d, %d, %d) vertices \n", idx, idy, idz);
701 ierr = LOG_CELL_VERTICES(cell, rank); CHKERRQ(ierr);
702
704 PetscFunctionReturn(0);
705}
PetscErrorCode LOG_CELL_VERTICES(const Cell *cell, PetscMPIInt rank)
Prints the coordinates of a cell's vertices.
Definition logging.c:205
PetscErrorCode GetCellVerticesFromGrid(Cmpnts ***coor, PetscInt idx, PetscInt idy, PetscInt idz, Cell *cell)
Implementation of GetCellVerticesFromGrid().
Here is the call graph for this function:
Here is the caller graph for this function:

◆ EvaluateParticlePosition()

PetscErrorCode EvaluateParticlePosition ( const Cell cell,
PetscReal *  d,
const Cmpnts  p,
PetscInt *  position,
const PetscReal  threshold 
)

Implementation of EvaluateParticlePosition().

Determines the spatial relationship of a particle relative to a cubic cell.

Full API contract (arguments, ownership, side effects) is documented with the header declaration in include/walkingsearch.h.

See also
EvaluateParticlePosition()

Definition at line 715 of file walkingsearch.c.

716{
717 PetscErrorCode ierr;
718 PetscReal cellSize;
719 PetscReal cellThreshold;
720 PetscFunctionBeginUser;
722
723 // Validate input pointers to ensure they are not NULL, preventing potential segmentation faults.
724 if (cell == NULL || position == NULL) {
725 LOG_ALLOW(LOCAL,LOG_ERROR, "One or more input pointers are NULL.\n");
726 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "One or more input pointers are NULL.");
727 }
728
729 // Compute a local cell size
730 ierr = GetCellCharacteristicSize(cell, &cellSize); CHKERRQ(ierr);
731
732 // scale base threshold by cell size (ex: if threshold = 1e-6 and cell size = 0.01, cellThreshold = 1e-8)
733 cellThreshold = threshold*cellSize;
734
735 // Invoke the function to calculate signed distances from the particle to each face of the cell.
736 // The distances are stored in the array pointed to by 'd'.
737 ierr = CalculateDistancesToCellFaces(p, cell, d, cellThreshold); CHKERRQ(ierr);
738 CHKERRQ(ierr); // Check for errors in distance calculation.
739
740
741 // Catch degenerate-plane error manually:
742 if (ierr == PETSC_ERR_USER) {
744 "Skipping cell due to degenerate face.\n");
745 // We can set *position = -1 here
746 *position = -1; // treat as outside
747 return 0; // not a fatal error, just skip
748 } else {
749 CHKERRQ(ierr);
750 }
751
752
753 // Debugging output: Print the computed distances to each face for verification purposes.
754 LOG_ALLOW(LOCAL,LOG_DEBUG, "Face Distances:\n");
756 CHKERRQ(ierr); // Check for errors in printing distances.
757
758 // Determine the particle's position relative to the cell based on the computed distances.
759 // The function sets the value pointed to by 'position' accordingly:
760 // 0 for inside, 1 (2 or 3) for on the boundary, and -1 for outside.
761 ierr = DeterminePointPosition(d,position);
762 CHKERRQ(ierr); // Check for errors in position determination.
763
765 PetscFunctionReturn(0); // Indicate successful execution of the function.
766}
PetscBool is_function_allowed(const char *functionName)
Checks if a given function is in the allow-list.
Definition logging.c:183
LogLevel get_log_level()
Retrieves the current logging level from the environment variable LOG_LEVEL.
Definition logging.c:84
PetscErrorCode LOG_FACE_DISTANCES(PetscReal *d)
Prints the signed distances to each face of the cell.
Definition logging.c:230
PetscErrorCode GetCellCharacteristicSize(const Cell *cell, PetscReal *cellSize)
Implementation of GetCellCharacteristicSize().
PetscErrorCode CalculateDistancesToCellFaces(const Cmpnts p, const Cell *cell, PetscReal *d, const PetscReal threshold)
Implementation of CalculateDistancesToCellFaces().
PetscErrorCode DeterminePointPosition(PetscReal *d, PetscInt *result)
Implementation of DeterminePointPosition().
Here is the call graph for this function:
Here is the caller graph for this function:

◆ UpdateCellIndicesBasedOnDistances()

PetscErrorCode UpdateCellIndicesBasedOnDistances ( PetscReal  d[NUM_FACES],
PetscInt *  idx,
PetscInt *  idy,
PetscInt *  idz 
)

Internal helper implementation: UpdateCellIndicesBasedOnDistances().

Updates the cell indices based on the signed distances to each face.

Local to this translation unit.

Definition at line 774 of file walkingsearch.c.

775{
776
777 PetscFunctionBeginUser;
779 /*
780 PetscInt cxm,cxs; // maximum & minimum cell ID in x
781 PetscInt cym,cys; // maximum & minimum cell ID in y
782 PetscInt czm,czs; // maximum & minimum cell ID in z
783
784 cxs = info->xs; cxm = cxs + info->xm - 2;
785 cys = info->ys; cym = cys + info->ym - 2;
786 czs = info->zs; czm = czs + info->zm - 2;
787 */
788
789 // LOG_ALLOW(LOCAL, LOG_DEBUG, "Received d: "
790 // "d[LEFT=%d]=%.3e, d[RIGHT=%d]=%.3e, d[BOTTOM=%d]=%.3e, "
791 // "d[TOP=%d]=%.3e, d[FRONT=%d]=%.3e, d[BACK=%d]=%.3e\n",
792 // LEFT, d[LEFT], RIGHT, d[RIGHT], BOTTOM, d[BOTTOM],
793 // TOP, d[TOP], FRONT, d[FRONT], BACK, d[BACK]);
794 // LOG_ALLOW(LOCAL, LOG_DEBUG, "Raw d: "
795 // "[%.3e, %.3e, %.3e, %.3e, %.3e, %.3e]\n",
796 // d[0], d[1], d[2], d[3], d[4], d[5]);
797
798 // Validate input pointers
799 if (d == NULL) {
800 LOG_ALLOW(LOCAL,LOG_ERROR, "'d' is NULL.\n");
801 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Input array 'd' is NULL.");
802 }
803 if (idx == NULL || idy == NULL || idz == NULL) {
804 LOG_ALLOW(LOCAL,LOG_ERROR, "One or more index pointers are NULL.\n");
805 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "One or more index pointers are NULL.");
806 }
807
808 // Debug: Print current face distances
809 LOG_ALLOW(LOCAL,LOG_DEBUG, "Current Face Distances:\n");
811
812 // Update k-direction based on FRONT and BACK distances
813 if (d[FRONT] < 0.0) {
814 LOG_ALLOW(LOCAL,LOG_VERBOSE, "Condition met: d[FRONT] < 0.0, incrementing idz.\n");
815 (*idz) += 1;
816 }
817 else if(d[BACK] < 0.0){
818 LOG_ALLOW(LOCAL,LOG_VERBOSE, "Condition met: d[BACK] < 0.0, decrementing idz.\n");
819 (*idz) -= 1;
820 }
821
822 // Update i-direction based on LEFT and RIGHT distances
823 if (d[LEFT] < 0.0) {
824 LOG_ALLOW(LOCAL,LOG_VERBOSE, "Condition met: d[LEFT] < 0.0, decrementing idx.\n");
825 (*idx) -= 1;
826 }
827 else if (d[RIGHT] < 0.0) {
828 LOG_ALLOW(LOCAL,LOG_VERBOSE, "Condition met: d[RIGHT] < 0.0, incrementing idx.\n");
829 (*idx) += 1;
830 }
831
832 // Update j-direction based on BOTTOM and TOP distances
833 if (d[BOTTOM] < 0.0) {
834 LOG_ALLOW(LOCAL,LOG_VERBOSE, "Condition met: d[BOTTOM] < 0.0, decrementing idy.\n");
835 (*idy) -= 1;
836 }
837 else if (d[TOP] < 0.0) {
838 LOG_ALLOW(LOCAL,LOG_VERBOSE, "Condition met: d[TOP] < 0.0, incrementing idy.\n");
839 (*idy) += 1;
840 }
841
842 /*
843 // The 'cell' corners you can reference go from [xs .. xs+xm-1], but
844 // to form a valid cell in x, you need (idx+1) in range, so max is (xs+xm-2).
845 *idx = PetscMax(cxs, PetscMin(*idx, cxm));
846 *idy = PetscMax(cys, PetscMin(*idy, cym));
847 *idz = PetscMax(czs, PetscMin(*idz, czm));
848 */
849
850 LOG_ALLOW(LOCAL,LOG_DEBUG, "Updated Indices - idx, idy, idz: %d, %d, %d\n", *idx, *idy, *idz);
851
853 PetscFunctionReturn(0); // Indicate successful execution
854}
Here is the call graph for this function:
Here is the caller graph for this function:

◆ FinalizeTraversal()

PetscErrorCode FinalizeTraversal ( UserCtx user,
Particle particle,
PetscInt  traversal_steps,
PetscBool  cell_found,
PetscInt  idx,
PetscInt  idy,
PetscInt  idz 
)

Implementation of FinalizeTraversal().

Finalizes the traversal by reporting the results.

Full API contract (arguments, ownership, side effects) is documented with the header declaration in include/walkingsearch.h.

See also
FinalizeTraversal()

Definition at line 865 of file walkingsearch.c.

866{
867 PetscFunctionBeginUser;
868 // Validate input pointers
869 if (user == NULL || particle == NULL) {
870 LOG_ALLOW(LOCAL,LOG_ERROR, "One or more input pointers are NULL.\n");
871 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "One or more input pointers are NULL.");
872 }
873
874 if (cell_found) {
875 LOG_ALLOW(LOCAL,LOG_INFO, "Particle located in cell (%d, %d, %d) after %d traversal steps.\n",
876 idx, idy, idz, traversal_steps);
877 }
878 else {
879 LOG_ALLOW(LOCAL,LOG_WARNING, "Particle could not be located within the grid after %d traversal steps.\n", (PetscInt)traversal_steps);
880 particle->cell[0] = -1;
881 particle->cell[1] = -1;
882 particle->cell[2] = -1;
883 }
884
885 LOG_ALLOW(LOCAL, LOG_INFO, "Completed final traversal sync across all ranks.\n");
886
887
889 PetscFunctionReturn(0);
890}

◆ FindOwnerOfCell()

PetscErrorCode FindOwnerOfCell ( UserCtx user,
PetscInt  i,
PetscInt  j,
PetscInt  k,
PetscMPIInt *  owner_rank 
)

Internal helper implementation: FindOwnerOfCell().

Finds the MPI rank that owns a given global cell index.

Local to this translation unit.

Definition at line 898 of file walkingsearch.c.

899{
900 PetscErrorCode ierr;
901 PetscMPIInt size;
902
903 PetscFunctionBeginUser;
904
906
907 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size); CHKERRQ(ierr);
908
909 // --- 1. Input Validation ---
910 if (!user || !user->RankCellInfoMap) {
911 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "UserCtx or RankCellInfoMap is not initialized in FindOwnerOfCell.");
912 }
913 if (!owner_rank) {
914 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Output pointer owner_rank is NULL in FindOwnerOfCell.");
915 }
916
917 // --- 2. Linear Search through the Decomposition Map ---
918 // Initialize to a "not found" state.
919 *owner_rank = -1;
920
921 // Loop through the map, which contains the ownership info for every rank 'r'.
922 for (PetscMPIInt r = 0; r < size; ++r) {
923 const RankCellInfo *info = &user->RankCellInfoMap[r];
924
925 // A rank owns a cell if the cell's index is within its start (inclusive)
926 // and end (exclusive) range for all three dimensions.
927 if ((i >= info->xs_cell && i < info->xs_cell + info->xm_cell) &&
928 (j >= info->ys_cell && j < info->ys_cell + info->ym_cell) &&
929 (k >= info->zs_cell && k < info->zs_cell + info->zm_cell))
930 {
931 *owner_rank = r; // We found the owner.
932 break; // The search is over, exit the loop.
933 }
934 }
935
936 // --- 3. Logging for Diagnostics ---
937 // This is extremely useful for debugging particle migration issues.
938 if (*owner_rank == -1) {
939 LOG_ALLOW(LOCAL, LOG_DEBUG, "No owner found for global cell (%d, %d, %d). It is likely outside the domain.\n", i, j, k);
940 } else {
941 LOG_ALLOW(LOCAL, LOG_DEBUG, "Owner of cell (%d, %d, %d) is Rank %d.\n", i, j, k, *owner_rank);
942 }
943
945 PetscFunctionReturn(0);
946}
PetscInt ys_cell
Definition variables.h:187
PetscInt xs_cell
Definition variables.h:187
PetscInt zm_cell
Definition variables.h:188
PetscInt zs_cell
Definition variables.h:187
PetscInt xm_cell
Definition variables.h:188
RankCellInfo * RankCellInfoMap
Definition variables.h:881
PetscInt ym_cell
Definition variables.h:188
A lean struct to hold the global cell ownership range for a single MPI rank.
Definition variables.h:186
Here is the caller graph for this function:

◆ LocateParticleOrFindMigrationTarget()

PetscErrorCode LocateParticleOrFindMigrationTarget ( UserCtx user,
Particle particle,
ParticleLocationStatus status_out 
)

Implementation of LocateParticleOrFindMigrationTarget().

Locates a particle's host cell or identifies its migration target using a robust walk search.

Full API contract (arguments, ownership, side effects) is documented with the header declaration in include/walkingsearch.h.

See also
LocateParticleOrFindMigrationTarget()

Definition at line 957 of file walkingsearch.c.

960{
961 PetscErrorCode ierr;
962 PetscMPIInt rank;
963
964 // --- Search State Variables ---
965 PetscInt idx, idy, idz; // Current search cell global indices
966 PetscInt traversal_steps; // Counter to prevent infinite loops
967 PetscBool search_concluded = PETSC_FALSE; // Flag to terminate the main while loop
968
969 // --- Oscillation/Stuck Loop Detection Variables ---
970 PetscInt repeatedIndexCount = 0;
971 PetscInt prevIdx = PETSC_MIN_INT, prevIdy = PETSC_MIN_INT, prevIdz = PETSC_MIN_INT;
972 PetscInt last_position_result = -999;
973
974 PetscFunctionBeginUser;
976 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
977 if (user && user->simCtx) {
981 }
982 }
983
984 // --- 1. Initialize the Search ---
985 ierr = InitializeTraversalParameters(user, particle, &idx, &idy, &idz, &traversal_steps); CHKERRQ(ierr);
986
987 LOG_ALLOW(LOCAL,LOG_VERBOSE, " The Threshold for considering a particle to be at a face is %.16f.\n",DISTANCE_THRESHOLD);
988
989 LOG_ALLOW(LOCAL,LOG_TRACE," [PID %lld]Traversal Initiated at : i = %d, j = %d, k = %d.\n",(long long)particle->PID,idx,idy,idz);
990
991 // --- 2. Main Walking Search Loop ---
992 while (!search_concluded && traversal_steps < MAX_TRAVERSAL) {
993 traversal_steps++;
994
995 // --- 2a. GLOBAL Domain Boundary Check ---
996 // IMPROVED: Instead of failing immediately when hitting a boundary, clamp the indices
997 // to the boundary and continue searching. This allows the search to explore other
998 // directions that might lead to the particle.
999 PetscBool hit_boundary = PETSC_FALSE;
1000 if (idx < 0) {
1001 idx = 0;
1002 hit_boundary = PETSC_TRUE;
1004 LOG_ALLOW(LOCAL, LOG_DEBUG, "[PID %lld]: Hit global i-min boundary, clamped to idx=0.\n", (long long)particle->PID);
1005 }
1006 if (idx >= (user->IM - 1)) {
1007 idx = user->IM - 2;
1008 hit_boundary = PETSC_TRUE;
1010 LOG_ALLOW(LOCAL, LOG_DEBUG, "[PID %lld]: Hit global i-max boundary, clamped to idx=%d.\n", (long long)particle->PID, idx);
1011 }
1012 if (idy < 0) {
1013 idy = 0;
1014 hit_boundary = PETSC_TRUE;
1016 LOG_ALLOW(LOCAL, LOG_DEBUG, "[PID %lld]: Hit global j-min boundary, clamped to idy=0.\n", (long long)particle->PID);
1017 }
1018 if (idy >= (user->JM - 1)) {
1019 idy = user->JM - 2;
1020 hit_boundary = PETSC_TRUE;
1022 LOG_ALLOW(LOCAL, LOG_DEBUG, "[PID %lld]: Hit global j-max boundary, clamped to idy=%d.\n", (long long)particle->PID, idy);
1023 }
1024 if (idz < 0) {
1025 idz = 0;
1026 hit_boundary = PETSC_TRUE;
1028 LOG_ALLOW(LOCAL, LOG_DEBUG, "[PID %lld]: Hit global k-min boundary, clamped to idz=0.\n", (long long)particle->PID);
1029 }
1030 if (idz >= (user->KM - 1)) {
1031 idz = user->KM - 2;
1032 hit_boundary = PETSC_TRUE;
1034 LOG_ALLOW(LOCAL, LOG_DEBUG, "[PID %lld]: Hit global k-max boundary, clamped to idz=%d.\n", (long long)particle->PID, idz);
1035 }
1036
1037 if (hit_boundary) {
1038 LOG_ALLOW(LOCAL, LOG_INFO, "[PID %lld]: Hit GLOBAL domain boundary. Clamped to boundary cell (%d,%d,%d) and continuing search.\n",
1039 (long long)particle->PID, idx, idy, idz);
1040 }
1041
1042 // --- 2b. LOCAL GHOST REGION CHECK (PREVENTS SEGV) ---
1043 // Before trying to access cell data, check if we have it.
1044 PetscBool is_cell_local;
1045 ierr = CheckCellWithinLocalGrid(user, idx, idy, idz, &is_cell_local); CHKERRQ(ierr);
1046 if (!is_cell_local) {
1047 // We have walked outside the local rank's ghost region.
1048 // This definitively means the particle belongs to another rank.
1049 // Conclude the search here; the current (idx,idy,idz) is the handoff cell.
1050 LOG_ALLOW(LOCAL, LOG_INFO, "[PID %lld]: Walked outside local ghost region to cell (%d,%d,%d). Concluding search for handoff.\n",
1051 (long long)particle->PID, idx, idy, idz);
1052 search_concluded = PETSC_TRUE;
1053 continue; // Skip the rest of the loop; proceed to finalization.
1054 }
1055
1056 // --- 2c. Stuck Loop Detection & Enhanced Tie-Breaker ---
1057 if (idx == prevIdx && idy == prevIdy && idz == prevIdz) {
1058 repeatedIndexCount++;
1059 if (repeatedIndexCount > REPEAT_COUNT_THRESHOLD) {
1060 // Only apply tie-breaker if we are stuck for the right reason (on a boundary)
1061 if (last_position_result >= 1) {
1063 LOG_ALLOW(LOCAL, LOG_WARNING, "[PID %lld]: Stuck on boundary of cell (%d,%d,%d) for %d steps. Applying enhanced tie-breaker.\n",
1064 (long long)particle->PID, idx, idy, idz, repeatedIndexCount);
1065
1066 // Re-evaluate at the stuck cell to get definitive weights
1067 Cell final_cell;
1068 PetscReal final_d[NUM_FACES];
1069 PetscInt final_position; // Dummy variable
1070
1071 ierr = RetrieveCurrentCell(user, idx, idy, idz, &final_cell); CHKERRQ(ierr);
1072 ierr = EvaluateParticlePosition(&final_cell, final_d, particle->loc, &final_position, DISTANCE_THRESHOLD);
1073
1074 if (ierr == 0) { // If evaluation succeeded
1075 ierr = UpdateParticleWeights(final_d, particle); CHKERRQ(ierr);
1076 search_concluded = PETSC_TRUE; // Conclude search, accepting this cell.
1077 } else { // Evaluation failed (e.g., degenerate cell)
1078 LOG_ALLOW(LOCAL, LOG_ERROR, "[PID %lld]: Tie-breaker failed during final evaluation at cell (%d,%d,%d). Search fails.\n",
1079 (long long)particle->PID, idx, idy, idz);
1080 idx = -1; // Invalidate result
1081 search_concluded = PETSC_TRUE;
1082 }
1083 } else { // Stuck for the wrong reason (not on a boundary)
1084 LOG_ALLOW(LOCAL, LOG_ERROR, "[PID %lld]: Search is stuck at cell (%d,%d,%d) but not on a boundary. This indicates a logic error. Failing search.\n",
1085 (long long)particle->PID, idx, idy, idz);
1086 idx = -1; // Invalidate result
1087 search_concluded = PETSC_TRUE;
1088 }
1089 if(search_concluded) continue;
1090 }
1091 } else {
1092 repeatedIndexCount = 0;
1093 }
1094 prevIdx = idx; prevIdy = idy; prevIdz = idz;
1095
1096 // --- 2d. Geometric Evaluation ---
1097 Cell current_cell;
1098 PetscReal distances[NUM_FACES];
1099 PetscInt position_in_cell;
1100
1101 ierr = RetrieveCurrentCell(user, idx, idy, idz, &current_cell); CHKERRQ(ierr);
1102 ierr = EvaluateParticlePosition(&current_cell, distances, particle->loc, &position_in_cell, DISTANCE_THRESHOLD); CHKERRQ(ierr);
1103 last_position_result = position_in_cell;
1104
1105 // --- 2e. Decision Making ---
1106 if (position_in_cell >= 0) { // Particle is INSIDE or ON THE BOUNDARY
1107 search_concluded = PETSC_TRUE;
1108 ierr = UpdateParticleWeights(distances, particle); CHKERRQ(ierr);
1109 } else { // Particle is OUTSIDE
1110 ierr = UpdateCellIndicesBasedOnDistances(distances, &idx, &idy, &idz); CHKERRQ(ierr);
1111 }
1112 }
1113
1114 user->simCtx->searchMetrics.traversalStepsSum += traversal_steps;
1116 traversal_steps);
1117
1118 // --- 3. Finalize and Determine Actionable Status ---
1119 if (idx == -1 || (!search_concluded && traversal_steps >= MAX_TRAVERSAL)) {
1120 if (idx != -1) {
1122 LOG_ALLOW(LOCAL, LOG_ERROR, "[PID %lld]: Search FAILED, exceeded MAX_TRAVERSAL limit of %d.\n",
1123 (long long)particle->PID, MAX_TRAVERSAL);
1124 }
1125 *status_out = LOST;
1126 particle->cell[0] = -1; particle->cell[1] = -1; particle->cell[2] = -1;
1127 } else {
1128 // Search succeeded in finding a candidate cell, now determine its owner.
1129 PetscMPIInt owner_rank;
1130 ierr = FindOwnerOfCell(user, idx, idy, idz, &owner_rank); CHKERRQ(ierr);
1131
1132 LOG_ALLOW(LOCAL,LOG_TRACE," [PID %ld] Owner rank : %d.\n",particle->PID,owner_rank);
1133
1134 // Always update the particle's cell index. It's a good guess for the receiving rank.
1135 particle->cell[0] = idx; particle->cell[1] = idy; particle->cell[2] = idz;
1136
1137 if (owner_rank == rank) {
1138 *status_out = ACTIVE_AND_LOCATED;
1139 } else if (owner_rank != -1) {
1140 // Particle belongs to another rank. Return the direct, actionable status.
1141 *status_out = MIGRATING_OUT;
1142 particle->destination_rank = owner_rank;
1143 } else { // Found a valid index, but no owner in the map.
1144 *status_out = LOST;
1145 particle->cell[0] = -1; particle->cell[1] = -1; particle->cell[2] = -1;
1146 }
1147 }
1148
1149 LOG_ALLOW(LOCAL,LOG_DEBUG,"[Rank %d][PID %ld] Search complete.\n",rank,particle->PID);
1150
1151 // --- 4. Report the Final Outcome ---
1152 ierr = ReportSearchOutcome(particle, *status_out, traversal_steps); CHKERRQ(ierr);
1153
1155 PetscFunctionReturn(0);
1156}
PetscErrorCode UpdateParticleWeights(PetscReal *d, Particle *particle)
Updates a particle's interpolation weights based on distances to cell faces.
@ LOG_TRACE
Very fine-grained tracing information for in-depth debugging.
Definition logging.h:32
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:814
@ LOST
Definition variables.h:139
@ ACTIVE_AND_LOCATED
Definition variables.h:137
@ MIGRATING_OUT
Definition variables.h:138
PetscInt64 boundaryClampCount
Definition variables.h:231
PetscInt KM
Definition variables.h:820
PetscInt64 traversalStepsSum
Definition variables.h:226
PetscInt currentSettlementPass
Definition variables.h:235
PetscInt64 reSearchCount
Definition variables.h:227
PetscMPIInt destination_rank
Definition variables.h:172
PetscInt64 maxTraversalSteps
Definition variables.h:228
PetscInt JM
Definition variables.h:820
SearchMetricsState searchMetrics
Definition variables.h:752
PetscInt IM
Definition variables.h:820
PetscInt64 searchAttempts
Definition variables.h:222
PetscInt64 tieBreakCount
Definition variables.h:230
PetscInt64 maxTraversalFailCount
Definition variables.h:229
Defines the vertices of a single hexahedral grid cell.
Definition variables.h:160
#define DISTANCE_THRESHOLD
#define MAX_TRAVERSAL
PetscErrorCode EvaluateParticlePosition(const Cell *cell, PetscReal *d, const Cmpnts p, PetscInt *position, const PetscReal threshold)
Implementation of EvaluateParticlePosition().
PetscErrorCode ReportSearchOutcome(const Particle *particle, ParticleLocationStatus status, PetscInt traversal_steps)
Implementation of ReportSearchOutcome().
#define REPEAT_COUNT_THRESHOLD
PetscErrorCode UpdateCellIndicesBasedOnDistances(PetscReal d[NUM_FACES], PetscInt *idx, PetscInt *idy, PetscInt *idz)
Internal helper implementation: UpdateCellIndicesBasedOnDistances().
PetscErrorCode FindOwnerOfCell(UserCtx *user, PetscInt i, PetscInt j, PetscInt k, PetscMPIInt *owner_rank)
Internal helper implementation: FindOwnerOfCell().
PetscErrorCode RetrieveCurrentCell(UserCtx *user, PetscInt idx, PetscInt idy, PetscInt idz, Cell *cell)
Implementation of RetrieveCurrentCell().
PetscErrorCode InitializeTraversalParameters(UserCtx *user, Particle *particle, PetscInt *idx, PetscInt *idy, PetscInt *idz, PetscInt *traversal_steps)
Implementation of InitializeTraversalParameters().
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ReportSearchOutcome()

PetscErrorCode ReportSearchOutcome ( const Particle particle,
ParticleLocationStatus  status,
PetscInt  traversal_steps 
)

Implementation of ReportSearchOutcome().

Logs the final outcome of the particle location search.

Full API contract (arguments, ownership, side effects) is documented with the header declaration in include/walkingsearch.h.

See also
ReportSearchOutcome()

Definition at line 1167 of file walkingsearch.c.

1170{
1171 PetscFunctionBeginUser;
1173 switch (status) {
1174 case ACTIVE_AND_LOCATED:
1175 LOG_ALLOW(LOCAL, LOG_INFO, "Search SUCCESS [PID %lld]: Located in global cell (%d, %d, %d) after %d steps.\n",
1176 (long long)particle->PID, particle->cell[0], particle->cell[1], particle->cell[2], traversal_steps);
1177 break;
1178 case MIGRATING_OUT:
1179 LOG_ALLOW(LOCAL, LOG_INFO, "Search SUCCESS [PID %lld]: Identified for migration to Rank %d. Handoff cell is (%d, %d, %d) after %d steps.\n",
1180 (long long)particle->PID, particle->destination_rank, particle->cell[0], particle->cell[1], particle->cell[2], traversal_steps);
1181 break;
1182 case LOST:
1183 LOG_ALLOW(LOCAL, LOG_WARNING, "Search FAILED [PID %lld]: Particle is LOST after %d steps.\n",
1184 (long long)particle->PID, traversal_steps);
1185 break;
1186 default:
1187 LOG_ALLOW(LOCAL, LOG_WARNING, "Search ended with unexpected status %d for PID %lld.\n", status, (long long)particle->PID);
1188 break;
1189 }
1191 PetscFunctionReturn(0);
1192}
Here is the caller graph for this function: