PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
ParticlePhysics.c File Reference
#include "ParticlePhysics.h"
Include dependency graph for ParticlePhysics.c:

Go to the source code of this file.

Macros

#define ERROR_MSG_BUFFER_SIZE   256
 
#define __FUNCT__   "UpdateParticleField"
 
#define __FUNCT__   "UpdateFieldForAllParticles"
 
#define __FUNCT__   "UpdateAllParticleFields"
 

Functions

PetscErrorCode UpdateParticleField (const char *fieldName, PetscReal dt, PetscReal *psi_io, PetscReal diffusivity, PetscReal mean_val, PetscReal cell_vol, PetscReal C_model)
 Updates a single particle's field based on its state and physics model.
 
PetscErrorCode UpdateFieldForAllParticles (UserCtx *user, const char *fieldName)
 Loops over all local particles and updates a specified field.
 
PetscErrorCode UpdateAllParticleFields (UserCtx *user)
 Orchestrates the update of all physical properties for particles.
 

Macro Definition Documentation

◆ ERROR_MSG_BUFFER_SIZE

#define ERROR_MSG_BUFFER_SIZE   256

Definition at line 4 of file ParticlePhysics.c.

◆ __FUNCT__ [1/3]

#define __FUNCT__   "UpdateParticleField"

Definition at line 8 of file ParticlePhysics.c.

◆ __FUNCT__ [2/3]

#define __FUNCT__   "UpdateFieldForAllParticles"

Definition at line 8 of file ParticlePhysics.c.

◆ __FUNCT__ [3/3]

#define __FUNCT__   "UpdateAllParticleFields"

Definition at line 8 of file ParticlePhysics.c.

Function Documentation

◆ UpdateParticleField()

PetscErrorCode UpdateParticleField ( const char *  fieldName,
PetscReal  dt,
PetscReal *  psi_io,
PetscReal  diffusivity,
PetscReal  mean_val,
PetscReal  cell_vol,
PetscReal  C_model 
)

Updates a single particle's field based on its state and physics model.

Updates a single particle's field based on its state and the specified field name.

Implements the IEM (Interaction by Exchange with the Mean) model for scalar mixing.

Physics: dPsi/dt = -Omega * (Psi - <Psi>) Solution: Psi_new = <Psi> + (Psi_old - <Psi>) * exp(-Omega * dt)

Parameters
[in]fieldNameName of the field (e.g., "Psi").
[in]dtTime step size.
[in,out]psi_ioPointer to the particle's scalar value (Psi).
[in]diffusivityParticle diffusivity (Gamma + Gamma_t).
[in]mean_valLocal Eulerian mean value (<Psi>).
[in]cell_volVolume of the host cell (1/Jacobian).
[in]C_modelModel constant (C_IEM).

Definition at line 25 of file ParticlePhysics.c.

32{
33 PetscFunctionBeginUser;
34
35 if (strcmp(fieldName, "Psi") == 0) {
36 // --- IEM Mixing Model ---
37
38 // 1. Calculate Characteristic Length Scale (Delta)
39 // Delta^2 approx Volume^(2/3)
40 // Safety: Ensure volume is positive
41 if (cell_vol < 1.0e-14) cell_vol = 1.0e-14;
42
43 PetscReal delta2 = PetscPowReal(cell_vol, 0.6666667);
44
45 // 2. Calculate Mixing Frequency (Omega)
46 // Omega = C_phi * (Gamma_eff) / Delta^2
47 PetscReal omega = C_model * diffusivity / delta2;
48
49 // 3. Analytical Integration
50 // Exact solution for linear decay ODE
51 PetscReal decay = PetscExpReal(-omega * dt);
52
53 // Update Value IN-PLACE
54 PetscReal psi_old = *psi_io;
55 *psi_io = mean_val + (psi_old - mean_val) * decay;
56 }
57 else {
58 // Logic for other fields can be added here
59 // e.g. Temperature might use a different C_model or source term
60 }
61
62 PetscFunctionReturn(0);
63}
Here is the caller graph for this function:

◆ UpdateFieldForAllParticles()

PetscErrorCode UpdateFieldForAllParticles ( UserCtx user,
const char *  fieldName 
)

Loops over all local particles and updates a specified field.

Prepares necessary Eulerian and Lagrangian data structures before looping.

Parameters
[in,out]userPointer to the UserCtx.
[in]fieldNameThe name of the field to update (e.g., "Psi").
Returns
PetscErrorCode 0 on success.

Definition at line 76 of file ParticlePhysics.c.

77{
78 PetscErrorCode ierr;
79 DM swarm = user->swarm;
80 DM da = user->da; // Scalar DM for accessing lPsi and lAj
81 PetscInt n_local;
82 PetscReal dt = user->simCtx->dt;
83
84 // Model Constant (Default to 2.0 if not set)
85 PetscReal C_IEM = 2.0;//(user->simCtx->C_IEM > 1.0e-6) ? user->simCtx->C_IEM : 2.0;
86
87 // Pointers for Swarm Data
88 PetscReal *psi_arr = NULL;
89 PetscReal *diff_arr = NULL;
90 PetscInt *cell_arr = NULL;
91
92 // Pointers for Eulerian Data (Read-Only)
93 PetscReal ***grid_mean = NULL;
94 PetscReal ***grid_aj = NULL;
95
96 // Flags to track what we accessed (for cleanup)
97 PetscBool accessed_eulerian = PETSC_FALSE;
98
99 PetscFunctionBeginUser;
101
102 ierr = DMSwarmGetLocalSize(swarm, &n_local); CHKERRQ(ierr);
103 if (n_local == 0) {
104 PetscFunctionReturn(0);
105 }
106
107 // 1. Access Target Field (Write Access)
108 ierr = DMSwarmGetField(swarm, fieldName, NULL, NULL, (void**)&psi_arr); CHKERRQ(ierr);
109
110 // 2. Access Auxiliary Data based on Physics Requirements
111 // For "Psi", we need Diffusivity, CellID, Grid Mean, and Grid Metrics.
112 if (strcmp(fieldName, "Psi") == 0) {
113
114 // --- Lagrangian Data ---
115 ierr = DMSwarmGetField(swarm, "Diffusivity", NULL, NULL, (void**)&diff_arr); CHKERRQ(ierr);
116 ierr = DMSwarmGetField(swarm, "DMSwarm_CellID", NULL, NULL, (void**)&cell_arr); CHKERRQ(ierr);
117
118 // --- Eulerian Data ---
119 // We need Local vectors (lPsi, lAj) to handle ghost indexing safely
120 if (!user->lPsi || !user->lAj) {
121 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "UserCtx lPsi or lAj not initialized.");
122 }
123
124 ierr = DMDAVecGetArrayRead(da, user->lPsi, &grid_mean); CHKERRQ(ierr);
125 ierr = DMDAVecGetArrayRead(da, user->lAj, &grid_aj); CHKERRQ(ierr);
126 accessed_eulerian = PETSC_TRUE;
127 }
128
129 // 3. Loop over particles
130 for (PetscInt p = 0; p < n_local; ++p) {
131
132 // Default values for generic fields
133 PetscReal p_diff = 0.0;
134 PetscReal p_mean = 0.0;
135 PetscReal p_vol = 1.0;
136
137 if (strcmp(fieldName, "Psi") == 0) {
138 // Unpack indices
139 PetscInt i = cell_arr[3*p + 0];
140 PetscInt j = cell_arr[3*p + 1];
141 PetscInt k = cell_arr[3*p + 2];
142
143 // Fetch Data
144 p_diff = diff_arr[p];
145 p_mean = grid_mean[k][j][i]; // Direct Lookup
146
147 // Calculate Volume = 1.0 / Jacobian
148 PetscReal jac = grid_aj[k][j][i];
149 p_vol = (jac > 1.0e-14) ? (1.0 / jac) : 1.0e-14;
150 }
151
152 // Call Kernel
153 ierr = UpdateParticleField(fieldName, dt, &psi_arr[p], p_diff, p_mean, p_vol, C_IEM);
154 CHKERRQ(ierr);
155 }
156
157 // 4. Restore Arrays
158 ierr = DMSwarmRestoreField(swarm, fieldName, NULL, NULL, (void**)&psi_arr); CHKERRQ(ierr);
159
160 if (strcmp(fieldName, "Psi") == 0) {
161 ierr = DMSwarmRestoreField(swarm, "Diffusivity", NULL, NULL, (void**)&diff_arr); CHKERRQ(ierr);
162 ierr = DMSwarmRestoreField(swarm, "DMSwarm_CellID", NULL, NULL, (void**)&cell_arr); CHKERRQ(ierr);
163 }
164
165 if (accessed_eulerian) {
166 ierr = DMDAVecRestoreArrayRead(da, user->lPsi, &grid_mean); CHKERRQ(ierr);
167 ierr = DMDAVecRestoreArrayRead(da, user->lAj, &grid_aj); CHKERRQ(ierr);
168 }
169
170 LOG_ALLOW(GLOBAL, LOG_INFO, "Updated particle physics for field '%s'.\n", fieldName);
171
173 PetscFunctionReturn(0);
174}
PetscErrorCode UpdateParticleField(const char *fieldName, PetscReal dt, PetscReal *psi_io, PetscReal diffusivity, PetscReal mean_val, PetscReal cell_vol, PetscReal C_model)
Updates a single particle's field based on its state and physics model.
#define GLOBAL
Scope for global logging across all processes.
Definition logging.h:46
#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:200
#define PROFILE_FUNCTION_END
Marks the end of a profiled code block.
Definition logging.h:746
@ LOG_INFO
Informational messages about program execution.
Definition logging.h:31
#define PROFILE_FUNCTION_BEGIN
Marks the beginning of a profiled code block (typically a function).
Definition logging.h:737
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:731
PetscReal dt
Definition variables.h:599
Vec lPsi
Definition variables.h:801
Vec lAj
Definition variables.h:776
Here is the call graph for this function:
Here is the caller graph for this function:

◆ UpdateAllParticleFields()

PetscErrorCode UpdateAllParticleFields ( UserCtx user)

Orchestrates the update of all physical properties for particles.

This function serves as the top-level entry point for updating particle-specific physical quantities after their position and the surrounding fluid velocity are known. It calls a sequence of more specific update routines for each property.

For example, it can be configured to update:

Parameters
[in,out]userPointer to the UserCtx containing the swarm and simulation context.
Returns
PetscErrorCode 0 on success.

Definition at line 178 of file ParticlePhysics.c.

179{
180 PetscErrorCode ierr;
181 PetscFunctionBeginUser;
183
184 LOG_ALLOW(GLOBAL, LOG_INFO, "Updating all particle physical properties...\n");
185
186 // Update Scalar Mixing (IEM)
187 // IMPORTANT: Ensure ScatterParticleFieldToEulerField("Psi") has been called
188 // immediately before this to populate user->lPsi with the current mean.
189 ierr = UpdateFieldForAllParticles(user, "Psi"); CHKERRQ(ierr);
190
191 LOG_ALLOW(GLOBAL, LOG_INFO, "All particle physical properties updated.\n");
192
194 PetscFunctionReturn(0);
195}
PetscErrorCode UpdateFieldForAllParticles(UserCtx *user, const char *fieldName)
Loops over all local particles and updates a specified field.
Here is the call graph for this function:
Here is the caller graph for this function: