PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
ParticlePhysics.c
Go to the documentation of this file.
1#include "ParticlePhysics.h"
3
4#ifndef ERROR_MSG_BUFFER_SIZE
5#define ERROR_MSG_BUFFER_SIZE 256
6#endif
7
8#undef __FUNCT__
9#define __FUNCT__ "UpdateParticleField"
10/**
11 * @brief Internal helper implementation: `UpdateParticleField()`.
12 * @details Local to this translation unit.
13 */
14PetscErrorCode UpdateParticleField(const char *fieldName,
15 PetscReal dt,
16 PetscReal *psi_io,
17 PetscReal diffusivity,
18 PetscReal mean_val,
19 PetscReal cell_vol,
20 PetscReal C_model)
21{
22 PetscFunctionBeginUser;
23
24 if (strcmp(fieldName, "Psi") == 0) {
25 if (cell_vol < 1.0e-14) cell_vol = 1.0e-14;
26
27 PetscReal delta2 = PetscPowReal(cell_vol, 0.6666667);
28 PetscReal omega = C_model * diffusivity / delta2;
29 PetscReal decay = PetscExpReal(-omega * dt);
30
31 PetscReal psi_old = *psi_io;
32 *psi_io = mean_val + (psi_old - mean_val) * decay;
33 }
34
35 PetscFunctionReturn(0);
36}
37
38#undef __FUNCT__
39#define __FUNCT__ "UpdateFieldForAllParticles"
40/**
41 * @brief Internal helper implementation: `UpdateFieldForAllParticles()`.
42 * @details Local to this translation unit.
43 */
44PetscErrorCode UpdateFieldForAllParticles(UserCtx *user, const char *fieldName)
45{
46 PetscErrorCode ierr;
47 DM swarm = user->swarm;
48 DM da = user->da;
49 PetscInt n_local;
50 PetscReal dt = user->simCtx->dt;
51 PetscReal C_IEM = 2.0;
52
53 PetscReal *psi_arr = NULL;
54 PetscReal *diff_arr = NULL;
55 PetscInt *cell_arr = NULL;
56
57 PetscReal ***grid_mean = NULL;
58 PetscReal ***grid_aj = NULL;
59
60 PetscBool accessed_eulerian = PETSC_FALSE;
61
62 PetscFunctionBeginUser;
64
65 ierr = DMSwarmGetLocalSize(swarm, &n_local); CHKERRQ(ierr);
66 if (n_local == 0) {
67 PetscFunctionReturn(0);
68 }
69
70 ierr = DMSwarmGetField(swarm, fieldName, NULL, NULL, (void**)&psi_arr); CHKERRQ(ierr);
71
72 if (strcmp(fieldName, "Psi") == 0) {
73 ierr = DMSwarmGetField(swarm, "Diffusivity", NULL, NULL, (void**)&diff_arr); CHKERRQ(ierr);
74 ierr = DMSwarmGetField(swarm, "DMSwarm_CellID", NULL, NULL, (void**)&cell_arr); CHKERRQ(ierr);
75
76 if (!user->lPsi || !user->lAj) {
77 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "UserCtx lPsi or lAj not initialized.");
78 }
79
80 ierr = DMDAVecGetArrayRead(da, user->lPsi, &grid_mean); CHKERRQ(ierr);
81 ierr = DMDAVecGetArrayRead(da, user->lAj, &grid_aj); CHKERRQ(ierr);
82 accessed_eulerian = PETSC_TRUE;
83 }
84
85 for (PetscInt p = 0; p < n_local; ++p) {
86 PetscReal p_diff = 0.0;
87 PetscReal p_mean = 0.0;
88 PetscReal p_vol = 1.0;
89
90 if (strcmp(fieldName, "Psi") == 0) {
91 PetscInt i = cell_arr[3*p + 0];
92 PetscInt j = cell_arr[3*p + 1];
93 PetscInt k = cell_arr[3*p + 2];
94
95 p_diff = diff_arr[p];
96 p_mean = grid_mean[k][j][i];
97
98 PetscReal jac = grid_aj[k][j][i];
99 p_vol = (jac > 1.0e-14) ? (1.0 / jac) : 1.0e-14;
100 }
101
102 ierr = UpdateParticleField(fieldName, dt, &psi_arr[p], p_diff, p_mean, p_vol, C_IEM);
103 CHKERRQ(ierr);
104 }
105
106 ierr = DMSwarmRestoreField(swarm, fieldName, NULL, NULL, (void**)&psi_arr); CHKERRQ(ierr);
107
108 if (strcmp(fieldName, "Psi") == 0) {
109 ierr = DMSwarmRestoreField(swarm, "Diffusivity", NULL, NULL, (void**)&diff_arr); CHKERRQ(ierr);
110 ierr = DMSwarmRestoreField(swarm, "DMSwarm_CellID", NULL, NULL, (void**)&cell_arr); CHKERRQ(ierr);
111 }
112
113 if (accessed_eulerian) {
114 ierr = DMDAVecRestoreArrayRead(da, user->lPsi, &grid_mean); CHKERRQ(ierr);
115 ierr = DMDAVecRestoreArrayRead(da, user->lAj, &grid_aj); CHKERRQ(ierr);
116 }
117
118 LOG_ALLOW(GLOBAL, LOG_INFO, "Updated particle physics for field '%s'.\n", fieldName);
119
121 PetscFunctionReturn(0);
122}
123
124#undef __FUNCT__
125#define __FUNCT__ "UpdateAllParticleFields"
126/**
127 * @brief Implementation of \ref UpdateAllParticleFields().
128 * @details Full API contract (arguments, ownership, side effects) is documented with
129 * the header declaration in `include/ParticlePhysics.h`.
130 * @see UpdateAllParticleFields()
131 */
132PetscErrorCode UpdateAllParticleFields(UserCtx *user)
133{
134 PetscErrorCode ierr;
135 PetscFunctionBeginUser;
137
138 LOG_ALLOW(GLOBAL, LOG_INFO, "Updating all particle physical properties...\n");
139
141 ierr = ApplyVerificationScalarOverrideToParticles(user); CHKERRQ(ierr);
142 LOG_ALLOW(GLOBAL, LOG_INFO, "Verification scalar override active; skipped model-driven Psi update.\n");
144 PetscFunctionReturn(0);
145 }
146
147 ierr = UpdateFieldForAllParticles(user, "Psi"); CHKERRQ(ierr);
148
149 LOG_ALLOW(GLOBAL, LOG_INFO, "All particle physical properties updated.\n");
150
152 PetscFunctionReturn(0);
153}
PetscErrorCode UpdateParticleField(const char *fieldName, PetscReal dt, PetscReal *psi_io, PetscReal diffusivity, PetscReal mean_val, PetscReal cell_vol, PetscReal C_model)
Internal helper implementation: UpdateParticleField().
PetscErrorCode UpdateFieldForAllParticles(UserCtx *user, const char *fieldName)
Internal helper implementation: UpdateFieldForAllParticles().
PetscErrorCode UpdateAllParticleFields(UserCtx *user)
Implementation of UpdateAllParticleFields().
Header file for Particle related physics modules.
#define GLOBAL
Scope for global logging across all processes.
Definition logging.h:45
#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
#define PROFILE_FUNCTION_END
Marks the end of a profiled code block.
Definition logging.h:779
@ LOG_INFO
Informational messages about program execution.
Definition logging.h:30
#define PROFILE_FUNCTION_BEGIN
Marks the beginning of a profiled code block (typically a function).
Definition logging.h:770
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:814
PetscReal dt
Definition variables.h:658
Vec lPsi
Definition variables.h:883
Vec lAj
Definition variables.h:858
User-defined context containing data specific to a single computational grid level.
Definition variables.h:811
PetscErrorCode ApplyVerificationScalarOverrideToParticles(UserCtx *user)
Populates the particle Psi field from a verification-only source override.
PetscBool VerificationScalarOverrideActive(const SimCtx *simCtx)
Reports whether a verification-only scalar override is active.