Class TwoFluidPipe
- All Implemented Interfaces:
Serializable, Runnable, PipeLineInterface, ProcessEquipmentInterface, TwoPortInterface, SimulationInterface, NamedInterface
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:
-
Nested Class Summary
Nested ClassesModifier and TypeClassDescriptionstatic enumBoundary condition type.static enumInsulation type presets with typical U-values. -
Field Summary
FieldsModifier and TypeFieldDescriptionprivate LiquidAccumulationTrackerLiquid accumulation tracker.private doubleCFL number for time stepping (0 < CFL < 1).Track which slugs have already been counted at outlet (by slug ID).private intCurrent step count.private doublePipe inner diameter (m).private doubleSpatial step size (m).private double[]Elevation profile at each section (m).private booleanEnable heat transfer from surroundings.private booleanEnable Joule-Thomson effect.private booleanEnable slug tracking.private TwoFluidConservationEquationsConservation equations solver.private FlowRegimeDetectorFlow regime detector.private double[]Gas velocity profile (m/s).private doubleHeat transfer coefficient (W/(m²·K)).private double[]Heat transfer coefficient profile along pipe (W/(m²·K)).private doubleHydrate formation temperature (K).private boolean[]Sections flagged for hydrate risk.private booleanInclude energy equation.private booleanInclude mass transfer (flash/condensation).private TwoFluidPipe.BoundaryConditionInlet boundary condition type.private TwoFluidPipe.InsulationTypeCurrent insulation type.private doubleprivate doubleTotal pipe length (m).private double[]Liquid holdup profile.private double[]Liquid velocity profile (m/s).private static final org.apache.logging.log4j.Loggerprivate doubleMaximum simulation time (s).private doubleprivate doubleprivate intNumber of computational cells.private TwoFluidPipe.BoundaryConditionOutlet boundary condition type.private doubleOutlet pressure (Pa).private booleanFlag indicating if outlet pressure was explicitly set.private intOutlet slug statistics.private double[]Pressure profile (Pa).private SystemInterfaceReference fluid for flash calculations.private doublePipe wall roughness (m).private TwoFluidSection[]Pipe sections with state.private static final longprivate doubleCurrent simulation time (s).private SlugTrackerSlug tracker.private doubleSoil/burial thermal resistance (m²·K/W).private doubleSurface temperature for heat transfer (K).private double[]Surface temperature profile along pipe (K).private double[]Temperature profile (K).private intUpdate thermodynamics every N steps.private TimeIntegratorTime integrator.private doubleprivate doublePipe wall density (kg/m³) - steel default.private doublePipe wall specific heat capacity (J/(kg·K)) - steel default.private double[]Pipe wall temperature profile (K).private doublePipe wall thickness (m).private doubleWax appearance temperature (K).private boolean[]Sections flagged for wax risk.Fields inherited from class Pipeline
equilibriumHeatTransfer, equilibriumMassTransfer, fileName, flowPattern, legHeights, legPositions, numberOfLegs, numberOfNodesInLeg, outerHeatTransferCoeffs, outerTemperature, pipe, pipeDiameters, pipelineMechanicalDesign, pipeWallRoughness, system, times, wallHeatTransferCoeffsFields inherited from class TwoPortEquipment
inStream, outStreamFields inherited from class ProcessEquipmentBaseClass
conditionAnalysisMessage, energyStream, hasController, isSolved, properties, reportFields inherited from class SimulationBaseClass
calcIdentifier, calculateSteadyState, timeFields inherited from class NamedBaseClass
name -
Constructor Summary
ConstructorsConstructorDescriptionTwoFluidPipe(String name) Constructor with name only.TwoFluidPipe(String name, StreamInterface inStream) Constructor with inlet stream. -
Method Summary
Modifier and TypeMethodDescriptionprivate voidApply boundary conditions.private doubleCalculate stable time step using CFL condition.private double[]calculateLocalHoldup(TwoFluidSection sec, TwoFluidSection prev, double mDotGas, double mDotLiq, double area) Calculate local liquid holdup using drift-flux model with terrain effects.private doubleEstimate pressure gradient for steady-state initialization.Get accumulation tracker for detailed analysis.doubleGet pipe diameter.doubleGet distance to first hydrate risk location.intGet first section index with hydrate risk.Get flow regime at each section.double[]Get gas velocity profile.doubleGet the heat transfer coefficient.double[]Get the heat transfer coefficient profile.doubleGet hydrate formation temperature.intGet number of sections with hydrate risk.boolean[]Get sections with hydrate formation risk.Get the current insulation type.doubleGet time of last slug arrival at outlet.doubleGet pipe length.double[]Get liquid holdup profile.doublegetLiquidInventory(String unit) Get total liquid inventory in the pipe.double[]Get liquid velocity profile.doubleGet maximum simulation time.doubleGet maximum slug length observed at outlet.doubleGet maximum slug volume observed at outlet.intGet number of sections.double[]Get oil holdup profile along the pipeline.double[]Get oil velocity profile along the pipeline.double[]Get oil-water velocity slip profile along the pipeline.intGet number of slugs that have arrived at outlet.double[]Get position array for plotting.double[]Get pressure profile.doubleGet pipe wall roughness.doubleGet current simulation time.Get slug statistics summary string.Get slug tracker for slug statistics.doubleGet soil thermal resistance.doubleGet the surface temperature used for heat transfer calculations.double[]Get the surface temperature profile.double[]Get temperature profile.double[]getTemperatureProfile(String unit) Get temperature profile with specified unit.doubleGet total liquid volume delivered by slugs at outlet.double[]Get the pipe wall temperature profile.doubleGet pipe wall thickness.double[]Get water cut profile along the pipeline.double[]Get water holdup profile along the pipeline.double[]Get water velocity profile along the pipeline.doubleGet wax appearance temperature.boolean[]Get sections with wax deposition risk.booleanCheck if any section has hydrate risk.booleanCheck if any section has wax risk.private voidInitialize pipe sections with inlet conditions.private voidInitialize sub-models.booleanCheck if heat transfer is enabled.booleanCheck if Joule-Thomson effect is enabled.booleanCheck if water-oil velocity slip modeling is enabled.voidIn this method all thermodynamic and unit operations will be calculated in a steady state calculation.private voidRun steady-state initialization.voidrunTransient(double dt, UUID id) Run transient simulation for specified time step.voidsetCflNumber(double cfl) Set CFL number for time stepping.voidsetDiameter(double diameter) Set pipe diameter.voidsetElevationProfile(double[] elevations) Set elevation profile.voidsetEnableJouleThomson(boolean enable) Enable or disable Joule-Thomson effect.voidsetEnableSlugTracking(boolean enable) Enable/disable slug tracking.voidsetEnableWaterOilSlip(boolean enable) Enable water-oil velocity slip modeling.voidsetHeatTransferCoefficient(double heatTransferCoefficient) Set heat transfer coefficient for convective heat transfer.voidsetHeatTransferProfile(double[] profile) Set heat transfer coefficient profile along the pipe.voidsetHydrateFormationTemperature(double temperature, String unit) Set hydrate formation temperature for risk monitoring.voidsetIncludeEnergyEquation(boolean include) Enable/disable energy equation.voidsetIncludeMassTransfer(boolean include) Enable/disable mass transfer (flashing/condensation).voidSet insulation type using predefined U-values.voidsetLength(double length) Set pipe length.voidsetMaxSimulationTime(double time) Set maximum simulation time for transient calculations.voidsetNumberOfSections(int numberOfSections) Set number of computational sections.voidsetOutletPressure(double pressure) Set outlet pressure.voidsetOutletPressure(double pressure, String unit) Set outlet pressure with unit.voidsetRoughness(double roughness) Set pipe wall roughness.voidsetSoilThermalResistance(double resistance) Set soil/burial thermal resistance.voidsetSurfaceTemperature(double temperature, String unit) Set surface temperature for heat transfer calculations.voidsetSurfaceTemperatureProfile(double[] profile) Set surface temperature profile along the pipe.voidsetThermodynamicUpdateInterval(int interval) Set thermodynamic update interval.voidsetWallProperties(double thickness, double density, double heatCapacity) Set pipe wall properties for transient thermal calculations.voidsetWaxAppearanceTemperature(double temperature, String unit) Set wax appearance temperature for risk monitoring.private voidTrack slugs arriving at outlet and collect statistics.private voidUpdate outlet stream with current outlet conditions.private voidUpdate result arrays from section states.private voidupdateTemperatureProfile(double massFlow, double area) Update temperature profile along the pipe accounting for heat transfer.private voidUpdate thermodynamic properties using flash calculations.private voidupdateTransientTemperature(double massFlow, double area, double dt) Update temperature profile for transient simulation including pipe wall thermal mass.private voidupdateWaterOilHoldups(TwoFluidSection sec, TwoFluidSection prev, double alphaL, double area) Update water and oil holdups for three-phase flow with terrain effects.private voidvalidateAndCorrectState(double[][] U_new, double[][] U_prev) Validate and correct state vector to prevent numerical instabilities.private voidValidate section states and fix any numerical issues.Methods inherited from class Pipeline
displayResult, getCapacityDuty, getCapacityMax, getEntropyProduction, getMechanicalDesign, getOutletPressure, getPipe, getSuperficialVelocity, getTimes, initMechanicalDesign, setEquilibriumHeatTransfer, setEquilibriumMassTransfer, setHeightProfile, setInitialFlowPattern, setLegPositions, setNumberOfLegs, setNumberOfNodesInLeg, setOuterTemperatures, setOutputFileName, setPipeDiameters, setPipeOuterHeatTransferCoefficients, setPipeWallHeatTransferCoefficients, setPipeWallRoughness, setTimeSeries, toJson, toJsonMethods inherited from class TwoPortEquipment
getInletPressure, getInletStream, getInletTemperature, getMassBalance, getOutletPressure, getOutletStream, getOutletTemperature, setInletPressure, setInletStream, setInletTemperature, setOutletStream, setOutletTemperature, validateSetupMethods inherited from class ProcessEquipmentBaseClass
copy, equals, getConditionAnalysisMessage, getController, getEnergyStream, getExergyChange, getMassBalance, getMinimumFlow, getPressure, getPressure, getProperty, getReport_json, getResultTable, getSpecification, getTemperature, getTemperature, getThermoSystem, hashCode, isActive, isActive, isSetEnergyStream, reportResults, run_step, runConditionAnalysis, setController, setEnergyStream, setEnergyStream, setFlowValveController, setMinimumFlow, setPressure, setRegulatorOutSignal, setSpecification, setTemperature, solvedMethods inherited from class SimulationBaseClass
getCalculateSteadyState, getCalculationIdentifier, getTime, increaseTime, isRunInSteps, setCalculateSteadyState, setCalculationIdentifier, setRunInSteps, setTimeMethods inherited from class NamedBaseClass
getName, getTagName, setName, setTagNameMethods inherited from class Object
clone, finalize, getClass, notify, notifyAll, toString, wait, wait, waitMethods inherited from interface NamedInterface
getName, getTagName, setName, setTagNameMethods inherited from interface ProcessEquipmentInterface
getExergyChange, getFluid, getRestCapacity, needRecalculationMethods inherited from interface SimulationInterface
getCalculateSteadyState, getCalculationIdentifier, getTime, increaseTime, isRunInSteps, run, run_step, run_step, runTransient, setCalculateSteadyState, setCalculationIdentifier, setRunInSteps, setTime, solvedMethods inherited from interface TwoPortInterface
getInletPressure, getInletStream, getInletTemperature, getInStream, getOutletPressure, getOutletStream, getOutletTemperature, getOutStream, setInletPressure, setInletStream, setInletTemperature, setOutletStream, setOutletTemperature
-
Field Details
-
serialVersionUID
private static final long serialVersionUID- See Also:
-
logger
private static final org.apache.logging.log4j.Logger logger -
length
private double lengthTotal pipe length (m). -
diameter
private double diameterPipe inner diameter (m). -
roughness
private double roughnessPipe wall roughness (m). -
numberOfSections
private int numberOfSectionsNumber of computational cells. -
elevationProfile
private double[] elevationProfileElevation profile at each section (m). -
sections
Pipe sections with state. -
dx
private double dxSpatial step size (m). -
simulationTime
private double simulationTimeCurrent simulation time (s). -
maxSimulationTime
private double maxSimulationTimeMaximum simulation time (s). -
cflNumber
private double cflNumberCFL number for time stepping (0 < CFL < 1). -
equations
Conservation equations solver. -
timeIntegrator
Time integrator. -
flowRegimeDetector
Flow regime detector. -
accumulationTracker
Liquid accumulation tracker. -
slugTracker
Slug tracker. -
inletBCType
Inlet boundary condition type. -
outletBCType
Outlet boundary condition type. -
outletPressure
private double outletPressureOutlet pressure (Pa). -
outletPressureSet
private boolean outletPressureSetFlag indicating if outlet pressure was explicitly set. -
includeEnergyEquation
private boolean includeEnergyEquationInclude energy equation. -
includeMassTransfer
private boolean includeMassTransferInclude mass transfer (flash/condensation). -
enableHeatTransfer
private boolean enableHeatTransferEnable heat transfer from surroundings. -
surfaceTemperature
private double surfaceTemperatureSurface temperature for heat transfer (K). -
heatTransferCoefficient
private double heatTransferCoefficientHeat transfer coefficient (W/(m²·K)). -
heatTransferProfile
private double[] heatTransferProfileHeat transfer coefficient profile along pipe (W/(m²·K)). -
surfaceTemperatureProfile
private double[] surfaceTemperatureProfileSurface temperature profile along pipe (K). -
wallThickness
private double wallThicknessPipe wall thickness (m). -
wallDensity
private double wallDensityPipe wall density (kg/m³) - steel default. -
wallHeatCapacity
private double wallHeatCapacityPipe wall specific heat capacity (J/(kg·K)) - steel default. -
wallTemperatureProfile
private double[] wallTemperatureProfilePipe wall temperature profile (K). -
soilThermalResistance
private double soilThermalResistanceSoil/burial thermal resistance (m²·K/W). -
enableJouleThomson
private boolean enableJouleThomsonEnable Joule-Thomson effect. -
hydrateFormationTemperature
private double hydrateFormationTemperatureHydrate formation temperature (K). -
waxAppearanceTemperature
private double waxAppearanceTemperatureWax appearance temperature (K). -
hydrateRiskSections
private boolean[] hydrateRiskSectionsSections flagged for hydrate risk. -
waxRiskSections
private boolean[] waxRiskSectionsSections flagged for wax risk. -
insulationType
Current insulation type. -
enableSlugTracking
private boolean enableSlugTrackingEnable slug tracking. -
outletSlugCount
private int outletSlugCountOutlet slug statistics. -
totalSlugVolumeAtOutlet
private double totalSlugVolumeAtOutlet -
lastSlugArrivalTime
private double lastSlugArrivalTime -
maxSlugLengthAtOutlet
private double maxSlugLengthAtOutlet -
maxSlugVolumeAtOutlet
private double maxSlugVolumeAtOutlet -
countedOutletSlugs
-
thermodynamicUpdateInterval
private int thermodynamicUpdateIntervalUpdate thermodynamics every N steps. -
currentStep
private int currentStepCurrent step count. -
pressureProfile
private double[] pressureProfilePressure profile (Pa). -
temperatureProfile
private double[] temperatureProfileTemperature profile (K). -
liquidHoldupProfile
private double[] liquidHoldupProfileLiquid holdup profile. -
gasVelocityProfile
private double[] gasVelocityProfileGas velocity profile (m/s). -
liquidVelocityProfile
private double[] liquidVelocityProfileLiquid velocity profile (m/s). -
referenceFluid
Reference fluid for flash calculations.
-
-
Constructor Details
-
TwoFluidPipe
-
TwoFluidPipe
Constructor with inlet stream.- Parameters:
name- Equipment nameinStream- 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)
- Parameters:
sec- Current sectionprev- 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
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 sectionprev- Previous section (upstream)alphaL- Total liquid holduparea- Pipe cross-sectional area
-
run
-
runTransient
Run transient simulation for specified time step.- Specified by:
runTransientin interfaceSimulationInterface- Overrides:
runTransientin classPipeline- 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 validateU_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
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
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
Get accumulation tracker for detailed analysis.- Returns:
- Accumulation tracker
-
getSlugTracker
-
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
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:
setOutletPressurein interfaceTwoPortInterface- Overrides:
setOutletPressurein classTwoPortEquipment- Parameters:
pressure- Pressure (Pa)
-
setOutletPressure
Set outlet pressure with unit.- Parameters:
pressure- Pressure valueunit- 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
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 unitunit- 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
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
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
Set hydrate formation temperature for risk monitoring.- Parameters:
temperature- Hydrate formation temperatureunit- Temperature unit ("K" or "C")
-
getHydrateFormationTemperature
public double getHydrateFormationTemperature()Get hydrate formation temperature.- Returns:
- Hydrate formation temperature [K], or 0 if not set
-
setWaxAppearanceTemperature
Set wax appearance temperature for risk monitoring.- Parameters:
temperature- Wax appearance temperatureunit- 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
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
-