Class TwoFluidPipe

All Implemented Interfaces:
Serializable, Runnable, PipeLineInterface, ProcessEquipmentInterface, TwoPortInterface, SimulationInterface, NamedInterface

public class TwoFluidPipe extends Pipeline
Two-fluid transient multiphase pipe model.

Implements a full two-fluid model for 1D transient multiphase pipeline flow. Unlike the drift-flux based TransientPipe, this model solves separate momentum equations for each phase, providing more accurate predictions for:

  • Countercurrent flow
  • Slug flow dynamics
  • Terrain-induced liquid accumulation
  • Transient pressure waves

Conservation Equations

  • Gas Mass: ∂/∂t(α_g·ρ_g·A) + ∂/∂x(α_g·ρ_g·v_g·A) = Γ_g
  • Liquid Mass: ∂/∂t(α_L·ρ_L·A) + ∂/∂x(α_L·ρ_L·v_L·A) = Γ_L
  • Gas Momentum: ∂/∂t(α_g·ρ_g·v_g·A) + ∂/∂x(α_g·ρ_g·v_g²·A + α_g·P·A) = S_g
  • Liquid Momentum: ∂/∂t(α_L·ρ_L·v_L·A) + ∂/∂x(α_L·ρ_L·v_L²·A + α_L·P·A) = S_L
  • Mixture Energy: (optional)

Usage Example

// Create two-phase fluid
SystemInterface fluid = new SystemSrkEos(300, 50);
fluid.addComponent("methane", 0.85);
fluid.addComponent("n-pentane", 0.15);
fluid.setMixingRule("classic");
fluid.setMultiPhaseCheck(true);

// Create inlet stream
Stream inlet = new Stream("inlet", fluid);
inlet.setFlowRate(10, "kg/sec");
inlet.run();

// Create two-fluid pipe
TwoFluidPipe pipe = new TwoFluidPipe("Pipeline", inlet);
pipe.setLength(5000); // 5 km
pipe.setDiameter(0.3); // 300 mm
pipe.setNumberOfSections(100);

// Set terrain profile
double[] elevations = new double[100];
for (int i = 0; i < 100; i++) {
  elevations[i] = 50.0 * Math.sin(i * Math.PI / 50); // Undulating terrain
}
pipe.setElevationProfile(elevations);

// Initialize steady state
pipe.run();

// Transient simulation
UUID id = UUID.randomUUID();
for (int step = 0; step < 1000; step++) {
  pipe.runTransient(0.1, id); // 0.1 second steps
}

// Get results
double[] pressures = pipe.getPressureProfile();
double[] holdups = pipe.getLiquidHoldupProfile();
double liquidInventory = pipe.getLiquidInventory("m3");

References

  • Bendiksen, K.H. et al. (1991) - The Dynamic Two-Fluid Model OLGA
  • Taitel, Y. and Dukler, A.E. (1976) - Flow regime transitions
  • Issa, R.I. and Kempf, M.H.W. (2003) - Simulation of slug flow
Version:
1.0
Author:
Even Solbraa
See Also:
  • Field Details

    • serialVersionUID

      private static final long serialVersionUID
      See Also:
    • logger

      private static final org.apache.logging.log4j.Logger logger
    • length

      private double length
      Total pipe length (m).
    • diameter

      private double diameter
      Pipe inner diameter (m).
    • roughness

      private double roughness
      Pipe wall roughness (m).
    • numberOfSections

      private int numberOfSections
      Number of computational cells.
    • elevationProfile

      private double[] elevationProfile
      Elevation profile at each section (m).
    • sections

      private TwoFluidSection[] sections
      Pipe sections with state.
    • dx

      private double dx
      Spatial step size (m).
    • simulationTime

      private double simulationTime
      Current simulation time (s).
    • maxSimulationTime

      private double maxSimulationTime
      Maximum simulation time (s).
    • cflNumber

      private double cflNumber
      CFL number for time stepping (0 < CFL < 1).
    • equations

      private TwoFluidConservationEquations equations
      Conservation equations solver.
    • timeIntegrator

      private TimeIntegrator timeIntegrator
      Time integrator.
    • flowRegimeDetector

      private FlowRegimeDetector flowRegimeDetector
      Flow regime detector.
    • accumulationTracker

      private LiquidAccumulationTracker accumulationTracker
      Liquid accumulation tracker.
    • slugTracker

      private SlugTracker slugTracker
      Slug tracker.
    • inletBCType

      private TwoFluidPipe.BoundaryCondition inletBCType
      Inlet boundary condition type.
    • outletBCType

      private TwoFluidPipe.BoundaryCondition outletBCType
      Outlet boundary condition type.
    • outletPressure

      private double outletPressure
      Outlet pressure (Pa).
    • outletPressureSet

      private boolean outletPressureSet
      Flag indicating if outlet pressure was explicitly set.
    • includeEnergyEquation

      private boolean includeEnergyEquation
      Include energy equation.
    • includeMassTransfer

      private boolean includeMassTransfer
      Include mass transfer (flash/condensation).
    • enableHeatTransfer

      private boolean enableHeatTransfer
      Enable heat transfer from surroundings.
    • surfaceTemperature

      private double surfaceTemperature
      Surface temperature for heat transfer (K).
    • heatTransferCoefficient

      private double heatTransferCoefficient
      Heat transfer coefficient (W/(m²·K)).
    • heatTransferProfile

      private double[] heatTransferProfile
      Heat transfer coefficient profile along pipe (W/(m²·K)).
    • surfaceTemperatureProfile

      private double[] surfaceTemperatureProfile
      Surface temperature profile along pipe (K).
    • wallThickness

      private double wallThickness
      Pipe wall thickness (m).
    • wallDensity

      private double wallDensity
      Pipe wall density (kg/m³) - steel default.
    • wallHeatCapacity

      private double wallHeatCapacity
      Pipe wall specific heat capacity (J/(kg·K)) - steel default.
    • wallTemperatureProfile

      private double[] wallTemperatureProfile
      Pipe wall temperature profile (K).
    • soilThermalResistance

      private double soilThermalResistance
      Soil/burial thermal resistance (m²·K/W).
    • enableJouleThomson

      private boolean enableJouleThomson
      Enable Joule-Thomson effect.
    • hydrateFormationTemperature

      private double hydrateFormationTemperature
      Hydrate formation temperature (K).
    • waxAppearanceTemperature

      private double waxAppearanceTemperature
      Wax appearance temperature (K).
    • hydrateRiskSections

      private boolean[] hydrateRiskSections
      Sections flagged for hydrate risk.
    • waxRiskSections

      private boolean[] waxRiskSections
      Sections flagged for wax risk.
    • insulationType

      private TwoFluidPipe.InsulationType insulationType
      Current insulation type.
    • enableSlugTracking

      private boolean enableSlugTracking
      Enable slug tracking.
    • outletSlugCount

      private int outletSlugCount
      Outlet slug statistics.
    • totalSlugVolumeAtOutlet

      private double totalSlugVolumeAtOutlet
    • lastSlugArrivalTime

      private double lastSlugArrivalTime
    • maxSlugLengthAtOutlet

      private double maxSlugLengthAtOutlet
    • maxSlugVolumeAtOutlet

      private double maxSlugVolumeAtOutlet
    • countedOutletSlugs

      private Set<Integer> countedOutletSlugs
      Track which slugs have already been counted at outlet (by slug ID).
    • thermodynamicUpdateInterval

      private int thermodynamicUpdateInterval
      Update thermodynamics every N steps.
    • currentStep

      private int currentStep
      Current step count.
    • pressureProfile

      private double[] pressureProfile
      Pressure profile (Pa).
    • temperatureProfile

      private double[] temperatureProfile
      Temperature profile (K).
    • liquidHoldupProfile

      private double[] liquidHoldupProfile
      Liquid holdup profile.
    • gasVelocityProfile

      private double[] gasVelocityProfile
      Gas velocity profile (m/s).
    • liquidVelocityProfile

      private double[] liquidVelocityProfile
      Liquid velocity profile (m/s).
    • referenceFluid

      private SystemInterface referenceFluid
      Reference fluid for flash calculations.
  • Constructor Details

    • TwoFluidPipe

      public TwoFluidPipe(String name)
      Constructor with name only.
      Parameters:
      name - Equipment name
    • TwoFluidPipe

      public TwoFluidPipe(String name, StreamInterface inStream)
      Constructor with inlet stream.
      Parameters:
      name - Equipment name
      inStream - Inlet stream
  • Method Details

    • initSubModels

      private void initSubModels()
      Initialize sub-models.
    • initializeSections

      private void initializeSections()
      Initialize pipe sections with inlet conditions.
    • runSteadyState

      private void runSteadyState()
      Run steady-state initialization.
    • updateTemperatureProfile

      private void updateTemperatureProfile(double massFlow, double area)
      Update temperature profile along the pipe accounting for heat transfer.

      Steady-state energy balance: m_dot * Cp * dT/dx = -h * π * D * (T - T_surface) - μ_JT * dP/dx

      Parameters:
      massFlow - Total mass flow rate [kg/s]
      area - Pipe cross-sectional area [m²]
    • updateTransientTemperature

      private void updateTransientTemperature(double massFlow, double area, double dt)
      Update temperature profile for transient simulation including pipe wall thermal mass.

      Solves coupled fluid-wall energy equations:

      • Fluid: ρ_f * Cp_f * A * dT_f/dt = -m_dot * Cp_f * dT_f/dx - h_i * π * D * (T_f - T_w)
      • Wall: ρ_w * Cp_w * A_w * dT_w/dt = h_i * π * D * (T_f - T_w) - h_o * π * D_o * (T_w - T_amb)
      Parameters:
      massFlow - Total mass flow rate [kg/s]
      area - Pipe cross-sectional area [m²]
      dt - Time step [s]
    • calculateLocalHoldup

      private double[] calculateLocalHoldup(TwoFluidSection sec, TwoFluidSection prev, double mDotGas, double mDotLiq, double area)
      Calculate local liquid holdup using drift-flux model with terrain effects.

      Uses the drift-flux model: v_G = C_0 * v_m + v_gj where:

      • C_0 = distribution coefficient (~1.0-1.2 for stratified flow)
      • v_gj = drift velocity (gas rises relative to mixture)
      Terrain effects modify slip velocity based on inclination.
      Parameters:
      sec - Current section
      prev - Previous section (upstream)
      mDotGas - Gas mass flow rate [kg/s]
      mDotLiq - Liquid mass flow rate [kg/s]
      area - Pipe cross-sectional area [m²]
      Returns:
      Array with [liquidHoldup, gasHoldup]
    • estimatePressureGradient

      private double estimatePressureGradient(TwoFluidSection sec)
      Estimate pressure gradient for steady-state initialization.

      Uses Haaland equation (explicit approximation to Colebrook-White) for friction factor, consistent with AdiabaticPipe and PipeBeggsAndBrills.

      Parameters:
      sec - Current pipe section
      Returns:
      Pressure gradient estimate (Pa/m)
    • updateWaterOilHoldups

      private void updateWaterOilHoldups(TwoFluidSection sec, TwoFluidSection prev, double alphaL, double area)
      Update water and oil holdups for three-phase flow with terrain effects.

      Water is denser than oil, so it tends to accumulate in valleys (low spots) more than oil. This method calculates the local water cut which can vary along the pipe based on:

      • Gravity segregation: water settles faster in low-velocity regions
      • Terrain effects: water accumulates more in valleys
      • Slip between oil and water phases
      Parameters:
      sec - Current section
      prev - Previous section (upstream)
      alphaL - Total liquid holdup
      area - Pipe cross-sectional area
    • run

      public void run(UUID id)
      Description copied from class: Pipeline

      In this method all thermodynamic and unit operations will be calculated in a steady state calculation.

      Specified by:
      run in interface SimulationInterface
      Overrides:
      run in class Pipeline
      Parameters:
      id - UUID
    • runTransient

      public void runTransient(double dt, UUID id)
      Run transient simulation for specified time step.
      Specified by:
      runTransient in interface SimulationInterface
      Overrides:
      runTransient in class Pipeline
      Parameters:
      dt - Requested time step (s)
      id - Calculation identifier
    • validateAndCorrectState

      private void validateAndCorrectState(double[][] U_new, double[][] U_prev)
      Validate and correct state vector to prevent numerical instabilities.
      Parameters:
      U_new - New state to validate
      U_prev - Previous state for fallback
    • validateSectionStates

      private void validateSectionStates()
      Validate section states and fix any numerical issues.
    • calcStableTimeStep

      private double calcStableTimeStep()
      Calculate stable time step using CFL condition.
      Returns:
      stable time step [s]
    • updateThermodynamics

      private void updateThermodynamics()
      Update thermodynamic properties using flash calculations.
    • applyBoundaryConditions

      private void applyBoundaryConditions()
      Apply boundary conditions.
    • updateOutletStream

      private void updateOutletStream()
      Update outlet stream with current outlet conditions.
    • updateResultArrays

      private void updateResultArrays()
      Update result arrays from section states.
    • getLiquidInventory

      public double getLiquidInventory(String unit)
      Get total liquid inventory in the pipe.

      Calculates inventory from conservative mass per length, converted to volume using local liquid density. This ensures consistency with the solver's mass tracking.

      Parameters:
      unit - Volume unit ("m3", "bbl", "L")
      Returns:
      Liquid volume
    • getPressureProfile

      public double[] getPressureProfile()
      Get pressure profile.
      Returns:
      Pressure at each section (Pa)
    • getTemperatureProfile

      public double[] getTemperatureProfile()
      Get temperature profile.
      Returns:
      Temperature at each section (K)
    • getLiquidHoldupProfile

      public double[] getLiquidHoldupProfile()
      Get liquid holdup profile.

      For consistency with oil and water holdups, the liquid holdup is calculated as the sum of oil and water holdups.

      Returns:
      Holdup at each section (0-1)
    • getWaterCutProfile

      public double[] getWaterCutProfile()
      Get water cut profile along the pipeline.

      For three-phase flow, water cut can vary along the pipeline as water accumulates in low spots (valleys) due to its higher density compared to oil.

      Returns:
      Water cut at each section (0-1, fraction of liquid that is water)
    • getWaterHoldupProfile

      public double[] getWaterHoldupProfile()
      Get water holdup profile along the pipeline.
      Returns:
      Water holdup at each section (0-1, fraction of pipe area occupied by water)
    • getOilHoldupProfile

      public double[] getOilHoldupProfile()
      Get oil holdup profile along the pipeline.
      Returns:
      Oil holdup at each section (0-1, fraction of pipe area occupied by oil)
    • getGasVelocityProfile

      public double[] getGasVelocityProfile()
      Get gas velocity profile.
      Returns:
      Gas velocity at each section (m/s)
    • getLiquidVelocityProfile

      public double[] getLiquidVelocityProfile()
      Get liquid velocity profile.
      Returns:
      Liquid velocity at each section (m/s)
    • getOilVelocityProfile

      public double[] getOilVelocityProfile()
      Get oil velocity profile along the pipeline.

      When water-oil slip is enabled, this returns the independent oil velocity. Otherwise, it returns the combined liquid velocity.

      Returns:
      Oil velocity at each section (m/s)
    • getWaterVelocityProfile

      public double[] getWaterVelocityProfile()
      Get water velocity profile along the pipeline.

      When water-oil slip is enabled, this returns the independent water velocity. Otherwise, it returns the combined liquid velocity.

      Returns:
      Water velocity at each section (m/s)
    • getOilWaterSlipProfile

      public double[] getOilWaterSlipProfile()
      Get oil-water velocity slip profile along the pipeline.

      Returns the difference between oil and water velocities (vOil - vWater). Positive values indicate oil is flowing faster than water.

      Returns:
      Oil-water slip velocity at each section (m/s)
    • getFlowRegimeProfile

      public PipeSection.FlowRegime[] getFlowRegimeProfile()
      Get flow regime at each section.
      Returns:
      Array of flow regimes
    • getPositionProfile

      public double[] getPositionProfile()
      Get position array for plotting.
      Returns:
      Position along pipe (m)
    • getSimulationTime

      public double getSimulationTime()
      Get current simulation time.
      Returns:
      Time (s)
    • getAccumulationTracker

      public LiquidAccumulationTracker getAccumulationTracker()
      Get accumulation tracker for detailed analysis.
      Returns:
      Accumulation tracker
    • getSlugTracker

      public SlugTracker getSlugTracker()
      Get slug tracker for slug statistics.
      Returns:
      Slug tracker
    • trackOutletSlugs

      private void trackOutletSlugs()
      Track slugs arriving at outlet and collect statistics. Each slug is only counted once when it first reaches the outlet region.
    • getOutletSlugCount

      public int getOutletSlugCount()
      Get number of slugs that have arrived at outlet.
      Returns:
      Outlet slug count
    • getTotalSlugVolumeAtOutlet

      public double getTotalSlugVolumeAtOutlet()
      Get total liquid volume delivered by slugs at outlet.
      Returns:
      Total slug volume (m³)
    • getLastSlugArrivalTime

      public double getLastSlugArrivalTime()
      Get time of last slug arrival at outlet.
      Returns:
      Time (s)
    • getMaxSlugLengthAtOutlet

      public double getMaxSlugLengthAtOutlet()
      Get maximum slug length observed at outlet.
      Returns:
      Max length (m)
    • getMaxSlugVolumeAtOutlet

      public double getMaxSlugVolumeAtOutlet()
      Get maximum slug volume observed at outlet.
      Returns:
      Max volume (m³)
    • getSlugStatisticsSummary

      public String getSlugStatisticsSummary()
      Get slug statistics summary string.
      Returns:
      Statistics summary
    • setLength

      public void setLength(double length)
      Set pipe length.
      Parameters:
      length - Length (m)
    • getLength

      public double getLength()
      Get pipe length.
      Returns:
      Length (m)
    • setDiameter

      public void setDiameter(double diameter)
      Set pipe diameter.
      Parameters:
      diameter - Diameter (m)
    • getDiameter

      public double getDiameter()
      Get pipe diameter.
      Returns:
      Diameter (m)
    • setRoughness

      public void setRoughness(double roughness)
      Set pipe wall roughness.
      Parameters:
      roughness - Roughness (m)
    • getRoughness

      public double getRoughness()
      Get pipe wall roughness.
      Returns:
      Roughness (m)
    • setNumberOfSections

      public void setNumberOfSections(int numberOfSections)
      Set number of computational sections.
      Parameters:
      numberOfSections - Number of sections
    • getNumberOfSections

      public int getNumberOfSections()
      Get number of sections.
      Returns:
      Number of sections
    • setElevationProfile

      public void setElevationProfile(double[] elevations)
      Set elevation profile.
      Parameters:
      elevations - Elevation at each section (m)
    • setOutletPressure

      public void setOutletPressure(double pressure)
      Set outlet pressure.
      Specified by:
      setOutletPressure in interface TwoPortInterface
      Overrides:
      setOutletPressure in class TwoPortEquipment
      Parameters:
      pressure - Pressure (Pa)
    • setOutletPressure

      public void setOutletPressure(double pressure, String unit)
      Set outlet pressure with unit.
      Parameters:
      pressure - Pressure value
      unit - Pressure unit ("Pa", "bara", "barg", "psia")
    • setCflNumber

      public void setCflNumber(double cfl)
      Set CFL number for time stepping.
      Parameters:
      cfl - CFL number (0 < cfl < 1)
    • setMaxSimulationTime

      public void setMaxSimulationTime(double time)
      Set maximum simulation time for transient calculations.
      Parameters:
      time - Maximum simulation time in seconds
    • getMaxSimulationTime

      public double getMaxSimulationTime()
      Get maximum simulation time.
      Returns:
      Maximum simulation time in seconds
    • setIncludeEnergyEquation

      public void setIncludeEnergyEquation(boolean include)
      Enable/disable energy equation.
      Parameters:
      include - true to include energy equation
    • setSurfaceTemperature

      public void setSurfaceTemperature(double temperature, String unit)
      Set surface temperature for heat transfer calculations.

      Enables heat transfer modeling. The pipe loses/gains heat to reach this temperature.

      Parameters:
      temperature - Surface temperature in the specified unit
      unit - Temperature unit ("K" or "C")
    • setHeatTransferCoefficient

      public void setHeatTransferCoefficient(double heatTransferCoefficient)
      Set heat transfer coefficient for convective heat transfer.

      Heat transfer rate: Q = h * A * (T_pipe - T_surface)
      where h = heat transfer coefficient (W/(m²·K))
      A = pipe surface area (m²)
      T_pipe = bulk fluid temperature (K)
      T_surface = surrounding surface temperature (K)

      Typical values:

      • Insulated subsea pipe: 5-15 W/(m²·K)
      • Uninsulated subsea pipe: 20-30 W/(m²·K)
      • Exposed/above-ground pipe: 50-100 W/(m²·K)
      Parameters:
      heatTransferCoefficient - Heat transfer coefficient in W/(m²·K)
    • getSurfaceTemperature

      public double getSurfaceTemperature()
      Get the surface temperature used for heat transfer calculations.
      Returns:
      Surface temperature in Kelvin
    • getHeatTransferCoefficient

      public double getHeatTransferCoefficient()
      Get the heat transfer coefficient.
      Returns:
      Heat transfer coefficient in W/(m²·K)
    • isHeatTransferEnabled

      public boolean isHeatTransferEnabled()
      Check if heat transfer is enabled.
      Returns:
      true if heat transfer modeling is active
    • setIncludeMassTransfer

      public void setIncludeMassTransfer(boolean include)
      Enable/disable mass transfer (flashing/condensation).
      Parameters:
      include - true to include mass transfer
    • setEnableSlugTracking

      public void setEnableSlugTracking(boolean enable)
      Enable/disable slug tracking.
      Parameters:
      enable - true to enable slug tracking
    • setEnableWaterOilSlip

      public void setEnableWaterOilSlip(boolean enable)
      Enable water-oil velocity slip modeling.

      When enabled, uses 7-equation model with separate oil and water momentum equations, allowing water and oil phases to flow at different velocities. This is important for:

      • Uphill flow: water slips back relative to oil due to higher density
      • Downhill flow: water accelerates relative to oil
      • Stratified oil-water flow with shear at interface
      Parameters:
      enable - true to enable 7-equation model with oil-water slip
    • isWaterOilSlipEnabled

      public boolean isWaterOilSlipEnabled()
      Check if water-oil velocity slip modeling is enabled.
      Returns:
      true if 7-equation model is enabled
    • setThermodynamicUpdateInterval

      public void setThermodynamicUpdateInterval(int interval)
      Set thermodynamic update interval.
      Parameters:
      interval - Update every N time steps
    • setInsulationType

      public void setInsulationType(TwoFluidPipe.InsulationType type)
      Set insulation type using predefined U-values.

      This is a convenience method that sets appropriate heat transfer coefficient based on insulation type. Automatically enables heat transfer modeling.

      Parameters:
      type - Insulation type preset
    • getInsulationType

      public TwoFluidPipe.InsulationType getInsulationType()
      Get the current insulation type.
      Returns:
      Current insulation type
    • setHeatTransferProfile

      public void setHeatTransferProfile(double[] profile)
      Set heat transfer coefficient profile along the pipe.

      Allows different U-values at different positions (e.g., buried vs exposed sections).

      Parameters:
      profile - Array of U-values [W/(m²·K)], one per section
    • getHeatTransferProfile

      public double[] getHeatTransferProfile()
      Get the heat transfer coefficient profile.
      Returns:
      Array of U-values or null if constant
    • setSurfaceTemperatureProfile

      public void setSurfaceTemperatureProfile(double[] profile)
      Set surface temperature profile along the pipe.

      Allows different ambient temperatures at different positions (e.g., varying seabed depth).

      Parameters:
      profile - Array of surface temperatures [K], one per section
    • getSurfaceTemperatureProfile

      public double[] getSurfaceTemperatureProfile()
      Get the surface temperature profile.
      Returns:
      Array of surface temperatures or null if constant
    • setWallProperties

      public void setWallProperties(double thickness, double density, double heatCapacity)
      Set pipe wall properties for transient thermal calculations.
      Parameters:
      thickness - Wall thickness [m]
      density - Wall material density [kg/m³]
      heatCapacity - Wall specific heat capacity [J/(kg·K)]
    • getWallThickness

      public double getWallThickness()
      Get pipe wall thickness.
      Returns:
      Wall thickness [m]
    • setSoilThermalResistance

      public void setSoilThermalResistance(double resistance)
      Set soil/burial thermal resistance.

      For buried pipelines, this adds thermal resistance between pipe outer wall and ambient. The effective U-value becomes: U_eff = 1 / (1/U + R_soil)

      Parameters:
      resistance - Soil thermal resistance [m²·K/W]
    • getSoilThermalResistance

      public double getSoilThermalResistance()
      Get soil thermal resistance.
      Returns:
      Soil thermal resistance [m²·K/W]
    • setEnableJouleThomson

      public void setEnableJouleThomson(boolean enable)
      Enable or disable Joule-Thomson effect.

      When enabled, temperature drops due to pressure reduction are calculated. This is important for gas pipelines with significant pressure drop.

      Parameters:
      enable - true to enable J-T effect
    • isJouleThomsonEnabled

      public boolean isJouleThomsonEnabled()
      Check if Joule-Thomson effect is enabled.
      Returns:
      true if J-T effect is modeled
    • setHydrateFormationTemperature

      public void setHydrateFormationTemperature(double temperature, String unit)
      Set hydrate formation temperature for risk monitoring.
      Parameters:
      temperature - Hydrate formation temperature
      unit - Temperature unit ("K" or "C")
    • getHydrateFormationTemperature

      public double getHydrateFormationTemperature()
      Get hydrate formation temperature.
      Returns:
      Hydrate formation temperature [K], or 0 if not set
    • setWaxAppearanceTemperature

      public void setWaxAppearanceTemperature(double temperature, String unit)
      Set wax appearance temperature for risk monitoring.
      Parameters:
      temperature - Wax appearance temperature
      unit - Temperature unit ("K" or "C")
    • getWaxAppearanceTemperature

      public double getWaxAppearanceTemperature()
      Get wax appearance temperature.
      Returns:
      Wax appearance temperature [K], or 0 if not set
    • getHydrateRiskSections

      public boolean[] getHydrateRiskSections()
      Get sections with hydrate formation risk.
      Returns:
      Array of booleans, true where temperature is below hydrate formation temperature
    • getWaxRiskSections

      public boolean[] getWaxRiskSections()
      Get sections with wax deposition risk.
      Returns:
      Array of booleans, true where temperature is below wax appearance temperature
    • hasHydrateRisk

      public boolean hasHydrateRisk()
      Check if any section has hydrate risk.
      Returns:
      true if any section temperature is below hydrate formation temperature
    • hasWaxRisk

      public boolean hasWaxRisk()
      Check if any section has wax risk.
      Returns:
      true if any section temperature is below wax appearance temperature
    • getTemperatureProfile

      public double[] getTemperatureProfile(String unit)
      Get temperature profile with specified unit.
      Parameters:
      unit - Temperature unit ("K", "C", or "F")
      Returns:
      Temperature profile in the specified unit
    • getWallTemperatureProfile

      public double[] getWallTemperatureProfile()
      Get the pipe wall temperature profile.
      Returns:
      Wall temperature profile [K], or null if not calculated
    • getHydrateRiskSectionCount

      public int getHydrateRiskSectionCount()
      Get number of sections with hydrate risk.
      Returns:
      Count of sections below hydrate formation temperature
    • getFirstHydrateRiskSection

      public int getFirstHydrateRiskSection()
      Get first section index with hydrate risk.
      Returns:
      Section index where hydrate risk first occurs, or -1 if no risk
    • getDistanceToHydrateRisk

      public double getDistanceToHydrateRisk()
      Get distance to first hydrate risk location.
      Returns:
      Distance [m] from inlet to first hydrate risk, or -1 if no risk