Skip to content

Commit 8f2acba

Browse files
committed
Add a minimal free-surface gradient interface to WaveBase
Introduces GetElevationGradientXY(pos, t)
1 parent 51232a3 commit 8f2acba

7 files changed

Lines changed: 111 additions & 0 deletions

File tree

include/hydroc/waves/irregular_wave.h

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -53,6 +53,7 @@ class IrregularWaves : public WaveBase {
5353
double GetElevation(const Eigen::Vector3d& position, double time) override;
5454
Eigen::Vector3d GetVelocity(const Eigen::Vector3d& position, double time) override;
5555
Eigen::Vector3d GetAcceleration(const Eigen::Vector3d& position, double time) override;
56+
Eigen::Vector2d GetElevationGradientXY(const Eigen::Vector3d& position, double time) const override;
5657

5758
private:
5859
IrregularWaveParams params_;

include/hydroc/waves/regular_wave.h

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -33,6 +33,7 @@ class RegularWave : public WaveBase {
3333
double GetElevation(const Eigen::Vector3d& position, double time) override;
3434
Eigen::Vector3d GetVelocity(const Eigen::Vector3d& position, double time) override;
3535
Eigen::Vector3d GetAcceleration(const Eigen::Vector3d& position, double time) override;
36+
Eigen::Vector2d GetElevationGradientXY(const Eigen::Vector3d& position, double time) const override;
3637

3738
private:
3839
unsigned int num_bodies_ = 0;

include/hydroc/waves/wave_base.h

Lines changed: 35 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -14,6 +14,22 @@ enum class WaveMode {
1414
irregular = 2
1515
};
1616

17+
/**
18+
* @brief Abstract base class for all wave models.
19+
*
20+
* Coordinate conventions:
21+
* - Waves propagate in the +X direction (unidirectional).
22+
* - Z is vertical (positive upward), with z = mwl_ at mean water level.
23+
* - Y is horizontal, perpendicular to wave propagation.
24+
*
25+
* Units:
26+
* - Positions: meters [m]
27+
* - Time: seconds [s]
28+
* - Elevation η: meters [m]
29+
* - Gradients ∂η/∂x, ∂η/∂y: dimensionless [m/m]
30+
* - Velocities: meters per second [m/s]
31+
* - Accelerations: meters per second squared [m/s²]
32+
*/
1733
class WaveBase {
1834
public:
1935
virtual ~WaveBase() = default;
@@ -25,6 +41,25 @@ class WaveBase {
2541
virtual Eigen::Vector3d GetVelocity(const Eigen::Vector3d& position, double time) = 0;
2642
virtual Eigen::Vector3d GetAcceleration(const Eigen::Vector3d& position, double time) = 0;
2743

44+
/**
45+
* @brief Returns the free-surface gradient (∂η/∂x, ∂η/∂y) at the given position and time.
46+
*
47+
* Used for computing surface normals in visualization. The normal vector can be
48+
* constructed as: n = normalize(−∂η/∂x, −∂η/∂y, 1).
49+
*
50+
* @param position World coordinates (x, y, z) in meters. Only x, y are used.
51+
* @param time Simulation time in seconds.
52+
* @return Eigen::Vector2d containing (∂η/∂x, ∂η/∂y), dimensionless.
53+
*
54+
* @note Default implementation returns (0, 0) for flat surface (NoWave).
55+
* @note Current wave models are unidirectional (+X), so ∂η/∂y = 0.
56+
*/
57+
virtual Eigen::Vector2d GetElevationGradientXY(const Eigen::Vector3d& position, double time) const {
58+
(void)position;
59+
(void)time;
60+
return Eigen::Vector2d::Zero();
61+
}
62+
2863
double mwl_ = 0.0;
2964
double g_ = 9.81;
3065
double water_depth_ = 0.0;

src/hydro/waves/irregular_wave.cpp

Lines changed: 7 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -189,6 +189,13 @@ double IrregularWaves::GetElevation(const Eigen::Vector3d& position, double time
189189
return GetEtaIrregular(position, time, spectrum_frequencies_, spectral_densities_, spectral_widths_, wave_phases_, wavenumbers_);
190190
}
191191

192+
Eigen::Vector2d IrregularWaves::GetElevationGradientXY(const Eigen::Vector3d& position, double time) const {
193+
// Irregular waves propagate in +X direction only; ∂η/∂y = 0.
194+
double deta_dx = GetEtaGradientXIrregular(position, time, spectrum_frequencies_, spectral_densities_,
195+
spectral_widths_, wave_phases_, wavenumbers_);
196+
return Eigen::Vector2d(deta_dx, 0.0);
197+
}
198+
192199
Eigen::VectorXd IrregularWaves::GetForceAtTime(double t) {
193200
unsigned int total_dofs = params_.num_bodies_ * 6;
194201
Eigen::VectorXd f(total_dofs);

src/hydro/waves/regular_wave.cpp

Lines changed: 7 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -77,6 +77,13 @@ double RegularWave::GetElevation(const Eigen::Vector3d& position, double time) {
7777
return GetEta(position, time, regular_wave_omega_, regular_wave_amplitude_, regular_wave_phase_, wavenumber_);
7878
}
7979

80+
Eigen::Vector2d RegularWave::GetElevationGradientXY(const Eigen::Vector3d& position, double time) const {
81+
// Regular waves propagate in +X direction only; ∂η/∂y = 0.
82+
double deta_dx = GetEtaGradientX(position, time, regular_wave_omega_, regular_wave_amplitude_,
83+
regular_wave_phase_, wavenumber_);
84+
return Eigen::Vector2d(deta_dx, 0.0);
85+
}
86+
8087
Eigen::VectorXd RegularWave::GetForceAtTime(double t) {
8188
unsigned int dof = num_bodies_ * 6;
8289
Eigen::VectorXd f(dof);

src/hydro/waves/wave_utilities.cpp

Lines changed: 29 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -62,6 +62,18 @@ double GetEta(const Eigen::Vector3d& position,
6262
return amplitude * cos(wavenumber * x_pos - omega * time + phase);
6363
}
6464

65+
double GetEtaGradientX(const Eigen::Vector3d& position,
66+
double time,
67+
double omega,
68+
double amplitude,
69+
double phase,
70+
double wavenumber) {
71+
// Analytic derivative of η = A·cos(k·x − ω·t + φ):
72+
// ∂η/∂x = −A·k·sin(k·x − ω·t + φ)
73+
double x_pos = position.x();
74+
return -amplitude * wavenumber * sin(wavenumber * x_pos - omega * time + phase);
75+
}
76+
6577
double GetEtaIrregular(const Eigen::Vector3d& position,
6678
double time,
6779
const Eigen::VectorXd& freqs_hz,
@@ -78,6 +90,23 @@ double GetEtaIrregular(const Eigen::Vector3d& position,
7890
return eta;
7991
}
8092

93+
double GetEtaGradientXIrregular(const Eigen::Vector3d& position,
94+
double time,
95+
const Eigen::VectorXd& freqs_hz,
96+
const Eigen::VectorXd& spectral_densities,
97+
const Eigen::VectorXd& spectral_widths,
98+
const Eigen::VectorXd& wave_phases,
99+
const Eigen::VectorXd& wavenumbers) {
100+
// Sum of component gradients: ∂η/∂x = Σ[−Aᵢ·kᵢ·sin(kᵢ·x − ωᵢ·t + φᵢ)]
101+
double deta_dx = 0.0;
102+
for (Eigen::Index i = 0; i < freqs_hz.size(); ++i) {
103+
double amplitude = std::sqrt(2.0 * spectral_densities[i] * spectral_widths[i]);
104+
double omega = 2.0 * M_PI * freqs_hz[i];
105+
deta_dx += GetEtaGradientX(position, time, omega, amplitude, wave_phases[i], wavenumbers[i]);
106+
}
107+
return deta_dx;
108+
}
109+
81110
std::vector<double> GetEtaIrregularTimeSeries(const Eigen::Vector3d& position,
82111
const std::vector<double> time_index,
83112
const Eigen::VectorXd& freqs_hz,

src/hydro/waves/wave_utilities.h

Lines changed: 31 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -44,6 +44,21 @@ double GetEta(const Eigen::Vector3d& position,
4444
double phase,
4545
double wavenumber);
4646

47+
/**
48+
* @brief Returns ∂η/∂x for a single regular wave component.
49+
*
50+
* Analytic derivative of GetEta(): ∂η/∂x = −A·k·sin(k·x − ω·t + φ).
51+
* Waves propagate in +X direction; ∂η/∂y = 0 by construction.
52+
*
53+
* @return Slope in x-direction (dimensionless).
54+
*/
55+
double GetEtaGradientX(const Eigen::Vector3d& position,
56+
double time,
57+
double omega,
58+
double amplitude,
59+
double phase,
60+
double wavenumber);
61+
4762
double GetEtaIrregular(const Eigen::Vector3d& position,
4863
double time,
4964
const Eigen::VectorXd& freqs_hz,
@@ -52,6 +67,22 @@ double GetEtaIrregular(const Eigen::Vector3d& position,
5267
const Eigen::VectorXd& wave_phases,
5368
const Eigen::VectorXd& wavenumbers);
5469

70+
/**
71+
* @brief Returns ∂η/∂x for irregular (spectral) waves by summing component gradients.
72+
*
73+
* Analytic derivative: ∂η/∂x = Σ[−Aᵢ·kᵢ·sin(kᵢ·x − ωᵢ·t + φᵢ)].
74+
* Waves propagate in +X direction; ∂η/∂y = 0 by construction.
75+
*
76+
* @return Slope in x-direction (dimensionless).
77+
*/
78+
double GetEtaGradientXIrregular(const Eigen::Vector3d& position,
79+
double time,
80+
const Eigen::VectorXd& freqs_hz,
81+
const Eigen::VectorXd& spectral_densities,
82+
const Eigen::VectorXd& spectral_widths,
83+
const Eigen::VectorXd& wave_phases,
84+
const Eigen::VectorXd& wavenumbers);
85+
5586
std::vector<double> GetEtaIrregularTimeSeries(const Eigen::Vector3d& position,
5687
const std::vector<double> time_index,
5788
const Eigen::VectorXd& freqs_hz,

0 commit comments

Comments
 (0)