Interface PhaseInterface

All Superinterfaces:
Cloneable, Serializable, ThermodynamicConstantsInterface
All Known Subinterfaces:
PhaseCPAInterface, PhaseEosInterface
All Known Implementing Classes:
Phase, PhaseAmmoniaEos, PhaseBNS, PhaseBWRSEos, PhaseCSPsrkEos, PhaseDefault, PhaseDesmukhMather, PhaseDuanSun, PhaseElectrolyteCPA, PhaseElectrolyteCPAMM, PhaseElectrolyteCPAOld, PhaseElectrolyteCPAstatoil, PhaseEos, PhaseEOSCGEos, PhaseGE, PhaseGENRTL, PhaseGENRTLmodifiedHV, PhaseGENRTLmodifiedWS, PhaseGERG2004Eos, PhaseGERG2008Eos, PhaseGEUnifac, PhaseGEUnifacPSRK, PhaseGEUnifacUMRPRU, PhaseGEUniquac, PhaseGEUniquacmodifiedHV, PhaseGEWilson, PhaseHydrate, PhaseIdealGas, PhaseKentEisenberg, PhaseLeachmanEos, PhaseModifiedFurstElectrolyteEos, PhaseModifiedFurstElectrolyteEosMod2004, PhasePCSAFT, PhasePCSAFTa, PhasePCSAFTRahmat, PhasePitzer, PhasePrCPA, PhasePrEos, PhasePrEosvolcor, PhasePureComponentSolid, PhaseRK, PhaseSolid, PhaseSolidComplex, PhaseSoreideWhitson, PhaseSpanWagnerEos, PhaseSrkCPA, PhaseSrkCPAs, PhaseSrkEos, PhaseSrkEosvolcor, PhaseSrkPenelouxEos, PhaseTSTEos, PhaseUMRCPA, PhaseVegaEos, PhaseWaterIAPWS, PhaseWax

public interface PhaseInterface extends ThermodynamicConstantsInterface, Cloneable

PhaseInterface interface.

Version:
$Id: $Id
Author:
Even Solbraa
  • Method Details

    • addComponent

      void addComponent(String name, double moles, double molesInPhase, int compIndex)

      Add component to component array and update moles variables.

      Parameters:
      name - Name of component.
      moles - Total number of moles of component.
      molesInPhase - Number of moles in phase.
      compIndex - Index number of component in phase object component array.
    • setMoleFractions

      void setMoleFractions(double[] x)

      Set x and normalize for all Components in phase.

      Parameters:
      x - Mole fractions of component in a phase.
    • getPhaseFraction

      default double getPhaseFraction()

      getPhaseFraction. Alias for getBeta()

      Returns:
      Beta value
    • getComposition

      double[] getComposition(String unit)

      Returns the composition vector in unit molefraction/wtfraction/molespersec/volumefraction.

      Parameters:
      unit - Supported units are molefraction, wtfraction, molespersec, volumefraction
      Returns:
      composition array with unit
    • getCp0

      double getCp0()

      getCp0.

      Returns:
      a double
    • getDensity_AGA8

      double getDensity_AGA8()
      Get density of a phase using the AGA8-Detail EoS.
      Returns:
      density with unit kg/m3
    • getJouleThomsonCoefficient

      double getJouleThomsonCoefficient()
      Get the Joule Thomson Coefficient of a phase.
      Returns:
      Joule Thomson coefficient in K/bar
    • getJouleThomsonCoefficient

      double getJouleThomsonCoefficient(String unit)
      Get the Joule Thomson Coefficient of a phase.
      Parameters:
      unit - Supported units are K/bar, C/bar
      Returns:
      Joule Thomson coefficient in specified unit
    • getMolarComposition

      double[] getMolarComposition()
      Returns the mole composition vector in unit mole fraction.
      Returns:
      an array of type double
    • getVolume

      double getVolume()
      method to return phase volume note: without Peneloux volume correction.
      Returns:
      volume in unit m3*1e5
    • getVolume

      double getVolume(String unit)
      method to return fluid volume.
      Parameters:
      unit - Supported units are m3, litre
      Returns:
      volume in specified unit
    • getGamma

      double getGamma()
      method to return heat capacity ratio/adiabatic index/Poisson constant. The method calculates it as Cp (real) /Cv (real).
      Returns:
      gamma
    • getGamma2

      default double getGamma2()
      method to return heat capacity ratio calculated as Cp/(Cp-R*getNumberOfMolesInPhase).
      Returns:
      kappa
    • getComponentWithIndex

      ComponentInterface getComponentWithIndex(int index)

      getComponentWithIndex.

      Parameters:
      index - a int
      Returns:
      a ComponentInterface object
    • getWtFractionOfWaxFormingComponents

      double getWtFractionOfWaxFormingComponents()

      getWtFractionOfWaxFormingComponents.

      Returns:
      a double
    • getCompressibilityX

      double getCompressibilityX()

      getCompressibilityX.

      Returns:
      a double
    • getCompressibilityY

      double getCompressibilityY()

      getCompressibilityY.

      Returns:
      a double
    • getIsothermalCompressibility

      double getIsothermalCompressibility()

      getIsothermalCompressibility.

      Returns:
      a double
    • getIsobaricThermalExpansivity

      double getIsobaricThermalExpansivity()

      getIsobaricThermalExpansivity.

      Returns:
      a double
    • getdrhodN

      double getdrhodN()

      getdrhodN.

      Returns:
      a double
    • setInitType

      void setInitType(int initType)

      setInitType.

      Parameters:
      initType - a int
    • init

      default void init()

      Init using current phase properties.

    • init

      default void init(double totalNumberOfMoles, int numberOfComponents, int initType, double beta)

      init. Uses existing phase type.

      Parameters:
      totalNumberOfMoles - Total number of moles in all phases of Stream.
      numberOfComponents - Number of components in system.
      initType - a int. Use 0 to init, and 1 to reset.
      beta - Mole fraction of this phase in system.
    • init

      void init(double totalNumberOfMoles, int numberOfComponents, int initType, PhaseType pt, double beta)

      init.

      Parameters:
      totalNumberOfMoles - Total number of moles in all phases of Stream.
      numberOfComponents - Number of components in system.
      initType - a int. Use 0 to init, and 1 to reset.
      pt - Type of phase.
      beta - Mole fraction of this phase in system.
    • initPhysicalProperties

      void initPhysicalProperties()

      Init / calculate all physical properties of phase.

    • initPhysicalProperties

      void initPhysicalProperties(PhysicalPropertyType ppt)

      Init / calculate specified physical property of phase.

      Parameters:
      ppt - PhysicalPropertyType enum object.
    • initPhysicalProperties

      default void initPhysicalProperties(String name)
      Initialize / calculate a specified physical properties for all phases and interfaces.
      Parameters:
      name - Name of physical property.
    • getdrhodP

      double getdrhodP()

      getdrhodP.

      Returns:
      a double
    • getdrhodT

      double getdrhodT()

      getdrhodT.

      Returns:
      a double
    • getEnthalpydP

      double getEnthalpydP()

      getEnthalpydP.

      Returns:
      a double
    • getEnthalpydT

      double getEnthalpydT()

      getEnthalpydT.

      Returns:
      a double
    • getEntropydP

      double getEntropydP()

      getEntropydP.

      Returns:
      a double
    • getEntropydT

      double getEntropydT()

      getEntropydT.

      Returns:
      a double
    • getMoleFraction

      double getMoleFraction()

      getMoleFraction.

      Returns:
      a double
    • getcomponentArray

      ComponentInterface[] getcomponentArray()

      Get component array of Phase.

      Returns:
      an array of ComponentInterface objects
    • getComponentNames

      String[] getComponentNames()
      Get normalized names of components in phase.
      Returns:
      Array of names of components in phase.
    • getMass

      double getMass()

      getMass.

      Returns:
      a double
    • getWtFraction

      double getWtFraction(SystemInterface system)

      getWtFraction.

      Parameters:
      system - a SystemInterface object
      Returns:
      a double
    • getMolarVolume

      double getMolarVolume()
      method to return molar volume of the phase note: without Peneloux volume correction.
      Returns:
      molar volume volume in unit m3/mol*1e5
    • getMolarVolume

      double getMolarVolume(String unit)
      method to return molar volume of the fluid: eventual volume correction included.
      Parameters:
      unit - Supported units are m3/mol, litre/mol
      Returns:
      molar volume volume in unit
    • getFlowRate

      double getFlowRate(String flowunit)
      method to return flow rate of a phase.
      Parameters:
      flowunit - Supported units are kg/sec, kg/min, kg/hr, m3/sec, m3/min, m3/hr, ft3/sec, mole/sec, mole/min, mole/hr
      Returns:
      flow rate in specified unit
    • setComponentArray

      void setComponentArray(ComponentInterface[] components)

      setComponentArray.

      Parameters:
      components - an array of ComponentInterface objects
    • getDensity_GERG2008

      double getDensity_GERG2008()
      Get density of a phase using the GERG-2008 EoS.
      Returns:
      density with unit kg/m3
    • getProperties_GERG2008

      double[] getProperties_GERG2008()

      Get GERG properties of a phase using the GERG-2008 EoS.

      Returns:
      an array of type double
    • getAlpha0_GERG2008

      org.netlib.util.doubleW[] getAlpha0_GERG2008()
      Overloaded method to get the Leachman a0matrix with default hydrogen type ('normal').
      Returns:
      matrix of the reduced ideal helmholtz free energy and its derivatives
    • getAlphares_GERG2008

      org.netlib.util.doubleW[][] getAlphares_GERG2008()
      Overloaded method to get the Leachman armatrix with default hydrogen type ('normal').
      Returns:
      matrix of the reduced residual helmholtz free energy and its derivatives
    • getDensity_EOSCG

      double getDensity_EOSCG()
      Get density of a phase using the EOS-CG EoS.
      Returns:
      density with unit kg/m3
    • getProperties_EOSCG

      double[] getProperties_EOSCG()
      Get EOS-CG properties of a phase using the EOS-CG model.
      Returns:
      an array of type double
    • getAlpha0_EOSCG

      org.netlib.util.doubleW[] getAlpha0_EOSCG()
      Get EOS-CG ideal Helmholtz contribution and derivatives.
      Returns:
      matrix of the reduced ideal helmholtz free energy and its derivatives
    • getAlphares_EOSCG

      org.netlib.util.doubleW[][] getAlphares_EOSCG()
      Get EOS-CG residual Helmholtz contribution and derivatives.
      Returns:
      matrix of the reduced residual helmholtz free energy and its derivatives
    • getDensity_Leachman

      double getDensity_Leachman(String hydrogenType)
      method to get Leachman density of a phase using the Leachman EoS.
      Parameters:
      hydrogenType - Supported types are 'normal', 'para', 'ortho'
      Returns:
      density with unit kg/m3
    • getDensity_Leachman

      double getDensity_Leachman()
      Overloaded method to get the Leachman density with default hydrogen type ('normal').
      Returns:
      density with unit kg/m3
    • getProperties_Leachman

      double[] getProperties_Leachman(String hydrogenType)

      method to get Leachman properties of a phase using the Leachman EoS.

      Parameters:
      hydrogenType - a String object
      Returns:
      an array of type double
    • getProperties_Leachman

      double[] getProperties_Leachman()
      Overloaded method to get the Leachman properties with default hydrogen type ('normal').
      Returns:
      density with unit kg/m3
    • getAlpha0_Leachman

      org.netlib.util.doubleW[] getAlpha0_Leachman()
      Overloaded method to get the Leachman a0matrix with default hydrogen type ('normal').
      Returns:
      matrix of the reduced ideal helmholtz free energy and its derivatives
    • getAlpha0_Leachman

      org.netlib.util.doubleW[] getAlpha0_Leachman(String hydrogenType)

      method to get Leachman a0matrix of a phase using the Leachman EoS.

      Parameters:
      hydrogenType - Supported types are 'normal', 'para', 'ortho'
      Returns:
      matrix of the reduced ideal helmholtz free energy and its derivatives
    • getAlphares_Leachman

      org.netlib.util.doubleW[][] getAlphares_Leachman(String hydrogenType)

      method to get Leachman armatrix of a phase using the Leachman EoS.

      Parameters:
      hydrogenType - Supported types are 'normal', 'para', 'ortho'
      Returns:
      matrix of the reduced residual helmholtz free energy and its derivatives
    • getAlphares_Leachman

      org.netlib.util.doubleW[][] getAlphares_Leachman()
      Overloaded method to get the Leachman armatrix with default hydrogen type ('normal').
      Returns:
      matrix of the reduced residual helmholtz free energy and its derivatives
    • getDensity_Vega

      double getDensity_Vega()
      method to get helium density of a phase using the Vega EoS.
      Returns:
      density with unit kg/m3
    • getProperties_Vega

      double[] getProperties_Vega()

      method to get helium properties of a phase using the Vega EoS.

      Returns:
      an array of type double
    • getAlpha0_Vega

      org.netlib.util.doubleW[] getAlpha0_Vega()
      Overloaded method to get the Leachman a0matrix with default hydrogen type ('normal').
      Returns:
      matrix of the reduced ideal helmholtz free energy and its derivatives
    • getAlphares_Vega

      org.netlib.util.doubleW[][] getAlphares_Vega()
      Overloaded method to get the Leachman armatrix with default hydrogen type ('normal').
      Returns:
      matrix of the reduced residual helmholtz free energy and its derivatives
    • getDensity

      double getDensity()
      Get density of a phase note: does not use Peneloux volume correction.
      Returns:
      density with unit kg/m3
    • getDensity

      double getDensity(String unit)
      Get density of a fluid note: with Peneloux volume correction.
      Parameters:
      unit - Supported units are kg/m3, mol/m3
      Returns:
      density in specified unit
    • getWaterDensity

      double getWaterDensity(String unit)
      method to get density of a water phase using specific water method.
      Parameters:
      unit - Supported units are kg/m3, mol/m3
      Returns:
      density in specified unit
    • removeComponent

      void removeComponent(String name, double moles, double molesInPhase)

      Remove component from Phase.

      Parameters:
      name - Name of component.
      moles - Total number of moles of component.
      molesInPhase - Number of moles in phase.
    • getFugacity

      double getFugacity(int compNumb)

      getFugacity.

      Parameters:
      compNumb - a int
      Returns:
      a double
    • getFugacity

      double getFugacity(String compName)

      getFugacity.

      Parameters:
      compName - a String object
      Returns:
      a double
    • getTotalVolume

      double getTotalVolume()
      method to return phase volume note: without Peneloux volume correction.
      Returns:
      volume in unit m3*1e5
    • getCorrectedVolume

      double getCorrectedVolume()
      method to return phase volume with Peneloux volume correction need to call initPhysicalProperties() before this method is called.
      Returns:
      volume in unit m3
    • hasTBPFraction

      boolean hasTBPFraction()

      hasTBPFraction.

      Returns:
      a boolean
    • getMolalMeanIonicActivity

      double getMolalMeanIonicActivity(int comp1, int comp2)
      Calculates the mean ionic activity coefficient on the molality scale for an electrolyte defined by two ionic species. The conversion from mole fraction scale (γ±,x) to molality scale (γ±,m) follows: γ±,m = γ±,x * x_water

      Reference: Robinson, R.A. and Stokes, R.H. "Electrolyte Solutions", 2nd ed., Butterworths, London, 1965.

      Parameters:
      comp1 - component index of the first ion (e.g., cation)
      comp2 - component index of the second ion (e.g., anion)
      Returns:
      mean ionic activity coefficient on the molality scale
    • getGibbsEnergy

      double getGibbsEnergy()

      getGibbsEnergy.

      Returns:
      a double
    • getMixGibbsEnergy

      double getMixGibbsEnergy()

      getMixGibbsEnergy.

      Returns:
      a double
    • getExcessGibbsEnergy

      double getExcessGibbsEnergy()

      getExcessGibbsEnergy.

      Returns:
      a double
    • getExcessGibbsEnergySymetric

      double getExcessGibbsEnergySymetric()

      getExcessGibbsEnergySymetric.

      Returns:
      a double
    • hasPlusFraction

      boolean hasPlusFraction()

      hasPlusFraction.

      Returns:
      a boolean
    • getSresTP

      double getSresTP()

      getSresTP.

      Returns:
      a double
    • setProperties

      void setProperties(PhaseInterface phase)

      setProperties. Transfer properties from another phase object.

      Parameters:
      phase - a PhaseInterface object
    • useVolumeCorrection

      void useVolumeCorrection(boolean volcor)

      useVolumeCorrection.

      Parameters:
      volcor - a boolean
    • useVolumeCorrection

      boolean useVolumeCorrection()

      useVolumeCorrection.

      Returns:
      a boolean
    • getBeta

      double getBeta()

      Getter for property beta. Beta is the mole fraction of a phase of all the moles of a system.

      Returns:
      Beta value
    • setBeta

      void setBeta(double b)

      Setter for property beta. Beta is the mole fraction of a phase of all the moles of a system.

      Parameters:
      b - Beta value to set.
    • getWtFrac

      double getWtFrac(int component)

      getWtFrac.

      Parameters:
      component - a int
      Returns:
      a double
    • getWtFrac

      double getWtFrac(String componentName)

      getWtFrac.

      Parameters:
      componentName - a String object
      Returns:
      a double
    • setMixingRuleGEModel

      void setMixingRuleGEModel(String name)

      setMixingRuleGEModel.

      Parameters:
      name - a String object
    • getComponent

      ComponentInterface getComponent(String name)

      Get Component by name.

      Parameters:
      name - Name of component
      Returns:
      a ComponentInterface object
    • getComponent

      ComponentInterface getComponent(int i)

      Get Component by index.

      Parameters:
      i - Component index
      Returns:
      a ComponentInterface object
    • getActivityCoefficient

      double getActivityCoefficient(int k)

      getActivityCoefficient.

      Parameters:
      k - a int
      Returns:
      a double
    • getActivityCoefficient

      double getActivityCoefficient(int k, int p)

      getActivityCoefficient.

      Parameters:
      k - a int
      p - a int
      Returns:
      a double
    • getActivityCoefficient

      double getActivityCoefficient(int k, String scale)
      Get activity coefficient on a specified concentration scale.
      Parameters:
      k - component index
      scale - concentration scale: "molefraction" (default), "molality" (mol/kg solvent), or "molarity" (mol/L solution)
      Returns:
      activity coefficient on the specified scale
    • setPressure

      void setPressure(double pres)

      Set the pressure in bara (absolute pressure in bar).

      Parameters:
      pres - a double
    • getpH

      double getpH()
      Calculate pH of the aqueous phase using the IUPAC standard definition.

      Uses activity-based pH calculation since NeqSim's chemical reaction equilibrium constants are defined on the mole fraction scale: pH = -log10(gamma_x * x_H3O+) where gamma_x is the activity coefficient and x_H3O+ is the mole fraction of H3O+.

      Returns:
      pH value
    • getpH

      double getpH(String method)
      Calculate pH of the phase using specified method.

      Available methods:

      • activity (default): pH = -log10(gamma_x * x_H3O+) - consistent with mole fraction-based equilibrium constants
      • molality (IUPAC standard): pH = -log10(gamma_m * m_H3O+) - correct for all concentrations
      • molarity: pH = -log10([H3O+]) where [H3O+] is in mol/L - ignores activity coefficient
      Parameters:
      method - The calculation method: "activity" (default), "molality", or "molarity"
      Returns:
      pH value
    • normalize

      void normalize()
      Normalize property x.

      Property x is the mole fraction of a component in a specific phase. Normalizing, means that the sum of x for all Components in a phase equal 1.0.

    • getLogPureComponentFugacity

      double getLogPureComponentFugacity(int k)

      getLogPureComponentFugacity.

      Parameters:
      k - a int
      Returns:
      a double
    • getPureComponentFugacity

      double getPureComponentFugacity(int k)

      getPureComponentFugacity.

      Parameters:
      k - a int
      Returns:
      a double
    • getPureComponentFugacity

      double getPureComponentFugacity(int k, boolean pure)

      getPureComponentFugacity.

      Parameters:
      k - a int
      pure - a boolean
      Returns:
      a double
    • addMoles

      default void addMoles(int component, double dn)

      Change the number of moles of component of phase,i.e., numberOfMolesInPhase but do not change the total number of moles of component in system. NB! Phase fraction beta is not updated by this method. Must be done separately to keep consistency between phase and component calculation of total number of moles in system.

      Parameters:
      component - Component number to change
      dn - Number of moles of component added to phase
    • addMolesChemReac

      default void addMolesChemReac(int component, double dn)

      Change the number of moles of component of phase, i.e., numberOfMolesInPhase, and total number of moles of component in system, i.e., numberOfMoles with the same amount. NB! Phase fraction beta is not updated by this method. Must be done separately to keep consistency between phase and component calculation of total number of moles in system.

      Parameters:
      component - Component number to change
      dn - Number of moles of component added to phase and system
    • addMolesChemReac

      void addMolesChemReac(int component, double dn, double totdn)

      Change the number of moles of component of phase, i.e., numberOfMolesInPhase and Component properties for the number of moles of component of phase, i.e., numberOfMolesInPhase, and total number of moles of component in system, i.e., numberOfMoles with separate amounts. NB! Phase fraction beta is not updated by this method. Must be done separately to keep consistency between phase and component calculation of total number of moles in system.

      Parameters:
      component - Component number to change
      dn - Number of moles of component to add to phase
      totdn - Number of moles of component to add to system
    • setEmptyFluid

      void setEmptyFluid()

      Set the flow rate (moles) of all components to zero.

    • getPhysicalProperties

      PhysicalProperties getPhysicalProperties()

      getPhysicalProperties.

      Returns:
      a PhysicalProperties object
    • resetPhysicalProperties

      void resetPhysicalProperties()

      resetPhysicalProperties.

    • setPhysicalProperties

      default void setPhysicalProperties()

      Specify the type of the physical properties.

    • setPhysicalProperties

      void setPhysicalProperties(PhysicalPropertyModel ppm)

      Specify the type of the physical properties.

      Parameters:
      ppm - PhysicalPropertyModel enum object
    • getPhysicalPropertyModel

      PhysicalPropertyModel getPhysicalPropertyModel()

      Getter for property ppm.

      Returns:
      a PhysicalPropertyModel object
    • setPpm

      void setPpm(PhysicalPropertyModel ppm)

      Setter for property ppm.

      Parameters:
      ppm - PhysicalPropertyModel enum object
    • setPhysicalPropertyModel

      void setPhysicalPropertyModel(PhysicalPropertyModel ppm)

      setPhysicalPropertyModel.

      Parameters:
      ppm - a PhysicalPropertyModel object
    • getInitType

      int getInitType()

      getInitType.

      Returns:
      a int
    • setAttractiveTerm

      void setAttractiveTerm(int i)

      setAttractiveTerm.

      Parameters:
      i - a int
    • setMixingRule

      void setMixingRule(MixingRuleTypeInterface mr)

      setMixingRule.

      Parameters:
      mr - a MixingRuleTypeInterface
    • setMixingRule

      default void setMixingRule(int mr)

      setMixingRule.

      Parameters:
      mr - a int
    • resetMixingRule

      void resetMixingRule(MixingRuleTypeInterface mr)

      resetMixingRule.

      Parameters:
      mr - a int
    • setTemperature

      void setTemperature(double temperature)
      Set the temperature of the phase.
      Parameters:
      temperature - in unit Kelvin
    • molarVolume

      double molarVolume(double pressure, double temperature, double A, double B, PhaseType pt) throws IsNaNException, TooManyIterationsException

      molarVolume.

      Parameters:
      pressure - a double
      temperature - a double
      A - a double
      B - a double
      pt - the PhaseType of the phase
      Returns:
      a double
      Throws:
      IsNaNException - if any.
      TooManyIterationsException - if any.
    • geta

      double geta(PhaseInterface phase, double temperature, double pressure, int numbcomp)

      geta.

      Parameters:
      phase - a PhaseInterface object
      temperature - a double
      pressure - a double
      numbcomp - a int
      Returns:
      a double
    • getb

      double getb(PhaseInterface phase, double temperature, double pressure, int numbcomp)

      getb.

      Parameters:
      phase - a PhaseInterface object
      temperature - a double
      pressure - a double
      numbcomp - a int
      Returns:
      a double
    • getAntoineVaporPressure

      double getAntoineVaporPressure(double temp)

      getAntoineVaporPressure.

      Parameters:
      temp - a double
      Returns:
      a double
    • calcA

      double calcA(PhaseInterface phase, double temperature, double pressure, int numbcomp)

      calcA.

      Parameters:
      phase - a PhaseInterface object
      temperature - a double
      pressure - a double
      numbcomp - a int
      Returns:
      a double
    • calcB

      double calcB(PhaseInterface phase, double temperature, double pressure, int numbcomp)

      calcB.

      Parameters:
      phase - a PhaseInterface object
      temperature - a double
      pressure - a double
      numbcomp - a int
      Returns:
      a double
    • calcAi

      double calcAi(int compNumb, PhaseInterface phase, double temperature, double pressure, int numbcomp)

      calcAi.

      Parameters:
      compNumb - a int
      phase - a PhaseInterface object
      temperature - a double
      pressure - a double
      numbcomp - a int
      Returns:
      a double
    • calcAiT

      double calcAiT(int compNumb, PhaseInterface phase, double temperature, double pressure, int numbcomp)

      calcAiT.

      Parameters:
      compNumb - a int
      phase - a PhaseInterface object
      temperature - a double
      pressure - a double
      numbcomp - a int
      Returns:
      a double
    • calcAij

      double calcAij(int compNumb, int j, PhaseInterface phase, double temperature, double pressure, int numbcomp)

      calcAij.

      Parameters:
      compNumb - a int
      j - a int
      phase - a PhaseInterface object
      temperature - a double
      pressure - a double
      numbcomp - a int
      Returns:
      a double
    • calcBij

      double calcBij(int compNumb, int j, PhaseInterface phase, double temperature, double pressure, int numbcomp)

      calcBij.

      Parameters:
      compNumb - a int
      j - a int
      phase - a PhaseInterface object
      temperature - a double
      pressure - a double
      numbcomp - a int
      Returns:
      a double
    • calcAT

      double calcAT(int compNumb, PhaseInterface phase, double temperature, double pressure, int numbcomp)

      calcAT.

      Parameters:
      compNumb - a int
      phase - a PhaseInterface object
      temperature - a double
      pressure - a double
      numbcomp - a int
      Returns:
      a double
    • calcBi

      double calcBi(int compNumb, PhaseInterface phase, double temperature, double pressure, int numbcomp)

      calcBi.

      Parameters:
      compNumb - a int
      phase - a PhaseInterface object
      temperature - a double
      pressure - a double
      numbcomp - a int
      Returns:
      a double
    • calcR

      double calcR()

      calcR.

      Returns:
      a double
    • getg

      double getg()

      getg.

      Returns:
      a double
    • getEnthalpy

      double getEnthalpy()
      method to return enthalpy of a phase in unit Joule.
      Returns:
      a double
    • getEnthalpy

      double getEnthalpy(String unit)
      method to return phase enthalpy in a specified unit.
      Parameters:
      unit - Supported units are J, J/mol, kJ/kmol, J/kg and kJ/kg
      Returns:
      enthalpy in specified unit
    • getEntropy

      double getEntropy()
      method to return entropy of the phase.
      Returns:
      a double
    • getEntropy

      double getEntropy(String unit)
      method to return entropy of the phase.
      Parameters:
      unit - Supported units are J/K, J/moleK, J/kgK and kJ/kgK
      Returns:
      entropy in specified unit
    • getViscosity

      double getViscosity()
      method to return viscosity of the phase.
      Returns:
      viscosity in unit kg/msec
    • getViscosity

      double getViscosity(String unit)
      method to return viscosity of the phase in a specified unit.
      Parameters:
      unit - Supported units are kg/msec, Pas, cP (centipoise)
      Returns:
      viscosity in specified unit
    • getThermalConductivity

      double getThermalConductivity()
      method to return conductivity of a phase.
      Returns:
      conductivity in unit W/m*K
    • getThermalConductivity

      double getThermalConductivity(String unit)
      method to return conductivity in a specified unit.
      Parameters:
      unit - Supported units are W/mK, W/cmK
      Returns:
      conductivity in specified unit
    • getCp

      double getCp()
      method to return specific heat capacity (Cp).
      Returns:
      Cp in unit J/K
    • getCp

      double getCp(String unit)
      method to return specific heat capacity (Cp) in a specified unit.
      Parameters:
      unit - Supported units are J/K, J/molK, J/kgK and kJ/kgK
      Returns:
      Cp in specified unit
    • getHresTP

      double getHresTP()

      getHresTP.

      Returns:
      a double
    • getGresTP

      double getGresTP()

      getGresTP.

      Returns:
      a double
    • getCv

      double getCv()
      method to return specific heat capacity (Cv).
      Returns:
      Cv in unit J/K
    • getCv

      double getCv(String unit)
      method to return specific heat capacity (Cv) in a specified unit.
      Parameters:
      unit - Supported units are J/K, J/molK, J/kgK and kJ/kgK
      Returns:
      Cv in specified unit
    • getKappa

      double getKappa()
      method to return real gas isentropic exponent (kappa = - Cp/Cv*(v/p)*dp/dv method to return heat capacity ratio/adiabatic index/Poisson constant.
      Returns:
      kappa
    • getCpres

      double getCpres()

      getCpres.

      Returns:
      a double
    • getZ

      double getZ()

      getZ.

      Returns:
      a double
    • getPseudoCriticalPressure

      double getPseudoCriticalPressure()

      getPseudoCriticalPressure.

      Returns:
      a double
    • getPseudoCriticalTemperature

      double getPseudoCriticalTemperature()

      getPseudoCriticalTemperature.

      Returns:
      a double
    • getPhase

      PhaseInterface getPhase()

      getPhase.

      Returns:
      a PhaseInterface object
    • getNumberOfComponents

      int getNumberOfComponents()

      Get number of components added to Phase.

      Returns:
      the number of components in Phase.
    • setNumberOfComponents

      void setNumberOfComponents(int k)

      setNumberOfComponents.

      Parameters:
      k - a int
    • getComponents

      ComponentInterface[] getComponents()

      getComponents.

      Returns:
      an array of ComponentInterface objects
    • getNumberOfMolesInPhase

      double getNumberOfMolesInPhase()

      Get the number of moles the phase contains.

      Returns:
      The number of moles in the phase.
    • calcMolarVolume

      void calcMolarVolume(boolean test)

      calcMolarVolume.

      Parameters:
      test - a boolean
    • setTotalVolume

      void setTotalVolume(double volume)

      setTotalVolume.

      Parameters:
      volume - a double
    • setMolarVolume

      void setMolarVolume(double molarVolume)

      setMolarVolume.

      Parameters:
      molarVolume - a double
    • getInfiniteDiluteFugacity

      double getInfiniteDiluteFugacity(int k, int p)

      getInfiniteDiluteFugacity.

      Parameters:
      k - a int
      p - a int
      Returns:
      a double
    • getHelmholtzEnergy

      double getHelmholtzEnergy()

      getHelmholtzEnergy.

      Returns:
      a double
    • getNumberOfMolecularComponents

      int getNumberOfMolecularComponents()

      getNumberOfMolecularComponents.

      Returns:
      a int
    • getNumberOfIonicComponents

      int getNumberOfIonicComponents()

      getNumberOfIonicComponents.

      Returns:
      a int
    • getA

      double getA()

      getA.

      Returns:
      a double
    • getB

      double getB()

      getB.

      Returns:
      a double
    • getAT

      double getAT()

      getAT.

      Returns:
      a double
    • getATT

      double getATT()

      getATT.

      Returns:
      a double
    • clone

      clone.

      Returns:
      a PhaseInterface object
    • getTemperature

      double getTemperature()
      Get temperature of phase.
      Returns:
      temperature in unit K
    • getTemperature

      double getTemperature(String unit)
      Get temperature in a specified unit.
      Parameters:
      unit - Supported units are K, C, R
      Returns:
      temperature in specified unit
    • getPressure

      double getPressure()
      Get pressure of phase.
      Returns:
      pressure in unit bara
    • getPressure

      double getPressure(String unit)
      Get pressure of phase in a specified unit.
      Parameters:
      unit - Supported units are bara, barg, Pa, MPa, psi, psia, psig
      Returns:
      pressure in specified unit
    • getMolarMass

      double getMolarMass()
      Get molar mass of phase.
      Returns:
      molar mass in unit kg/mol
    • getMolarMass

      double getMolarMass(String unit)
      Get molar mass of a fluid phase.
      Parameters:
      unit - Supported units are kg/mol, gr/mol
      Returns:
      molar mass in specified unit
    • getInternalEnergy

      double getInternalEnergy()

      getInternalEnergy.

      Returns:
      a double
    • getInternalEnergy

      double getInternalEnergy(String unit)

      getInternalEnergy.

      Parameters:
      unit - a String object
      Returns:
      a double
    • getdPdrho

      double getdPdrho()

      getdPdrho.

      Returns:
      a double
    • getdPdTVn

      double getdPdTVn()

      getdPdTVn.

      Returns:
      a double
    • getdPdVTn

      double getdPdVTn()

      getdPdVTn.

      Returns:
      a double
    • Fn

      double Fn()

      Fn.

      Returns:
      a double
    • FT

      double FT()

      FT.

      Returns:
      a double
    • FV

      double FV()

      FV.

      Returns:
      a double
    • FD

      double FD()

      FD.

      Returns:
      a double
    • FB

      double FB()

      FB.

      Returns:
      a double
    • gb

      double gb()

      gb.

      Returns:
      a double
    • fb

      double fb()

      fb.

      Returns:
      a double
    • gV

      double gV()

      gV.

      Returns:
      a double
    • fv

      double fv()

      fv.

      Returns:
      a double
    • FnV

      double FnV()

      FnV.

      Returns:
      a double
    • FnB

      double FnB()

      FnB.

      Returns:
      a double
    • FTT

      double FTT()

      FTT.

      Returns:
      a double
    • FBT

      double FBT()

      FBT.

      Returns:
      a double
    • FDT

      double FDT()

      FDT.

      Returns:
      a double
    • FBV

      double FBV()

      FBV.

      Returns:
      a double
    • FBB

      double FBB()

      FBB.

      Returns:
      a double
    • FDV

      double FDV()

      FDV.

      Returns:
      a double
    • FBD

      double FBD()

      FBD.

      Returns:
      a double
    • FTV

      double FTV()

      FTV.

      Returns:
      a double
    • FVV

      double FVV()

      FVV.

      Returns:
      a double
    • gVV

      double gVV()

      gVV.

      Returns:
      a double
    • gBV

      double gBV()

      gBV.

      Returns:
      a double
    • gBB

      double gBB()

      gBB.

      Returns:
      a double
    • fVV

      double fVV()

      fVV.

      Returns:
      a double
    • fBV

      double fBV()

      fBV.

      Returns:
      a double
    • fBB

      double fBB()

      fBB.

      Returns:
      a double
    • dFdT

      double dFdT()

      Calculate derivative of F per Temperature, i.e., dF/dT.

      Returns:
      a double
    • dFdV

      double dFdV()

      Calculate derivative of F per Volume, i.e., dF/dV.

      Returns:
      a double
    • dFdTdV

      double dFdTdV()

      Calculate derivative of F per Temperature and Volume, i.e., dF/dT * 1/dV.

      Returns:
      a double
    • dFdVdV

      double dFdVdV()

      dFdVdV.

      Returns:
      a double
    • dFdTdT

      double dFdTdT()

      dFdTdT.

      Returns:
      a double
    • getOsmoticCoefficientOfWater

      double getOsmoticCoefficientOfWater()

      getOsmoticCoefficientOfWater.

      Returns:
      a double
    • getOsmoticCoefficient

      double getOsmoticCoefficient(int watNumb)

      getOsmoticCoefficient.

      Parameters:
      watNumb - a int
      Returns:
      a double
    • getOsmoticCoefficientOfWaterMolality

      double getOsmoticCoefficientOfWaterMolality()
      Get the osmotic coefficient of water on the molality scale. This is the definition used by Robinson and Stokes (1965):
      φ = -ln(a_w) / (M_w * Σm_i)
      
      where:
      • a_w = water activity
      • M_w = molar mass of water (kg/mol)
      • Σm_i = sum of ion molalities (mol/kg solvent)
      Returns:
      osmotic coefficient on molality scale
    • getMeanIonicActivity

      double getMeanIonicActivity(int comp1, int comp2)

      getMeanIonicActivity.

      Parameters:
      comp1 - a int
      comp2 - a int
      Returns:
      a double
    • getLogInfiniteDiluteFugacity

      double getLogInfiniteDiluteFugacity(int k, int p)

      getLogInfiniteDiluteFugacity.

      Parameters:
      k - a int
      p - a int
      Returns:
      a double
    • getLogInfiniteDiluteFugacity

      double getLogInfiniteDiluteFugacity(int k)

      getLogInfiniteDiluteFugacity.

      Parameters:
      k - a int
      Returns:
      a double
    • getMixingRule

      MixingRulesInterface getMixingRule()
      Get mixing rule.
      Returns:
      a MixingRulesInterface
    • getMixingRuleType

      MixingRuleTypeInterface getMixingRuleType()
      Get mixing rule type.
      Returns:
      a MixingRuleTypeInterface
    • initRefPhases

      void initRefPhases(boolean onlyPure)

      initRefPhases.

      Parameters:
      onlyPure - a boolean
    • getRefPhase

      PhaseInterface getRefPhase(int index)

      Indexed getter for property refPhase.

      Parameters:
      index - a int
      Returns:
      a PhaseInterface object
    • getRefPhase

      PhaseInterface[] getRefPhase()

      Getter for property refPhase.

      Returns:
      an array of PhaseInterface objects
    • setRefPhase

      void setRefPhase(int index, PhaseInterface refPhase)

      Indexed setter for property refPhase.

      Parameters:
      index - a int
      refPhase - a PhaseInterface object
    • setRefPhase

      void setRefPhase(PhaseInterface[] refPhase)

      Setter for property refPhase.

      Parameters:
      refPhase - an array of PhaseInterface objects
    • setParams

      void setParams(PhaseInterface phase, double[][] alpha, double[][] Dij, double[][] DijT, String[][] mixRule, double[][] intparam)

      setParams.

      Parameters:
      phase - a PhaseInterface object
      alpha - an array of type double
      Dij - an array of type double
      DijT - an array of type double
      mixRule - an array of String objects
      intparam - an array of type double
    • getType

      PhaseType getType()
      Getter for property pt.
      Returns:
      PhaseType enum object.
    • setType

      void setType(PhaseType pt)
      Setter for property pt.
      Parameters:
      pt - PhaseType to set.
    • getPhaseTypeName

      default String getPhaseTypeName()

      Getter for property phaseTypeName.

      Returns:
      a String object
    • setPhaseTypeName

      default void setPhaseTypeName(String phaseTypeName)

      Setter for property phaseTypeName.

      Parameters:
      phaseTypeName - a String object
    • isMixingRuleDefined

      boolean isMixingRuleDefined()

      Check if mixing rule is defined.

      Returns:
      Returns true if MixingRule is defined and false if not.
    • getActivityCoefficientSymetric

      double getActivityCoefficientSymetric(int k)

      getActivityCoefficientSymetric.

      Parameters:
      k - a int
      Returns:
      a double
    • getActivityCoefficientUnSymetric

      double getActivityCoefficientUnSymetric(int k)

      getActivityCoefficientUnSymetric.

      Parameters:
      k - a int
      Returns:
      a double
    • hasComponent

      default boolean hasComponent(String name)
      Verify if phase has a component.
      Parameters:
      name - Name of component to look for. NB! Converts name to normalized name.
      Returns:
      True if component is found.
    • hasComponent

      boolean hasComponent(String name, boolean normalized)
      Verify if phase has a component.
      Parameters:
      name - Name of component to look for.
      normalized - Set true to convert input name to normalized component name.
      Returns:
      True if component is found.
    • isAsphalteneRich

      default boolean isAsphalteneRich()
      Check if this phase is rich in asphaltene components.

      This method detects asphaltene-rich phases regardless of whether the phase is modeled as solid (PhaseType.ASPHALTENE) or liquid (PhaseType.LIQUID_ASPHALTENE, Pedersen's liquid-liquid approach). A phase is considered asphaltene-rich if:

      • The phase type is ASPHALTENE or LIQUID_ASPHALTENE, or
      • The total mole fraction of asphaltene components exceeds 0.5

      Asphaltene components are identified by name containing "asphaltene" (case-insensitive).

      Returns:
      true if the phase is asphaltene-rich
    • getLogActivityCoefficient

      double getLogActivityCoefficient(int k, int p)

      getLogActivityCoefficient.

      Parameters:
      k - a int
      p - a int
      Returns:
      a double
    • isConstantPhaseVolume

      boolean isConstantPhaseVolume()

      isConstantPhaseVolume.

      Returns:
      a boolean
    • setConstantPhaseVolume

      void setConstantPhaseVolume(boolean constantPhaseVolume)

      setConstantPhaseVolume.

      Parameters:
      constantPhaseVolume - a boolean
    • getSoundSpeed

      double getSoundSpeed()
      Get the speed of sound of a phase note: implemented in phaseEos.
      Returns:
      speed of sound in m/s
    • getSoundSpeed

      double getSoundSpeed(String unit)
      Get the speed of sound of a system. The sound speed is implemented based on a molar average over the phases
      Parameters:
      unit - Supported units are m/s, km/h
      Returns:
      speed of sound in m/s
    • getModelName

      String getModelName()
      method to return name of thermodynamic model.
      Returns:
      String model name
    • getZvolcorr

      double getZvolcorr()
      method to return Z volume corrected gas compressibility.
      Returns:
      double Z volume corrected