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__   "UpdateAllParticlePhysics"
 

Functions

PetscErrorCode UpdateParticleField (const char *fieldName, PetscReal t, Cmpnts pos, Cmpnts vel, PetscReal *psi_io)
 Updates a single particle's field based on its state and the specified field name.
 
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 7 of file ParticlePhysics.c.

◆ __FUNCT__ [1/3]

#define __FUNCT__   "UpdateParticleField"

Definition at line 11 of file ParticlePhysics.c.

◆ __FUNCT__ [2/3]

#define __FUNCT__   "UpdateFieldForAllParticles"

Definition at line 11 of file ParticlePhysics.c.

◆ __FUNCT__ [3/3]

#define __FUNCT__   "UpdateAllParticlePhysics"

Definition at line 11 of file ParticlePhysics.c.

Function Documentation

◆ UpdateParticleField()

PetscErrorCode UpdateParticleField ( const char *  fieldName,
PetscReal  t,
Cmpnts  pos,
Cmpnts  vel,
PetscReal *  psi_io 
)

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

This function serves as a switchboard for various particle property calculations. Given a particle's data, it computes a new value for the specified field.

Currently supported fields:

  • "Psi": Calculates the particle's kinetic energy (0.5 * |v|^2).
Parameters
[in]fieldNameThe name of the field to update (e.g., "Psi").
[in]tThe current simulation time.
[in]posThe particle's physical position.
[in]velThe particle's velocity.
[in,out]psi_ioA pointer to the particle's current Psi value (used for accumulation if needed).
Returns
PetscErrorCode 0 on success.

Definition at line 29 of file ParticlePhysics.c.

34{
35 PetscFunctionBeginUser;
36
37 if (strcmp(fieldName, "Psi") == 0) {
38 // --- Calculate Kinetic Energy ---
39 // Formula: Psi = 0.5 * (vx^2 + vy^2 + vz^2)
40 *psi_io = 0.5 * (vel.x * vel.x + vel.y * vel.y + vel.z * vel.z);
41
42 }
43 // else if (strcmp(fieldName, "SomeOtherField") == 0) {
44 // // ... logic for another field ...
45 // }
46 else {
47 // If the field name is not recognized, we can choose to do nothing or raise an error.
48 // For now, we'll just log a warning and do nothing.
49 LOG_ALLOW(GLOBAL, LOG_WARNING, "UpdateParticleField: Field name '%s' not recognized. No update performed.\n", fieldName);
50 }
51
52 PetscFunctionReturn(0);
53}
#define GLOBAL
Scope for global logging across all processes.
Definition logging.h:47
#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:201
@ LOG_WARNING
Non-critical issues that warrant attention.
Definition logging.h:30
PetscScalar x
Definition variables.h:101
PetscScalar z
Definition variables.h:101
PetscScalar y
Definition variables.h:101
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.

This function orchestrates the update of a single particle field across the entire local swarm. It gets access to the necessary particle data arrays and calls the UpdateParticleField kernel for each particle.

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

Definition at line 69 of file ParticlePhysics.c.

70{
71 PetscErrorCode ierr;
72 DM swarm = user->swarm;
73 PetscInt n_local;
74 PetscReal current_time = user->simCtx->ti;
75
76 // Pointers to swarm data arrays
77 const PetscReal *pos_arr;
78 const PetscReal *vel_arr;
79 PetscReal *psi_arr; // This will be our target field array
80
81 PetscFunctionBeginUser;
83
84 ierr = DMSwarmGetLocalSize(swarm, &n_local); CHKERRQ(ierr);
85 if (n_local == 0) {
86 PetscFunctionReturn(0); // Nothing to do on this rank
87 }
88
89 // Get read access to position and velocity
90 ierr = DMSwarmGetField(swarm, "position", NULL, NULL, (void**)&pos_arr); CHKERRQ(ierr);
91 ierr = DMSwarmGetField(swarm, "velocity", NULL, NULL, (void**)&vel_arr); CHKERRQ(ierr);
92
93 // Get write access to the target field
94 ierr = DMSwarmGetField(swarm, fieldName, NULL, NULL, (void**)&psi_arr); CHKERRQ(ierr);
95
96 // Loop over each local particle
97 for (PetscInt p = 0; p < n_local; ++p) {
98 // Unpack data for the current particle
99 Cmpnts current_pos = {pos_arr[3*p + 0], pos_arr[3*p + 1], pos_arr[3*p + 2]};
100 Cmpnts current_vel = {vel_arr[3*p + 0], vel_arr[3*p + 1], vel_arr[3*p + 2]};
101
102 // Call the kernel to calculate and update the value in the array directly.
103 // We pass the address of the specific particle's entry in the psi_arr.
104 ierr = UpdateParticleField(fieldName, current_time, current_pos, current_vel, &psi_arr[p]); CHKERRQ(ierr);
105 }
106
107 // Restore all field arrays
108 ierr = DMSwarmRestoreField(swarm, "position", NULL, NULL, (void**)&pos_arr); CHKERRQ(ierr);
109 ierr = DMSwarmRestoreField(swarm, "velocity", NULL, NULL, (void**)&vel_arr); CHKERRQ(ierr);
110 ierr = DMSwarmRestoreField(swarm, fieldName, NULL, NULL, (void**)&psi_arr); CHKERRQ(ierr);
111
112 LOG_ALLOW(GLOBAL, LOG_INFO, "Updated particle field '%s' for all local particles.\n", fieldName);
113
115 PetscFunctionReturn(0);
116}
PetscErrorCode UpdateParticleField(const char *fieldName, PetscReal t, Cmpnts pos, Cmpnts vel, PetscReal *psi_io)
Updates a single particle's field based on its state and the specified field name.
#define PROFILE_FUNCTION_END
Marks the end of a profiled code block.
Definition logging.h:740
@ LOG_INFO
Informational messages about program execution.
Definition logging.h:32
#define PROFILE_FUNCTION_BEGIN
Marks the beginning of a profiled code block (typically a function).
Definition logging.h:731
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:664
PetscReal ti
Definition variables.h:547
A 3D point or vector with PetscScalar components.
Definition variables.h:100
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 137 of file ParticlePhysics.c.

138{
139 PetscErrorCode ierr;
140
141 PetscFunctionBeginUser;
143
144 LOG_ALLOW(GLOBAL, LOG_INFO, "Updating all particle physical properties...\n");
145
146 // --- Call updater for 'Psi' field (currently calculating Kinetic Energy) ---
147 ierr = UpdateFieldForAllParticles(user, "Psi"); CHKERRQ(ierr);
148
149 // --- FUTURE EXPANSION: Add calls for other fields here ---
150 /*
151 if (simCtx->enable_temperature_solve) {
152 ierr = UpdateFieldForAllParticles(user, "Temperature"); CHKERRQ(ierr);
153 }
154 if (simCtx->enable_species_transport) {
155 ierr = UpdateFieldForAllParticles(user, "SpeciesMassFraction"); CHKERRQ(ierr);
156 }
157 */
158
159 LOG_ALLOW(GLOBAL, LOG_INFO, "All particle physical properties updated.\n");
160
162 PetscFunctionReturn(0);
163}
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: