Skip to content

Commit 403d334

Browse files
committed
unify RadiationComponent construction across legacy and HydroSystem paths
Introduce CreateRadiationComponent() in TestHydro as the single source of truth for constructing RadiationComponent with the current BEM radiation data, convolution mode, tapered options, and diagnostics settings.
1 parent d615169 commit 403d334

2 files changed

Lines changed: 47 additions & 42 deletions

File tree

include/hydroc/hydro_forces.h

Lines changed: 6 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -471,6 +471,12 @@ class TestHydro {
471471
// to ensure consistent construction.
472472
std::unique_ptr<hydrochrono::hydro::HydrostaticsComponent> CreateHydrostaticsComponent() const;
473473

474+
// Factory: creates RadiationComponent with current BEM data and convolution settings.
475+
// Used by both EnsureRadiationComponent() and EnsureHydroSystemAndCoupler()
476+
// to ensure consistent construction.
477+
// Note: Each instance owns its own velocity history (they are independent).
478+
std::unique_ptr<hydrochrono::hydro::RadiationComponent> CreateRadiationComponent() const;
479+
474480
// Internal helpers for HydroSystem + ChronoHydroCoupler path
475481
// Constructs hydro_system_ and chrono_coupler_ once; subsequent calls are no-ops.
476482
void EnsureHydroSystemAndCoupler();

src/hydro_forces.cpp

Lines changed: 41 additions & 42 deletions
Original file line numberDiff line numberDiff line change
@@ -450,29 +450,7 @@ void TestHydro::EnsureRadiationComponent() {
450450
return; // Already created
451451
}
452452

453-
const int rirf_steps = file_info_.GetRIRFDims(2);
454-
455-
// Convert TestHydro::RadiationConvolutionMode to hydrochrono::hydro::RadiationConvolutionMode
456-
hydrochrono::hydro::RadiationConvolutionMode component_mode;
457-
if (convolution_mode_ == RadiationConvolutionMode::TaperedDirect) {
458-
component_mode = hydrochrono::hydro::RadiationConvolutionMode::TaperedDirect;
459-
} else {
460-
component_mode = hydrochrono::hydro::RadiationConvolutionMode::Baseline;
461-
}
462-
463-
// Convert TaperedDirectOptions to hydrochrono::hydro::TaperedDirectOptions
464-
hydrochrono::hydro::TaperedDirectOptions component_opts;
465-
component_opts.smoothing = tapered_opts_.smoothing;
466-
component_opts.window_length = tapered_opts_.window_length;
467-
component_opts.rirf_end_time = tapered_opts_.rirf_end_time;
468-
component_opts.taper_start_percent = tapered_opts_.taper_start_percent;
469-
component_opts.taper_end_percent = tapered_opts_.taper_end_percent;
470-
component_opts.taper_final_amplitude = tapered_opts_.taper_final_amplitude;
471-
component_opts.export_plot_csv = tapered_opts_.export_plot_csv;
472-
473-
radiation_component_ = std::make_unique<hydrochrono::hydro::RadiationComponent>(
474-
file_info_, num_bodies_, rirf_steps, rirf_time_vector, rirf_width_vector,
475-
component_mode, component_opts, diagnostics_output_dir_);
453+
radiation_component_ = CreateRadiationComponent();
476454
}
477455

478456
void TestHydro::EnsureExcitationComponent() {
@@ -502,6 +480,44 @@ std::unique_ptr<hydrochrono::hydro::HydrostaticsComponent> TestHydro::CreateHydr
502480
file_info_, num_bodies_, equilibrium_, cb_minus_cg_, gravitational_acceleration);
503481
}
504482

483+
std::unique_ptr<hydrochrono::hydro::RadiationComponent> TestHydro::CreateRadiationComponent() const {
484+
// Single source of truth for RadiationComponent construction.
485+
// Used by EnsureRadiationComponent() and EnsureHydroSystemAndCoupler().
486+
// Each RadiationComponent instance owns its own velocity history.
487+
//
488+
// Configuration inputs:
489+
// - file_info_: BEM data including RIRF kernels
490+
// - num_bodies_: number of bodies in system
491+
// - rirf_time_vector, rirf_width_vector: time discretization for RIRF
492+
// - convolution_mode_: Baseline or TaperedDirect
493+
// - tapered_opts_: preprocessing options for TaperedDirect mode
494+
// - diagnostics_output_dir_: where to write debug CSVs
495+
496+
const int rirf_steps = file_info_.GetRIRFDims(2);
497+
498+
// Convert TestHydro::RadiationConvolutionMode to hydrochrono::hydro::RadiationConvolutionMode
499+
hydrochrono::hydro::RadiationConvolutionMode component_mode;
500+
if (convolution_mode_ == RadiationConvolutionMode::TaperedDirect) {
501+
component_mode = hydrochrono::hydro::RadiationConvolutionMode::TaperedDirect;
502+
} else {
503+
component_mode = hydrochrono::hydro::RadiationConvolutionMode::Baseline;
504+
}
505+
506+
// Convert TestHydro::TaperedDirectOptions to hydrochrono::hydro::TaperedDirectOptions
507+
hydrochrono::hydro::TaperedDirectOptions component_opts;
508+
component_opts.smoothing = tapered_opts_.smoothing;
509+
component_opts.window_length = tapered_opts_.window_length;
510+
component_opts.rirf_end_time = tapered_opts_.rirf_end_time;
511+
component_opts.taper_start_percent = tapered_opts_.taper_start_percent;
512+
component_opts.taper_end_percent = tapered_opts_.taper_end_percent;
513+
component_opts.taper_final_amplitude = tapered_opts_.taper_final_amplitude;
514+
component_opts.export_plot_csv = tapered_opts_.export_plot_csv;
515+
516+
return std::make_unique<hydrochrono::hydro::RadiationComponent>(
517+
file_info_, num_bodies_, rirf_steps, rirf_time_vector, rirf_width_vector,
518+
component_mode, component_opts, diagnostics_output_dir_);
519+
}
520+
505521
// ------------------------------------------------------------
506522
// Internal helpers for HydroSystem + ChronoHydroCoupler path
507523
// ------------------------------------------------------------
@@ -521,25 +537,8 @@ void TestHydro::EnsureHydroSystemAndCoupler() {
521537
// Hydrostatics component (uses shared factory for consistent construction)
522538
components.push_back(CreateHydrostaticsComponent());
523539

524-
// Radiation component (owns its own velocity history)
525-
const int rirf_steps = file_info_.GetRIRFDims(2);
526-
hydrochrono::hydro::RadiationConvolutionMode component_mode;
527-
if (convolution_mode_ == RadiationConvolutionMode::TaperedDirect) {
528-
component_mode = hydrochrono::hydro::RadiationConvolutionMode::TaperedDirect;
529-
} else {
530-
component_mode = hydrochrono::hydro::RadiationConvolutionMode::Baseline;
531-
}
532-
hydrochrono::hydro::TaperedDirectOptions component_opts;
533-
component_opts.smoothing = tapered_opts_.smoothing;
534-
component_opts.window_length = tapered_opts_.window_length;
535-
component_opts.rirf_end_time = tapered_opts_.rirf_end_time;
536-
component_opts.taper_start_percent = tapered_opts_.taper_start_percent;
537-
component_opts.taper_end_percent = tapered_opts_.taper_end_percent;
538-
component_opts.taper_final_amplitude = tapered_opts_.taper_final_amplitude;
539-
component_opts.export_plot_csv = tapered_opts_.export_plot_csv;
540-
components.push_back(std::make_unique<hydrochrono::hydro::RadiationComponent>(
541-
file_info_, num_bodies_, rirf_steps, rirf_time_vector, rirf_width_vector,
542-
component_mode, component_opts, diagnostics_output_dir_));
540+
// Radiation component (uses shared factory for consistent construction)
541+
components.push_back(CreateRadiationComponent());
543542

544543
// Excitation component (uses shared factory for consistent construction)
545544
components.push_back(CreateExcitationComponent());

0 commit comments

Comments
 (0)