TRIOMA API Documentation
In this section the documentation for all modules in TRIOMA:
Correlations
- src.TRIOMA.tools.correlations.Nu_DittusBoelter(Re, Pr)
Calculate the Nusselt number using the Dittus-Boelter correlation. Parameters: - Re (float): Reynolds number - Pr (float): Prandtl number Returns: - Nu (float): Nusselt number
- src.TRIOMA.tools.correlations.Nu_Gnielinsky(Re, Pr, f)
Calculate the Nusselt number using the Gnielinsky correlation. Parameters: - Re (float): Reynolds number - Pr (float): Prandtl number - f (float): friction factor Returns: - Nu (float): Nusselt number
- src.TRIOMA.tools.correlations.Nu_SiederTate(Re, Pr, mu_w, mu_c)
Calculate the Nusselt number using the Sieder-Tate correlation.
Parameters: - Re (float): Reynolds number - Pr (float): Prandtl number - mu_w (float): viscosity of the wall fluid - mu_c (float): viscosity of the core fluid
Returns: - Nu (float): Nusselt number
- src.TRIOMA.tools.correlations.Pr(c_p, mu, k)
Calculate the Prandtl number. Parameters: - c_p (float): specific heat capacity - mu (float): viscosity - k (float): thermal conductivity Returns: - Pr (float): Prandtl number
- src.TRIOMA.tools.correlations.Re(rho, u, L, mu)
Calculate the Reynolds number. Parameters: - rho (float): density - u (float): velocity - L (float): length - mu (float): viscosity Returns: - Re (float): Reynolds number
- src.TRIOMA.tools.correlations.Schmidt(D, mu, rho)
Calculate the Schmidt number. Parameters: - D (float): diffusivity - mu (float): viscosity - rho (float): density Returns: - Sc (float): Schmidt number
- src.TRIOMA.tools.correlations.Sherwood(Sc, Re)
Calculate the Sherwood number. Parameters: - Sc (float): Schmidt number - Re (float): Reynolds number Returns: - Sh (float): Sherwood number
- src.TRIOMA.tools.correlations.Sherwood_HT_analogy(Re, Sc)
Calculate the Sherwood number using the heat transfer analogy. Parameters: - Re (float): Reynolds number - Sc (float): Schmidt number Returns: - Sh (float): Sherwood number
- src.TRIOMA.tools.correlations.Sherwood_bubbles(Sc, Re)
Calculate the Sherwood number for bubbles. Parameters: - Sc (float): Schmidt number - Re (float): Reynolds number Returns: - Sh (float): Sherwood number
- src.TRIOMA.tools.correlations.f_Haaland(Re, e_D)
Calculate the friction factor using the Haaland correlation. Parameters: - Re (float): Reynolds number - e_D (float): relative roughness Returns: - f (float): friction factor
- src.TRIOMA.tools.correlations.f_Pethukov(Re, Pr)
Calculate the friction factor using the Pethukov correlation. Parameters: - Re (float): Reynolds number - Pr (float): Prandtl number Returns: - f (float): friction factor
- src.TRIOMA.tools.correlations.get_deltaTML(T_in_hot, T_out_hot, T_in_cold, T_out_cold)
Calculate the log mean temperature difference. Parameters: - T_in_hot (float): hot inlet temperature - T_out_hot (float): hot outlet temperature - T_in_cold (float): cold inlet temperature - T_out_cold (float): cold outlet temperature Returns: - deltaTML (float): log mean temperature difference
- src.TRIOMA.tools.correlations.get_h_from_Nu(Nu, k, D)
Calculate the heat transfer coefficient from the Nusselt number. Parameters: - Nu (float): Nusselt number - k (float): thermal conductivity - D (float): diameter Returns: - h (float): heat transfer coefficient
- src.TRIOMA.tools.correlations.get_k_from_Sh(Sh, L, D)
Calculate the mass transfer coefficient from the Sherwood number. Parameters: - Sh (float): Sherwood number - L (float): length - D (float): diameter Returns: - k (float): mass transfer coefficient
- src.TRIOMA.tools.correlations.get_length_HX(deltaTML, d_hyd, U, Q)
Calculate the length of the heat exchanger. Parameters: - deltaTML (float): log mean temperature difference - d_hyd (float): hydraulic diameter - U (float): overall heat transfer coefficient - Q (float): heat transfer rate Returns: - L (float): length of the heat exchanger
Component Class
- class src.TRIOMA.tools.component_tools.Component(geometry: Geometry = None, c_in: float = None, eff: float = None, fluid: Fluid = None, membrane: Membrane = None, name: str = None, p_out: float = 1e-15, loss: bool = False, inv: float = None, delta_p: float = None, pumping_power: float = None, V: float = None)
Represents a component in a plant to make a high level T transport analysis.
- Args:
Geometry (Geometry): The geometry of the component. c_in (float): The concentration of the component at the inlet. fluid (Fluid): The fluid associated with the component. Defaults to None. membrane (Membrane): The membrane associated with the component. Defaults to None.
- T_leak() float
Calculates the leakage of the component.
- Returns:
float: The leakage of the component.
- analytical_efficiency(p_out=1e-15)
Calculate the analytical efficiency of a component.
Parameters: - p_out: Pressure of H isotope at the outlet of the component. Defaults to 1E-15.
Returns: - eff_an: Analytical efficiency of the component (from Humrickhouse papers)
This function calculates the analytical efficiency of a component based on the given length (L) of the component. It uses various properties of the component, such as membrane properties, fluid properties, and adimensionals.
If the fluid is a molten salt (MS=True), the analytical efficiency is calculated using the following formula: eff_an = 1 - xi * (lambertw(z=np.exp(beta - tau - 1), tol=1e-10) ** 2 + 2 * lambertw(z=np.exp(beta - tau - 1), tol=1e-10)).
If the fluid is a liquid metal (MS= False), the analytical efficiency is calculated using the following formula: eff_an = 1 - np.exp(-tau * zeta / (1 + zeta))* (1-p_in/p_out)^0.5.
The output of the function is the analytical efficiency of the component as Component.eff_an.
- connect_to_component(component2: Component | BreedingBlanket = None)
sets the inlet conc of the object component equal to the outlet of self
- converge_split_HX(tol: float = 0.001, T_in_hot: float = None, T_out_hot: float = None, T_in_cold: float = None, T_out_cold: float = None, R_sec: float = None, Q: float = None, plotvar: bool = False) Circuit
Splits the component into N components to better discretize Temperature effects Tries to find the optimal number of components to split the component into
- define_component_volumes()
Calculates the volumes of the component.
- friction_factor(Re: float) float
Calculates the friction factor for the component.
- Args:
Re (float): Reynolds number.
- Returns:
float: The friction factor.
- get_adimensionals()
Calculates the adimensional parameters H and W.
Updates the H and W attributes of the Component object.
- get_efficiency(plotvar: bool = False, c_guess: float = None, p_out=1e-15)
Calculates the efficiency of the component.
- get_flux(c: float = None, c_guess: float = 1e-09, p_out=1e-15)
Calculates the Tritium flux of the component. It can make some approximations based on W and H to make the solver faster
- Args:
c (float): The concentration in the component fluid.
- Returns:
float: The permeation flux.
- get_global_HX_coeff(R_conv_sec: float = 0)
Calculates the global heat exchange coefficient of the component. It can take the secondary resistance to convection as input. defaults to no resistance
- Returns:
float: The global heat exchange coefficient of the component.
- get_pipe_flowrate()
Calculates the volumetric flow rate of the component [m^3/s].
- Returns:
float: The flow rate of the component.
- get_pressure_drop() float
Calculates the pressure drop across the component.
- Returns:
float: The pressure drop across the component.
- get_pumping_power() float
Calculates the pumping power required for the component.
- Returns:
float: The pumping power required for the component in W.
- get_regime(print_var: bool = False)
Gets the regime of the component.
- Returns:
str: The regime of the component.
- get_total_flowrate()
Calculates the total flow rate of the component.
- outlet_c_comp() float
Calculates the concentration of the component at the outlet.
- Returns:
float: The concentration of the component at the outlet.
- split_HX(N: int = 25, T_in_hot: int = None, T_out_hot: int = None, T_in_cold: int = None, T_out_cold: int = None, R_sec: int = None, Q: int = None, plotvar: bool = False) Circuit
Splits the component into N components to better discretize Temperature effects
- use_analytical_efficiency(p_out=1e-15)
Evaluates the analytical efficiency and substitutes it in the efficiency attribute of the component.
- Args:
L (float): the length of the pipe component
- Returns:
None
Circuit Class
- class src.TRIOMA.tools.component_tools.Circuit(components: list = None, closed: bool = False)
Represent a circuit of components connected in series
This class represents a circuit consisting of multiple components connected in series. It provides methods to update attributes, add components, calculate circuit efficiency, and plot the circuit.
- Attributes:
components (list): A list of components in the circuit. closed (bool): A boolean indicating whether the circuit is a closed loop or not.
- Methods:
- update_attribute(attr_name, new_value):
Updates the value of the specified attribute.
- add_component(component):
Adds a component to the circuit.
- get_eff_circuit():
Calculates the efficiency of the circuit.
- get_gains_and_losses():
Calculates the gains and losses of the circuit.
- plot_circuit():
Plots the circuit using matplotlib.
- Example usage:
circuit = Circuit() circuit.add_component(component1) circuit.add_component(component2) circuit.get_eff_circuit() circuit.plot_circuit()
- add_component(component: Component | BreedingBlanket | Circuit | GLC)
Adds a component to the circuit.
- Args:
component (Component): The component to add.
- get_circuit_pumping_power()
Calculates the pumping power required for the circuit.
- Returns:
float: The pumping power required for the circuit in W.
- get_eff_circuit()
Calculates the efficiency of the circuit based on the components present.
- Returns:
circtuit.eff (float): The efficiency of the circuit.
- Raises:
None
- Example Usage:
circuit.get_eff_circuit()
- get_gains_and_losses()
Calculates the gains and losses of the components in the circuit.
- Returns:
circuit.extraction_perc (float): The extraction percentage of the circuit. circuit.loss_perc (float): The loss percentage of the circuit.
- get_inventory(flag_an=True)
Calculates the inventory (in mol) of the circuit based on the components present.
- Returns:
circuit.inv (float): The inventory of the circuit.
- Raises:
None
- Example Usage:
circuit.get_circuit_inventory()
- inspect_circuit(name=None)
Inspects the circuit components.
- Parameters:
name (str, optional): The name of the component to inspect. If not provided, all components will be inspected.
- Returns:
None
- plot_circuit()
Plot the circuit diagram for the components in the circuit.
This function uses matplotlib to create a circuit diagram for the components in the circuit. Each component is represented by a rectangle with arrows indicating the flow direction. The color of the rectangle represents the position of the component in the circuit.
- Returns:
None
- solve_circuit(tol=1e-06)
Solve the circuit by calculating the concentration of the components at the outlet. If the circuit is a closed loop, the concentration of the first component is set to the concentration of the last component until the stationary regime is reached.
Fluid Class
- class src.TRIOMA.tools.component_tools.Fluid(T: float = None, D: float = None, D_0: float = None, E_d: float = None, Solubility: float = None, Solubility_0: float = None, E_s: float = None, MS: bool = True, d_Hyd: float = None, k_t: float = None, mu: float = None, rho: float = None, U0: float = None, k: float = None, cp: float = None, inv: float = None, recirculation: float = 0, V: float = None)
Represents a fluid in a component for Tritium transport analysis
- Args:
T (float): Temperature of the fluid. D (float): Tritium Diffusivity of the fluid. D_0 (float): Preexponential Diffusivity of the fluid. E_d (float): Activation energy of the fluid diffusivity. Solubility (float): Solubility of the fluid. Solubility_0 (float): Preexponential solubility of the fluid. E_s (float): Activation energy of the fluid solubility. MS (bool): Indicates whether the fluid is a molten salt or a liquid metal. d_Hyd (float, optional): Hydraulic diameter of the fluid. Defaults to None. k_t (float, optional): Mass transport coefficient of the fluid. Defaults to None. mu (float, optional): Viscosity of the fluid. Defaults to None. rho (float, optional): Density of the fluid. Defaults to None. U0 (float, optional): Velocity of the fluid. Defaults to None. inv (float, optional): Inventory of the fluid. Defaults to None. recirculation(float,optional): fraction of recirculated flowrate. Defaults to 0. 1 = 100%, 0.5=50%.
- get_kt(turbulator=None)
Calculates the mass transport coefficient (k_t) for the fluid.
If the hydraulic diameter (d_Hyd) is defined, the mass transport coefficient is calculated using correlations. Otherwise, an error message is printed.
- Returns:
None
- set_properties_from_fluid_material(fluid_material: FluidMaterial = None)
Sets the properties of the fluid from a FluidMaterial object.
- Args:
fluid_material (FluidMaterial): The FluidMaterial object to set the properties from.
Membrane Class
- class src.TRIOMA.tools.component_tools.Membrane(T: float = None, D: float = None, thick: float = None, K_S: float = None, k_d: float = None, k_r: float = None, k: float = None, D_0: float = None, E_d: float = None, K_S_0: float = None, E_S: float = None, inv: float = None, V: float = None)
Represents a metallic membrane of a component for H transport.
- Attributes:
T (float): Temperature of the membrane. D (float): Diffusion coefficient of the membrane. thick (float): Thickness of the membrane. K_S (float): Solubility coefficient of the membrane. k_d (float, optional): Dissociation rate constant of the membrane. Defaults to None. k_r (float, optional): Recombination rate constant of the membrane. Defaults to None. k (float, optional): Thermal conductivity of the membrane. Defaults to None. D_0 (float, optional): Pre-exponential factor of the membrane. Defaults to None.Overwrites D if defined E_d (float, optional): Activation energy of the diffusivity in the membrane in eV. Defaults to None. Overwrites D if defined K_S_0 (float, optional): Pre-exponential factor of the solubility in the membrane. Defaults to None.Overwrites K_S if defined E_S (float, optional): Activation energy of the solubility in the membrane in eV. Defaults to None. Overwrites K_S if defined inv (float, optional): Inventory of the membrane in mol. Defaults to None.
- set_properties_from_solid_material(solid_material: SolidMaterial = None)
Sets the properties of the membrane from a SolidMaterial object.
- Args:
solid_material (SolidMaterial): The SolidMaterial object to set the properties from.
Solid Material Class
- class src.TRIOMA.tools.component_tools.SolidMaterial(T: float = None, D: float = None, K_S: float = None, k: float = None)
Represents a solid material used in a component.
- Attributes:
D (float): The Diffusivity of the solid material. K_S (float): The Sievert constant of the solid material.
Fluid Material Class
- class src.TRIOMA.tools.component_tools.FluidMaterial(T: float = None, D: float = None, Solubility: float = None, MS: bool = None, mu: float = None, rho: float = None, k: float = None, cp: float = None)
Represents a fluid material with various properties.
- Attributes:
T (float): Temperature of the fluid material. D (float): Density of the fluid material. Solubility (float): Solubility of the fluid material. MS (float): Molecular weight of the fluid material. mu (float): Viscosity of the fluid material. rho (float): Density of the fluid material. k (float): Thermal conductivity of the fluid material. cp (float): Specific heat capacity of the fluid material.
Breeding Blanket Class
- class src.TRIOMA.tools.component_tools.BreedingBlanket(c_in: float = None, Q: float = None, TBR: float = None, T_out: float = None, T_in: float = None, fluid: Fluid = None, name: str = None)
Represents a breeding blanket component in a fuel cycle system.
- Attributes:
c_in(float): Inlet concentration of tritium in the breeding blanket component. Q (float): Heat generated by the breeding blanket component. TBR (float): Tritium breeding ratio of the breeding blanket component. T_out (float): Outlet temperature of the breeding blanket component. T_in (float): Inlet temperature of the breeding blanket component. fluid (Fluid): Fluid used in the breeding blanket component.
- get_cout(print_var: bool = False)
Calculates the outlet concentration of tritium in the breeding blanket component.
- Args:
print_var (bool): If True, prints the intermediate variables.
- Returns:
None
- get_flowrate()
Calculates the flow rate of the coolant in the breeding blanket component.