PICurv 0.1.0
A Parallel Particle-In-Cell Solver for Curvilinear LES
Loading...
Searching...
No Matches
wallfunction.c
Go to the documentation of this file.
1/**
2 * @file wallfunction.c
3 * @brief Wall function implementations for near-wall turbulence modeling
4 *
5 * This file contains various wall function models that bridge the gap between
6 * the wall and the first computational grid point. Wall functions allow simulations
7 * to avoid resolving the viscous sublayer, reducing computational cost while
8 * maintaining reasonable accuracy for turbulent wall-bounded flows.
9 *
10 * PHYSICAL BACKGROUND:
11 * The near-wall region in turbulent flows is characterized by three layers:
12 * 1. Viscous sublayer (y+ < 5): Dominated by viscous effects, u+ = y+
13 * 2. Buffer layer (5 < y+ < 30): Transition region
14 * 3. Log-law region (y+ > 30): Inertial layer, u+ = (1/κ) ln(y+) + B
15 * where:
16 * u+ = u / u_τ (normalized velocity)
17 * y+ = y * u_τ / ν (normalized wall distance)
18 * u_τ = sqrt(τ_w / ρ) (friction velocity)
19 * κ ≈ 0.41 (von Karman constant)
20 * B ≈ 5.5 (log-law intercept for smooth walls)
21 *
22 * IMPLEMENTED MODELS:
23 * - No-slip: Linear interpolation for stationary walls
24 * - Free-slip: Normal interpolation, tangential extrapolation
25 * - Werner-Wengle: Algebraic wall function
26 * - Log-law: Standard equilibrium wall function with roughness
27 * - Cabot: Non-equilibrium wall function with pressure gradient effects
28 *
29 * @author Original implementation adapted from legacy code
30 * @date Enhanced with comprehensive documentation
31 */
32
33#include "wallfunction.h"
34#include <math.h>
35#include <stdlib.h>
36#include <stdio.h>
37
38// ============================================================================
39// PHYSICAL CONSTANTS
40// ============================================================================
41
42/** @brief von Karman constant (universal turbulence constant) */
43#define KAPPA 0.41
44
45/** @brief Log-law intercept constant B for smooth walls */
46#define LOGLAW_B 5.5
47
48/** @brief Viscous sublayer thickness y+ threshold */
49#define VISCOUS_SUBLAYER_YPLUS 11.81
50
51/** @brief Smooth-to-rough transition y+ threshold */
52#define ROUGHNESS_TRANSITION_YPLUS 2.25
53
54/** @brief Fully rough regime y+ threshold */
55#define FULLY_ROUGH_YPLUS 90.0
56
57/** @brief Eddy viscosity damping coefficient (van Driest damping) */
58#define DAMPING_COEFFICIENT 19.0
59
60// ============================================================================
61// BASIC BOUNDARY CONDITION FUNCTIONS
62// ============================================================================
63
64/**
65 * @brief Internal helper implementation: `noslip()`.
66 * @details Local to this translation unit.
67 */
68void noslip(UserCtx *user,
69 double distance_reference, double distance_boundary,
70 Cmpnts velocity_wall, Cmpnts velocity_reference,
71 Cmpnts *velocity_boundary,
72 double normal_x, double normal_y, double normal_z)
73{
74 (void)user;
75 (void)normal_x;
76 (void)normal_y;
77 (void)normal_z;
78 // Compute velocity difference between reference point and wall
79 double delta_u = velocity_reference.x - velocity_wall.x;
80 double delta_v = velocity_reference.y - velocity_wall.y;
81 double delta_w = velocity_reference.z - velocity_wall.z;
82
83 // Linear interpolation factor
84 double interpolation_factor = distance_boundary / distance_reference;
85
86 // Apply linear interpolation
87 (*velocity_boundary).x = interpolation_factor * delta_u;
88 (*velocity_boundary).y = interpolation_factor * delta_v;
89 (*velocity_boundary).z = interpolation_factor * delta_w;
90
91 // Add wall velocity (shift to absolute frame)
92 (*velocity_boundary).x += velocity_wall.x;
93 (*velocity_boundary).y += velocity_wall.y;
94 (*velocity_boundary).z += velocity_wall.z;
95}
96
97/**
98 * @brief Internal helper implementation: `freeslip()`.
99 * @details Local to this translation unit.
100 */
101void freeslip(UserCtx *user,
102 double distance_reference, double distance_boundary,
103 Cmpnts velocity_wall, Cmpnts velocity_reference,
104 Cmpnts *velocity_boundary,
105 double normal_x, double normal_y, double normal_z)
106{
107 (void)user;
108 // Extract normal components of velocities
109 double wall_normal_velocity = velocity_wall.x * normal_x +
110 velocity_wall.y * normal_y +
111 velocity_wall.z * normal_z;
112
113 double reference_normal_velocity = velocity_reference.x * normal_x +
114 velocity_reference.y * normal_y +
115 velocity_reference.z * normal_z;
116
117 // Interpolate normal velocity component
118 double boundary_normal_velocity = wall_normal_velocity +
119 (reference_normal_velocity - wall_normal_velocity) * (distance_boundary / distance_reference);
120
121 // Extract tangential velocity components (extrapolated from reference point)
122 double tangential_u = velocity_reference.x - reference_normal_velocity * normal_x;
123 double tangential_v = velocity_reference.y - reference_normal_velocity * normal_y;
124 double tangential_w = velocity_reference.z - reference_normal_velocity * normal_z;
125
126 // Reconstruct total velocity: U = U_t + U_n * n
127 (*velocity_boundary).x = tangential_u + boundary_normal_velocity * normal_x;
128 (*velocity_boundary).y = tangential_v + boundary_normal_velocity * normal_y;
129 (*velocity_boundary).z = tangential_w + boundary_normal_velocity * normal_z;
130}
131
132// ============================================================================
133// ROUGHNESS-CORRECTED LOG-LAW FUNCTIONS
134// ============================================================================
135
136/**
137 * @brief Internal helper implementation: `E_coeff()`.
138 * @details Local to this translation unit.
139 */
140double E_coeff(double friction_velocity, double roughness_height, double kinematic_viscosity)
141{
142 // Compute roughness Reynolds number
143 double roughness_reynolds = friction_velocity * roughness_height / kinematic_viscosity;
144
145 double roughness_correction;
146
147 if (roughness_reynolds <= ROUGHNESS_TRANSITION_YPLUS) {
148 // Hydraulically smooth regime: no correction needed
149 roughness_correction = 0.0;
150 }
151 else if (roughness_reynolds < FULLY_ROUGH_YPLUS) {
152 // Transitional regime: smooth interpolation using sine function
153 // This provides a gradual transition between smooth and rough wall behavior
154 double max_correction = LOGLAW_B - 8.5 + (1.0 / KAPPA) * log(roughness_reynolds);
155
156 roughness_correction = max_correction *
157 sin(0.4258 * (log(roughness_reynolds) - 0.811));
158 }
159 else {
160 // Fully rough regime: maximum correction
161 roughness_correction = LOGLAW_B - 8.5 + (1.0 / KAPPA) * log(roughness_reynolds);
162 }
163
164 // Return E = exp[κ(B - ΔB)]
165 return exp(KAPPA * (LOGLAW_B - roughness_correction));
166}
167
168/**
169 * @brief Internal helper implementation: `u_hydset_roughness()`.
170 * @details Local to this translation unit.
171 */
172double u_hydset_roughness(double kinematic_viscosity, double wall_distance,
173 double friction_velocity, double roughness_height)
174{
175 // Normalized wall distance
176 double yplus = friction_velocity * wall_distance / kinematic_viscosity;
177
178 // Threshold values for regime transitions
179 const double viscous_limit = VISCOUS_SUBLAYER_YPLUS;
180 const double loglayer_limit = 300.0;
181
182 double tangential_velocity;
183
184 if (yplus <= viscous_limit) {
185 // Viscous sublayer: linear velocity profile
186 tangential_velocity = friction_velocity * yplus;
187 }
188 else if (yplus <= loglayer_limit) {
189 // Log-layer: logarithmic velocity profile with roughness correction
190 double E = E_coeff(friction_velocity, roughness_height, kinematic_viscosity);
191 tangential_velocity = (friction_velocity / KAPPA) * log(E * yplus);
192 }
193 else {
194 // Outside valid range (wake region or beyond)
195 tangential_velocity = -1.0;
196 }
197
198 return tangential_velocity;
199}
200
201// ============================================================================
202// NEWTON-RAPHSON SOLVER FOR FRICTION VELOCITY (LOG-LAW)
203// ============================================================================
204
205/**
206 * @brief Internal helper implementation: `f_hydset()`.
207 * @details Local to this translation unit.
208 */
209double f_hydset(double kinematic_viscosity, double known_velocity,
210 double wall_distance, double friction_velocity_guess,
211 double roughness_height)
212{
213 double yplus = friction_velocity_guess * wall_distance / kinematic_viscosity;
214 double residual;
215
216 if (yplus <= VISCOUS_SUBLAYER_YPLUS) {
217 // Viscous sublayer: u = u_τ * y+
218 residual = friction_velocity_guess * yplus - known_velocity;
219 }
220 else {
221 // Log-layer: u = (u_τ / κ) * ln(E * y+)
222 double E = E_coeff(friction_velocity_guess, roughness_height, kinematic_viscosity);
223 residual = friction_velocity_guess * (1.0 / KAPPA * log(E * yplus)) - known_velocity;
224 }
225
226 return residual;
227}
228
229/**
230 * @brief Implementation of \ref df_hydset().
231 * @details Full API contract (arguments, ownership, side effects) is documented with
232 * the header declaration in `include/wallfunction.h`.
233 * @see df_hydset()
234 */
235double df_hydset(double kinematic_viscosity, double known_velocity,
236 double wall_distance, double friction_velocity_guess,
237 double roughness_height)
238{
239 const double perturbation = 1.e-7;
240
241 double f_plus = f_hydset(kinematic_viscosity, known_velocity, wall_distance,
242 friction_velocity_guess + perturbation, roughness_height);
243
244 double f_minus = f_hydset(kinematic_viscosity, known_velocity, wall_distance,
245 friction_velocity_guess - perturbation, roughness_height);
246
247 return (f_plus - f_minus) / (2.0 * perturbation);
248}
249
250/**
251 * @brief Implementation of \ref find_utau_hydset().
252 * @details Full API contract (arguments, ownership, side effects) is documented with
253 * the header declaration in `include/wallfunction.h`.
254 * @see find_utau_hydset()
255 */
256double find_utau_hydset(double kinematic_viscosity, double known_velocity,
257 double wall_distance, double initial_guess,
258 double roughness_height)
259{
260 double friction_velocity = initial_guess;
261 double friction_velocity_old = initial_guess;
262
263 const int max_iterations = 30;
264 const double convergence_tolerance = 1.e-7;
265
266 int iteration;
267 for (iteration = 0; iteration < max_iterations; iteration++) {
268 // Newton-Raphson update
269 double residual = f_hydset(kinematic_viscosity, known_velocity, wall_distance,
270 friction_velocity_old, roughness_height);
271 double derivative = df_hydset(kinematic_viscosity, known_velocity, wall_distance,
272 friction_velocity_old, roughness_height);
273
274 friction_velocity = friction_velocity_old - residual / derivative;
275
276 // Check convergence
277 if (fabs(friction_velocity - friction_velocity_old) < convergence_tolerance) {
278 break;
279 }
280
281 friction_velocity_old = friction_velocity;
282 }
283
284 // Warn if convergence not achieved
285 if (iteration == max_iterations) {
287 "WARNING: u_tau iteration did not converge after %d iterations. "
288 "Final difference: %le\n",
289 iteration, fabs(friction_velocity - friction_velocity_old));
290 }
291
292 return friction_velocity;
293}
294
295// ============================================================================
296// TURBULENT EDDY VISCOSITY MODELS
297// ============================================================================
298
299/**
300 * @brief Internal helper implementation: `nu_t()`.
301 * @details Local to this translation unit.
302 */
303double nu_t(double yplus)
304{
305 // van Driest damping function: [1 - exp(-y+ / A+)]²
306 double damping_function = 1.0 - exp(-yplus / DAMPING_COEFFICIENT);
307
308 // Mixing length model: ν_t / ν = κ * y+ * D²
309 return KAPPA * yplus * damping_function * damping_function;
310}
311
312/**
313 * @brief Implementation of \ref integrate_1().
314 * @details Full API contract (arguments, ownership, side effects) is documented with
315 * the header declaration in `include/wallfunction.h`.
316 * @see integrate_1()
317 */
318double integrate_1(double kinematic_viscosity, double wall_distance,
319 double friction_velocity, int integration_mode)
320{
321 const int num_points = 30;
322 double step_size = wall_distance / (num_points - 1);
323 double integral_values[num_points];
324
325 // Compute integrand at each point
326 for (int i = 0; i < num_points; i++) {
327 double current_distance = i * step_size;
328 double yplus = current_distance * friction_velocity / kinematic_viscosity;
329 double eddy_viscosity_ratio = nu_t(yplus);
330
331 if (integration_mode == 0) {
332 // Mode 0: integrand = 1 / [(1 + ν_t/ν) * ν]
333 integral_values[i] = 1.0 / ((1.0 + eddy_viscosity_ratio) * kinematic_viscosity);
334 }
335 else {
336 // Mode 1: integrand = y / [(1 + ν_t/ν) * ν]
337 integral_values[i] = current_distance /
338 ((1.0 + eddy_viscosity_ratio) * kinematic_viscosity);
339 }
340 }
341
342 // Trapezoidal rule integration with corrected endpoints
343 double integral_sum = 0.0;
344
345 // Interior points: weight = (f[i-1] + 2*f[i] + f[i+1]) * dy/4
346 for (int i = 1; i < num_points - 1; i++) {
347 integral_sum += (integral_values[i-1] + 2.0 * integral_values[i] +
348 integral_values[i+1]) * 0.25 * step_size;
349 }
350
351 // Endpoint correction for higher accuracy
352 integral_sum += 0.1667 * step_size *
353 (2.0 * integral_values[0] + integral_values[1] +
354 2.0 * integral_values[num_points-1] + integral_values[num_points-2]);
355
356 return integral_sum;
357}
358
359// ============================================================================
360// CABOT NON-EQUILIBRIUM WALL FUNCTION
361// ============================================================================
362
363/**
364 * @brief Implementation of \ref taw().
365 * @details Full API contract (arguments, ownership, side effects) is documented with
366 * the header declaration in `include/wallfunction.h`.
367 * @see taw()
368 */
369double taw(double kinematic_viscosity, double friction_velocity,
370 double wall_distance, double velocity, double pressure_gradient_tangent)
371{
372 // Compute integrated viscosity profiles
373 double F1 = integrate_1(kinematic_viscosity, wall_distance, friction_velocity, 0);
374 double Fy = integrate_1(kinematic_viscosity, wall_distance, friction_velocity, 1);
375
376 // Apply momentum balance: τ_w * F1 = u - (dp/dx) * Fy
377 return (1.0 / F1) * (velocity - pressure_gradient_tangent * Fy);
378}
379
380/**
381 * @brief Implementation of \ref u_Cabot().
382 * @details Full API contract (arguments, ownership, side effects) is documented with
383 * the header declaration in `include/wallfunction.h`.
384 * @see u_Cabot()
385 */
386double u_Cabot(double kinematic_viscosity, double wall_distance,
387 double friction_velocity, double pressure_gradient_tangent,
388 double wall_shear_stress)
389{
390 double F1 = integrate_1(kinematic_viscosity, wall_distance, friction_velocity, 0);
391 double Fy = integrate_1(kinematic_viscosity, wall_distance, friction_velocity, 1);
392
393 return wall_shear_stress * F1 + pressure_gradient_tangent * Fy;
394}
395
396/**
397 * @brief Internal helper implementation: `f_Cabot()`.
398 * @details Local to this translation unit.
399 */
400double f_Cabot(double kinematic_viscosity, double velocity, double wall_distance,
401 double friction_velocity_guess, double pressure_gradient_tangent,
402 double pressure_gradient_normal)
403{
404 (void)pressure_gradient_normal;
405 double wall_shear = taw(kinematic_viscosity, friction_velocity_guess,
406 wall_distance, velocity, pressure_gradient_tangent);
407
408 return friction_velocity_guess - sqrt(fabs(wall_shear));
409}
410
411/**
412 * @brief Implementation of \ref df_Cabot().
413 * @details Full API contract (arguments, ownership, side effects) is documented with
414 * the header declaration in `include/wallfunction.h`.
415 * @see df_Cabot()
416 */
417double df_Cabot(double kinematic_viscosity, double velocity, double wall_distance,
418 double friction_velocity_guess, double pressure_gradient_tangent,
419 double pressure_gradient_normal)
420{
421 const double perturbation = 1.e-7;
422
423 double f_plus = f_Cabot(kinematic_viscosity, velocity, wall_distance,
424 friction_velocity_guess + perturbation,
425 pressure_gradient_tangent, pressure_gradient_normal);
426
427 double f_minus = f_Cabot(kinematic_viscosity, velocity, wall_distance,
428 friction_velocity_guess - perturbation,
429 pressure_gradient_tangent, pressure_gradient_normal);
430
431 return (f_plus - f_minus) / (2.0 * perturbation);
432}
433
434/**
435 * @brief Implementation of \ref find_utau_Cabot().
436 * @details Full API contract (arguments, ownership, side effects) is documented with
437 * the header declaration in `include/wallfunction.h`.
438 * @see find_utau_Cabot()
439 */
440void find_utau_Cabot(double kinematic_viscosity, double velocity, double wall_distance,
441 double initial_guess, double pressure_gradient_tangent,
442 double pressure_gradient_normal, double *friction_velocity,
443 double *wall_shear_velocity, double *wall_shear_normal)
444{
445 double current_guess = initial_guess;
446 double new_guess;
447
448 const int max_iterations = 30;
449 const double convergence_tolerance = 1.e-10;
450
451 int iteration;
452 for (iteration = 0; iteration < max_iterations; iteration++) {
453 double residual = f_Cabot(kinematic_viscosity, velocity, wall_distance,
454 current_guess, pressure_gradient_tangent,
455 pressure_gradient_normal);
456
457 double derivative = df_Cabot(kinematic_viscosity, velocity, wall_distance,
458 current_guess, pressure_gradient_tangent,
459 pressure_gradient_normal);
460
461 new_guess = fabs(current_guess - residual / derivative);
462
463 if (fabs(current_guess - new_guess) < convergence_tolerance) {
464 break;
465 }
466
467 current_guess = new_guess;
468 }
469
470 if (fabs(current_guess - new_guess) > 1.e-5 && iteration >= 29) {
471 printf("\n!!!!!!!! Cabot Iteration Failed !!!!!!!!\n");
472 }
473
474 *friction_velocity = new_guess;
475 *wall_shear_velocity = taw(kinematic_viscosity, new_guess, wall_distance,
476 velocity, pressure_gradient_tangent);
477 *wall_shear_normal = taw(kinematic_viscosity, new_guess, wall_distance,
478 0.0, pressure_gradient_normal);
479}
480
481// ============================================================================
482// WERNER-WENGLE ALGEBRAIC WALL FUNCTION
483// ============================================================================
484
485/**
486 * @brief Internal helper implementation: `u_Werner()`.
487 * @details Local to this translation unit.
488 */
489double u_Werner(double kinematic_viscosity, double wall_distance,
490 double friction_velocity)
491{
492 double yplus = friction_velocity * wall_distance / kinematic_viscosity;
493
494 // Werner-Wengle constants
495 const double power_law_coefficient = 8.3;
496 const double power_law_exponent = 1.0 / 7.0;
497
498 double velocity;
499
500 if (yplus <= VISCOUS_SUBLAYER_YPLUS) {
501 // Viscous sublayer: u+ = y+
502 velocity = yplus * friction_velocity;
503 }
504 else {
505 // Power-law region: u+ = A * (y+)^B
506 velocity = power_law_coefficient * pow(yplus, power_law_exponent) * friction_velocity;
507 }
508
509 return velocity;
510}
511
512/**
513 * @brief Internal helper implementation: `f_Werner()`.
514 * @details Local to this translation unit.
515 */
516double f_Werner(double kinematic_viscosity, double velocity,
517 double wall_distance, double friction_velocity)
518{
519 const double power_law_coefficient = 8.3;
520 const double power_law_exponent = 1.0 / 7.0;
521
522 // Transition point between viscous and power-law regions
523 double transition_distance = VISCOUS_SUBLAYER_YPLUS * kinematic_viscosity / friction_velocity;
524
525 // Transition velocity
526 double transition_velocity = kinematic_viscosity / (2.0 * transition_distance) *
527 pow(power_law_coefficient, 2.0 / (1.0 - power_law_exponent));
528
529 double residual;
530
531 if (fabs(velocity) <= transition_velocity) {
532 // Viscous sublayer regime
533 residual = friction_velocity * friction_velocity -
534 velocity / wall_distance * kinematic_viscosity;
535 }
536 else {
537 // Power-law regime (more complex inversion formula)
538 double term1 = 0.5 * (1.0 - power_law_exponent) *
539 pow(power_law_coefficient, (1.0 + power_law_exponent) /
540 (1.0 - power_law_exponent)) *
541 pow(kinematic_viscosity / wall_distance, 1.0 + power_law_exponent);
542
543 double term2 = (1.0 + power_law_exponent) / power_law_coefficient *
544 pow(kinematic_viscosity / wall_distance, power_law_exponent) *
545 fabs(velocity);
546
547 residual = friction_velocity * friction_velocity -
548 velocity / fabs(velocity) * pow(term1 + term2, 2.0 / (1.0 + power_law_exponent));
549 }
550
551 return residual;
552}
553
554/**
555 * @brief Implementation of \ref df_Werner().
556 * @details Full API contract (arguments, ownership, side effects) is documented with
557 * the header declaration in `include/wallfunction.h`.
558 * @see df_Werner()
559 */
560double df_Werner(double kinematic_viscosity, double velocity,
561 double wall_distance, double friction_velocity)
562{
563 const double perturbation = 1.e-7;
564
565 double f_plus = f_Werner(kinematic_viscosity, velocity, wall_distance,
566 friction_velocity + perturbation);
567
568 double f_minus = f_Werner(kinematic_viscosity, velocity, wall_distance,
569 friction_velocity - perturbation);
570
571 return (f_plus - f_minus) / (2.0 * perturbation);
572}
573
574/**
575 * @brief Implementation of \ref find_utau_Werner().
576 * @details Full API contract (arguments, ownership, side effects) is documented with
577 * the header declaration in `include/wallfunction.h`.
578 * @see find_utau_Werner()
579 */
580double find_utau_Werner(double kinematic_viscosity, double velocity,
581 double wall_distance, double initial_guess)
582{
583 double current_guess = initial_guess;
584 double new_guess;
585
586 const int max_iterations = 20;
587 const double convergence_tolerance = 1.e-7;
588
589 int iteration;
590 for (iteration = 0; iteration < max_iterations; iteration++) {
591 double residual = f_Werner(kinematic_viscosity, velocity, wall_distance, current_guess);
592 double derivative = df_Werner(kinematic_viscosity, velocity, wall_distance, current_guess);
593
594 new_guess = current_guess - residual / derivative;
595
596 if (fabs(current_guess - new_guess) < convergence_tolerance) {
597 break;
598 }
599
600 current_guess = new_guess;
601 }
602
603 if (fabs(current_guess - new_guess) > 1.e-5 && iteration >= 19) {
604 printf("\n!!!!!!!! Werner-Wengle Iteration Failed !!!!!!!!\n");
605 }
606
607 return new_guess;
608}
609
610// ============================================================================
611// SIMPLE LOG-LAW WALL FUNCTION (NO ROUGHNESS)
612// ============================================================================
613
614/**
615 * @brief Implementation of \ref u_loglaw().
616 * @details Full API contract (arguments, ownership, side effects) is documented with
617 * the header declaration in `include/wallfunction.h`.
618 * @see u_loglaw()
619 */
620double u_loglaw(double wall_distance, double friction_velocity, double roughness_length)
621{
622 return friction_velocity * (1.0 / KAPPA) * log((roughness_length + wall_distance) / roughness_length);
623}
624
625/**
626 * @brief Internal helper implementation: `find_utau_loglaw()`.
627 * @details Local to this translation unit.
628 */
629double find_utau_loglaw(double velocity, double wall_distance, double roughness_length)
630{
631 return KAPPA * velocity / log((wall_distance + roughness_length) / roughness_length);
632}
633
634// ============================================================================
635// UTILITY FUNCTIONS
636// ============================================================================
637
638// ============================================================================
639// HIGH-LEVEL WALL FUNCTION INTERFACE FUNCTIONS
640// ============================================================================
641
642/**
643 * @brief Implementation of \ref wall_function().
644 * @details Full API contract (arguments, ownership, side effects) is documented with
645 * the header declaration in `include/wallfunction.h`.
646 * @see wall_function()
647 */
648void wall_function(UserCtx *user, double distance_reference, double distance_boundary,
649 Cmpnts velocity_wall, Cmpnts velocity_reference,
650 Cmpnts *velocity_boundary, PetscReal *friction_velocity,
651 double normal_x, double normal_y, double normal_z)
652{
653 SimCtx *simCtx = user->simCtx;
654 double kinematic_viscosity = 1.0 / simCtx->ren;
655
656 // Decompose velocity into components
657 double delta_u = velocity_reference.x - velocity_wall.x;
658 double delta_v = velocity_reference.y - velocity_wall.y;
659 double delta_w = velocity_reference.z - velocity_wall.z;
660
661 double normal_velocity = delta_u * normal_x + delta_v * normal_y + delta_w * normal_z;
662
663 double tangential_u = delta_u - normal_velocity * normal_x;
664 double tangential_v = delta_v - normal_velocity * normal_y;
665 double tangential_w = delta_w - normal_velocity * normal_z;
666
667 double tangential_magnitude = sqrt(tangential_u * tangential_u +
668 tangential_v * tangential_v +
669 tangential_w * tangential_w);
670
671 // Apply Werner-Wengle wall function
672 double utau_guess = 0.05;
673 double tangential_modeled = u_Werner(kinematic_viscosity, distance_boundary, utau_guess);
674 if (friction_velocity) {
675 // Werner-Wengle path currently uses a fixed u_tau estimate.
676 *friction_velocity = (PetscReal)utau_guess;
677 }
678
679 // Scale tangential components
680 if (tangential_magnitude > 1.e-10) {
681 tangential_u *= tangential_modeled / tangential_magnitude;
682 tangential_v *= tangential_modeled / tangential_magnitude;
683 tangential_w *= tangential_modeled / tangential_magnitude;
684 }
685 else {
686 tangential_u = tangential_v = tangential_w = 0.0;
687 }
688
689 // Reconstruct velocity: U = U_t + U_n * n
690 (*velocity_boundary).x = tangential_u + (distance_boundary / distance_reference) * normal_velocity * normal_x;
691 (*velocity_boundary).y = tangential_v + (distance_boundary / distance_reference) * normal_velocity * normal_y;
692 (*velocity_boundary).z = tangential_w + (distance_boundary / distance_reference) * normal_velocity * normal_z;
693
694 // Add wall velocity
695 (*velocity_boundary).x += velocity_wall.x;
696 (*velocity_boundary).y += velocity_wall.y;
697 (*velocity_boundary).z += velocity_wall.z;
698}
699
700/**
701 * @brief Internal helper implementation: `wall_function_loglaw()`.
702 * @details Local to this translation unit.
703 */
704void wall_function_loglaw(UserCtx *user, double roughness_height,
705 double distance_reference, double distance_boundary,
706 Cmpnts velocity_wall, Cmpnts velocity_reference,
707 Cmpnts *velocity_boundary, PetscReal *friction_velocity,
708 double normal_x, double normal_y, double normal_z)
709{
710 SimCtx *simCtx = user->simCtx;
711 double kinematic_viscosity = 1.0 / simCtx->ren;
712
713 // Decompose velocity
714 double delta_u = velocity_reference.x - velocity_wall.x;
715 double delta_v = velocity_reference.y - velocity_wall.y;
716 double delta_w = velocity_reference.z - velocity_wall.z;
717
718 double normal_velocity = delta_u * normal_x + delta_v * normal_y + delta_w * normal_z;
719
720 double tangential_u = delta_u - normal_velocity * normal_x;
721 double tangential_v = delta_v - normal_velocity * normal_y;
722 double tangential_w = delta_w - normal_velocity * normal_z;
723
724 double tangential_magnitude = sqrt(tangential_u * tangential_u +
725 tangential_v * tangential_v +
726 tangential_w * tangential_w);
727
728 // Solve for friction velocity
729 *friction_velocity = find_utau_hydset(kinematic_viscosity, tangential_magnitude,
730 distance_reference, 0.001, roughness_height);
731
732 // Compute wall function velocity
733 double tangential_modeled = u_hydset_roughness(kinematic_viscosity, distance_boundary,
734 *friction_velocity, roughness_height);
735
736 // Check if outside valid range (y+ > 300)
737 if (tangential_modeled < 0.0) {
738 tangential_modeled = tangential_magnitude;
739 }
740
741 // Scale tangential components
742 if (tangential_magnitude > 1.e-10) {
743 tangential_u *= tangential_modeled / tangential_magnitude;
744 tangential_v *= tangential_modeled / tangential_magnitude;
745 tangential_w *= tangential_modeled / tangential_magnitude;
746 }
747 else {
748 tangential_u = tangential_v = tangential_w = 0.0;
749 }
750
751 // Reconstruct total velocity
752 (*velocity_boundary).x = tangential_u + (distance_boundary / distance_reference) * normal_velocity * normal_x;
753 (*velocity_boundary).y = tangential_v + (distance_boundary / distance_reference) * normal_velocity * normal_y;
754 (*velocity_boundary).z = tangential_w + (distance_boundary / distance_reference) * normal_velocity * normal_z;
755
756 // Add wall velocity
757 (*velocity_boundary).x += velocity_wall.x;
758 (*velocity_boundary).y += velocity_wall.y;
759 (*velocity_boundary).z += velocity_wall.z;
760}
761
762/**
763 * @brief Internal helper implementation: `wall_function_Cabot()`.
764 * @details Local to this translation unit.
765 */
766void wall_function_Cabot(UserCtx *user, double roughness_height,
767 double distance_reference, double distance_boundary,
768 Cmpnts velocity_wall, Cmpnts velocity_reference,
769 Cmpnts *velocity_boundary, PetscReal *friction_velocity,
770 double normal_x, double normal_y, double normal_z,
771 double pressure_gradient_x, double pressure_gradient_y,
772 double pressure_gradient_z, int iteration_count)
773{
774 (void)roughness_height;
775 SimCtx *simCtx = user->simCtx;
776 double kinematic_viscosity = 1.0 / simCtx->ren;
777
778 // Decompose velocity
779 double delta_u = velocity_reference.x - velocity_wall.x;
780 double delta_v = velocity_reference.y - velocity_wall.y;
781 double delta_w = velocity_reference.z - velocity_wall.z;
782
783 double normal_velocity = delta_u * normal_x + delta_v * normal_y + delta_w * normal_z;
784
785 double tangential_u = delta_u - normal_velocity * normal_x;
786 double tangential_v = delta_v - normal_velocity * normal_y;
787 double tangential_w = delta_w - normal_velocity * normal_z;
788
789 double tangential_magnitude = sqrt(tangential_u * tangential_u +
790 tangential_v * tangential_v +
791 tangential_w * tangential_w);
792
793 // Compute tangential pressure gradient
794 double pressure_gradient_tangent = 0.0;
795 if (tangential_magnitude > 1.e-10) {
796 pressure_gradient_tangent = (pressure_gradient_x * tangential_u +
797 pressure_gradient_y * tangential_v +
798 pressure_gradient_z * tangential_w) / tangential_magnitude;
799 }
800
801 // Normal pressure gradient (currently set to zero)
802 double pressure_gradient_normal = 0.0;
803
804 // Solve for friction velocity (only on first iteration or periodically)
805 if (iteration_count == 0 || iteration_count > 4) {
806 double utau, wall_shear1, wall_shear2;
807 find_utau_Cabot(kinematic_viscosity, tangential_magnitude, distance_reference,
808 0.01, pressure_gradient_tangent, pressure_gradient_normal,
809 &utau, &wall_shear1, &wall_shear2);
810 *friction_velocity = utau;
811 }
812
813 // Compute wall function velocity with pressure gradient correction
814 double tangential_modeled = u_Cabot(kinematic_viscosity, distance_boundary,
815 *friction_velocity, pressure_gradient_tangent,
816 taw(kinematic_viscosity, *friction_velocity,
817 distance_reference, tangential_magnitude,
818 pressure_gradient_tangent));
819
820 // Scale tangential components
821 if (tangential_magnitude > 1.e-10) {
822 tangential_u *= tangential_modeled / tangential_magnitude;
823 tangential_v *= tangential_modeled / tangential_magnitude;
824 tangential_w *= tangential_modeled / tangential_magnitude;
825 }
826 else {
827 tangential_u = tangential_v = tangential_w = 0.0;
828 }
829
830 // Reconstruct total velocity
831 (*velocity_boundary).x = tangential_u + (distance_boundary / distance_reference) * normal_velocity * normal_x;
832 (*velocity_boundary).y = tangential_v + (distance_boundary / distance_reference) * normal_velocity * normal_y;
833 (*velocity_boundary).z = tangential_w + (distance_boundary / distance_reference) * normal_velocity * normal_z;
834
835 // Add wall velocity
836 (*velocity_boundary).x += velocity_wall.x;
837 (*velocity_boundary).y += velocity_wall.y;
838 (*velocity_boundary).z += velocity_wall.z;
839}
#define LOCAL
Logging scope definitions for controlling message output.
Definition logging.h:44
#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:199
@ LOG_DEBUG
Detailed debugging information.
Definition logging.h:31
SimCtx * simCtx
Back-pointer to the master simulation context.
Definition variables.h:814
PetscReal ren
Definition variables.h:691
PetscScalar x
Definition variables.h:101
PetscScalar z
Definition variables.h:101
PetscScalar y
Definition variables.h:101
A 3D point or vector with PetscScalar components.
Definition variables.h:100
The master context for the entire simulation.
Definition variables.h:643
User-defined context containing data specific to a single computational grid level.
Definition variables.h:811
double df_Cabot(double kinematic_viscosity, double velocity, double wall_distance, double friction_velocity_guess, double pressure_gradient_tangent, double pressure_gradient_normal)
Implementation of df_Cabot().
double f_Cabot(double kinematic_viscosity, double velocity, double wall_distance, double friction_velocity_guess, double pressure_gradient_tangent, double pressure_gradient_normal)
Internal helper implementation: f_Cabot().
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)
Internal helper implementation: wall_function_loglaw().
double find_utau_loglaw(double velocity, double wall_distance, double roughness_length)
Internal helper implementation: find_utau_loglaw().
double u_Werner(double kinematic_viscosity, double wall_distance, double friction_velocity)
Internal helper implementation: u_Werner().
double u_hydset_roughness(double kinematic_viscosity, double wall_distance, double friction_velocity, double roughness_height)
Internal helper implementation: u_hydset_roughness().
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)
Implementation of wall_function().
#define DAMPING_COEFFICIENT
Eddy viscosity damping coefficient (van Driest damping)
double taw(double kinematic_viscosity, double friction_velocity, double wall_distance, double velocity, double pressure_gradient_tangent)
Implementation of taw().
double find_utau_hydset(double kinematic_viscosity, double known_velocity, double wall_distance, double initial_guess, double roughness_height)
Implementation of find_utau_hydset().
double u_Cabot(double kinematic_viscosity, double wall_distance, double friction_velocity, double pressure_gradient_tangent, double wall_shear_stress)
Implementation of u_Cabot().
#define ROUGHNESS_TRANSITION_YPLUS
Smooth-to-rough transition y+ threshold.
double df_Werner(double kinematic_viscosity, double velocity, double wall_distance, double friction_velocity)
Implementation of df_Werner().
double f_Werner(double kinematic_viscosity, double velocity, double wall_distance, double friction_velocity)
Internal helper implementation: f_Werner().
#define VISCOUS_SUBLAYER_YPLUS
Viscous sublayer thickness y+ threshold.
double u_loglaw(double wall_distance, double friction_velocity, double roughness_length)
Implementation of u_loglaw().
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)
Internal helper implementation: wall_function_Cabot().
#define LOGLAW_B
Log-law intercept constant B for smooth walls.
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)
Internal helper implementation: noslip().
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)
Implementation of find_utau_Cabot().
double integrate_1(double kinematic_viscosity, double wall_distance, double friction_velocity, int integration_mode)
Implementation of integrate_1().
#define FULLY_ROUGH_YPLUS
Fully rough regime y+ threshold.
double f_hydset(double kinematic_viscosity, double known_velocity, double wall_distance, double friction_velocity_guess, double roughness_height)
Internal helper implementation: f_hydset().
double E_coeff(double friction_velocity, double roughness_height, double kinematic_viscosity)
Internal helper implementation: E_coeff().
#define KAPPA
von Karman constant (universal turbulence constant)
double nu_t(double yplus)
Internal helper implementation: nu_t().
double df_hydset(double kinematic_viscosity, double known_velocity, double wall_distance, double friction_velocity_guess, double roughness_height)
Implementation of df_hydset().
double find_utau_Werner(double kinematic_viscosity, double velocity, double wall_distance, double initial_guess)
Implementation of find_utau_Werner().
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)
Internal helper implementation: freeslip().