PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
wallfunction.h
Go to the documentation of this file.
1#ifndef WALLFUNCTION_H
2#define WALLFUNCTION_H
3
4#include <petscpf.h>
5#include <petscdmswarm.h>
6#include <stdlib.h>
7#include <time.h>
8#include <math.h>
9#include <petsctime.h>
10#include <petscsys.h>
11#include <petscdmcomposite.h>
12#include <petscsystypes.h>
13
14// Include additional headers
15#include "variables.h" // Shared type definitions
16#include "ParticleSwarm.h" // Particle swarm functions
17#include "walkingsearch.h" // Particle location functions
18#include "grid.h" // Grid functions
19#include "logging.h" // Logging macros
20#include "io.h" // Data Input and Output functions
21#include "interpolation.h" // Interpolation routines
22#include "ParticleMotion.h" // Functions related to motion of particles
23#include "Boundaries.h" // Functions related to Boundary condition
24
25// ===========================================================================
26// Function Declarations
27// ===========================================================================
28
29/**
30 * @brief Applies no-slip wall boundary condition with linear interpolation
31 *
32 * This function enforces a no-slip boundary condition (zero velocity at the wall)
33 * by linearly interpolating between the wall velocity (typically zero) and the
34 * velocity at a reference point in the flow.
35 *
36 * MATHEMATICAL FORMULATION:
37 * For a point at distance sb from the wall, with a reference velocity Uc
38 * at distance sc from the wall:
39 * U_boundary = U_wall + (U_reference - U_wall) * (sb / sc)
40 *
41 * PHYSICAL INTERPRETATION:
42 * This provides a first-order approximation assuming linear velocity variation
43 * in the near-wall region, which is valid in the viscous sublayer.
44 *
45 * @param[in] user Simulation context (unused but required for interface)
46 * @param[in] distance_reference Wall-normal distance to reference point (sc)
47 * @param[in] distance_boundary Wall-normal distance to boundary point (sb)
48 * @param[in] velocity_wall Velocity at the wall (Ua), typically zero
49 * @param[in] velocity_reference Velocity at reference point (Uc)
50 * @param[out] velocity_boundary Computed velocity at boundary point (Ub)
51 * @param[in] normal_x X-component of wall normal vector
52 * @param[in] normal_y Y-component of wall normal vector
53 * @param[in] normal_z Z-component of wall normal vector
54 *
55 * @note For moving walls, velocity_wall would be non-zero
56 * @note The normal vector should point INTO the fluid domain
57 */
58void noslip(UserCtx *user,
59 double distance_reference, double distance_boundary,
60 Cmpnts velocity_wall, Cmpnts velocity_reference,
61 Cmpnts *velocity_boundary,
62 double normal_x, double normal_y, double normal_z);
63
64/**
65 * @brief Applies free-slip wall boundary condition
66 *
67 * Free-slip conditions allow tangential flow but enforce zero normal velocity.
68 * This is appropriate for inviscid walls or symmetry planes where there is
69 * no shear stress but flow cannot penetrate the boundary.
70 *
71 * DECOMPOSITION:
72 * Velocity is decomposed into normal and tangential components:
73 * U = U_n * n + U_t
74 * where U_n = U · n (normal component)
75 * U_t = U - U_n * n (tangential component)
76 *
77 * BOUNDARY CONDITIONS:
78 * - Normal component: Interpolated from interior (∂U_n/∂n ≠ 0)
79 * - Tangential component: Extrapolated from interior (∂U_t/∂n = 0)
80 *
81 * @param[in] user Simulation context
82 * @param[in] distance_reference Wall-normal distance to reference point
83 * @param[in] distance_boundary Wall-normal distance to boundary point
84 * @param[in] velocity_wall Velocity at the wall
85 * @param[in] velocity_reference Velocity at reference point
86 * @param[out] velocity_boundary Computed velocity at boundary point
87 * @param[in] normal_x X-component of wall normal vector
88 * @param[in] normal_y Y-component of wall normal vector
89 * @param[in] normal_z Z-component of wall normal vector
90 */
91void freeslip(UserCtx *user,
92 double distance_reference, double distance_boundary,
93 Cmpnts velocity_wall, Cmpnts velocity_reference,
94 Cmpnts *velocity_boundary,
95 double normal_x, double normal_y, double normal_z);
96
97/**
98 * @brief Computes roughness-modified log-law coefficient E
99 *
100 * The coefficient E accounts for wall roughness effects on the log-law:
101 * u+ = (1/κ) ln(E * y+)
102 *
103 * ROUGHNESS REGIMES:
104 * 1. Hydraulically smooth (ks+ < 2.25): E = exp(κ*B), no roughness effect
105 * 2. Transitional (2.25 < ks+ < 90): Smooth interpolation between regimes
106 * 3. Fully rough (ks+ > 90): E modified by roughness height
107 * where ks+ = u_τ * ks / ν (roughness Reynolds number)
108 *
109 * @param[in] friction_velocity Friction velocity u_τ
110 * @param[in] roughness_height Equivalent sand grain roughness height ks
111 * @param[in] kinematic_viscosity Kinematic viscosity ν
112 * @return Roughness-modified log-law coefficient E
113 *
114 * @note This follows the Nikuradse sand grain roughness correlation
115 */
116double E_coeff(double friction_velocity, double roughness_height, double kinematic_viscosity);
117
118/**
119 * @brief Computes velocity from log-law for rough walls
120 *
121 * Calculates the tangential velocity at a given wall distance using
122 * the roughness-modified log-law of the wall.
123 *
124 * APPLICABLE REGIMES:
125 * - Viscous sublayer (y+ < 11.53): u+ = y+
126 * - Log-layer (11.53 < y+ < 300): u+ = (1/κ) ln(E * y+)
127 * - Outer layer (y+ > 300): Returns -1 (invalid, wake region)
128 *
129 * @param[in] kinematic_viscosity Kinematic viscosity ν
130 * @param[in] wall_distance Distance from wall y
131 * @param[in] friction_velocity Friction velocity u_τ
132 * @param[in] roughness_height Equivalent sand grain roughness ks
133 * @return Tangential velocity u_t, or -1 if outside valid range
134 *
135 * @note The upper limit y+ < 300 ensures we stay within the log-layer
136 * @note For y+ > 300, wake effects become important and log-law is invalid
137 */
138double u_hydset_roughness(double kinematic_viscosity, double wall_distance,
139 double friction_velocity, double roughness_height);
140
141/**
142 * @brief Residual function for friction velocity equation (log-law with roughness)
143 *
144 * This function computes the residual for Newton-Raphson iteration:
145 * f(u_τ) = u_predicted(u_τ) - u_known
146 * where u_predicted comes from the log-law or linear law depending on y+.
147 *
148 * @param[in] kinematic_viscosity Kinematic viscosity
149 * @param[in] known_velocity Known velocity at distance y
150 * @param[in] wall_distance Distance from wall
151 * @param[in] friction_velocity_guess Current guess for u_τ
152 * @param[in] roughness_height Wall roughness height
153 * @return Residual value (should be zero at solution)
154 */
155double f_hydset(double kinematic_viscosity, double known_velocity,
156 double wall_distance, double friction_velocity_guess,
157 double roughness_height);
158
159/**
160 * @brief Numerical derivative of residual function
161 *
162 * Computes df/du_τ using finite differences for Newton-Raphson iteration.
163 *
164 * @param[in] kinematic_viscosity Kinematic viscosity
165 * @param[in] known_velocity Known velocity
166 * @param[in] wall_distance Distance from wall
167 * @param[in] friction_velocity_guess Current guess for u_τ
168 * @param[in] roughness_height Wall roughness
169 * @return Derivative df/du_τ
170 */
171double df_hydset(double kinematic_viscosity, double known_velocity,
172 double wall_distance, double friction_velocity_guess,
173 double roughness_height);
174
175/**
176 * @brief Solves for friction velocity using Newton-Raphson iteration
177 *
178 * Given a known velocity at a known distance from the wall, this function
179 * iteratively solves for the friction velocity u_τ that satisfies the
180 * roughness-modified log-law or linear law.
181 *
182 * ALGORITHM:
183 * Newton-Raphson: u_τ^(n+1) = u_τ^n - f(u_τ^n) / f'(u_τ^n)
184 * Convergence criterion: |u_τ^(n+1) - u_τ^n| < 1e-7
185 *
186 * @param[in] kinematic_viscosity Kinematic viscosity
187 * @param[in] known_velocity Velocity at reference point
188 * @param[in] wall_distance Distance from wall to reference point
189 * @param[in] initial_guess Initial guess for u_τ
190 * @param[in] roughness_height Wall roughness height
191 * @return Converged friction velocity u_τ
192 *
193 * @note Maximum 30 iterations; warns if convergence not achieved
194 * @note Typical convergence in 3-5 iterations for good initial guess
195 */
196double find_utau_hydset(double kinematic_viscosity, double known_velocity,
197 double wall_distance, double initial_guess,
198 double roughness_height);
199
200
201/**
202 * @brief Computes turbulent eddy viscosity ratio (ν_t / ν)
203 *
204 * Uses the mixing length model with van Driest damping:
205 * ν_t / ν = κ * y+ * [1 - exp(-y+ / A+)]²
206 * where A+ ≈ 19 is the damping coefficient.
207 *
208 * PHYSICAL INTERPRETATION:
209 * - Near wall (y+ → 0): ν_t → 0 (viscous effects dominate)
210 * - Far from wall (y+ >> A+): ν_t ~ κ * y+ (mixing length theory)
211 * - Damping function smoothly transitions between these limits
212 *
213 * @param[in] yplus Normalized wall distance y+ = y * u_τ / ν
214 * @return Eddy viscosity ratio ν_t / ν
215 *
216 * @note This is the standard mixing length model with van Driest damping
217 * @note Valid in the inner layer (y+ < 50-100)
218 */
219double nu_t(double yplus);
220
221/**
222 * @brief Integrates eddy viscosity profile from wall to distance y
223 *
224 * Computes integrals needed for non-equilibrium wall functions:
225 * If mode=0: ∫[0 to y] dy / (ν + ν_t)
226 * If mode=1: ∫[0 to y] y dy / (ν + ν_t)
227 *
228 * These integrals appear in the solution of the momentum equation
229 * with pressure gradients in the near-wall region.
230 *
231 * NUMERICAL METHOD:
232 * - Trapezoidal rule with 30 integration points
233 * - Higher-order correction at endpoints for accuracy
234 *
235 * @param[in] kinematic_viscosity Kinematic viscosity ν
236 * @param[in] wall_distance Distance from wall y
237 * @param[in] friction_velocity Friction velocity u_τ
238 * @param[in] integration_mode 0 for F integral, 1 for Fy integral
239 * @return Value of the integral
240 *
241 * @note Used in Cabot wall function for pressure gradient effects
242 */
243double integrate_1(double kinematic_viscosity, double wall_distance,
244 double friction_velocity, int integration_mode);
245
246/**
247 * @brief Computes wall shear stress with pressure gradient effects
248 *
249 * Solves the integrated momentum equation in the near-wall region:
250 * τ_w = (u - dp/dx * F_y) / F_1
251 * where F_1 and F_y are integrals of the effective viscosity profile.
252 *
253 * This accounts for non-equilibrium effects due to streamwise pressure gradients.
254 *
255 * @param[in] kinematic_viscosity Kinematic viscosity
256 * @param[in] friction_velocity Friction velocity u_τ
257 * @param[in] wall_distance Distance from wall
258 * @param[in] velocity Velocity at wall_distance
259 * @param[in] pressure_gradient_tangent Tangential pressure gradient dp/ds
260 * @return Wall shear stress τ_w
261 *
262 * @note This is more accurate than equilibrium wall functions in separated flows
263 */
264double taw(double kinematic_viscosity, double friction_velocity,
265 double wall_distance, double velocity, double pressure_gradient_tangent);
266
267/**
268 * @brief Computes velocity using Cabot wall function
269 *
270 * Reconstructs velocity from wall shear stress and pressure gradient:
271 * u = τ_w * F1 + (dp/dx) * Fy
272 *
273 * @param[in] kinematic_viscosity Kinematic viscosity
274 * @param[in] wall_distance Distance from wall
275 * @param[in] friction_velocity Friction velocity
276 * @param[in] pressure_gradient_tangent Tangential pressure gradient
277 * @param[in] wall_shear_stress Wall shear stress
278 * @return Velocity at wall_distance
279 */
280double u_Cabot(double kinematic_viscosity, double wall_distance,
281 double friction_velocity, double pressure_gradient_tangent,
282 double wall_shear_stress);
283
284/**
285 * @brief Residual function for Cabot wall function
286 *
287 * Computes: f(u_τ) = u_τ - sqrt(|τ_w(u_τ)|)
288 * This enforces consistency between u_τ and τ_w.
289 *
290 * @param[in] kinematic_viscosity Kinematic viscosity
291 * @param[in] velocity Velocity
292 * @param[in] wall_distance Distance from wall
293 * @param[in] friction_velocity_guess Guess for u_τ
294 * @param[in] pressure_gradient_tangent Tangential pressure gradient
295 * @param[in] pressure_gradient_normal Normal pressure gradient (currently unused)
296 * @return Residual value
297 */
298double f_Cabot(double kinematic_viscosity, double velocity, double wall_distance,
299 double friction_velocity_guess, double pressure_gradient_tangent,
300 double pressure_gradient_normal);
301
302/**
303 * @brief Numerical derivative for Cabot wall function
304 */
305double df_Cabot(double kinematic_viscosity, double velocity, double wall_distance,
306 double friction_velocity_guess, double pressure_gradient_tangent,
307 double pressure_gradient_normal);
308
309/**
310 * @brief Solves for friction velocity using Cabot wall function
311 *
312 * This non-equilibrium wall function accounts for pressure gradient effects,
313 * making it more accurate in separated or strongly accelerating/decelerating flows.
314 *
315 * @param[in] kinematic_viscosity Kinematic viscosity
316 * @param[in] velocity Velocity at reference point
317 * @param[in] wall_distance Distance from wall
318 * @param[in] initial_guess Initial guess for u_τ
319 * @param[in] pressure_gradient_tangent Tangential pressure gradient
320 * @param[in] pressure_gradient_normal Normal pressure gradient
321 * @param[out] friction_velocity Converged friction velocity
322 * @param[out] wall_shear_velocity Wall shear stress for velocity
323 * @param[out] wall_shear_normal Wall shear stress for normal pressure gradient
324 */
325void find_utau_Cabot(double kinematic_viscosity, double velocity, double wall_distance,
326 double initial_guess, double pressure_gradient_tangent,
327 double pressure_gradient_normal, double *friction_velocity,
328 double *wall_shear_velocity, double *wall_shear_normal);
329
330/**
331 * @brief Computes velocity using Werner-Wengle wall function
332 *
333 * Algebraic wall function that provides explicit relation:
334 * u+ = y+ for y+ < 11.81 (viscous sublayer)
335 * u+ = A * (y+)^B for y+ ≥ 11.81 (power law)
336 * where A = 8.3, B = 1/7 are empirical constants.
337 *
338 * ADVANTAGES:
339 * - Explicit (no iteration required)
340 * - Computationally cheap
341 * - Reasonable accuracy for equilibrium boundary layers
342 *
343 * LIMITATIONS:
344 * - Less accurate than log-law for high Reynolds numbers
345 * - No pressure gradient effects
346 *
347 * @param[in] kinematic_viscosity Kinematic viscosity
348 * @param[in] wall_distance Distance from wall
349 * @param[in] friction_velocity Friction velocity
350 * @return Tangential velocity
351 */
352double u_Werner(double kinematic_viscosity, double wall_distance,
353 double friction_velocity);
354
355
356/**
357 * @brief Residual function for Werner-Wengle iteration
358 *
359 * Computes residual: f(u_τ) = u_τ² - g(u, y, ν)
360 * where g is derived from the velocity profile inversion.
361 *
362 * @param[in] kinematic_viscosity Kinematic viscosity
363 * @param[in] velocity Known velocity
364 * @param[in] wall_distance Distance from wall
365 * @param[in] friction_velocity Guess for friction velocity
366 * @return Residual value
367 */
368double f_Werner(double kinematic_viscosity, double velocity,
369 double wall_distance, double friction_velocity);
370
371/**
372 * @brief Numerical derivative for Werner-Wengle iteration
373 */
374double df_Werner(double kinematic_viscosity, double velocity,
375 double wall_distance, double friction_velocity);
376
377/**
378 * @brief Solves for friction velocity using Werner-Wengle wall function
379 *
380 * @param[in] kinematic_viscosity Kinematic viscosity
381 * @param[in] velocity Velocity at reference point
382 * @param[in] wall_distance Distance from wall
383 * @param[in] initial_guess Initial guess for u_τ
384 * @return Converged friction velocity
385 */
386double find_utau_Werner(double kinematic_viscosity, double velocity,
387 double wall_distance, double initial_guess);
388
389/**
390 * @brief Computes velocity using simple log-law (smooth wall with roughness offset)
391 *
392 * Simple logarithmic profile:
393 * u = (u_τ / κ) * ln((y + y0) / y0)
394 * where y0 is a roughness length scale.
395 *
396 * @param[in] wall_distance Distance from wall
397 * @param[in] friction_velocity Friction velocity
398 * @param[in] roughness_length Roughness length scale y0
399 * @return Velocity
400 *
401 * @note This is a simplified model; use u_hydset_roughness for better accuracy
402 */
403double u_loglaw(double wall_distance, double friction_velocity, double roughness_length);
404
405/**
406 * @brief Solves for friction velocity using simple log-law (explicit formula)
407 *
408 * Explicit inversion: u_τ = κ * u / ln((y + y0) / y0)
409 *
410 * @param[in] velocity Known velocity
411 * @param[in] wall_distance Distance from wall
412 * @param[in] roughness_length Roughness length scale
413 * @return Friction velocity
414 */
415double find_utau_loglaw(double velocity, double wall_distance, double roughness_length);
416
417/**
418 * @brief Applies standard wall function with Werner-Wengle model
419 *
420 * This is a high-level interface that:
421 * 1. Decomposes velocity into normal/tangential components
422 * 2. Applies wall function to tangential velocity
423 * 3. Interpolates normal velocity
424 * 4. Reconstructs total velocity
425 *
426 * @param[in] user Simulation context
427 * @param[in] distance_reference Distance to reference point
428 * @param[in] distance_boundary Distance to boundary point
429 * @param[in] velocity_wall Wall velocity
430 * @param[in] velocity_reference Reference velocity
431 * @param[out] velocity_boundary Output boundary velocity
432 * @param[out] friction_velocity Output friction velocity
433 * @param[in] normal_x X-component of wall normal
434 * @param[in] normal_y Y-component of wall normal
435 * @param[in] normal_z Z-component of wall normal
436 */
437void wall_function(UserCtx *user, double distance_reference, double distance_boundary,
438 Cmpnts velocity_wall, Cmpnts velocity_reference,
439 Cmpnts *velocity_boundary, PetscReal *friction_velocity,
440 double normal_x, double normal_y, double normal_z);
441
442/**
443 * @brief Applies log-law wall function with roughness correction
444 *
445 * This is the recommended wall function interface for most applications.
446 * It uses the roughness-corrected log-law which is accurate for both
447 * smooth and rough walls.
448 *
449 * @param[in] user Simulation context
450 * @param[in] roughness_height Wall roughness height ks
451 * @param[in] distance_reference Distance to reference point
452 * @param[in] distance_boundary Distance to boundary point
453 * @param[in] velocity_wall Wall velocity (typically zero)
454 * @param[in] velocity_reference Reference velocity
455 * @param[out] velocity_boundary Output boundary velocity
456 * @param[out] friction_velocity Output friction velocity u_τ
457 * @param[in] normal_x X-component of wall normal
458 * @param[in] normal_y Y-component of wall normal
459 * @param[in] normal_z Z-component of wall normal
460 *
461 * @note Sets velocity to reference value if y+ > 300 (outside log-layer)
462 */
463void wall_function_loglaw(UserCtx *user, double roughness_height,
464 double distance_reference, double distance_boundary,
465 Cmpnts velocity_wall, Cmpnts velocity_reference,
466 Cmpnts *velocity_boundary, PetscReal *friction_velocity,
467 double normal_x, double normal_y, double normal_z);
468
469/**
470 * @brief Applies Cabot non-equilibrium wall function with pressure gradients
471 *
472 * This wall function accounts for pressure gradient effects and is most
473 * accurate in separated or strongly accelerating/decelerating flows.
474 *
475 * @param[in] user Simulation context
476 * @param[in] roughness_height Wall roughness (currently unused in this version)
477 * @param[in] distance_reference Distance to reference point
478 * @param[in] distance_boundary Distance to boundary point
479 * @param[in] velocity_wall Wall velocity
480 * @param[in] velocity_reference Reference velocity
481 * @param[out] velocity_boundary Output boundary velocity
482 * @param[out] friction_velocity Output friction velocity
483 * @param[in] normal_x X-component of wall normal
484 * @param[in] normal_y Y-component of wall normal
485 * @param[in] normal_z Z-component of wall normal
486 * @param[in] pressure_gradient_x X-component of pressure gradient
487 * @param[in] pressure_gradient_y Y-component of pressure gradient
488 * @param[in] pressure_gradient_z Z-component of pressure gradient
489 * @param[in] iteration_count Current iteration count
490 *
491 * @note Only updates friction velocity every 5 iterations for stability
492 * @note Normal pressure gradient currently set to zero
493 */
494void wall_function_Cabot(UserCtx *user, double roughness_height,
495 double distance_reference, double distance_boundary,
496 Cmpnts velocity_wall, Cmpnts velocity_reference,
497 Cmpnts *velocity_boundary, PetscReal *friction_velocity,
498 double normal_x, double normal_y, double normal_z,
499 double pressure_gradient_x, double pressure_gradient_y,
500 double pressure_gradient_z, int iteration_count);
501
502#endif // WALLFUNCTION_H
Header file for Particle Motion and migration related functions.
Header file for Particle Swarm management functions.
Public interface for grid, solver, and metric setup routines.
Public interface for data input/output routines.
Logging utilities and macros for PETSc-based applications.
Main header file for a complex fluid dynamics solver.
A 3D point or vector with PetscScalar components.
Definition variables.h:100
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.
double df_Cabot(double kinematic_viscosity, double velocity, double wall_distance, double friction_velocity_guess, double pressure_gradient_tangent, double pressure_gradient_normal)
Numerical derivative for Cabot wall function.
double f_Cabot(double kinematic_viscosity, double velocity, double wall_distance, double friction_velocity_guess, double pressure_gradient_tangent, double pressure_gradient_normal)
Residual function for Cabot wall function.
void wall_function_loglaw(UserCtx *user, double roughness_height, double distance_reference, double distance_boundary, Cmpnts velocity_wall, Cmpnts velocity_reference, Cmpnts *velocity_boundary, PetscReal *friction_velocity, double normal_x, double normal_y, double normal_z)
Applies log-law wall function with roughness correction.
double find_utau_loglaw(double velocity, double wall_distance, double roughness_length)
Solves for friction velocity using simple log-law (explicit formula)
double u_Werner(double kinematic_viscosity, double wall_distance, double friction_velocity)
Computes velocity using Werner-Wengle wall function.
double u_hydset_roughness(double kinematic_viscosity, double wall_distance, double friction_velocity, double roughness_height)
Computes velocity from log-law for rough walls.
void wall_function(UserCtx *user, double distance_reference, double distance_boundary, Cmpnts velocity_wall, Cmpnts velocity_reference, Cmpnts *velocity_boundary, PetscReal *friction_velocity, double normal_x, double normal_y, double normal_z)
Applies standard wall function with Werner-Wengle model.
double taw(double kinematic_viscosity, double friction_velocity, double wall_distance, double velocity, double pressure_gradient_tangent)
Computes wall shear stress with pressure gradient effects.
double find_utau_hydset(double kinematic_viscosity, double known_velocity, double wall_distance, double initial_guess, double roughness_height)
Solves for friction velocity using Newton-Raphson iteration.
double u_Cabot(double kinematic_viscosity, double wall_distance, double friction_velocity, double pressure_gradient_tangent, double wall_shear_stress)
Computes velocity using Cabot wall function.
double df_Werner(double kinematic_viscosity, double velocity, double wall_distance, double friction_velocity)
Numerical derivative for Werner-Wengle iteration.
double f_Werner(double kinematic_viscosity, double velocity, double wall_distance, double friction_velocity)
Residual function for Werner-Wengle iteration.
double u_loglaw(double wall_distance, double friction_velocity, double roughness_length)
Computes velocity using simple log-law (smooth wall with roughness offset)
void wall_function_Cabot(UserCtx *user, double roughness_height, double distance_reference, double distance_boundary, Cmpnts velocity_wall, Cmpnts velocity_reference, Cmpnts *velocity_boundary, PetscReal *friction_velocity, double normal_x, double normal_y, double normal_z, double pressure_gradient_x, double pressure_gradient_y, double pressure_gradient_z, int iteration_count)
Applies Cabot non-equilibrium wall function with pressure gradients.
void noslip(UserCtx *user, double distance_reference, double distance_boundary, Cmpnts velocity_wall, Cmpnts velocity_reference, Cmpnts *velocity_boundary, double normal_x, double normal_y, double normal_z)
Applies no-slip wall boundary condition with linear interpolation.
void find_utau_Cabot(double kinematic_viscosity, double velocity, double wall_distance, double initial_guess, double pressure_gradient_tangent, double pressure_gradient_normal, double *friction_velocity, double *wall_shear_velocity, double *wall_shear_normal)
Solves for friction velocity using Cabot wall function.
double integrate_1(double kinematic_viscosity, double wall_distance, double friction_velocity, int integration_mode)
Integrates eddy viscosity profile from wall to distance y.
double f_hydset(double kinematic_viscosity, double known_velocity, double wall_distance, double friction_velocity_guess, double roughness_height)
Residual function for friction velocity equation (log-law with roughness)
double E_coeff(double friction_velocity, double roughness_height, double kinematic_viscosity)
Computes roughness-modified log-law coefficient E.
double nu_t(double yplus)
Computes turbulent eddy viscosity ratio (ν_t / ν)
double df_hydset(double kinematic_viscosity, double known_velocity, double wall_distance, double friction_velocity_guess, double roughness_height)
Numerical derivative of residual function.
double find_utau_Werner(double kinematic_viscosity, double velocity, double wall_distance, double initial_guess)
Solves for friction velocity using Werner-Wengle wall function.
void freeslip(UserCtx *user, double distance_reference, double distance_boundary, Cmpnts velocity_wall, Cmpnts velocity_reference, Cmpnts *velocity_boundary, double normal_x, double normal_y, double normal_z)
Applies free-slip wall boundary condition.