6#ifndef ERROR_MSG_BUFFER_SIZE
7#define ERROR_MSG_BUFFER_SIZE 256
12#define __FUNCT__ "UpdateParticlePosition"
24 PetscFunctionBeginUser;
30 position->
x += velocity->
x * dt;
31 position->
y += velocity->
y * dt;
32 position->
z += velocity->
z * dt;
36 PetscFunctionReturn(0);
40#define __FUNCT__ "UpdateAllParticlePositions"
51 DM swarm = user->
swarm;
53 PetscReal *pos = NULL;
54 PetscReal *vel = NULL;
58 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
60 PetscFunctionBeginUser;
65 ierr = DMSwarmGetLocalSize(swarm, &nLocal); CHKERRQ(ierr);
68 PetscFunctionReturn(0);
71 ierr = DMSwarmGetField(swarm,
"position", NULL, NULL, (
void**)&pos); CHKERRQ(ierr);
72 ierr = DMSwarmGetField(swarm,
"velocity", NULL, NULL, (
void**)&vel); CHKERRQ(ierr);
77 for (p = 0; p < nLocal; p++) {
79 temp_pos.
x = pos[3*p];
80 temp_pos.
y = pos[3*p + 1];
81 temp_pos.
z = pos[3*p + 2];
84 temp_vel.
x = vel[3*p];
85 temp_vel.
y = vel[3*p + 1];
86 temp_vel.
z = vel[3*p + 2];
91 pos[3*p] = temp_pos.
x;
92 pos[3*p + 1] = temp_pos.
y;
93 pos[3*p + 2] = temp_pos.
z;
95 vel[3*p] = temp_vel.
x;
96 vel[3*p + 1] = temp_vel.
y;
97 vel[3*p + 2] = temp_vel.
z;
101 ierr = DMSwarmRestoreField(swarm,
"position", NULL, NULL, (
void**)&pos); CHKERRQ(ierr);
102 ierr = DMSwarmRestoreField(swarm,
"velocity", NULL, NULL, (
void**)&vel); CHKERRQ(ierr);
108 PetscFunctionReturn(0);
128#define __FUNCT__ "CheckAndRemoveOutOfBoundsParticles"
166 PetscInt *removedCountLocal,
167 PetscInt *removedCountGlobal,
171 DM swarm = user->
swarm;
172 PetscInt nLocalInitial;
173 PetscReal *pos_p = NULL;
174 PetscInt64 *pid_p = NULL;
175 PetscInt local_removed_count = 0;
176 PetscMPIInt global_removed_count_mpi = 0;
177 PetscMPIInt rank, size;
179 PetscFunctionBeginUser;
180 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
181 ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size); CHKERRQ(ierr);
185 *removedCountLocal = 0;
186 if (removedCountGlobal) *removedCountGlobal = 0;
188 ierr = DMSwarmGetLocalSize(swarm, &nLocalInitial); CHKERRQ(ierr);
192 if (nLocalInitial > 0) {
194 ierr = DMSwarmGetField(swarm,
"position", NULL, NULL, (
void **)&pos_p); CHKERRQ(ierr);
195 ierr = DMSwarmGetField(swarm,
"DMSwarm_pid", NULL, NULL, (
void **)&pid_p); CHKERRQ(ierr);
198 for (PetscInt p = nLocalInitial - 1; p >= 0; p--) {
199 PetscBool isInsideAnyBox = PETSC_FALSE;
200 Cmpnts current_pos = {pos_p[3*p + 0], pos_p[3*p + 1], pos_p[3*p + 2]};
203 for (PetscMPIInt proc = 0; proc < size; proc++) {
205 isInsideAnyBox = PETSC_TRUE;
210 if (!isInsideAnyBox) {
211 LOG_ALLOW(
LOCAL,
LOG_DEBUG,
"Rank %d: Removing out-of-bounds particle [PID %lld] at local index %d. Pos: (%g, %g, %g)\n",
212 rank, (
long long)pid_p[p], p, current_pos.
x, current_pos.
y, current_pos.
z);
219 ierr = DMSwarmRestoreField(swarm,
"position", NULL, NULL, (
void **)&pos_p); CHKERRQ(ierr);
220 ierr = DMSwarmRestoreField(swarm,
"DMSwarm_pid", NULL, NULL, (
void **)&pid_p); CHKERRQ(ierr);
223 ierr = DMSwarmRemovePointAtIndex(swarm, p); CHKERRQ(ierr);
224 local_removed_count++;
227 PetscInt nLocalCurrent;
228 ierr = DMSwarmGetLocalSize(swarm, &nLocalCurrent); CHKERRQ(ierr);
230 if (nLocalCurrent > 0 && p > 0) {
231 ierr = DMSwarmGetField(swarm,
"position", NULL, NULL, (
void **)&pos_p); CHKERRQ(ierr);
232 ierr = DMSwarmGetField(swarm,
"DMSwarm_pid", NULL, NULL, (
void **)&pid_p); CHKERRQ(ierr);
247 if (pos_p) { ierr = DMSwarmRestoreField(swarm,
"position", NULL, NULL, (
void **)&pos_p); CHKERRQ(ierr); }
248 if (pid_p) { ierr = DMSwarmRestoreField(swarm,
"DMSwarm_pid", NULL, NULL, (
void **)&pid_p); CHKERRQ(ierr); }
251 PetscInt nLocalFinal;
252 ierr = DMSwarmGetLocalSize(swarm, &nLocalFinal); CHKERRQ(ierr);
253 LOG_ALLOW(
LOCAL,
LOG_INFO,
"[Rank %d] Finished removing %d out-of-bounds particles. Final local size: %d.\n", rank, local_removed_count, nLocalFinal);
256 *removedCountLocal = local_removed_count;
257 if (removedCountGlobal) {
258 ierr = MPI_Allreduce(&local_removed_count, &global_removed_count_mpi, 1, MPI_INT, MPI_SUM, PetscObjectComm((PetscObject)swarm)); CHKERRQ(ierr);
259 *removedCountGlobal = global_removed_count_mpi;
264 PetscFunctionReturn(0);
268#define __FUNCT__ "CheckAndRemoveLostParticles"
300 PetscInt *removedCountLocal,
301 PetscInt *removedCountGlobal)
304 DM swarm = user->
swarm;
305 PetscInt nLocalInitial;
306 PetscInt *status_p = NULL;
307 PetscInt64 *pid_p = NULL;
308 PetscReal *pos_p = NULL;
309 PetscInt local_removed_count = 0;
310 PetscMPIInt global_removed_count_mpi = 0;
313 PetscFunctionBeginUser;
315 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
319 *removedCountLocal = 0;
320 if (removedCountGlobal) *removedCountGlobal = 0;
322 ierr = DMSwarmGetLocalSize(swarm, &nLocalInitial); CHKERRQ(ierr);
326 if (nLocalInitial > 0) {
328 ierr = DMSwarmGetField(swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status_p); CHKERRQ(ierr);
329 ierr = DMSwarmGetField(swarm,
"DMSwarm_pid", NULL, NULL, (
void **)&pid_p); CHKERRQ(ierr);
330 ierr = DMSwarmGetField(swarm,
"position", NULL, NULL, (
void **)&pos_p); CHKERRQ(ierr);
333 for (PetscInt p = nLocalInitial - 1; p >= 0; p--) {
334 if (status_p[p] ==
LOST) {
335 LOG_ALLOW(
LOCAL,
LOG_DEBUG,
"Rank %d: Removing LOST particle [PID %lld] at local index %d. Position: (%.4f, %.4f, %.4f).\n",
336 rank, (
long long)pid_p[p], p, pos_p[3*p], pos_p[3*p+1], pos_p[3*p+2]);
343 ierr = DMSwarmRestoreField(swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status_p); CHKERRQ(ierr);
344 ierr = DMSwarmRestoreField(swarm,
"DMSwarm_pid", NULL, NULL, (
void **)&pid_p); CHKERRQ(ierr);
345 ierr = DMSwarmRestoreField(swarm,
"position", NULL, NULL, (
void **)&pos_p); CHKERRQ(ierr);
348 ierr = DMSwarmRemovePointAtIndex(swarm, p); CHKERRQ(ierr);
349 local_removed_count++;
352 PetscInt nLocalCurrent;
353 ierr = DMSwarmGetLocalSize(swarm, &nLocalCurrent); CHKERRQ(ierr);
355 if (nLocalCurrent > 0 && p > 0) {
356 ierr = DMSwarmGetField(swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status_p); CHKERRQ(ierr);
357 ierr = DMSwarmGetField(swarm,
"DMSwarm_pid", NULL, NULL, (
void **)&pid_p); CHKERRQ(ierr);
358 ierr = DMSwarmGetField(swarm,
"position", NULL, NULL, (
void **)&pos_p); CHKERRQ(ierr);
374 if (status_p) { ierr = DMSwarmRestoreField(swarm,
"DMSwarm_location_status", NULL, NULL, (
void **)&status_p); CHKERRQ(ierr); }
375 if (pid_p) { ierr = DMSwarmRestoreField(swarm,
"DMSwarm_pid", NULL, NULL, (
void **)&pid_p); CHKERRQ(ierr); }
376 if (pos_p) { ierr = DMSwarmRestoreField(swarm,
"position", NULL, NULL, (
void **)&pos_p); CHKERRQ(ierr); }
379 PetscInt nLocalFinal;
380 ierr = DMSwarmGetLocalSize(swarm, &nLocalFinal); CHKERRQ(ierr);
381 LOG_ALLOW(
LOCAL,
LOG_INFO,
"Rank %d: Finished removing %d LOST particles. Final local size: %d.\n", rank, local_removed_count, nLocalFinal);
384 *removedCountLocal = local_removed_count;
385 if (removedCountGlobal) {
386 ierr = MPI_Allreduce(&local_removed_count, &global_removed_count_mpi, 1, MPI_INT, MPI_SUM, PetscObjectComm((PetscObject)swarm)); CHKERRQ(ierr);
387 *removedCountGlobal = global_removed_count_mpi;
393 PetscFunctionReturn(0);
398#define __FUNCT__ "SetMigrationRanks"
411 DM swarm = user->
swarm;
413 PetscInt *rankField = NULL;
415 PetscFunctionBeginUser;
419 ierr = DMSwarmGetField(swarm,
"DMSwarm_rank", NULL, NULL, (
void **)&rankField); CHKERRQ(ierr);
422 for(p_idx = 0; p_idx < migrationCount; ++p_idx) {
423 rankField[migrationList[p_idx].
local_index] = migrationList[p_idx].target_rank;
426 ierr = DMSwarmRestoreField(swarm,
"DMSwarm_rank", NULL, NULL, (
void **)&rankField); CHKERRQ(ierr);
429 PetscFunctionReturn(0);
433#define __FUNCT__ "PerformMigration"
448 DM swarm = user->
swarm;
451 PetscFunctionBeginUser;
453 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
458 ierr = DMSwarmMigrate(swarm, PETSC_TRUE); CHKERRQ(ierr);
462 PetscFunctionReturn(0);
470#define __FUNCT__ "CalculateParticleCountPerCell"
488 DM swarm = user->
swarm;
491 PetscInt *global_cell_id_arr;
492 PetscScalar ***count_arr_3d;
496 PetscInt particles_counted_locally = 0;
498 PetscFunctionBeginUser;
500 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
503 if (!da) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL,
"UserCtx->da is NULL.");
504 if (!swarm) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL,
"UserCtx->swarm is NULL.");
505 if (!countVec) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL,
"UserCtx->ParticleCount is NULL.");
508 ierr = DMDAGetInfo(da, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &count_dof, NULL, NULL, NULL, NULL, NULL); CHKERRQ(ierr);
509 if (count_dof != 1) { PetscSNPrintf(msg,
sizeof(msg),
"countDM must have DOF=1, got %d.", count_dof); SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, msg); }
512 ierr = VecSet(countVec, 0.0); CHKERRQ(ierr);
516 ierr = DMSwarmGetLocalSize(swarm, &nlocal); CHKERRQ(ierr);
517 ierr = DMSwarmGetField(swarm,
"DMSwarm_CellID", NULL, NULL, (
void **)&global_cell_id_arr); CHKERRQ(ierr);
518 ierr = DMSwarmGetField(swarm,
"DMSwarm_pid",NULL,NULL,(
void **)&PID_arr);CHKERRQ(ierr);
522 ierr = DMDAVecGetArray(da, countVec, &count_arr_3d); CHKERRQ(ierr);
525 PetscInt xs, ys, zs, xm, ym, zm;
526 ierr = DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm); CHKERRQ(ierr);
529 LOG_ALLOW(
LOCAL,
LOG_DEBUG,
"CalculateParticleCountPerCell (Rank %d): Processing %d local particles using GLOBAL CellIDs.\n",rank,nlocal);
530 for (p = 0; p < nlocal; p++) {
532 PetscInt64 i_geom = global_cell_id_arr[p * 3 + 0];
533 PetscInt64 j_geom = global_cell_id_arr[p * 3 + 1];
534 PetscInt64 k_geom = global_cell_id_arr[p * 3 + 2];
537 PetscInt i = (PetscInt)i_geom + 1;
538 PetscInt j = (PetscInt)j_geom + 1;
539 PetscInt k = (PetscInt)k_geom + 1;
546 LOG_LOOP_ALLOW(
LOCAL,
LOG_VERBOSE,p,100,
"[Rank %d] Read CellID for p=%d, PID = %ld: (%ld, %ld, %ld)\n", rank, p,PID_arr[p],i, j, k);
556 if (i >= xs && i < xs + xm &&
557 j >= ys && j < ys + ym &&
558 k >= zs && k < zs + zm )
563 count_arr_3d[k][j][i] += 1.0;
564 particles_counted_locally++;
567 LOG_ALLOW(
LOCAL,
LOG_DEBUG,
"CalculateParticleCountPerCell (Rank %d): Skipping particle %ld with global CellID (%ld, %ld, %ld) - likely outside local+ghost range.\n",rank, PID_arr[p] , i, j, k);
570 LOG_ALLOW(
LOCAL,
LOG_DEBUG,
"CalculateParticleCountPerCell (Rank %d): Local counting finished. Processed %d particles locally.\n", rank, particles_counted_locally);
573 ierr = DMDAVecRestoreArray(da, countVec, &count_arr_3d); CHKERRQ(ierr);
574 ierr = DMSwarmRestoreField(swarm,
"DMSwarm_CellID", NULL, NULL, (
void **)&global_cell_id_arr); CHKERRQ(ierr);
575 ierr = DMSwarmRestoreField(swarm,
"DMSwarm_pid",NULL,NULL,(
void **)&PID_arr);CHKERRQ(ierr);
579 ierr = VecAssemblyBegin(countVec); CHKERRQ(ierr);
580 ierr = VecAssemblyEnd(countVec); CHKERRQ(ierr);
583 PetscReal total_counted_particles = 0.0, max_count_in_cell = 0.0;
584 ierr = VecSum(countVec, &total_counted_particles); CHKERRQ(ierr);
585 PetscInt max_idx_global = -1;
586 ierr = VecMax(countVec, &max_idx_global, &max_count_in_cell); CHKERRQ(ierr);
588 total_counted_particles, max_count_in_cell);
591 if (max_idx_global >= 0) {
595 ierr = DMDAGetInfo(da, NULL, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL); CHKERRQ(ierr);
598 PetscInt Kmax = max_idx_global / (M * N);
599 PetscInt Jmax = (max_idx_global % (M * N)) / M;
600 PetscInt Imax = max_idx_global % M;
602 (
int)Imax, (
int)Jmax, (
int)Kmax, (
int)max_idx_global);
605 PetscScalar count_at_origin = 0.0;
606 PetscScalar ***count_arr_for_check;
607 ierr = DMDAVecGetArrayRead(da, countVec, &count_arr_for_check); CHKERRQ(ierr);
609 PetscInt xs, ys, zs, xm, ym, zm;
610 ierr = DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm); CHKERRQ(ierr);
611 if (0 >= xs && 0 < xs+xm && 0 >= ys && 0 < ys+ym && 0 >= zs && 0 < zs+zm) {
612 count_at_origin = count_arr_for_check[0][0][0];
615 count_at_origin = -999.0;
617 ierr = DMDAVecRestoreArrayRead(da, countVec, &count_arr_for_check); CHKERRQ(ierr);
629 PetscFunctionReturn(0);
635#define __FUNCT__ "ResizeSwarmGlobally"
641 PetscInt N_current, nlocal_current;
645 PetscFunctionBeginUser;
647 ierr = PetscObjectGetComm((PetscObject)swarm, &comm); CHKERRQ(ierr);
648 ierr = MPI_Comm_rank(comm, &rank); CHKERRQ(ierr);
649 ierr = DMSwarmGetSize(swarm, &N_current); CHKERRQ(ierr);
650 ierr = DMSwarmGetLocalSize(swarm, &nlocal_current); CHKERRQ(ierr);
652 PetscInt delta = N_target - N_current;
656 PetscFunctionReturn(0);
660 PetscInt num_to_remove_global = -delta;
661 LOG_ALLOW(
GLOBAL,
LOG_INFO,
"Current size %d > target size %d. Removing %d particles globally.\n", N_current, N_target, num_to_remove_global);
671 ierr = MPI_Exscan(&nlocal_current, &rstart, 1, MPIU_INT, MPI_SUM, comm); CHKERRMPI(ierr);
673 rend = rstart + nlocal_current;
677 PetscInt nlocal_remove_count = 0;
678 if (rend > N_target) {
679 PetscInt start_remove_local_idx = (N_target > rstart) ? (N_target - rstart) : 0;
680 nlocal_remove_count = nlocal_current - start_remove_local_idx;
683 if (nlocal_remove_count < 0) nlocal_remove_count = 0;
684 if (nlocal_remove_count > nlocal_current) nlocal_remove_count = nlocal_current;
686 LOG_ALLOW(
LOCAL,
LOG_DEBUG,
"Rank %d: Global range [%d, %d). Target size %d. Need to remove %d local particles (from end).\n", rank, rstart, rend, N_target, nlocal_remove_count);
689 PetscInt removal_ops_done = 0;
690 for (PetscInt p = nlocal_current - 1; p >= 0 && removal_ops_done < nlocal_remove_count; --p) {
691 ierr = DMSwarmRemovePointAtIndex(swarm, p); CHKERRQ(ierr);
695 if (removal_ops_done != nlocal_remove_count) {
696 SETERRQ(comm, PETSC_ERR_PLIB,
"Rank %d: Failed to remove the expected number of local particles (%d != %d)", rank, removal_ops_done, nlocal_remove_count);
701 ierr = MPI_Barrier(comm); CHKERRMPI(ierr);
704 PetscInt num_to_add_global = delta;
705 LOG_ALLOW(
GLOBAL,
LOG_INFO,
"Current size %d < target size %d. Adding %d particles globally.\n", N_current, N_target, num_to_add_global);
706 ierr = DMSwarmAddNPoints(swarm, num_to_add_global); CHKERRQ(ierr);
712 ierr = DMSwarmGetSize(swarm, &N_final); CHKERRQ(ierr);
713 if (N_final != N_target) {
714 SETERRQ(comm, PETSC_ERR_PLIB,
"Failed to resize swarm: expected %d particles, got %d", N_target, N_final);
718 PetscFunctionReturn(0);
722#define __FUNCT__ "PreCheckAndResizeSwarm"
742 char filename[PETSC_MAX_PATH_LEN];
744 PetscInt N_current = 0;
747 const char *refFieldName =
"position";
748 const PetscInt bs = 3;
756 const int placeholder_int = 0;
761 PetscFunctionBeginUser;
763 ierr = PetscObjectGetComm((PetscObject)user->
swarm, &comm); CHKERRQ(ierr);
764 ierr = MPI_Comm_rank(comm, &rank); CHKERRQ(ierr);
772 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONGSTATE,
"Invalid execution mode for reading simulation fields.");
782 ierr = PetscSNPrintf(filename,
sizeof(filename),
"%s/%s%05" PetscInt_FMT
"_%d.%s",
789 PetscBool fileExists = PETSC_FALSE;
790 ierr = PetscTestFile(filename,
'r', &fileExists); CHKERRQ(ierr);
801 ierr = VecCreate(PETSC_COMM_SELF, &tmpVec); CHKERRQ(ierr);
802 ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF, filename, FILE_MODE_READ, &viewer); CHKERRQ(ierr);
803 ierr = VecLoad(tmpVec, viewer); CHKERRQ(ierr);
804 ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
806 ierr = VecGetSize(tmpVec, &vecSize); CHKERRQ(ierr);
807 ierr = VecDestroy(&tmpVec); CHKERRQ(ierr);
809 if (vecSize % bs != 0) {
811 LOG_ALLOW(
GLOBAL,
LOG_ERROR,
"Rank 0: Vector size %d from file '%s' is not divisible by block size %d.\n", vecSize, filename, bs);
813 N_file = vecSize / bs;
820 ierr = MPI_Bcast(&N_file, 1, MPIU_INT, 0, comm); CHKERRMPI(ierr);
824 SETERRQ(comm, PETSC_ERR_FILE_OPEN,
"Mandatory reference file '%s' not found for timestep %d (as determined by Rank 0).", filename, ti);
827 SETERRQ(comm, PETSC_ERR_FILE_READ,
"Reference file '%s' has incorrect format (as determined by Rank 0).", filename);
830 SETERRQ(comm, PETSC_ERR_PLIB,
"Received invalid particle count %d from Rank 0.", N_file);
835 ierr = DMSwarmGetSize(user->
swarm, &N_current); CHKERRQ(ierr);
837 if (N_file != N_current) {
838 LOG_ALLOW(
GLOBAL,
LOG_INFO,
"Swarm size %d differs from file size %d. Resizing swarm globally.\n", N_current, N_file);
848 PetscFunctionReturn(0);
853#define __FUNCT__ "ReinitializeParticlesOnInletSurface"
872 DM swarm = user->
swarm;
873 PetscReal *positions_field = NULL;
874 PetscInt64 *particleIDs = NULL;
875 PetscInt *cell_ID_field = NULL;
876 const Cmpnts ***coor_nodes_local_array;
879 PetscInt xs_gnode_rank, ys_gnode_rank, zs_gnode_rank;
880 PetscInt IM_nodes_global, JM_nodes_global, KM_nodes_global;
882 PetscRandom rand_logic_reinit_i, rand_logic_reinit_j, rand_logic_reinit_k;
883 PetscInt nlocal_current;
884 PetscInt particles_actually_reinitialized_count = 0;
885 PetscBool can_this_rank_service_inlet = PETSC_FALSE;
887 PetscFunctionBeginUser;
893 PetscFunctionReturn(0);
896 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
897 ierr = DMSwarmGetLocalSize(swarm, &nlocal_current); CHKERRQ(ierr);
900 if (nlocal_current == 0) {
901 LOG_ALLOW(
LOCAL,
LOG_DEBUG,
"[T=%.4f, Step=%d] Rank %d has no local particles to re-initialize on inlet.\n", currentTime, step, rank);
902 PetscFunctionReturn(0);
906 ierr = DMDAGetLocalInfo(user->
da, &info); CHKERRQ(ierr);
907 ierr = DMDAGetInfo(user->
da, NULL, &IM_nodes_global, &JM_nodes_global, &KM_nodes_global, NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL); CHKERRQ(ierr);
908 ierr = DMDAGetGhostCorners(user->
da, &xs_gnode_rank, &ys_gnode_rank, &zs_gnode_rank, NULL, NULL, NULL); CHKERRQ(ierr);
911 IM_nodes_global -= 1; JM_nodes_global -= 1; KM_nodes_global -= 1;
913 const PetscInt IM_cells_global = IM_nodes_global > 0 ? IM_nodes_global - 1 : 0;
914 const PetscInt JM_cells_global = JM_nodes_global > 0 ? JM_nodes_global - 1 : 0;
915 const PetscInt KM_cells_global = KM_nodes_global > 0 ? KM_nodes_global - 1 : 0;
920 ierr =
CanRankServiceInletFace(user, &info, IM_nodes_global, JM_nodes_global, KM_nodes_global, &can_this_rank_service_inlet); CHKERRQ(ierr);
923 ierr = DMGetCoordinatesLocal(user->
da, &Coor_local); CHKERRQ(ierr);
924 ierr = DMDAVecGetArrayRead(user->
fda, Coor_local, (
void*)&coor_nodes_local_array); CHKERRQ(ierr);
925 ierr = DMSwarmGetField(swarm,
"position", NULL, NULL, (
void**)&positions_field); CHKERRQ(ierr);
926 ierr = DMSwarmGetField(swarm,
"DMSwarm_pid", NULL, NULL, (
void**)&particleIDs); CHKERRQ(ierr);
927 ierr = DMSwarmGetField(swarm,
"DMSwarm_CellID", NULL, NULL, (
void**)&cell_ID_field); CHKERRQ(ierr);
929 if (!can_this_rank_service_inlet) {
933 LOG_ALLOW(
LOCAL,
LOG_DEBUG,
"[T=%.4f, Step=%d] Rank %d is resetting %d local particles to inlet center (%.6f, %.6f, %.6f) for migration.\n", currentTime, step, rank, nlocal_current, user->
simCtx->
CMx_c, user->
simCtx->
CMy_c, user->
simCtx->
CMz_c);
935 for(PetscInt p = 0; p < nlocal_current; p++){
940 cell_ID_field[3*p+0] = -1;
941 cell_ID_field[3*p+1] = -1;
942 cell_ID_field[3*p+2] = -1;
946 ierr = DMDAVecRestoreArrayRead(user->
fda, Coor_local, (
void*)&coor_nodes_local_array); CHKERRQ(ierr);
947 ierr = DMSwarmRestoreField(swarm,
"position", NULL, NULL, (
void**)&positions_field); CHKERRQ(ierr);
948 ierr = DMSwarmRestoreField(swarm,
"DMSwarm_pid", NULL, NULL, (
void**)&particleIDs); CHKERRQ(ierr);
949 ierr = DMSwarmRestoreField(swarm,
"DMSwarm_CellID", NULL, NULL, (
void**)&cell_ID_field); CHKERRQ(ierr);
951 PetscFunctionReturn(0);
962 for (PetscInt p = 0; p < nlocal_current; p++) {
963 PetscInt ci_metric_lnode, cj_metric_lnode, ck_metric_lnode;
964 PetscReal xi_metric_logic, eta_metric_logic, zta_metric_logic;
970 IM_nodes_global, JM_nodes_global, KM_nodes_global,
971 &rand_logic_reinit_i, &rand_logic_reinit_j, &rand_logic_reinit_k,
972 &ci_metric_lnode, &cj_metric_lnode, &ck_metric_lnode,
973 &xi_metric_logic, &eta_metric_logic, &zta_metric_logic); CHKERRQ(ierr);
975 PetscBool garbage_flag = PETSC_FALSE;
977 IM_cells_global, JM_cells_global, KM_cells_global,
979 &ci_metric_lnode, &cj_metric_lnode, &ck_metric_lnode,
980 &xi_metric_logic, &eta_metric_logic, &zta_metric_logic,&garbage_flag); CHKERRQ(ierr);
982 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE,
"ReinitializeParticlesOnInletSurface only supports ParticleInitialization modes 0 and 3.");
986 ci_metric_lnode, cj_metric_lnode, ck_metric_lnode,
987 xi_metric_logic, eta_metric_logic, zta_metric_logic,
988 &phys_coords); CHKERRQ(ierr);
991 positions_field[3*p+0] = phys_coords.
x;
992 positions_field[3*p+1] = phys_coords.
y;
993 positions_field[3*p+2] = phys_coords.
z;
995 particles_actually_reinitialized_count++;
997 cell_ID_field[3*p+0] = -1;
998 cell_ID_field[3*p+1] = -1;
999 cell_ID_field[3*p+2] = -1;
1002 "Rank %d: PID %ld (idx %ld) RE-PLACED. CellOriginNode(locDAIdx):(%d,%d,%d). LogicCoords: (%.2e,%.2f,%.2f). PhysCoords: (%.6f,%.6f,%.6f).\n",
1003 rank, particleIDs[p], (
long)p,
1004 ci_metric_lnode, cj_metric_lnode, ck_metric_lnode,
1005 xi_metric_logic, eta_metric_logic, zta_metric_logic,
1006 phys_coords.
x, phys_coords.
y, phys_coords.
z);
1010 if (particles_actually_reinitialized_count > 0) {
1011 LOG_ALLOW(
GLOBAL,
LOG_INFO,
"[T=%.4f, Step=%d] Rank %d (on inlet face %d) successfully re-initialized %d of %d local particles.\n", currentTime, step, rank, user->
identifiedInletBCFace, particles_actually_reinitialized_count, nlocal_current);
1012 }
else if (nlocal_current > 0) {
1013 LOG_ALLOW(
GLOBAL,
LOG_WARNING,
"[T=%.4f, Step=%d] Rank %d claimed to service inlet face %d, but re-initialized 0 of %d local particles. This may indicate an issue if particles were expected to be re-placed.\n", currentTime, step, rank, user->
identifiedInletBCFace, nlocal_current);
1017 ierr = PetscRandomDestroy(&rand_logic_reinit_i); CHKERRQ(ierr);
1018 ierr = PetscRandomDestroy(&rand_logic_reinit_j); CHKERRQ(ierr);
1019 ierr = PetscRandomDestroy(&rand_logic_reinit_k); CHKERRQ(ierr);
1021 ierr = DMSwarmRestoreField(swarm,
"position", NULL, NULL, (
void**)&positions_field); CHKERRQ(ierr);
1022 ierr = DMSwarmRestoreField(swarm,
"DMSwarm_pid", NULL, NULL, (
void**)&particleIDs); CHKERRQ(ierr);
1023 ierr = DMSwarmRestoreField(swarm,
"DMSwarm_CellID", NULL, NULL, (
void**)&cell_ID_field); CHKERRQ(ierr);
1024 ierr = DMDAVecRestoreArrayRead(user->
fda, Coor_local, (
void*)&coor_nodes_local_array); CHKERRQ(ierr);
1028 PetscFunctionReturn(0);
1033 PetscInt64 val_a = *(
const PetscInt64*)a;
1034 PetscInt64 val_b = *(
const PetscInt64*)b;
1035 if (val_a < val_b)
return -1;
1036 if (val_a > val_b)
return 1;
1041#define __FUNCT__ "GetLocalPIDSnapshot"
1065 PetscInt64 **pids_snapshot_out)
1067 PetscErrorCode ierr;
1070 PetscFunctionBeginUser;
1075 if (!pids_snapshot_out) {
1076 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL,
"Output pointer pids_snapshot_out is NULL.");
1079 if (n_local > 0 && !pid_field) {
1080 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL,
"Input pid_field pointer is NULL for n_local > 0.");
1083 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
1088 *pids_snapshot_out = NULL;
1091 PetscFunctionReturn(0);
1095 ierr = PetscMalloc1(n_local, pids_snapshot_out); CHKERRQ(ierr);
1099 ierr = PetscMemcpy(*pids_snapshot_out, pid_field, n_local *
sizeof(PetscInt64)); CHKERRQ(ierr);
1104 ierr = PetscSortInt64(n_local, *pids_snapshot_out); CHKERRQ(ierr);
1109 PetscFunctionReturn(0);
1113#define __FUNCT__ "AddToMigrationList"
1136 PetscInt *capacity_p,
1138 PetscInt particle_local_idx,
1139 PetscMPIInt destination_rank)
1141 PetscErrorCode ierr;
1144 PetscFunctionBeginUser;
1149 if (!migration_list_p || !capacity_p || !count_p) {
1150 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL,
"Null pointer provided to AddToMigrationList for list management.");
1154 if (*count_p >= *capacity_p) {
1155 PetscInt old_capacity = *capacity_p;
1157 PetscInt new_capacity = (old_capacity == 0) ? 16 : old_capacity * 2;
1162 ierr = PetscRealloc(new_capacity *
sizeof(
MigrationInfo), migration_list_p); CHKERRQ(ierr);
1164 *capacity_p = new_capacity;
1166 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
1168 rank, old_capacity, new_capacity);
1175 list[*count_p].local_index = particle_local_idx;
1176 list[*count_p].target_rank = destination_rank;
1183 PetscFunctionReturn(0);
1188#define __FUNCT__ "FlagNewComersForLocation"
1218 PetscInt n_local_before,
1219 const PetscInt64 pids_before[])
1221 PetscErrorCode ierr;
1223 PetscInt n_local_after;
1224 PetscInt newcomer_count = 0;
1227 PetscInt64 *pid_field_after = NULL;
1228 PetscInt *status_field_after = NULL;
1229 PetscInt *cell_field_after = NULL;
1231 PetscFunctionBeginUser;
1237 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE,
"Input DMSwarm is NULL in FlagNewcomersForLocation.");
1240 if (n_local_before > 0 && !pids_before) {
1241 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL,
"Input pids_before array is NULL for n_local_before > 0.");
1244 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
1247 ierr = DMSwarmGetLocalSize(swarm, &n_local_after); CHKERRQ(ierr);
1250 rank, n_local_before, n_local_after);
1253 if (n_local_after == 0) {
1254 PetscFunctionReturn(0);
1259 ierr = DMSwarmGetField(swarm,
"DMSwarm_pid", NULL, NULL, (
void**)&pid_field_after); CHKERRQ(ierr);
1260 ierr = DMSwarmGetField(swarm,
"DMSwarm_location_status", NULL, NULL, (
void**)&status_field_after); CHKERRQ(ierr);
1261 ierr = DMSwarmGetField(swarm,
"DMSwarm_CellID", NULL, NULL, (
void**)&cell_field_after); CHKERRQ(ierr);
1262 if (!pid_field_after || !status_field_after) {
1263 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE,
"Failed to get required swarm fields in FlagNewcomersForLocation.");
1268 for (PetscInt p_idx = 0; p_idx < n_local_after; ++p_idx) {
1269 PetscInt64 current_pid = pid_field_after[p_idx];
1270 PetscBool is_found_in_before_list;
1273 ierr =
BinarySearchInt64(n_local_before, pids_before, current_pid, &is_found_in_before_list); CHKERRQ(ierr);
1276 if (!is_found_in_before_list) {
1285 rank, current_pid, p_idx);
1291 ierr = DMSwarmRestoreField(swarm,
"DMSwarm_pid", NULL, NULL, (
void**)&pid_field_after); CHKERRQ(ierr);
1292 ierr = DMSwarmRestoreField(swarm,
"DMSwarm_location_status", NULL, NULL, (
void**)&status_field_after); CHKERRQ(ierr);
1293 ierr = DMSwarmRestoreField(swarm,
"DMSwarm_CellID", NULL, NULL, (
void**)&cell_field_after); CHKERRQ(ierr);
1295 if (newcomer_count > 0) {
1301 PetscFunctionReturn(0);
1305#define __FUNCT__ "GuessParticleOwnerWithBBox"
1333 PetscMPIInt *guess_rank_out)
1335 PetscErrorCode ierr;
1336 PetscMPIInt rank, size;
1340 PetscFunctionBeginUser;
1345 if (!user || !particle || !guess_rank_out || !bboxlist) {
1346 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL,
"Null pointer provided to GuessParticleOwnerWithBBox.");
1348 if (!localBBox|| !neighbors) {
1349 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE,
"Required user->bboxl or user->neighbors is not initialized.");
1352 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
1353 ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size); CHKERRQ(ierr);
1355 *guess_rank_out = MPI_PROC_NULL;
1363 *guess_rank_out = rank;
1364 LOG_ALLOW(
LOCAL,
LOG_DEBUG,
"[PID %ld]: Fast path guess SUCCESS. Particle is within the local (Rank %d) bounding box.\n",
1365 particle->
PID, rank);
1368 PetscFunctionReturn(0);
1381 *guess_rank_out = neighbors->
rank_xm;
1383 *guess_rank_out = neighbors->
rank_xp;
1385 *guess_rank_out = neighbors->
rank_ym;
1387 *guess_rank_out = neighbors->
rank_yp;
1389 *guess_rank_out = neighbors->
rank_zm;
1391 *guess_rank_out = neighbors->
rank_zp;
1395 if (*guess_rank_out != MPI_PROC_NULL) {
1397 particle->
PID, *guess_rank_out);
1400 PetscFunctionReturn(0);
1408 for (PetscMPIInt r = 0; r < size; ++r) {
1409 if (r == rank)
continue;
1412 PetscBool is_in = PETSC_TRUE;
1414 LOG_ALLOW(
LOCAL,
LOG_VERBOSE,
"[Rank %d] Checking PID %lld at (%.4f, %.4f, %.4f) against Rank %d's box: [(%.4f, %.4f, %.4f) to (%.4f, %.4f, %.4f)]. Result: %s\n",
1415 (
int)rank, (
long long)particle->
PID,
1420 is_in ?
"INSIDE" :
"OUTSIDE");
1422 *guess_rank_out = r;
1424 particle->
PID, *guess_rank_out);
1427 PetscFunctionReturn(0);
1437 PetscFunctionReturn(0);
1441#define __FUNCT__ "LocateAllParticlesInGrid"
1476 PetscErrorCode ierr;
1477 PetscInt passes = 0;
1478 const PetscInt MAX_MIGRATION_PASSES = 10;
1479 PetscInt global_migrations_this_pass;
1482 PetscFunctionBeginUser;
1484 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank); CHKERRQ(ierr);
1485 LOG_ALLOW(
GLOBAL,
LOG_INFO,
"LocateAllParticlesInGrid (Orchestrator) - Beginning particle settlement process.\n");
1495 PetscInt local_migration_count = 0;
1496 PetscInt migrationListCapacity = 0;
1497 PetscInt nlocal_before;
1498 PetscInt64 *pids_before_snapshot = NULL;
1499 PetscInt local_lost_count = 0;
1501 ierr = DMSwarmGetLocalSize(user->
swarm, &nlocal_before); CHKERRQ(ierr);
1502 LOG_ALLOW(
LOCAL,
LOG_DEBUG,
"[Rank %d] Pass %d begins with %d local particles.\n", rank, passes, nlocal_before);
1506 if (nlocal_before > 0) {
1508 PetscReal *pos_p, *weights_p, *vel_p;
1509 PetscInt *cell_p, *status_p;
1511 ierr = DMSwarmGetField(user->
swarm,
"position", NULL, NULL, (
void**)&pos_p); CHKERRQ(ierr);
1512 ierr = DMSwarmGetField(user->
swarm,
"velocity", NULL, NULL, (
void**)&vel_p); CHKERRQ(ierr);
1513 ierr = DMSwarmGetField(user->
swarm,
"weight", NULL, NULL, (
void**)&weights_p); CHKERRQ(ierr);
1514 ierr = DMSwarmGetField(user->
swarm,
"DMSwarm_CellID", NULL, NULL, (
void**)&cell_p); CHKERRQ(ierr);
1515 ierr = DMSwarmGetField(user->
swarm,
"DMSwarm_pid", NULL, NULL, (
void**)&pid_p); CHKERRQ(ierr);
1516 ierr = DMSwarmGetField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void**)&status_p); CHKERRQ(ierr);
1522 for (PetscInt p_idx = 0; p_idx < nlocal_before; p_idx++) {
1527 "Local Particle idx=%d, PID=%ld, status=%s, cell=(%d, %d, %d)\n",
1536 LOG_ALLOW(
LOCAL,
LOG_VERBOSE,
" [rank %d][PID %ld] skipped in pass %d as it is already located at (%d,%d,%d).\n",rank,pid_p[p_idx],passes,cell_p[3*p_idx],cell_p[3*p_idx + 1],cell_p[3*p_idx + 2]);
1545 ierr =
UnpackSwarmFields(p_idx, pid_p, weights_p, pos_p, cell_p, vel_p, status_p, ¤t_particle); CHKERRQ(ierr);
1554 if (current_particle.
cell[0] >= 0) {
1597 PetscMPIInt guessed_owner_rank = MPI_PROC_NULL;
1601 if (guessed_owner_rank != MPI_PROC_NULL && guessed_owner_rank != rank) {
1602 LOG_ALLOW(
LOCAL,
LOG_VERBOSE,
"[PID %ld] Guess SUCCESS: Found migration target Rank %d. Finalizing.\n", current_particle.
PID, guessed_owner_rank);
1610 if (guessed_owner_rank == rank) {
1611 LOG_ALLOW(
LOCAL,
LOG_DEBUG,
"[PID %ld] Guess determined particle is local. Proceeding to robust walk to find cell.\n", current_particle.
PID);
1613 LOG_ALLOW(
LOCAL,
LOG_WARNING,
"[PID %ld] Guess FAILED to find an owner. Proceeding to robust walk for definitive search.\n", current_particle.
PID);
1622 ierr =
UpdateSwarmFields(p_idx, ¤t_particle, weights_p, cell_p, status_p); CHKERRQ(ierr);
1626 }
else if (final_status ==
LOST) {
1633 ierr = DMSwarmRestoreField(user->
swarm,
"position", NULL, NULL, (
void**)&pos_p); CHKERRQ(ierr);
1634 ierr = DMSwarmRestoreField(user->
swarm,
"velocity", NULL, NULL, (
void**)&vel_p); CHKERRQ(ierr);
1635 ierr = DMSwarmRestoreField(user->
swarm,
"weight", NULL, NULL, (
void**)&weights_p); CHKERRQ(ierr);
1636 ierr = DMSwarmRestoreField(user->
swarm,
"DMSwarm_CellID", NULL, NULL, (
void**)&cell_p); CHKERRQ(ierr);
1637 ierr = DMSwarmRestoreField(user->
swarm,
"DMSwarm_pid", NULL, NULL, (
void**)&pid_p); CHKERRQ(ierr);
1638 ierr = DMSwarmRestoreField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void**)&status_p); CHKERRQ(ierr);
1642 LOG_ALLOW(
LOCAL,
LOG_INFO,
"[Rank %d] Pass %d: Identified %d particles to migrate out.\n", rank, passes, local_migration_count);
1647 ierr = MPI_Allreduce(&local_migration_count, &global_migrations_this_pass, 1, MPIU_INT, MPI_SUM, PETSC_COMM_WORLD); CHKERRQ(ierr);
1649 if(global_migrations_this_pass > 0 ){
1651 LOG_ALLOW(
GLOBAL,
LOG_INFO,
"Pass %d: Migrating %d particles globally.\n", passes, global_migrations_this_pass);
1653 ierr =
SetMigrationRanks(user, migrationList, local_migration_count); CHKERRQ(ierr);
1663 ierr = PetscFree(pids_before_snapshot);
1664 ierr = PetscFree(migrationList);
1666 LOG_ALLOW(
GLOBAL,
LOG_INFO,
"End of pass %d. Total particles migrated globally: %d.\n", passes, global_migrations_this_pass);
1668 }
while (global_migrations_this_pass > 0 && passes < MAX_MIGRATION_PASSES);
1671 if (passes >= MAX_MIGRATION_PASSES) {
1672 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED,
"Particle migration failed to converge after %d passes. Check for particles oscillating between ranks.", MAX_MIGRATION_PASSES);
1678 PetscFunctionReturn(0);
1683#define __FUNCT__ "ResetAllParticleStatuses"
1698 PetscErrorCode ierr;
1702 PetscFunctionBeginUser;
1706 ierr = DMSwarmGetLocalSize(user->
swarm, &n_local); CHKERRQ(ierr);
1710 ierr = DMSwarmGetField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void**)&status_p); CHKERRQ(ierr);
1712 for (PetscInt p = 0; p < n_local; ++p) {
1721 ierr = DMSwarmRestoreField(user->
swarm,
"DMSwarm_location_status", NULL, NULL, (
void**)&status_p); CHKERRQ(ierr);
1726 PetscFunctionReturn(0);
PetscErrorCode GetRandomCellAndLogicalCoordsOnInletFace(UserCtx *user, const DMDALocalInfo *info, PetscInt xs_gnode_rank, PetscInt ys_gnode_rank, PetscInt zs_gnode_rank, PetscInt IM_nodes_global, PetscInt JM_nodes_global, PetscInt KM_nodes_global, PetscRandom *rand_logic_i_ptr, PetscRandom *rand_logic_j_ptr, PetscRandom *rand_logic_k_ptr, PetscInt *ci_metric_lnode_out, PetscInt *cj_metric_lnode_out, PetscInt *ck_metric_lnode_out, PetscReal *xi_metric_logic_out, PetscReal *eta_metric_logic_out, PetscReal *zta_metric_logic_out)
Assuming the current rank services the inlet face, this function selects a random cell (owned by this...
PetscErrorCode CanRankServiceInletFace(UserCtx *user, const DMDALocalInfo *info, PetscInt IM_nodes_global, PetscInt JM_nodes_global, PetscInt KM_nodes_global, PetscBool *can_service_inlet_out)
Determines if the current MPI rank owns any part of the globally defined inlet face,...
PetscErrorCode GetDeterministicFaceGridLocation(UserCtx *user, const DMDALocalInfo *info, PetscInt xs_gnode_rank, PetscInt ys_gnode_rank, PetscInt zs_gnode_rank, PetscInt IM_cells_global, PetscInt JM_cells_global, PetscInt KM_cells_global, PetscInt64 particle_global_id, PetscInt *ci_metric_lnode_out, PetscInt *cj_metric_lnode_out, PetscInt *ck_metric_lnode_out, PetscReal *xi_metric_logic_out, PetscReal *eta_metric_logic_out, PetscReal *zta_metric_logic_out, PetscBool *placement_successful_out)
Places particles in a deterministic grid/raster pattern on a specified domain face.
PetscErrorCode MetricLogicalToPhysical(UserCtx *user, const Cmpnts ***X, PetscInt i, PetscInt j, PetscInt k, PetscReal xi, PetscReal eta, PetscReal zta, Cmpnts *Xp)
Public wrapper: map (cell index, ξ,η,ζ) to (x,y,z).
PetscErrorCode UpdateParticlePosition(UserCtx *user, Cmpnts *position, const Cmpnts *velocity)
Updates a particle's position based on its velocity and the timestep dt (stored in user->dt).
static int compare_PetscInt64(const void *a, const void *b)
PetscErrorCode GuessParticleOwnerWithBBox(UserCtx *user, const Particle *particle, const BoundingBox *bboxlist, PetscMPIInt *guess_rank_out)
Provides a fast, heuristic-based guess for a particle's owner rank using bounding boxes.
PetscErrorCode ResizeSwarmGlobally(DM swarm, PetscInt N_target)
PetscErrorCode AddToMigrationList(MigrationInfo **migration_list_p, PetscInt *capacity_p, PetscInt *count_p, PetscInt particle_local_idx, PetscMPIInt destination_rank)
Safely adds a new migration task to a dynamically sized list.
PetscErrorCode SetMigrationRanks(UserCtx *user, const MigrationInfo *migrationList, PetscInt migrationCount)
Sets the target rank field (DMSwarmPICField_rank) for particles scheduled for migration.
PetscErrorCode CheckAndRemoveOutOfBoundsParticles(UserCtx *user, PetscInt *removedCountLocal, PetscInt *removedCountGlobal, const BoundingBox *bboxlist)
Checks for particles outside the global physical domain and removes them.
PetscErrorCode GetLocalPIDSnapshot(const PetscInt64 pid_field[], PetscInt n_local, PetscInt64 **pids_snapshot_out)
Creates a sorted snapshot of all Particle IDs (PIDs) from a raw data array.
PetscErrorCode UpdateAllParticlePositions(UserCtx *user)
Loops over all local particles in the DMSwarm, updating their positions based on velocity and the glo...
PetscErrorCode LocateAllParticlesInGrid(UserCtx *user, BoundingBox *bboxlist)
Orchestrates the complete particle location and migration process for one timestep.
#define ERROR_MSG_BUFFER_SIZE
PetscErrorCode ResetAllParticleStatuses(UserCtx *user)
This function is designed to be called at the end of a full timestep, after all particle-based calcul...
PetscErrorCode ReinitializeParticlesOnInletSurface(UserCtx *user, PetscReal currentTime, PetscInt step)
Re-initializes the positions of particles currently on this rank if this rank owns part of the design...
PetscErrorCode CheckAndRemoveLostParticles(UserCtx *user, PetscInt *removedCountLocal, PetscInt *removedCountGlobal)
Removes particles that have been definitively flagged as LOST by the location algorithm.
PetscErrorCode FlagNewcomersForLocation(DM swarm, PetscInt n_local_before, const PetscInt64 pids_before[])
Identifies newly arrived particles after migration and flags them for a location search.
static PetscBool IsParticleInBox(const BoundingBox *bbox, const Cmpnts *pos)
Checks if a particle position is within the bounds of a given bounding box.
PetscErrorCode PreCheckAndResizeSwarm(UserCtx *user, PetscInt ti, const char *ext)
Checks particle count from a saved file and resizes the swarm globally.
PetscErrorCode PerformMigration(UserCtx *user)
Performs particle migration based on the pre-populated DMSwarmPICField_rank field.
Header file for Particle Motion and migration related functions.
PetscErrorCode UnpackSwarmFields(PetscInt i, const PetscInt64 *PIDs, const PetscReal *weights, const PetscReal *positions, const PetscInt *cellIndices, PetscReal *velocities, PetscInt *LocStatus, Particle *particle)
Initializes a Particle struct with data from DMSwarm fields.
PetscErrorCode UpdateSwarmFields(PetscInt i, const Particle *particle, PetscReal *weights, PetscInt *cellIndices, PetscInt *status_field)
Updates DMSwarm fields with data from a Particle struct.
PetscErrorCode CalculateParticleCountPerCell(UserCtx *user)
Counts particles in each cell of the DMDA 'da' and stores the result in user->ParticleCount.
#define LOG_LOOP_ALLOW(scope, level, iterVar, interval, fmt,...)
Logs a message inside a loop, but only every interval iterations.
#define LOG_ALLOW_SYNC(scope, level, fmt,...)
----— DEBUG ---------------------------------------— #define LOG_ALLOW(scope, level,...
#define LOCAL
Logging scope definitions for controlling message output.
#define GLOBAL
Scope for global logging across all processes.
const char * BCFaceToString(BCFace face)
Helper function to convert BCFace enum to a string representation.
#define LOG_ALLOW(scope, level, fmt,...)
Logging macro that checks both the log level and whether the calling function is in the allowed-funct...
#define PROFILE_FUNCTION_END
Marks the end of a profiled code block.
const char * ParticleLocationStatusToString(ParticleLocationStatus level)
A function that outputs the name of the current level in the ParticleLocation enum.
@ LOG_ERROR
Critical errors that may halt the program.
@ LOG_INFO
Informational messages about program execution.
@ LOG_WARNING
Non-critical issues that warrant attention.
@ LOG_DEBUG
Detailed debugging information.
@ LOG_VERBOSE
Extremely detailed logs, typically for development use only.
#define PROFILE_FUNCTION_BEGIN
Marks the beginning of a profiled code block (typically a function).
PetscErrorCode InitializeLogicalSpaceRNGs(PetscRandom *rand_logic_i, PetscRandom *rand_logic_j, PetscRandom *rand_logic_k)
Initializes random number generators for logical space operations [0.0, 1.0).
PetscErrorCode BinarySearchInt64(PetscInt n, const PetscInt64 arr[], PetscInt64 key, PetscBool *found)
Performs a binary search for a key in a sorted array of PetscInt64.
PetscBool inletFaceDefined
BCFace identifiedInletBCFace
SimCtx * simCtx
Back-pointer to the master simulation context.
ParticleLocationStatus
Defines the state of a particle with respect to its location and migration status during the iterativ...
Cmpnts max_coords
Maximum x, y, z coordinates of the bounding box.
Cmpnts min_coords
Minimum x, y, z coordinates of the bounding box.
char * current_io_directory
PetscMPIInt destination_rank
ParticleLocationStatus location_status
char particle_subdir[PETSC_MAX_PATH_LEN]
char source_dir[PETSC_MAX_PATH_LEN]
@ EXEC_MODE_POSTPROCESSOR
char _io_context_buffer[PETSC_MAX_PATH_LEN]
PetscInt ParticleInitialization
char restart_dir[PETSC_MAX_PATH_LEN]
Defines a 3D axis-aligned bounding box.
A 3D point or vector with PetscScalar components.
Information needed to migrate a single particle between MPI ranks.
Defines a particle's core properties for Lagrangian tracking.
Stores the MPI ranks of neighboring subdomains.
The master context for the entire simulation.
User-defined context containing data specific to a single computational grid level.
Head of a generic C-style linked list.
PetscErrorCode LocateParticleOrFindMigrationTarget(UserCtx *user, Particle *particle, ParticleLocationStatus *status_out)
Locates a particle's host cell or identifies its migration target using a robust walk search.