Class TransientPipe

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

public class TransientPipe extends TwoPortEquipment implements PipeLineInterface
Transient multiphase pipe model using drift-flux formulation. Implements a 1D transient multiphase flow simulator for gas-liquid flow in pipelines. The model is suitable for analyzing terrain-induced slugging, liquid accumulation, and transient pressure behavior in production pipelines, risers, and flowlines. Three-Phase Flow Support: When both oil and aqueous (water) phases are present, the model automatically calculates volume-weighted average liquid properties (density, viscosity, enthalpy, sound speed) based on the individual phase volumes from the thermodynamic flash. This maintains the drift-flux framework while properly accounting for oil-water mixtures.

Features

  • Drift-flux model for gas-liquid slip (Zuber-Findlay formulation)
  • Mechanistic flow regime detection (Taitel-Dukler, Barnea)
  • Three-phase gas-oil-water flow with volume-weighted liquid property averaging
  • Liquid accumulation tracking at terrain low points
  • Lagrangian slug tracking for terrain-induced slugging
  • Integration with NeqSim thermodynamics (SRK, PR, CPA equations of state)
  • AUSM+ numerical flux scheme with adaptive CFL-based time stepping

Basic Usage

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

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

// Create and configure transient pipe
TransientPipe pipe = new TransientPipe("Pipeline", inlet);
pipe.setLength(1000); // 1000 m total length
pipe.setDiameter(0.2); // 200 mm inner diameter
pipe.setNumberOfSections(50); // 50 computational cells
pipe.setMaxSimulationTime(60); // 60 seconds simulation

// Run simulation
pipe.run();

// Access results
double[] pressures = pipe.getPressureProfile();
double[] holdups = pipe.getLiquidHoldupProfile();

Terrain Pipeline Example

TransientPipe pipe = new TransientPipe("TerrainPipe", inlet);
pipe.setLength(2000);
pipe.setDiameter(0.3);
pipe.setNumberOfSections(40);

// Define elevation profile with low point
double[] elevations = new double[40];
for (int i = 0; i < 40; i++) {
  double x = i * 50.0;
  if (x < 500)
    elevations[i] = 0;
  else if (x < 1000)
    elevations[i] = -20 * (x - 500) / 500;
  else if (x < 1500)
    elevations[i] = -20 + 20 * (x - 1000) / 500;
  else
    elevations[i] = 0;
}
pipe.setElevationProfile(elevations);

pipe.run();

// Check liquid accumulation
var accumTracker = pipe.getAccumulationTracker();
for (var zone : accumTracker.getAccumulationZones()) {
  System.out.println("Accumulation at: " + zone.getPosition() + " m");
}

Physical Model

The drift-flux model relates gas velocity to mixture velocity:
v_G = C₀ · v_m + v_d
where C₀ is the distribution coefficient (typically 1.0-1.2) and v_d is the drift velocity. Flow regime-dependent correlations from Bendiksen (1984) and Harmathy (1960) provide closure.

Numerical Method

The model solves conservation equations using an explicit finite volume scheme with AUSM+ flux splitting. Time stepping is adaptive based on CFL condition for numerical stability.

References

  • Taitel, Y. and Dukler, A.E. (1976) - AIChE Journal 22(1)
  • Barnea, D. (1987) - Int. J. Multiphase Flow 13(1)
  • Bendiksen, K.H. (1984) - Int. J. Multiphase Flow 10(4)
  • Zuber, N. and Findlay, J.A. (1965) - J. Heat Transfer 87(4)
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
    • GRAVITY

      private static final double GRAVITY
      See Also:
    • length

      private double length
    • diameter

      private double diameter
    • roughness

      private double roughness
    • numberOfSections

      private int numberOfSections
    • wallThickness

      private double wallThickness
    • pipeElasticity

      private double pipeElasticity
    • elevationProfile

      private double[] elevationProfile
    • inclinationProfile

      private double[] inclinationProfile
    • sections

      private PipeSection[] sections
    • dx

      private double dx
    • simulationTime

      private double simulationTime
    • maxSimulationTime

      private double maxSimulationTime
    • dt

      private double dt
    • cflNumber

      private double cflNumber
    • minTimeStep

      private double minTimeStep
    • maxTimeStep

      private double maxTimeStep
    • flowRegimeDetector

      private FlowRegimeDetector flowRegimeDetector
    • driftFluxModel

      private DriftFluxModel driftFluxModel
    • accumulationTracker

      private LiquidAccumulationTracker accumulationTracker
    • slugTracker

      private SlugTracker slugTracker
    • inletBCType

      private TransientPipe.BoundaryCondition inletBCType
    • outletBCType

      private TransientPipe.BoundaryCondition outletBCType
    • inletPressureValue

      private double inletPressureValue
    • outletPressureValue

      private double outletPressureValue
    • outletPressureExplicitlySet

      private boolean outletPressureExplicitlySet
    • inletMassFlow

      private double inletMassFlow
    • outletMassFlow

      private double outletMassFlow
    • includeHeatTransfer

      private boolean includeHeatTransfer
    • ambientTemperature

      private double ambientTemperature
    • overallHeatTransferCoeff

      private double overallHeatTransferCoeff
    • jouleThomsonCoeff

      private double jouleThomsonCoeff
    • referenceFluid

      private SystemInterface referenceFluid
    • updateThermodynamics

      private boolean updateThermodynamics
    • thermodynamicUpdateInterval

      private int thermodynamicUpdateInterval
    • currentStep

      private int currentStep
    • pressureProfile

      private double[] pressureProfile
    • temperatureProfile

      private double[] temperatureProfile
    • liquidHoldupProfile

      private double[] liquidHoldupProfile
    • gasVelocityProfile

      private double[] gasVelocityProfile
    • liquidVelocityProfile

      private double[] liquidVelocityProfile
    • pressureHistory

      private double[][] pressureHistory
    • historyInterval

      private int historyInterval
    • historyIndex

      private int historyIndex
    • massResidual

      private double massResidual
    • energyResidual

      private double energyResidual
    • isConverged

      private boolean isConverged
    • totalTimeSteps

      private int totalTimeSteps
  • Constructor Details

    • TransientPipe

      public TransientPipe()
      Default constructor. Creates a TransientPipe with default name "TransientPipe".
    • TransientPipe

      public TransientPipe(String name)
      Constructor with name.
      Parameters:
      name - Pipe name used for identification in process systems
    • TransientPipe

      public TransientPipe(String name, StreamInterface inletStream)
      Constructor with name and inlet stream. This is the recommended constructor for typical usage. The inlet stream provides initial conditions (composition, temperature, pressure, flow rate) for the simulation.
      Parameters:
      name - Pipe name used for identification in process systems
      inletStream - Inlet stream containing fluid properties and flow conditions
  • Method Details

    • initializeSubModels

      private void initializeSubModels()
      Initialize sub-models for flow regime detection, drift-flux, and slug tracking.
    • initializePipe

      public void initializePipe()
      Initialize the discretized pipe sections. Creates the computational mesh based on pipe length and number of sections. Also initializes elevation/inclination profiles and identifies low points for liquid accumulation tracking. This method is called automatically by run() if sections have not been initialized. It can also be called explicitly to inspect the mesh before running. Note: The number of sections determines spatial resolution. For accurate slug tracking, use at least 20 sections. Typical guideline: dx ≈ 10-50 pipe diameters.
    • initializeFromStream

      private void initializeFromStream()
      Initialize pipe state from inlet stream.
    • run

      public void run()
      Run the transient simulation.
      Specified by:
      run in interface Runnable
      Specified by:
      run in interface SimulationInterface
    • run

      public void run(UUID id)
      Run with calculation identifier.
      Specified by:
      run in interface SimulationInterface
      Parameters:
      id - Calculation identifier
    • runTransient

      public void runTransient(double dt, UUID id)
      Run transient simulation for a specified time step. This method advances the pipe simulation by the specified time step dt, making it suitable for use within a ProcessSystem.runTransient(double, java.util.UUID) loop. Unlike run(), which runs the complete simulation to maxSimulationTime, this method performs incremental time-stepping that can be coordinated with other equipment. On the first call, the pipe is initialized from the inlet stream. Subsequent calls advance the simulation state incrementally. The method uses adaptive sub-stepping internally to maintain numerical stability (CFL condition) while advancing by the requested dt. Usage Example:
      ProcessSystem process = new ProcessSystem();
      process.add(inlet);
      process.add(transientPipe);
      process.add(separator);
      process.run(); // Initial steady state
      
      // Transient loop
      for (int i = 0; i < 100; i++) {
        // Optionally change inlet conditions here
        process.runTransient(1.0); // Advance 1 second
        System.out.println("Outlet pressure: " + separator.getGasOutStream().getPressure());
      }
      
      Note: The inlet stream conditions are re-read at each call, allowing dynamic boundary conditions. If the inlet flow rate or composition changes, the pipe simulation will respond accordingly.
      Specified by:
      runTransient in interface SimulationInterface
      Parameters:
      dt - Time step to advance in seconds. The method will use internal sub-stepping if required for numerical stability.
      id - Calculation identifier for tracking
      See Also:
    • updateInletFromStream

      private void updateInletFromStream()
      Update inlet boundary conditions from the current inlet stream state. This method is called during runTransient(double, UUID) to capture any changes in the inlet stream conditions (flow rate, pressure, composition) that may have occurred since the last time step.
    • updateOutletFromStream

      private void updateOutletFromStream()
      Updates the outlet boundary conditions from the current outlet stream state. This method is called during runTransient(double, UUID) to capture any changes in the outlet stream conditions when using CONSTANT_FLOW outlet boundary condition. This allows external controllers or other process equipment to set the outlet flow rate on the outlet stream, which is then read back into the pipe model.
    • calculateTimeStep

      private double calculateTimeStep()
      Calculate adaptive time step based on CFL condition.

      Uses the full eigenvalue structure for hyperbolic systems. The characteristic wave speeds for multiphase flow are approximately |u ± c| where u is the mixture velocity and c is the mixture sound speed. We take the maximum of all four wave speeds: |u_G + c|, |u_G - c|, |u_L + c|, |u_L - c| to ensure stability for both forward and backward propagating waves.

      Returns:
      time step (s)
    • advanceTimeStep

      private void advanceTimeStep(double dt)
      Advance solution by one time step using finite volume method.
      Parameters:
      dt - time step (s)
    • getMassFlowForPressureProfile

      private double getMassFlowForPressureProfile()
      Get the controlling mass flow rate for pressure profile calculation. Returns the appropriate mass flow rate based on boundary conditions:
      • CONSTANT_FLOW inlet: use inlet mass flow
      • CONSTANT_FLOW outlet: use outlet mass flow
      • Both CONSTANT_PRESSURE: use inlet mass flow as fallback
      Returns:
      Mass flow rate in kg/s
    • updateBoundaryPressures

      private void updateBoundaryPressures()
      Update pressure at boundary sections only to maintain boundary condition consistency.

      Interior section pressures are computed from the momentum equation to preserve momentum conservation. Only boundary-adjacent sections are updated to ensure proper pressure gradient at boundaries.

    • updatePressureProfile

      private void updatePressureProfile()
      Update pressure profile based on friction and gravity pressure gradients.

      Legacy method: Marches from inlet to outlet (or outlet to inlet), computing pressure drop in each cell based on the local friction and gravity gradients. This ensures the pressure profile is consistent with the flow. For single-phase flow, uses the direct Darcy-Weisbach calculation.

      Note: This method is preserved for compatibility but is no longer called in transient mode to avoid overriding momentum-conserving pressure from the conservative solver.

    • calculateFluxes

      private double[][] calculateFluxes(PipeSection[] oldSections)
      Calculate fluxes at cell interfaces using AUSM+ scheme. For interior interfaces, the AUSM+ upwind scheme is used. For boundary interfaces, the fluxes are set according to the boundary conditions to ensure mass and momentum conservation:
      • CONSTANT_FLOW: Flux is set to the specified mass flow rate
      • CONSTANT_PRESSURE: Flux is computed from the AUSM scheme
      • CLOSED: Zero flux
      Parameters:
      oldSections - array of pipe sections from the previous time step
      Returns:
      2D array of fluxes [section][conserved variable index]
    • calculateAUSMFlux

      private double[] calculateAUSMFlux(PipeSection left, PipeSection right)
      Calculate AUSM+ numerical flux at interface.
      Parameters:
      left - left pipe section
      right - right pipe section
      Returns:
      flux array [mass_L, mass_G, momentum_L, momentum_G]
    • splitMachPlus

      private double splitMachPlus(double M)
    • splitMachMinus

      private double splitMachMinus(double M)
    • splitPressurePlus

      private double splitPressurePlus(double M, double p)
    • splitPressureMinus

      private double splitPressureMinus(double M, double p)
    • addSourceTerms

      private void addSourceTerms(double[] U, PipeSection section, double dt)
      Add source terms (gravity and friction) to the momentum equation. The momentum equation in conservative form is:
      ∂(ρu)/∂t + ∂(ρu² + p)/∂x = -ρg sin(θ) - τ_w/A
      
      Where the source terms are:
      • Gravity: S_g = -ρ_m · g · sin(θ) [kg/(m²·s²)]
      • Friction: S_f = -dP/dx_friction [kg/(m²·s²)]
      These are integrated as: U[2]_new = U[2]_old + S × dt
      Parameters:
      U - conserved variable array to be modified in-place
      section - pipe section for which to calculate source terms
      dt - time step size in seconds
    • updatePrimitiveVariables

      private void updatePrimitiveVariables(PipeSection section, PipeSection oldSection, double[] U, int sectionIndex, double dt, double momentumFluxLeft, double momentumFluxRight)
      Update primitive variables from conservative variables.

      This method converts from conservative variables (densities, momentum, energy) back to primitive variables (holdups, velocities, temperature, pressure). Interior cell pressures are computed from the momentum equation residual to properly capture acoustic wave propagation.

      Parameters:
      section - Current section to update
      oldSection - Previous time step section state
      U - Conservative variables [ρ_G·α_G, ρ_L·α_L, ρ_m·u, ρ_m·E]
      sectionIndex - Index of this section (0 to numberOfSections-1)
      dt - Time step (s)
      momentumFluxLeft - Momentum flux at left face (Pa)
      momentumFluxRight - Momentum flux at right face (Pa)
    • applyBoundaryConditions

      private void applyBoundaryConditions()
      Apply boundary conditions.
    • updateFlowRegimes

      private void updateFlowRegimes()
      Update flow regimes for all sections.
    • checkTerrainSlugging

      private void checkTerrainSlugging()
      Check for terrain-induced slugging.
    • updateThermodynamicProperties

      private void updateThermodynamicProperties()
      Update thermodynamic properties using NeqSim flash.
    • checkConvergence

      private void checkConvergence()
      Check for steady-state convergence.
    • storeResults

      private void storeResults()
      Store current results.
    • storeHistory

      private void storeHistory()
      Store pressure history for time series analysis.
    • updateOutletStream

      private void updateOutletStream()
      Update outlet stream with final conditions including flow rate. This method updates the outlet stream with:
      • Pressure from the outlet pipe section
      • Temperature from the outlet pipe section
      • Mass flow rate calculated from outlet section velocity and density
      The flow rate is computed directly from the conservation equations to maintain mass balance. During slug flow, the increased liquid holdup naturally increases the outlet density, which is reflected in the mass flow calculation.
    • getName

      public String getName()
      Description copied from class: NamedBaseClass

      Getter for the field name.

      Specified by:
      getName in interface NamedInterface
      Overrides:
      getName in class NamedBaseClass
      Returns:
      a String object
    • setName

      public void setName(String name)
      Description copied from class: NamedBaseClass

      Setter for the field name.

      Specified by:
      setName in interface NamedInterface
      Overrides:
      setName in class NamedBaseClass
      Parameters:
      name - a String object
    • setLength

      public void setLength(double length)
      Set total pipe length.
      Parameters:
      length - Pipe length in meters
    • getLength

      public double getLength()
      Get total pipe length.
      Returns:
      Pipe length in meters
    • setDiameter

      public void setDiameter(double diameter)
      Set pipe inner diameter.
      Parameters:
      diameter - Inner diameter in meters
    • getDiameter

      public double getDiameter()
      Get pipe inner diameter.
      Returns:
      Inner diameter in meters
    • getRoughness

      public double getRoughness()
      Get pipe wall roughness.
      Returns:
      Wall roughness in meters
    • getWallThickness

      public double getWallThickness()
      Get pipe wall thickness.
      Returns:
      Wall thickness in meters
    • getPipeElasticity

      public double getPipeElasticity()
      Get pipe material elasticity (Young's modulus).
      Returns:
      Elasticity in Pascals
    • isIncludeHeatTransfer

      public boolean isIncludeHeatTransfer()
      Check if heat transfer is enabled.
      Returns:
      true if heat transfer is included in calculations
    • setIncludeHeatTransfer

      public void setIncludeHeatTransfer(boolean include)
      Enable or disable heat transfer calculations.

      When enabled, the energy equation is solved to calculate temperature changes along the pipe due to:

      • Heat transfer to/from surroundings (ambient)
      • Joule-Thomson effect (gas expansion cooling)
      • Friction heating (viscous dissipation)
      • Elevation work
      Parameters:
      include - true to enable heat transfer calculations
    • getAmbientTemperature

      public double getAmbientTemperature()
      Get ambient temperature for heat transfer calculations.
      Returns:
      Ambient temperature in Kelvin
    • setAmbientTemperature

      public void setAmbientTemperature(double temperature)
      Set ambient temperature for heat transfer calculations.

      This is the external temperature used for heat loss/gain calculations. For subsea pipelines, this would be seawater temperature (~277-283 K). For buried onshore pipelines, this would be soil temperature.

      Parameters:
      temperature - Ambient temperature in Kelvin
    • getOverallHeatTransferCoeff

      public double getOverallHeatTransferCoeff()
      Get overall heat transfer coefficient.
      Returns:
      Heat transfer coefficient in W/(m²·K)
    • setOverallHeatTransferCoeff

      public void setOverallHeatTransferCoeff(double coeff)
      Set overall heat transfer coefficient.

      The overall heat transfer coefficient U accounts for all heat transfer resistances between the fluid and surroundings. Typical values:

      • Bare steel pipe in still air: 5-10 W/(m²·K)
      • Bare steel pipe in seawater: 200-500 W/(m²·K)
      • Insulated pipe (50mm): 1-3 W/(m²·K)
      • Buried pipe: 2-5 W/(m²·K)
      Parameters:
      coeff - Heat transfer coefficient in W/(m²·K)
    • getJouleThomsonCoeff

      public double getJouleThomsonCoeff()
      Get Joule-Thomson coefficient for temperature change during gas expansion.

      Returns the JT coefficient calculated from gas phase thermodynamics. This is automatically updated during simulation based on the fluid properties.

      Typical values (at ~300K, 50 bar):

      • Methane: ~4.5e-6 K/Pa (0.45 K/MPa)
      • Natural gas: 2-6e-6 K/Pa
      • CO2: ~1e-5 K/Pa
      • Hydrogen: ~-0.3e-6 K/Pa (warms on expansion)
      • Liquids: ~1e-8 K/Pa (very small)
      Returns:
      Joule-Thomson coefficient in K/Pa
    • getMassResidual

      public double getMassResidual()
      Get mass residual from last time step.
      Returns:
      Mass balance residual
    • getEnergyResidual

      public double getEnergyResidual()
      Get energy residual from last time step.
      Returns:
      Energy balance residual
    • isConverged

      public boolean isConverged()
      Check if simulation has converged.
      Returns:
      true if converged to steady-state
    • setRoughness

      public void setRoughness(double roughness)
      Set pipe wall roughness. Used for friction factor calculation. Typical values:
      • New steel pipe: 0.00004 m (40 μm)
      • Used steel pipe: 0.0001 m (100 μm)
      • Corroded pipe: 0.0003 m (300 μm)
      Parameters:
      roughness - Wall roughness in meters
    • setNumberOfSections

      public void setNumberOfSections(int n)
      Set number of computational cells (sections). More sections provide higher spatial resolution but require more computation. Guideline: dx ≈ 10-50 pipe diameters. For slug tracking, use at least 20 sections.
      Parameters:
      n - Number of sections
    • setElevationProfile

      public void setElevationProfile(double[] profile)
      Set elevation profile along the pipe. Array of elevations (meters) at each section node. The array length should match the number of sections. Positive values indicate upward elevation. Example:
      double[] elevations = new double[50];
      for (int i = 0; i < 50; i++) {
        elevations[i] = 10 * Math.sin(i * Math.PI / 49); // Terrain with low point
      }
      pipe.setElevationProfile(elevations);
      
      Parameters:
      profile - Array of elevations in meters
    • setInclinationProfile

      public void setInclinationProfile(double[] profile)
      Set inclination profile along the pipe. Array of inclination angles (radians) at each section. Positive values indicate upward inclination. Use this instead of elevation profile for constant-inclination sections.
      Parameters:
      profile - Array of inclination angles in radians
    • setMaxSimulationTime

      public void setMaxSimulationTime(double time)
      Set maximum simulation time. The simulation runs until this time is reached or steady-state is achieved.
      Parameters:
      time - Maximum simulation time in seconds
    • getSimulationTime

      public double getSimulationTime()
      Get current simulation time.
      Returns:
      Simulation time in seconds
    • setInletBoundaryCondition

      public void setInletBoundaryCondition(TransientPipe.BoundaryCondition bc)
      Set inlet boundary condition type.
      Parameters:
      bc - Boundary condition type
      See Also:
    • setOutletBoundaryCondition

      public void setOutletBoundaryCondition(TransientPipe.BoundaryCondition bc)
    • setinletPressureValue

      @Deprecated public void setinletPressureValue(double pressure)
      Deprecated.
      Set the inlet pressure value.
      Parameters:
      pressure - Inlet pressure in bara (bar absolute)
    • setInletPressure

      public void setInletPressure(double pressure)
      Set the inlet pressure value.
      Specified by:
      setInletPressure in interface TwoPortInterface
      Overrides:
      setInletPressure in class TwoPortEquipment
      Parameters:
      pressure - Inlet pressure in bara (bar absolute)
    • setoutletPressureValue

      @Deprecated public void setoutletPressureValue(double pressure)
      Deprecated.
      Set the outlet pressure value.
      Parameters:
      pressure - Outlet pressure in bara (bar absolute)
    • setOutletPressure

      public void setOutletPressure(double pressure)
      Set the outlet pressure value.
      Specified by:
      setOutletPressure in interface TwoPortInterface
      Overrides:
      setOutletPressure in class TwoPortEquipment
      Parameters:
      pressure - Outlet pressure in bara (bar absolute)
    • setInletMassFlow

      public void setInletMassFlow(double flow)
    • setOutletMassFlow

      public void setOutletMassFlow(double flow)
      Set the outlet mass flow rate. Use this when the outlet flow is controlled by downstream equipment (e.g., a valve). When using CONSTANT_FLOW outlet boundary condition, this value is used to calculate the pressure profile and outlet stream properties. This can be called before each runTransient() call to update the outlet flow from downstream valve Cv calculations.
      Parameters:
      flow - Outlet mass flow rate in kg/s
    • getOutletMassFlow

      public double getOutletMassFlow()
      Get the current outlet mass flow rate.
      Returns:
      Outlet mass flow rate in kg/s
    • setInletStream

      public void setInletStream(StreamInterface stream)
      Description copied from class: TwoPortEquipment
      Set inlet Stream of twoport.
      Specified by:
      setInletStream in interface TwoPortInterface
      Overrides:
      setInletStream in class TwoPortEquipment
      Parameters:
      stream - value to set
    • getOutletStream

      public StreamInterface getOutletStream()
      Description copied from class: TwoPortEquipment
      Get outlet Stream of twoport.
      Specified by:
      getOutletStream in interface TwoPortInterface
      Overrides:
      getOutletStream in class TwoPortEquipment
      Returns:
      outlet Stream of TwoPortEquipment
    • setOutletStream

      public void setOutletStream(StreamInterface stream)
      Description copied from class: TwoPortEquipment
      Set outlet Stream of twoport.
      Specified by:
      setOutletStream in interface TwoPortInterface
      Overrides:
      setOutletStream in class TwoPortEquipment
      Parameters:
      stream - value to set
    • getPressureProfile

      public double[] getPressureProfile()
    • getTemperatureProfile

      public double[] getTemperatureProfile()
      Get temperature profile at end of simulation.
      Returns:
      Array of temperatures (K) at each section, or null if not run
    • getLiquidHoldupProfile

      public double[] getLiquidHoldupProfile()
      Get liquid holdup profile at end of simulation. Liquid holdup (α_L) is the fraction of pipe cross-section occupied by liquid. Values range from 0 (all gas) to 1 (all liquid).
      Returns:
      Array of liquid holdups (dimensionless) at each section, or null if not run
    • getGasVelocityProfile

      public double[] getGasVelocityProfile()
      Get gas velocity profile at end of simulation.
      Returns:
      Array of gas velocities (m/s) at each section, or null if not run
    • getLiquidVelocityProfile

      public double[] getLiquidVelocityProfile()
      Get liquid velocity profile at end of simulation.
      Returns:
      Array of liquid velocities (m/s) at each section, or null if not run
    • getPressureHistory

      public double[][] getPressureHistory()
      Get pressure history over simulation time. The history is stored at intervals specified by historyInterval. The array dimensions are [time_index][position_index].
      Returns:
      2D array of pressures (Pa), or null if not run
    • getSections

      public PipeSection[] getSections()
      Get all pipe sections with their state variables. Each section contains detailed state information including pressure, temperature, holdups, velocities, flow regime, and other properties.
      Returns:
      Array of PipeSection objects
    • getSlugTracker

      public SlugTracker getSlugTracker()
      Get the slug tracker for detailed slug analysis. The slug tracker contains information about active slugs, slug statistics, and slug history. Example:
      SlugTracker tracker = pipe.getSlugTracker();
      int slugCount = tracker.getSlugCount();
      double avgLength = tracker.getAverageSlugLength();
      String stats = tracker.getStatisticsString();
      
      Returns:
      SlugTracker instance
    • getAccumulationTracker

      public LiquidAccumulationTracker getAccumulationTracker()
      Get the liquid accumulation tracker. Tracks liquid pooling at terrain low points. Useful for identifying potential liquid loading and slug initiation locations. Example:
      var tracker = pipe.getAccumulationTracker();
      for (var zone : tracker.getAccumulationZones()) {
        System.out.println("Low point at: " + zone.getPosition());
        System.out.println("Accumulated: " + zone.getAccumulatedVolume() + " m3");
      }
      
      Returns:
      LiquidAccumulationTracker instance
    • getTotalTimeSteps

      public int getTotalTimeSteps()
      Get total number of time steps taken during simulation.
      Returns:
      Number of time steps
    • setCflNumber

      public void setCflNumber(double cfl)
      Set CFL number for adaptive time stepping. The CFL (Courant-Friedrichs-Lewy) number controls the time step size relative to the grid spacing and wave speeds. Lower values (0.3-0.5) provide more stability but slower simulation. Higher values (0.7-0.9) are faster but may become unstable.
      Parameters:
      cfl - CFL number, clamped to range [0.1, 1.0]
    • setThermodynamicUpdateInterval

      public void setThermodynamicUpdateInterval(int interval)
      Set interval for thermodynamic property updates. Flash calculations are computationally expensive. This setting controls how often phase properties are recalculated. Higher values improve performance but may reduce accuracy for rapidly changing conditions.
      Parameters:
      interval - Update interval (number of time steps), minimum 1
    • setUpdateThermodynamics

      public void setUpdateThermodynamics(boolean update)
      Enable or disable thermodynamic updates during simulation. When disabled, phase properties remain constant at initial values. This is appropriate for isothermal simulations or when temperature/pressure changes are small.
      Parameters:
      update - true to enable updates, false to disable
    • getCalculationIdentifier

      public UUID getCalculationIdentifier()
      Description copied from class: SimulationBaseClass
      Getter for property calcIdentifier.
      Specified by:
      getCalculationIdentifier in interface SimulationInterface
      Overrides:
      getCalculationIdentifier in class SimulationBaseClass
      Returns:
      Value of calcIdentifier.
    • setCalculationIdentifier

      public void setCalculationIdentifier(UUID id)
      Description copied from class: SimulationBaseClass
      Setter for property calcIdentifier.
      Specified by:
      setCalculationIdentifier in interface SimulationInterface
      Overrides:
      setCalculationIdentifier in class SimulationBaseClass
      Parameters:
      id - Value to set.
    • getPipe

      public FlowSystemInterface getPipe()

      getPipe.

      This transient model uses internal discretization, not FlowSystemInterface.
      Specified by:
      getPipe in interface PipeLineInterface
      Returns:
      a FlowSystemInterface object
    • setNumberOfLegs

      public void setNumberOfLegs(int number)

      setNumberOfLegs.

      Specified by:
      setNumberOfLegs in interface PipeLineInterface
      Parameters:
      number - a int
    • setHeightProfile

      public void setHeightProfile(double[] heights)

      setHeightProfile.

      Specified by:
      setHeightProfile in interface PipeLineInterface
      Parameters:
      heights - an array of type double
    • setLegPositions

      public void setLegPositions(double[] positions)

      setLegPositions.

      Specified by:
      setLegPositions in interface PipeLineInterface
      Parameters:
      positions - an array of type double
    • setPipeDiameters

      public void setPipeDiameters(double[] diameters)

      setPipeDiameters.

      Specified by:
      setPipeDiameters in interface PipeLineInterface
      Parameters:
      diameters - an array of type double
    • setPipeWallRoughness

      public void setPipeWallRoughness(double[] rough)

      setPipeWallRoughness.

      Specified by:
      setPipeWallRoughness in interface PipeLineInterface
      Parameters:
      rough - an array of type double
    • setOuterTemperatures

      public void setOuterTemperatures(double[] outerTemp)

      setOuterTemperatures.

      Specified by:
      setOuterTemperatures in interface PipeLineInterface
      Parameters:
      outerTemp - an array of type double
    • setNumberOfNodesInLeg

      public void setNumberOfNodesInLeg(int number)

      setNumberOfNodesInLeg.

      Specified by:
      setNumberOfNodesInLeg in interface PipeLineInterface
      Parameters:
      number - a int
    • setOutputFileName

      public void setOutputFileName(String name)

      setOutputFileName.

      Specified by:
      setOutputFileName in interface PipeLineInterface
      Parameters:
      name - a String object
    • setInitialFlowPattern

      public void setInitialFlowPattern(String flowPattern)

      setInitialFlowPattern.

      Specified by:
      setInitialFlowPattern in interface PipeLineInterface
      Parameters:
      flowPattern - a String object
    • getInletStream

      public StreamInterface getInletStream()
      Get inlet Stream of twoport.
      Specified by:
      getInletStream in interface TwoPortInterface
      Overrides:
      getInletStream in class TwoPortEquipment
      Returns:
      inlet Stream of TwoPortEquipment