PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
rhs.h
Go to the documentation of this file.
1#ifndef RHS_H
2#define RHS_H
3
4#include "variables.h" // Provides definitions for UserCtx, SimCtx, IBMNodes, etc.
5#include "logging.h"
6#include "Metric.h"
7#include "BodyForces.h"
8
9/**
10 * @brief Computes the viscous contribution to the contravariant momentum RHS.
11 *
12 * This routine evaluates diffusive fluxes on the curvilinear grid and writes the
13 * resulting term into `Visc`. The caller is responsible for providing compatible
14 * vectors and for assembling any additional source terms afterwards.
15 *
16 * @param[in] user Block-level solver context containing metrics and model parameters.
17 * @param[in] Ucont Contravariant velocity field used by the discretization.
18 * @param[in] Ucat Cartesian velocity field used for derivative evaluation.
19 * @param[out] Visc Output vector receiving the viscous RHS contribution.
20 * @return PetscErrorCode 0 on success.
21 */
22PetscErrorCode Viscous(UserCtx *user, Vec Ucont, Vec Ucat, Vec Visc);
23
24/**
25 * @brief Computes the convective contribution to the contravariant momentum RHS.
26 *
27 * This routine evaluates the advection operator on the current velocity state and
28 * stores the contribution in `Conv` for subsequent combination with viscous and
29 * body-force terms.
30 *
31 * @param[in] user Block-level solver context containing metrics and numerics settings.
32 * @param[in] Ucont Contravariant velocity field used in face-normal flux construction.
33 * @param[in] Ucat Cartesian velocity field used by the convective stencil.
34 * @param[out] Conv Output vector receiving the convection RHS contribution.
35 * @return PetscErrorCode 0 on success.
36 */
37PetscErrorCode Convection(UserCtx *user, Vec Ucont, Vec Ucat, Vec Conv);
38
39/**
40 * @brief General dispatcher for applying all active body forces (momentum sources).
41 *
42 * This function serves as a central hub for adding momentum source terms to the
43 * contravariant right-hand-side (Rct) of the momentum equations. It is called once per RHS
44 * evaluation (e.g., once per Runge-Kutta stage).
45 *
46 * It introspects the simulation configuration to determine which, if any, body
47 * forces are active and calls their specific implementation functions.
48 *
49 * @param user The UserCtx containing the simulation state for a single block.
50 * @param Rct The PETSc Vec for the contravariant RHS, modified in place.
51 * @return PetscErrorCode 0 on success.
52 */
53PetscErrorCode ComputeBodyForces(UserCtx *user, Vec Rct);
54
55/**
56 * @brief Computes the Right-Hand Side (RHS) of the momentum equations.
57 *
58 * This function calculates the contribution of the convective and diffusive terms.
59 * It is called by the momentum solvers (e.g., RungeKutta).
60 *
61 * @param user The UserCtx for a single block.
62 * @param Rhs The PETSc Vec where the RHS result will be stored.
63 * @return PetscErrorCode 0 on success.
64 */
65extern PetscErrorCode ComputeRHS(UserCtx *user, Vec Rhs);
66/**
67 * @brief Computes the effective diffusivity scalar field (Gamma_eff) on the Eulerian grid.
68 *
69 * This function calculates the total diffusivity used to drive the stochastic
70 * motion of particles (Scalar FDF). It combines molecular diffusion and
71 * turbulent diffusion.
72 *
73 * Formula:
74 * Gamma_eff = nu/Sc + nu_t/Sc_t
75 *
76 * Where:
77 * - nu = 1/Re (kinematic viscosity)
78 * - nu_t (eddy viscosity from LES/RANS model)
79 * - Sc (molecular Schmidt number)
80 * - Sc_t (turbulent Schmidt number)
81 *
82 * @note If turbulence models are disabled, nu_t is assumed to be 0.
83 * @note This function updates the local ghost values of lDiffusivity at the end
84 * to ensure gradients can be computed correctly at subdomain boundaries.
85 *
86 * @param[in,out] user Pointer to the user context containing grid data and simulation parameters.
87 *
88 * @return PetscErrorCode 0 on success.
89 */
90PetscErrorCode ComputeEulerianDiffusivity(UserCtx *user);
91
92/**
93 * @brief Computes the Eulerian gradient of the effective diffusivity field.
94 *
95 * Reads the scalar diffusivity field and writes a vector gradient field used by
96 * particle stochastic transport updates.
97 *
98 * @param[in,out] user Pointer to user context containing diffusivity vectors.
99 * @return PetscErrorCode 0 on success.
100 */
101PetscErrorCode ComputeEulerianDiffusivityGradient(UserCtx *user);
102
103#endif // RHS_H
Logging utilities and macros for PETSc-based applications.
PetscErrorCode ComputeBodyForces(UserCtx *user, Vec Rct)
General dispatcher for applying all active body forces (momentum sources).
Definition rhs.c:1162
PetscErrorCode Viscous(UserCtx *user, Vec Ucont, Vec Ucat, Vec Visc)
Computes the viscous contribution to the contravariant momentum RHS.
Definition rhs.c:519
PetscErrorCode ComputeEulerianDiffusivity(UserCtx *user)
Computes the effective diffusivity scalar field (Gamma_eff) on the Eulerian grid.
Definition rhs.c:1959
PetscErrorCode ComputeEulerianDiffusivityGradient(UserCtx *user)
Computes the Eulerian gradient of the effective diffusivity field.
Definition rhs.c:2089
PetscErrorCode Convection(UserCtx *user, Vec Ucont, Vec Ucat, Vec Conv)
Computes the convective contribution to the contravariant momentum RHS.
Definition rhs.c:13
PetscErrorCode ComputeRHS(UserCtx *user, Vec Rhs)
Computes the Right-Hand Side (RHS) of the momentum equations.
Definition rhs.c:1185
Main header file for a complex fluid dynamics solver.
User-defined context containing data specific to a single computational grid level.
Definition variables.h:811