PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
ParticlePhysics.h
Go to the documentation of this file.
1/**
2 * @file ParticlePhysics.h
3 * @brief Header file for Particle related physics modules.
4 *
5 * This file contains declarations of functions responsible for solving for fields in the particle swarms within a simulation using PETSc's DMSwarm.
6 */
7
8 #ifndef PARTICLE_PHYSICS_H
9 #define PARTICLE_PHYSICS_H
10
11// Include necessary headers
12#include <petsc.h> // PETSc library header
13#include <petscdmswarm.h> // PETSc DMSwarm header
14#include <stdbool.h>
15#include <petscsys.h> // For PetscRealloc
16#include <math.h>
17#include "variables.h" // Common type definitions
18#include "logging.h" // Logging macros and definitions
19#include "walkingsearch.h" // Walking search function for particle migration
20
21/**
22 * @brief Updates a single particle's field based on its state and the specified field name.
23 *
24 * This function serves as a switchboard for various particle property calculations.
25 * Given a particle's data, it computes a new value for the specified field.
26 *
27 * Currently supported fields:
28 * - "Psi": Calculates the particle's kinetic energy (0.5 * |v|^2).
29 *
30 * @param[in] fieldName The name of the field to update (e.g., "Psi").
31 * @param[in] t The current simulation time.
32 * @param[in] pos The particle's physical position.
33 * @param[in] vel The particle's velocity.
34 * @param[in,out] psi_io A pointer to the particle's current Psi value (used for accumulation if needed).
35 *
36 * @return PetscErrorCode 0 on success.
37 */
38//PetscErrorCode UpdateParticleField(const char *fieldName,
39// PetscReal t,
40// Cmpnts pos,
41// Cmpnts vel,
42// PetscReal *psi_io);
43/**
44 * @brief Updates a single particle's field based on its state and physics model.
45 *
46 * Implements the IEM (Interaction by Exchange with the Mean) model for scalar mixing.
47 *
48 * Physics: dPsi/dt = -Omega * (Psi - <Psi>)
49 * Solution: Psi_new = <Psi> + (Psi_old - <Psi>) * exp(-Omega * dt)
50 *
51 * @param[in] fieldName Name of the field (e.g., "Psi").
52 * @param[in] dt Time step size.
53 * @param[in,out] psi_io Pointer to the particle's scalar value (Psi).
54 * @param[in] diffusivity Particle diffusivity (Gamma + Gamma_t).
55 * @param[in] mean_val Local Eulerian mean value (<Psi>).
56 * @param[in] cell_vol Volume of the host cell (1/Jacobian).
57 * @param[in] C_model Model constant (C_IEM).
58 */
59PetscErrorCode UpdateParticleField(const char *fieldName,
60 PetscReal dt,
61 PetscReal *psi_io,
62 PetscReal diffusivity,
63 PetscReal mean_val,
64 PetscReal cell_vol,
65 PetscReal C_model);
66
67
68/**
69 * @brief Loops over all local particles and updates a specified field.
70 *
71 * This function orchestrates the update of a single particle field across the entire
72 * local swarm. It gets access to the necessary particle data arrays and calls the
73 * `UpdateParticleField` kernel for each particle.
74 *
75 * @param[in,out] user Pointer to the UserCtx containing the swarm and simulation context.
76 * @param[in] fieldName The name of the field to update (e.g., "Psi").
77 *
78 * @return PetscErrorCode 0 on success.
79 */
80PetscErrorCode UpdateFieldForAllParticles(UserCtx *user, const char *fieldName);
81
82
83/**
84 * @brief Orchestrates the update of all physical properties for particles.
85 *
86 * This function serves as the top-level entry point for updating particle-specific
87 * physical quantities after their position and the surrounding fluid velocity are known.
88 * It calls a sequence of more specific update routines for each property.
89 *
90 * For example, it can be configured to update:
91 * - Particle kinetic energy (stored in "Psi")
92 * - Particle accumulated shear stress
93 * - Particle temperature
94 * - etc.
95 *
96 * @param[in,out] user Pointer to the UserCtx containing the swarm and simulation context.
97 *
98 * @return PetscErrorCode 0 on success.
99 */
100PetscErrorCode UpdateAllParticleFields(UserCtx *user);
101
102#endif // PARTICLE_PHYSICS_H
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 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.
Logging utilities and macros for PETSc-based applications.
Main header file for a complex fluid dynamics solver.
User-defined context containing data specific to a single computational grid level.
Definition variables.h:728
Header file for particle location functions using the walking search algorithm.