49#define VISCOUS_SUBLAYER_YPLUS 11.81
52#define ROUGHNESS_TRANSITION_YPLUS 2.25
55#define FULLY_ROUGH_YPLUS 90.0
58#define DAMPING_COEFFICIENT 19.0
69 double distance_reference,
double distance_boundary,
72 double normal_x,
double normal_y,
double normal_z)
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;
84 double interpolation_factor = distance_boundary / distance_reference;
87 (*velocity_boundary).x = interpolation_factor * delta_u;
88 (*velocity_boundary).y = interpolation_factor * delta_v;
89 (*velocity_boundary).z = interpolation_factor * delta_w;
92 (*velocity_boundary).x += velocity_wall.
x;
93 (*velocity_boundary).y += velocity_wall.
y;
94 (*velocity_boundary).z += velocity_wall.
z;
102 double distance_reference,
double distance_boundary,
104 Cmpnts *velocity_boundary,
105 double normal_x,
double normal_y,
double normal_z)
109 double wall_normal_velocity = velocity_wall.
x * normal_x +
110 velocity_wall.
y * normal_y +
111 velocity_wall.
z * normal_z;
113 double reference_normal_velocity = velocity_reference.
x * normal_x +
114 velocity_reference.
y * normal_y +
115 velocity_reference.
z * normal_z;
118 double boundary_normal_velocity = wall_normal_velocity +
119 (reference_normal_velocity - wall_normal_velocity) * (distance_boundary / distance_reference);
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;
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;
140double E_coeff(
double friction_velocity,
double roughness_height,
double kinematic_viscosity)
143 double roughness_reynolds = friction_velocity * roughness_height / kinematic_viscosity;
145 double roughness_correction;
149 roughness_correction = 0.0;
154 double max_correction =
LOGLAW_B - 8.5 + (1.0 /
KAPPA) * log(roughness_reynolds);
156 roughness_correction = max_correction *
157 sin(0.4258 * (log(roughness_reynolds) - 0.811));
161 roughness_correction =
LOGLAW_B - 8.5 + (1.0 /
KAPPA) * log(roughness_reynolds);
173 double friction_velocity,
double roughness_height)
176 double yplus = friction_velocity * wall_distance / kinematic_viscosity;
180 const double loglayer_limit = 300.0;
182 double tangential_velocity;
184 if (yplus <= viscous_limit) {
186 tangential_velocity = friction_velocity * yplus;
188 else if (yplus <= loglayer_limit) {
190 double E =
E_coeff(friction_velocity, roughness_height, kinematic_viscosity);
191 tangential_velocity = (friction_velocity /
KAPPA) * log(E * yplus);
195 tangential_velocity = -1.0;
198 return tangential_velocity;
209double f_hydset(
double kinematic_viscosity,
double known_velocity,
210 double wall_distance,
double friction_velocity_guess,
211 double roughness_height)
213 double yplus = friction_velocity_guess * wall_distance / kinematic_viscosity;
218 residual = friction_velocity_guess * yplus - known_velocity;
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;
235double df_hydset(
double kinematic_viscosity,
double known_velocity,
236 double wall_distance,
double friction_velocity_guess,
237 double roughness_height)
239 const double perturbation = 1.e-7;
241 double f_plus =
f_hydset(kinematic_viscosity, known_velocity, wall_distance,
242 friction_velocity_guess + perturbation, roughness_height);
244 double f_minus =
f_hydset(kinematic_viscosity, known_velocity, wall_distance,
245 friction_velocity_guess - perturbation, roughness_height);
247 return (f_plus - f_minus) / (2.0 * perturbation);
257 double wall_distance,
double initial_guess,
258 double roughness_height)
260 double friction_velocity = initial_guess;
261 double friction_velocity_old = initial_guess;
263 const int max_iterations = 30;
264 const double convergence_tolerance = 1.e-7;
267 for (iteration = 0; iteration < max_iterations; iteration++) {
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);
274 friction_velocity = friction_velocity_old - residual / derivative;
277 if (fabs(friction_velocity - friction_velocity_old) < convergence_tolerance) {
281 friction_velocity_old = friction_velocity;
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));
292 return friction_velocity;
309 return KAPPA * yplus * damping_function * damping_function;
318double integrate_1(
double kinematic_viscosity,
double wall_distance,
319 double friction_velocity,
int integration_mode)
321 const int num_points = 30;
322 double step_size = wall_distance / (num_points - 1);
323 double integral_values[num_points];
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);
331 if (integration_mode == 0) {
333 integral_values[i] = 1.0 / ((1.0 + eddy_viscosity_ratio) * kinematic_viscosity);
337 integral_values[i] = current_distance /
338 ((1.0 + eddy_viscosity_ratio) * kinematic_viscosity);
343 double integral_sum = 0.0;
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;
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]);
369double taw(
double kinematic_viscosity,
double friction_velocity,
370 double wall_distance,
double velocity,
double pressure_gradient_tangent)
373 double F1 =
integrate_1(kinematic_viscosity, wall_distance, friction_velocity, 0);
374 double Fy =
integrate_1(kinematic_viscosity, wall_distance, friction_velocity, 1);
377 return (1.0 / F1) * (velocity - pressure_gradient_tangent * Fy);
386double u_Cabot(
double kinematic_viscosity,
double wall_distance,
387 double friction_velocity,
double pressure_gradient_tangent,
388 double wall_shear_stress)
390 double F1 =
integrate_1(kinematic_viscosity, wall_distance, friction_velocity, 0);
391 double Fy =
integrate_1(kinematic_viscosity, wall_distance, friction_velocity, 1);
393 return wall_shear_stress * F1 + pressure_gradient_tangent * Fy;
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)
404 (void)pressure_gradient_normal;
405 double wall_shear =
taw(kinematic_viscosity, friction_velocity_guess,
406 wall_distance, velocity, pressure_gradient_tangent);
408 return friction_velocity_guess - sqrt(fabs(wall_shear));
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)
421 const double perturbation = 1.e-7;
423 double f_plus =
f_Cabot(kinematic_viscosity, velocity, wall_distance,
424 friction_velocity_guess + perturbation,
425 pressure_gradient_tangent, pressure_gradient_normal);
427 double f_minus =
f_Cabot(kinematic_viscosity, velocity, wall_distance,
428 friction_velocity_guess - perturbation,
429 pressure_gradient_tangent, pressure_gradient_normal);
431 return (f_plus - f_minus) / (2.0 * perturbation);
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)
445 double current_guess = initial_guess;
448 const int max_iterations = 30;
449 const double convergence_tolerance = 1.e-10;
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);
457 double derivative =
df_Cabot(kinematic_viscosity, velocity, wall_distance,
458 current_guess, pressure_gradient_tangent,
459 pressure_gradient_normal);
461 new_guess = fabs(current_guess - residual / derivative);
463 if (fabs(current_guess - new_guess) < convergence_tolerance) {
467 current_guess = new_guess;
470 if (fabs(current_guess - new_guess) > 1.e-5 && iteration >= 29) {
471 printf(
"\n!!!!!!!! Cabot Iteration Failed !!!!!!!!\n");
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);
489double u_Werner(
double kinematic_viscosity,
double wall_distance,
490 double friction_velocity)
492 double yplus = friction_velocity * wall_distance / kinematic_viscosity;
495 const double power_law_coefficient = 8.3;
496 const double power_law_exponent = 1.0 / 7.0;
502 velocity = yplus * friction_velocity;
506 velocity = power_law_coefficient * pow(yplus, power_law_exponent) * friction_velocity;
516double f_Werner(
double kinematic_viscosity,
double velocity,
517 double wall_distance,
double friction_velocity)
519 const double power_law_coefficient = 8.3;
520 const double power_law_exponent = 1.0 / 7.0;
526 double transition_velocity = kinematic_viscosity / (2.0 * transition_distance) *
527 pow(power_law_coefficient, 2.0 / (1.0 - power_law_exponent));
531 if (fabs(velocity) <= transition_velocity) {
533 residual = friction_velocity * friction_velocity -
534 velocity / wall_distance * kinematic_viscosity;
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);
543 double term2 = (1.0 + power_law_exponent) / power_law_coefficient *
544 pow(kinematic_viscosity / wall_distance, power_law_exponent) *
547 residual = friction_velocity * friction_velocity -
548 velocity / fabs(velocity) * pow(term1 + term2, 2.0 / (1.0 + power_law_exponent));
560double df_Werner(
double kinematic_viscosity,
double velocity,
561 double wall_distance,
double friction_velocity)
563 const double perturbation = 1.e-7;
565 double f_plus =
f_Werner(kinematic_viscosity, velocity, wall_distance,
566 friction_velocity + perturbation);
568 double f_minus =
f_Werner(kinematic_viscosity, velocity, wall_distance,
569 friction_velocity - perturbation);
571 return (f_plus - f_minus) / (2.0 * perturbation);
581 double wall_distance,
double initial_guess)
583 double current_guess = initial_guess;
586 const int max_iterations = 20;
587 const double convergence_tolerance = 1.e-7;
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);
594 new_guess = current_guess - residual / derivative;
596 if (fabs(current_guess - new_guess) < convergence_tolerance) {
600 current_guess = new_guess;
603 if (fabs(current_guess - new_guess) > 1.e-5 && iteration >= 19) {
604 printf(
"\n!!!!!!!! Werner-Wengle Iteration Failed !!!!!!!!\n");
620double u_loglaw(
double wall_distance,
double friction_velocity,
double roughness_length)
622 return friction_velocity * (1.0 /
KAPPA) * log((roughness_length + wall_distance) / roughness_length);
631 return KAPPA * velocity / log((wall_distance + roughness_length) / roughness_length);
650 Cmpnts *velocity_boundary, PetscReal *friction_velocity,
651 double normal_x,
double normal_y,
double normal_z)
654 double kinematic_viscosity = 1.0 / simCtx->
ren;
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;
661 double normal_velocity = delta_u * normal_x + delta_v * normal_y + delta_w * normal_z;
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;
667 double tangential_magnitude = sqrt(tangential_u * tangential_u +
668 tangential_v * tangential_v +
669 tangential_w * tangential_w);
672 double utau_guess = 0.05;
673 double tangential_modeled =
u_Werner(kinematic_viscosity, distance_boundary, utau_guess);
674 if (friction_velocity) {
676 *friction_velocity = (PetscReal)utau_guess;
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;
686 tangential_u = tangential_v = tangential_w = 0.0;
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;
695 (*velocity_boundary).x += velocity_wall.
x;
696 (*velocity_boundary).y += velocity_wall.
y;
697 (*velocity_boundary).z += velocity_wall.
z;
705 double distance_reference,
double distance_boundary,
707 Cmpnts *velocity_boundary, PetscReal *friction_velocity,
708 double normal_x,
double normal_y,
double normal_z)
711 double kinematic_viscosity = 1.0 / simCtx->
ren;
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;
718 double normal_velocity = delta_u * normal_x + delta_v * normal_y + delta_w * normal_z;
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;
724 double tangential_magnitude = sqrt(tangential_u * tangential_u +
725 tangential_v * tangential_v +
726 tangential_w * tangential_w);
729 *friction_velocity =
find_utau_hydset(kinematic_viscosity, tangential_magnitude,
730 distance_reference, 0.001, roughness_height);
733 double tangential_modeled =
u_hydset_roughness(kinematic_viscosity, distance_boundary,
734 *friction_velocity, roughness_height);
737 if (tangential_modeled < 0.0) {
738 tangential_modeled = tangential_magnitude;
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;
748 tangential_u = tangential_v = tangential_w = 0.0;
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;
757 (*velocity_boundary).x += velocity_wall.
x;
758 (*velocity_boundary).y += velocity_wall.
y;
759 (*velocity_boundary).z += velocity_wall.
z;
767 double distance_reference,
double distance_boundary,
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)
774 (void)roughness_height;
776 double kinematic_viscosity = 1.0 / simCtx->
ren;
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;
783 double normal_velocity = delta_u * normal_x + delta_v * normal_y + delta_w * normal_z;
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;
789 double tangential_magnitude = sqrt(tangential_u * tangential_u +
790 tangential_v * tangential_v +
791 tangential_w * tangential_w);
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;
802 double pressure_gradient_normal = 0.0;
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;
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));
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;
827 tangential_u = tangential_v = tangential_w = 0.0;
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;
836 (*velocity_boundary).x += velocity_wall.
x;
837 (*velocity_boundary).y += velocity_wall.
y;
838 (*velocity_boundary).z += velocity_wall.
z;
#define LOCAL
Logging scope definitions for controlling message output.
#define LOG_ALLOW(scope, level, fmt,...)
Logging macro that checks both the log level and whether the calling function is in the allowed-funct...
@ LOG_DEBUG
Detailed debugging information.
SimCtx * simCtx
Back-pointer to the master simulation context.
A 3D point or vector with PetscScalar components.
The master context for the entire simulation.
User-defined context containing data specific to a single computational grid level.
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().