API reference
This page provides a comprehensive API reference for the LineCableModels.jl package. It documents all public modules, types, functions, and constants, organized by functional area. Each section corresponds to a major module in the package, with detailed information about parameters, return values, and usage examples.
Contents
Data model
LineCableModels.DataModel — Module
LineCableModels.DataModelThe DataModel module provides data structures, constructors and utilities for modeling power cables within the LineCableModels.jl package. This module includes definitions for various cable components, and visualization tools for cable designs.
Overview
- Provides objects for detailed cable modeling with the
CableDesignand supporting types:CircStrands,Strip,Tubular,Semicon, andInsulator. - Includes objects for cable system modeling with the
LineCableSystemtype, and multiple formation patterns like trifoil and flat arrangements. - Contains functions for calculating the base electric properties of all elements within a
CableDesign, namely: resistance, inductance (via GMR), shunt capacitance, and shunt conductance (via loss factor). - Offers visualization tools for previewing cable cross-sections and system layouts.
- Provides a library system for storing and retrieving cable designs.
Dependencies
BaseColorsDataFramesDatesLineCableModels.CommonsLineCableModels.DataModel.BaseParamsLinearAlgebraMakieMeasurementsPlotsPrintfStatistics
Exports
LineCableModels.DataModel.CableComponent — Type
mutable struct CableComponent{T<:Union{Float64, Measurements.Measurement{Float64}}}Represents a CableComponent, i.e. a group of AbstractCablePart objects, with the equivalent geometric and material properties:
id::String: Cable component identification (e.g. core/sheath/armor).conductor_group::ConductorGroup: The conductor group containing all conductive parts.conductor_props::Material: Effective properties of the equivalent coaxial conductor.insulator_group::InsulatorGroup: The insulator group containing all insulating parts.insulator_props::Material: Effective properties of the equivalent coaxial insulator.
Cable components operate as containers for multiple cable parts, allowing the calculation of effective electromagnetic (EM) properties ($\sigma, \varepsilon, \mu$). This is performed by transforming the physical objects within the CableComponent into one equivalent coaxial homogeneous structure comprised of one conductor and one insulator, each one represented by effective Material types stored in conductor_props and insulator_props fields.
The effective properties approach is widely adopted in EMT-type simulations, and involves locking the internal and external radii of the conductor and insulator parts, respectively, and calculating the equivalent EM properties in order to match the previously determined values of R, L, C and G [4] [14].
In applications, the CableComponent type is mapped to the main cable structures described in manufacturer datasheets, e.g., core, sheath, armor and jacket.
LineCableModels.DataModel.CableComponent — Method
CableComponent(
id::String,
conductor_group::ConductorGroup,
insulator_group::InsulatorGroup
) -> CableComponent
Weakly-typed constructor that infers the scalar type T from the two groups, coerces them if necessary, and calls the strict kernel.
Arguments
id: Cable component identification.conductor_group: The conductor group (anyConductorGroup{S}).insulator_group: The insulator group (anyInsulatorGroup{R}).
Returns
- A
CableComponent{T}whereTis the resolved scalar type.
LineCableModels.DataModel.CableComponent — Method
Initializes a CableComponent object based on its constituent conductor and insulator groups. The constructor performs the following sequence of steps:
- Validate that the conductor and insulator groups have matching radii at their interface.
- Obtain the lumped-parameter values (R, L, C, G) from the conductor and insulator groups, which are computed within their respective constructors.
- Calculate the correction factors and equivalent electromagnetic properties of the conductor and insulator groups:
| Quantity | Symbol | Function |
|---|---|---|
| Resistivity (conductor) | $\rho_{con}$ | calc_equivalent_rho |
| Permeability (conductor) | $\mu_{con}$ | calc_equivalent_mu |
| Resistivity (insulator) | $\rho_{ins}$ | calc_sigma_lossfact |
| Permittivity (insulation) | $\varepsilon_{ins}$ | calc_equivalent_eps |
| Permeability (insulation) | $\mu_{ins}$ | calc_solenoid_correction |
Arguments
id: Cable component identification (e.g. core/sheath/armor).conductor_group: The conductor group containing all conductive parts.insulator_group: The insulator group containing all insulating parts.
Returns
A CableComponent instance with calculated equivalent properties:
id::String: Cable component identification.conductor_group::ConductorGroup{T}: The conductor group containing all conductive parts.conductor_props::Material{T}: Effective properties of the equivalent coaxial conductor.
* `rho`: Resistivity \[Ω·m\].
* `eps_r`: Relative permittivity \[dimensionless\].
* `mu_r`: Relative permeability \[dimensionless\].
* `T0`: Reference temperature \[°C\].
* `alpha`: Temperature coefficient of resistivity \[1/°C\].insulator_group::InsulatorGroup{T}: The insulator group containing all insulating parts.insulator_props::Material{T}: Effective properties of the equivalent coaxial insulator.
* `rho`: Resistivity \[Ω·m\].
* `eps_r`: Relative permittivity \[dimensionless\].
* `mu_r`: Relative permeability \[dimensionless\].
* `T0`: Reference temperature \[°C\].
* `alpha`: Temperature coefficient of resistivity \[1/°C\].Examples
conductor_group = ConductorGroup(...)
insulator_group = InsulatorGroup(...)
cable = CableComponent("component_id", conductor_group, insulator_group) # Create cable component with base parameters @ 50 HzSee also
LineCableModels.DataModel.CableDesign — Type
mutable struct CableDesign{T<:Union{Float64, Measurements.Measurement{Float64}}}Represents the design of a cable, including its unique identifier, nominal data, and components.
cable_id::String: Unique identifier for the cable design.nominal_data::Union{Nothing, NominalData{T}} where T<:Union{Float64, Measurements.Measurement{Float64}}: Informative reference data.components::Array{CableComponent{T}, 1} where T<:Union{Float64, Measurements.Measurement{Float64}}: Vector of cable components.
LineCableModels.DataModel.CableDesign — Method
CableDesign(
cable_id::String,
component::CableComponent;
nominal_data
) -> CableDesign
Weakly-typed constructor that infers the scalar type from the component (and nominal data if present), coerces values to that type, and calls the typed kernel.
LineCableModels.DataModel.CableDesign — Method
CableDesign(
cable_id::String,
conductor_group::ConductorGroup,
insulator_group::InsulatorGroup;
component_id,
nominal_data
) -> CableDesign
Constructs a CableDesign instance from conductor and insulator groups. Convenience wrapper that builds the component with reduced boilerplate.
LineCableModels.DataModel.CableDesign — Method
CableDesign(
cable_id::String,
component::CableComponent{T<:Union{Float64, Measurements.Measurement{Float64}}};
nominal_data
) -> CableDesign
Strict numeric kernel: constructs a CableDesign{T} from one component (typed) and optional nominal data (typed or nothing). Assumes all inputs are already at scalar type T.
Arguments
cable_id: Unique identifier for the cable design.component: InitialCableComponentfor the design.nominal_data: Reference data for the cable design. Default:NominalData().
Returns
- A
CableDesignobject with the specified properties.
Examples
conductor_group = ConductorGroup(central_conductor)
insulator_group = InsulatorGroup(main_insulator)
component = CableComponent(conductor_group, insulator_group)
design = CableDesign("example", component)See also
LineCableModels.DataModel.CablePosition — Type
CablePosition(
cable::Union{Nothing, CableDesign},
horz::Number,
vert::Number
) -> CablePosition
CablePosition(
cable::Union{Nothing, CableDesign},
horz::Number,
vert::Number,
conn::Union{Nothing, Dict{String, Int64}}
) -> CablePosition
Weakly-typed constructor that infers T from the cable and coordinates, builds/validates the phase mapping, coerces inputs to T, and calls the typed kernel.
LineCableModels.DataModel.CablePosition — Type
struct CablePosition{T<:Union{Float64, Measurements.Measurement{Float64}}}Represents a physically defined cable with position and phase mapping within a system.
design_data::CableDesign: TheCableDesignobject assigned to this cable position.horz::Union{Float64, Measurements.Measurement{Float64}}: Horizontal coordinate [m].vert::Union{Float64, Measurements.Measurement{Float64}}: Vertical coordinate [m].conn::Vector{Int64}: Phase mapping vector (aligned with design_data.components).
LineCableModels.DataModel.CablePosition — Method
Constructs a CablePosition instance with specified cable design, coordinates, and phase mapping.
Arguments
cable: ACableDesignobject defining the cable structure.horz: Horizontal coordinate [m].vert: Vertical coordinate [m].conn: A dictionary mapping component names to phase indices, ornothingfor default mapping.
Returns
- A
CablePositionobject with the assigned cable design, coordinates, and phase mapping.
The conn argument is a Dict that maps the cable components to their respective phases. The values (1, 2, 3) represent the phase numbers (A, B, C) in a three-phase system. Components mapped to phase 0 will be Kron-eliminated (grounded). Components set to the same phase will be bundled into an equivalent phase.
Examples
cable_design = CableDesign("example", nominal_data, components_dict)
xa, ya = 0.0, -1.0 # Coordinates in meters
# With explicit phase mapping
cablepos1 = CablePosition(cable_design, xa, ya, Dict("core" => 1))
# With default phase mapping (first component to phase 1, others to 0)
default_cablepos = CablePosition(cable_design, xa, ya)See also
LineCableModels.DataModel.CablesLibrary — Type
mutable struct CablesLibraryRepresents a library of cable designs stored as a dictionary.
data::Dict{String, CableDesign}: Dictionary mapping cable IDs to the respective CableDesign objects.
LineCableModels.DataModel.CablesLibrary — Method
CablesLibrary() -> CablesLibrary
Constructs an empty CablesLibrary instance.
Arguments
- None.
Returns
- A
CablesLibraryobject with an empty dictionary of cable designs.
Examples
# Create a new, empty library
library = CablesLibrary()See also
LineCableModels.DataModel.CircStrands — Type
struct CircStrands{T<:Union{Float64, Measurements.Measurement{Float64}}, U<:Int64} <: LineCableModels.DataModel.AbstractStrandsLayer{T<:Union{Float64, Measurements.Measurement{Float64}}}Represents an array of wires equally spaced around a circumference of arbitrary radius, with attributes:
r_in::Union{Float64, Measurements.Measurement{Float64}}: Internal radius of the wire array [m].r_ex::Union{Float64, Measurements.Measurement{Float64}}: External radius of the wire array [m].radius_wire::Union{Float64, Measurements.Measurement{Float64}}: Radius of each individual wire [m].num_wires::Int64: Number of wires in the array [dimensionless].lay_ratio::Union{Float64, Measurements.Measurement{Float64}}: Ratio defining the lay length of the wires (twisting factor) [dimensionless].mean_diameter::Union{Float64, Measurements.Measurement{Float64}}: Mean diameter of the wire array [m].pitch_length::Union{Float64, Measurements.Measurement{Float64}}: Pitch length of the wire array [m].lay_direction::Int64: Twisting direction of the strands (1 = unilay, -1 = contralay) [dimensionless].material_props::Material: Material object representing the physical properties of the wire material.temperature::Union{Float64, Measurements.Measurement{Float64}}: Temperature at which the properties are evaluated [°C].cross_section::Union{Float64, Measurements.Measurement{Float64}}: Cross-sectional area of all wires in the array [m²].resistance::Union{Float64, Measurements.Measurement{Float64}}: Electrical resistance per wire in the array [Ω/m].gmr::Union{Float64, Measurements.Measurement{Float64}}: Geometric mean radius of the wire array [m].
LineCableModels.DataModel.CircStrands — Method
CircStrands(
r_in::Union{Float64, Measurements.Measurement{Float64}},
radius_wire::Union{Float64, Measurements.Measurement{Float64}},
num_wires::Int64,
lay_ratio::Union{Float64, Measurements.Measurement{Float64}},
material_props::Material{T<:Union{Float64, Measurements.Measurement{Float64}}},
temperature::Union{Float64, Measurements.Measurement{Float64}},
lay_direction::Int64
) -> LineCableModels.DataModel.CircStrands{T, Int64} where T<:Union{Float64, Measurements.Measurement{Float64}}
Constructs a CircStrands instance based on specified geometric and material parameters.
Arguments
r_in: Internal radius of the wire array [m].radius_wire: Radius of each individual wire [m].num_wires: Number of wires in the array [dimensionless].lay_ratio: Ratio defining the lay length of the wires (twisting factor) [dimensionless].material_props: AMaterialobject representing the material properties.temperature: Temperature at which the properties are evaluated [°C].lay_direction: Twisting direction of the strands (1 = unilay, -1 = contralay) [dimensionless].
Returns
- A
CircStrandsobject with calculated geometric and electrical properties.
Examples
material_props = Material(1.7241e-8, 1.0, 0.999994, 20.0, 0.00393)
circstrands = CircStrands(0.01, Diameter(0.002), 7, 10, material_props, temperature=25)
println(circstrands.mean_diameter) # Outputs mean diameter in m
println(circstrands.resistance) # Outputs resistance in Ω/mSee also
LineCableModels.DataModel.ConductorGroup — Method
ConductorGroup(
component::CableComponent{T<:Union{Float64, Measurements.Measurement{Float64}}}
) -> ConductorGroup
Build a ConductorGroup equivalent for a CableComponent, preserving num_turns and num_wires from the original group.
Constructs a single-layer ConductorGroup{T} from the computed equivalent Tubular(component), but carries over bookkeeping fields needed by downstream corrections (e.g., solenoid correction using num_turns).
LineCableModels.DataModel.ConductorGroup — Method
Constructs a ConductorGroup instance initializing with the central conductor part.
Arguments
central_conductor: AnAbstractConductorPartobject located at the center of the conductor group.
Returns
- A
ConductorGroupobject initialized with geometric and electrical properties derived from the central conductor.
LineCableModels.DataModel.Diameter — Type
struct Diameter{T<:Real} <: LineCableModels.DataModel.AbstractRadiusRepresents the diameter of a cable component.
value::Real: Numerical value of the diameter [m].
LineCableModels.DataModel.Insulator — Type
struct Insulator{T<:Union{Float64, Measurements.Measurement{Float64}}} <: LineCableModels.DataModel.AbstractInsulatorPart{T<:Union{Float64, Measurements.Measurement{Float64}}}Represents an insulating layer with defined geometric, material, and electrical properties given by the attributes:
r_in::Union{Float64, Measurements.Measurement{Float64}}: Internal radius of the insulating layer [m].r_ex::Union{Float64, Measurements.Measurement{Float64}}: External radius of the insulating layer [m].material_props::Material: Material properties of the insulator.temperature::Union{Float64, Measurements.Measurement{Float64}}: Operating temperature of the insulator [°C].cross_section::Union{Float64, Measurements.Measurement{Float64}}: Cross-sectional area of the insulating layer [m²].resistance::Union{Float64, Measurements.Measurement{Float64}}: Electrical resistance of the insulating layer [Ω/m].gmr::Union{Float64, Measurements.Measurement{Float64}}: Geometric mean radius of the insulator [m].shunt_capacitance::Union{Float64, Measurements.Measurement{Float64}}: Shunt capacitance per unit length of the insulating layer [F/m].shunt_conductance::Union{Float64, Measurements.Measurement{Float64}}: Shunt conductance per unit length of the insulating layer [S·m].
LineCableModels.DataModel.Insulator — Method
Insulator(
component::CableComponent{T<:Union{Float64, Measurements.Measurement{Float64}}}
) -> Insulator
Constructs the equivalent coaxial insulation as an Insulator directly from a CableComponent, calling the strict positional constructor.
Arguments
component: TheCableComponentproviding geometry and material.
Returns
Insulator{T}with radii fromcomponent.insulator_groupand material fromcomponent.insulator_propsat the group temperature (fallback toT0).
LineCableModels.DataModel.Insulator — Method
Insulator(
r_in::Union{Float64, Measurements.Measurement{Float64}},
r_ex::Union{Float64, Measurements.Measurement{Float64}},
material_props::Material{T<:Union{Float64, Measurements.Measurement{Float64}}},
temperature::Union{Float64, Measurements.Measurement{Float64}}
) -> Insulator
Constructs an Insulator object with specified geometric and material parameters.
Arguments
r_in: Internal radius of the insulating layer [m].r_ex: External radius or thickness of the layer [m].material_props: Material properties of the insulating material.temperature: Operating temperature of the insulator [°C].
Returns
- An
Insulatorobject with calculated electrical properties.
Examples
material_props = Material(1e10, 3.0, 1.0, 20.0, 0.0)
insulator_layer = Insulator(0.01, 0.015, material_props, temperature=25)LineCableModels.DataModel.InsulatorGroup — Type
mutable struct InsulatorGroup{T<:Union{Float64, Measurements.Measurement{Float64}}} <: LineCableModels.DataModel.AbstractInsulatorPart{T<:Union{Float64, Measurements.Measurement{Float64}}}Represents a composite coaxial insulator group assembled from multiple insulating layers.
This structure serves as a container for different AbstractInsulatorPart elements (such as insulators and semiconductors) arranged in concentric layers. The InsulatorGroup aggregates these individual parts and provides equivalent electrical properties that represent the composite behavior of the entire assembly, stored in the attributes:
r_in::Union{Float64, Measurements.Measurement{Float64}}: Inner radius of the insulator group [m].r_ex::Union{Float64, Measurements.Measurement{Float64}}: Outer radius of the insulator group [m].cross_section::Union{Float64, Measurements.Measurement{Float64}}: Cross-sectional area of the entire insulator group [m²].shunt_capacitance::Union{Float64, Measurements.Measurement{Float64}}: Shunt capacitance per unit length of the insulator group [F/m].shunt_conductance::Union{Float64, Measurements.Measurement{Float64}}: Shunt conductance per unit length of the insulator group [S·m].layers::Array{LineCableModels.DataModel.AbstractInsulatorPart{T}, 1} where T<:Union{Float64, Measurements.Measurement{Float64}}: Vector of insulator layer components.
LineCableModels.DataModel.InsulatorGroup — Method
InsulatorGroup(
component::CableComponent{T<:Union{Float64, Measurements.Measurement{Float64}}}
) -> InsulatorGroup
Build an InsulatorGroup equivalent for a CableComponent while maintaining geometric coupling to the equivalent conductor group.
Stacks a single insulating layer of equivalent material and thickness over the new conductor group created from the same component.
LineCableModels.DataModel.InsulatorGroup — Method
Constructs an InsulatorGroup instance initializing with the initial insulator part.
Arguments
initial_insulator: AnAbstractInsulatorPartobject located at the innermost position of the insulator group.
Returns
- An
InsulatorGroupobject initialized with geometric and electrical properties derived from the initial insulator.
Examples
material_props = Material(1e10, 3.0, 1.0, 20.0, 0.0)
initial_insulator = Insulator(0.01, 0.015, material_props)
insulator_group = InsulatorGroup(initial_insulator)
println(insulator_group.layers) # Output: [initial_insulator]
println(insulator_group.shunt_capacitance) # Output: Capacitance in [F/m]LineCableModels.DataModel.LineCableSystem — Type
LineCableSystem(
system_id::String,
line_length::Number,
cable::CableDesign,
horz::Number,
vert::Number
) -> LineCableSystem
LineCableSystem(
system_id::String,
line_length::Number,
cable::CableDesign,
horz::Number,
vert::Number,
conn::Union{Nothing, Dict{String, Int64}}
) -> LineCableSystem
Weakly-typed convenience constructor. Builds a CablePosition from a CableDesign and coordinates, then constructs the system.
LineCableModels.DataModel.LineCableSystem — Type
mutable struct LineCableSystem{T<:Union{Float64, Measurements.Measurement{Float64}}}Represents a cable system configuration, defining the physical structure, cables, and their positions.
system_id::String: Unique identifier for the system.line_length::Union{Float64, Measurements.Measurement{Float64}}: Length of the cable system [m].num_cables::Int64: Number of cables in the system.num_phases::Int64: Number of actual phases in the system.cables::Array{CablePosition{T}, 1} where T<:Union{Float64, Measurements.Measurement{Float64}}: Cross-section cable positions.
LineCableModels.DataModel.LineCableSystem — Method
LineCableSystem(
system_id::String,
line_length::Number,
cable::CablePosition
) -> LineCableSystem
Weakly-typed constructor. Infers scalar type T from line_length and the cable (or its design), coerces as needed, and calls the strict kernel.
LineCableModels.DataModel.LineCableSystem — Method
Strict numeric kernel. Builds a typed LineCableSystem{T} from a vector of CablePosition{T}.
LineCableModels.DataModel.LineCableSystem — Method
Constructs a LineCableSystem with an initial cable position and system parameters.
Arguments
system_id: Identifier for the cable system.line_length: Length of the cable system [m].cable: InitialCablePositionobject defining a cable position and phase mapping.
Returns
- A
LineCableSystemobject initialized with a single cable position.
Examples
cable_design = CableDesign("example", nominal_data, components_dict)
cablepos1 = CablePosition(cable_design, 0.0, 0.0, Dict("core" => 1))
cable_system = LineCableSystem("test_case_1", 1000.0, cablepos1)
println(cable_system.num_phases) # Prints number of unique phase assignmentsSee also
LineCableModels.DataModel.NominalData — Type
struct NominalData{T<:Union{Float64, Measurements.Measurement{Float64}}}Stores nominal electrical and geometric parameters for a cable design.
designation_code::Union{Nothing, String}: Cable designation as per DIN VDE 0271/0276.U0::Union{Nothing, T} where T<:Union{Float64, Measurements.Measurement{Float64}}: Rated phase-to-earth voltage [kV].U::Union{Nothing, T} where T<:Union{Float64, Measurements.Measurement{Float64}}: Rated phase-to-phase voltage [kV].conductor_cross_section::Union{Nothing, T} where T<:Union{Float64, Measurements.Measurement{Float64}}: Cross-sectional area of the conductor [mm²].screen_cross_section::Union{Nothing, T} where T<:Union{Float64, Measurements.Measurement{Float64}}: Cross-sectional area of the screen [mm²].armor_cross_section::Union{Nothing, T} where T<:Union{Float64, Measurements.Measurement{Float64}}: Cross-sectional area of the armor [mm²].resistance::Union{Nothing, T} where T<:Union{Float64, Measurements.Measurement{Float64}}: Base (DC) resistance of the cable core [Ω/km].capacitance::Union{Nothing, T} where T<:Union{Float64, Measurements.Measurement{Float64}}: Capacitance of the main insulation [μF/km].inductance::Union{Nothing, T} where T<:Union{Float64, Measurements.Measurement{Float64}}: Inductance of the cable (trifoil formation) [mH/km].
LineCableModels.DataModel.NominalData — Method
NominalData(
;
designation_code,
U0,
U,
conductor_cross_section,
screen_cross_section,
armor_cross_section,
resistance,
capacitance,
inductance
) -> NominalData{Float64}
Weakly-typed constructor that infers the target scalar type T from the provided numeric kwargs (ignoring nothing and the string designation), coerces numerics to T, and calls the strict kernel.
If no numeric kwargs are provided, it defaults to Float64.
LineCableModels.DataModel.RectStrands — Type
struct RectStrands{T<:Union{Float64, Measurements.Measurement{Float64}}, S<:LineCableModels.DataModel.RectStrandsShape} <: LineCableModels.DataModel.AbstractStrandsLayer{T<:Union{Float64, Measurements.Measurement{Float64}}}Represents a concentric layer of rectangular strands with defined geometric, material, and electrical properties:
r_in::Union{Float64, Measurements.Measurement{Float64}}: Internal radial boundary [m].r_ex::Union{Float64, Measurements.Measurement{Float64}}: External radial boundary [m].material_props::Material: Material properties of the conductive strands.temperature::Union{Float64, Measurements.Measurement{Float64}}: Operating temperature of the layer [°C].resistance::Union{Float64, Measurements.Measurement{Float64}}: Equivalent electrical resistance of the layer [Ω/m].gmr::Union{Float64, Measurements.Measurement{Float64}}: Geometric mean radius (GMR) of the layer [m].shape::LineCableModels.DataModel.RectStrandsShape: Shape payload defining the internal geometric layout.
LineCableModels.DataModel.RectStrands — Method
RectStrands(
r_in::Union{Float64, Measurements.Measurement{Float64}},
thickness::Union{Float64, Measurements.Measurement{Float64}},
width::Union{Float64, Measurements.Measurement{Float64}},
num_wires::Int64,
lay_ratio::Union{Float64, Measurements.Measurement{Float64}},
material_props::Material{T<:Union{Float64, Measurements.Measurement{Float64}}},
temperature::Union{Float64, Measurements.Measurement{Float64}},
lay_direction::Int64
) -> LineCableModels.DataModel.RectStrands
Constructs a RectStrands instance.
Arguments
r_in: Internal radius of the layer [m].thickness: Radial thickness of the strands [m].width: Width of the individual rectangular strip [m].num_wires: Number of strands in the layer [dimensionless].lay_ratio: Ratio defining the lay length of the strands [dimensionless].material_props: AMaterialobject containing physical properties.temperature: Operating temperature [°C].lay_direction: Twisting direction (1 = unilay, -1 = contralay) [dimensionless].
Returns
- A
RectStrandsobject with calculated geometric and electrical properties.
Examples
material = Material(1.724e-8, 1.0, 1.0, 20.0, 0.00393)
layer = RectStrands(0.01, 0.002, 0.005, 10, 12.0, material, 25.0, 1)LineCableModels.DataModel.Sector — Type
struct Sector{T<:Union{Float64, Measurements.Measurement{Float64}}} <: LineCableModels.DataModel.AbstractConductorPart{T<:Union{Float64, Measurements.Measurement{Float64}}}Represents a single sector-shaped conductor with defined geometric and material properties.
r_in::Union{Float64, Measurements.Measurement{Float64}}: Inner radius (not applicable, typically 0 for the central point) [m].r_ex::Union{Float64, Measurements.Measurement{Float64}}: Outer radius (equivalent back radius) [m].params::SectorParams: Geometric parameters defining the sector's shape.rotation_angle_deg::Union{Float64, Measurements.Measurement{Float64}}: Rotation angle of this specific sector around the cable's center [degrees].material_props::Material: Material properties of the conductor.temperature::Union{Float64, Measurements.Measurement{Float64}}: Operating temperature of the conductor [°C].cross_section::Union{Float64, Measurements.Measurement{Float64}}: Cross-sectional area of the sector [m²].resistance::Union{Float64, Measurements.Measurement{Float64}}: Electrical resistance of the sector [Ω/m].gmr::Union{Float64, Measurements.Measurement{Float64}}: Geometric mean radius (GMR) of the sector (approximated) [m].vertices::Array{GeometryBasics.Point2{T}, 1} where T<:Union{Float64, Measurements.Measurement{Float64}}: Calculated vertices defining the polygon shape.centroid::GeometryBasics.Point2{T} where T<:Union{Float64, Measurements.Measurement{Float64}}: Geometric centroid of the sector shape.
LineCableModels.DataModel.SectorInsulator — Type
struct SectorInsulator{T<:Union{Float64, Measurements.Measurement{Float64}}} <: LineCableModels.DataModel.AbstractInsulatorPart{T<:Union{Float64, Measurements.Measurement{Float64}}}Represents an insulating layer surrounding a sector-shaped conductor.
r_in::Union{Float64, Measurements.Measurement{Float64}}: Inner radius (not applicable, defined by inner sector) [m].r_ex::Union{Float64, Measurements.Measurement{Float64}}: Outer radius (equivalent back radius of outer boundary) [m].inner_sector::Sector: The inner sector conductor that this insulator surrounds.thickness::Union{Float64, Measurements.Measurement{Float64}}: The thickness of the insulating layer [m].material_props::Material: Material properties of the insulator.temperature::Union{Float64, Measurements.Measurement{Float64}}: Operating temperature of the insulator [°C].cross_section::Union{Float64, Measurements.Measurement{Float64}}: Cross-sectional area of the insulating layer [m²].shunt_capacitance::Union{Float64, Measurements.Measurement{Float64}}: Shunt capacitance (approximated) [F/m].shunt_conductance::Union{Float64, Measurements.Measurement{Float64}}: Shunt conductance (approximated) [S·m].outer_vertices::Array{GeometryBasics.Point2{T}, 1} where T<:Union{Float64, Measurements.Measurement{Float64}}: Calculated vertices of the outer boundary polygon.
LineCableModels.DataModel.SectorParams — Type
struct SectorParams{T<:Union{Float64, Measurements.Measurement{Float64}}}Holds the geometric parameters that define the shape of a sector conductor.
n_sectors::Int64: Number of sectors in the full cable (e.g., 3 or 4).r_back::Union{Float64, Measurements.Measurement{Float64}}: Back radius of the sector (outermost curve) [m].d_sector::Union{Float64, Measurements.Measurement{Float64}}: Depth of the sector from the back to the base [m].r_corner::Union{Float64, Measurements.Measurement{Float64}}: Corner radius for rounding sharp edges [m].theta_cond_deg::Union{Float64, Measurements.Measurement{Float64}}: Angular width of the conductor's flat base/sides [degrees].d_insulation::Union{Float64, Measurements.Measurement{Float64}}: Insulation thickness [m].
LineCableModels.DataModel.Semicon — Type
struct Semicon{T<:Union{Float64, Measurements.Measurement{Float64}}} <: LineCableModels.DataModel.AbstractInsulatorPart{T<:Union{Float64, Measurements.Measurement{Float64}}}Represents a semiconducting layer with defined geometric, material, and electrical properties given by the attributes:
r_in::Union{Float64, Measurements.Measurement{Float64}}: Internal radius of the semiconducting layer [m].r_ex::Union{Float64, Measurements.Measurement{Float64}}: External radius of the semiconducting layer [m].material_props::Material: Material properties of the semiconductor.temperature::Union{Float64, Measurements.Measurement{Float64}}: Operating temperature of the semiconductor [°C].cross_section::Union{Float64, Measurements.Measurement{Float64}}: Cross-sectional area of the semiconducting layer [m²].resistance::Union{Float64, Measurements.Measurement{Float64}}: Electrical resistance of the semiconducting layer [Ω/m].gmr::Union{Float64, Measurements.Measurement{Float64}}: Geometric mean radius of the semiconducting layer [m].shunt_capacitance::Union{Float64, Measurements.Measurement{Float64}}: Shunt capacitance per unit length of the semiconducting layer [F/m].shunt_conductance::Union{Float64, Measurements.Measurement{Float64}}: Shunt conductance per unit length of the semiconducting layer [S·m].
LineCableModels.DataModel.Semicon — Method
Semicon(
r_in::Union{Float64, Measurements.Measurement{Float64}},
r_ex::Union{Float64, Measurements.Measurement{Float64}},
material_props::Material{T<:Union{Float64, Measurements.Measurement{Float64}}},
temperature::Union{Float64, Measurements.Measurement{Float64}}
) -> Semicon
Constructs a Semicon instance with calculated electrical and geometric properties.
Arguments
r_in: Internal radius of the semiconducting layer [m].r_ex: External radius or thickness of the layer [m].material_props: Material properties of the semiconducting material.temperature: Operating temperature of the layer [°C] (default: T₀).
Returns
- A
Semiconobject with initialized properties.
Examples
material_props = Material(1e6, 2.3, 1.0, 20.0, 0.00393)
semicon_layer = Semicon(0.01, Thickness(0.002), material_props, temperature=25)
println(semicon_layer.cross_section) # Expected output: ~6.28e-5 [m²]
println(semicon_layer.resistance) # Expected output: Resistance in [Ω/m]
println(semicon_layer.gmr) # Expected output: GMR in [m]
println(semicon_layer.shunt_capacitance) # Expected output: Capacitance in [F/m]
println(semicon_layer.shunt_conductance) # Expected output: Conductance in [S·m]LineCableModels.DataModel.Strip — Type
struct Strip{T<:Union{Float64, Measurements.Measurement{Float64}}} <: LineCableModels.DataModel.AbstractConductorPart{T<:Union{Float64, Measurements.Measurement{Float64}}}Represents a flat conductive strip with defined geometric and material properties given by the attributes:
r_in::Union{Float64, Measurements.Measurement{Float64}}: Internal radius of the strip [m].r_ex::Union{Float64, Measurements.Measurement{Float64}}: External radius of the strip [m].thickness::Union{Float64, Measurements.Measurement{Float64}}: Thickness of the strip [m].width::Union{Float64, Measurements.Measurement{Float64}}: Width of the strip [m].lay_ratio::Union{Float64, Measurements.Measurement{Float64}}: Ratio defining the lay length of the strip (twisting factor) [dimensionless].mean_diameter::Union{Float64, Measurements.Measurement{Float64}}: Mean diameter of the strip's helical path [m].pitch_length::Union{Float64, Measurements.Measurement{Float64}}: Pitch length of the strip's helical path [m].lay_direction::Int64: Twisting direction of the strip (1 = unilay, -1 = contralay) [dimensionless].material_props::Material: Material properties of the strip.temperature::Union{Float64, Measurements.Measurement{Float64}}: Temperature at which the properties are evaluated [°C].cross_section::Union{Float64, Measurements.Measurement{Float64}}: Cross-sectional area of the strip [m²].resistance::Union{Float64, Measurements.Measurement{Float64}}: Electrical resistance of the strip [Ω/m].gmr::Union{Float64, Measurements.Measurement{Float64}}: Geometric mean radius of the strip [m].
LineCableModels.DataModel.Strip — Method
Strip(
r_in::Union{Float64, Measurements.Measurement{Float64}},
r_ex::Union{Float64, Measurements.Measurement{Float64}},
width::Union{Float64, Measurements.Measurement{Float64}},
lay_ratio::Union{Float64, Measurements.Measurement{Float64}},
material_props::Material{T<:Union{Float64, Measurements.Measurement{Float64}}},
temperature::Union{Float64, Measurements.Measurement{Float64}},
lay_direction::Int64
) -> Strip
Constructs a Strip object with specified geometric and material parameters.
Arguments
r_in: Internal radius of the strip [m].r_ex: External radius or thickness of the strip [m].width: Width of the strip [m].lay_ratio: Ratio defining the lay length of the strip [dimensionless].material_props: Material properties of the strip.temperature: Temperature at which the properties are evaluated [°C]. Defaults toT₀.lay_direction: Twisting direction of the strip (1 = unilay, -1 = contralay) [dimensionless]. Defaults to 1.
Returns
- A
Stripobject with calculated geometric and electrical properties.
Examples
material_props = Material(1.7241e-8, 1.0, 0.999994, 20.0, 0.00393)
strip = Strip(0.01, Thickness(0.002), 0.05, 10, material_props, temperature=25)
println(strip.cross_section) # Output: 0.0001 [m²]
println(strip.resistance) # Output: Resistance value [Ω/m]See also
LineCableModels.DataModel.Thickness — Type
struct Thickness{T<:Real} <: LineCableModels.DataModel.AbstractRadiusRepresents the thickness of a cable component.
value::Real: Numerical value of the thickness [m].
LineCableModels.DataModel.Tubular — Type
struct Tubular{T<:Union{Float64, Measurements.Measurement{Float64}}} <: LineCableModels.DataModel.AbstractConductorPart{T<:Union{Float64, Measurements.Measurement{Float64}}}Represents a tubular or solid (r_in=0) conductor with geometric and material properties defined as:
r_in::Union{Float64, Measurements.Measurement{Float64}}: Internal radius of the tubular conductor [m].r_ex::Union{Float64, Measurements.Measurement{Float64}}: External radius of the tubular conductor [m].material_props::Material: AMaterialobject representing the physical properties of the conductor material.temperature::Union{Float64, Measurements.Measurement{Float64}}: Temperature at which the properties are evaluated [°C].cross_section::Union{Float64, Measurements.Measurement{Float64}}: Cross-sectional area of the tubular conductor [m²].resistance::Union{Float64, Measurements.Measurement{Float64}}: Electrical resistance (DC) of the tubular conductor [Ω/m].gmr::Union{Float64, Measurements.Measurement{Float64}}: Geometric mean radius of the tubular conductor [m].
LineCableModels.DataModel.Tubular — Method
Tubular(
component::CableComponent{T<:Union{Float64, Measurements.Measurement{Float64}}}
) -> Tubular
Constructs the equivalent coaxial conductor as a Tubular directly from a CableComponent, reusing the rigorously tested positional constructor.
Arguments
component: TheCableComponentproviding geometry and material.
Returns
Tubular{T}with radii fromcomponent.conductor_groupand material fromcomponent.conductor_propsat the group temperature (fallback toT0).
LineCableModels.DataModel.Tubular — Method
Tubular(
r_in::Union{Float64, Measurements.Measurement{Float64}},
r_ex::Union{Float64, Measurements.Measurement{Float64}},
material_props::Material{T<:Union{Float64, Measurements.Measurement{Float64}}},
temperature::Union{Float64, Measurements.Measurement{Float64}}
) -> Tubular
Initializes a Tubular object with specified geometric and material parameters.
Arguments
r_in: Internal radius of the tubular conductor [m].r_ex: External radius of the tubular conductor [m].material_props: AMaterialobject representing the physical properties of the conductor material.temperature: Temperature at which the properties are evaluated [°C]. Defaults toT₀.
Returns
- An instance of
Tubularinitialized with calculated geometric and electrical properties.
Examples
material_props = Material(1.7241e-8, 1.0, 0.999994, 20.0, 0.00393)
tubular = Tubular(0.01, 0.02, material_props, temperature=25)
println(tubular.cross_section) # Output: 0.000942 [m²]
println(tubular.resistance) # Output: Resistance value [Ω/m]See also
LineCableModels.DataModel.equivalent — Method
equivalent(
original_design::CableDesign;
new_id
) -> CableDesign
Builds a simplified CableDesign by replacing each component with a homogeneous equivalent and leverages shorthand constructors:
ConductorGroup(component::CableComponent{T}) = ConductorGroup(Tubular(component))InsulatorGroup(component::CableComponent{T}) = InsulatorGroup(Insulator(component))
The geometry is preserved from the original component, while materials are derived from the component's effective conductor and insulator properties.
LineCableModels.DataModel.flat_formation — Method
flat_formation(
xc::Union{Float64, Measurements.Measurement{Float64}},
yc::Union{Float64, Measurements.Measurement{Float64}},
s::Union{Float64, Measurements.Measurement{Float64}};
vertical
) -> NTuple{6, Union{Float64, Measurements.Measurement{Float64}}}
Calculates the coordinates of three conductors arranged in a flat (horizontal or vertical) formation.
Arguments
xc: X-coordinate of the reference point [m].yc: Y-coordinate of the reference point [m].s: Spacing between adjacent conductors [m].vertical: Boolean flag indicating whether the formation is vertical.
Returns
- A tuple containing:
xa,ya: Coordinates of the first conductor [m].xb,yb: Coordinates of the second conductor [m].xc,yc: Coordinates of the third conductor [m].
Examples
# Horizontal formation
xa, ya, xb, yb, xc, yc = flat_formation(0.0, 0.0, 0.1)
println((xa, ya)) # First conductor coordinates
println((xb, yb)) # Second conductor coordinates
println((xc, yc)) # Third conductor coordinates
# Vertical formation
xa, ya, xb, yb, xc, yc = flat_formation(0.0, 0.0, 0.1, vertical=true)LineCableModels.DataModel.trifoil_formation — Method
trifoil_formation(
x0::Union{Float64, Measurements.Measurement{Float64}},
y0::Union{Float64, Measurements.Measurement{Float64}},
r_ext::Union{Float64, Measurements.Measurement{Float64}}
) -> NTuple{6, Union{Float64, Measurements.Measurement{Float64}}}
Calculates the coordinates of three cables arranged in a trifoil pattern.
Arguments
x0: X-coordinate of the center point [m].y0: Y-coordinate of the center point [m].r_ext: External radius of the circular layout [m].
Returns
- A tuple containing:
xa,ya: Coordinates of the top cable [m].xb,yb: Coordinates of the bottom-left cable [m].xc,yc: Coordinates of the bottom-right cable [m].
Examples
xa, ya, xb, yb, xc, yc = trifoil_formation(0.0, 0.0, 0.035)
println((xa, ya)) # Coordinates of top cable
println((xb, yb)) # Coordinates of bottom-left cable
println((xc, yc)) # Coordinates of bottom-right cableBase parameters (R, L, C, G)
LineCableModels.DataModel.BaseParams — Module
LineCableModels.DataModel.BaseParamsThe BaseParams submodule provides fundamental functions for determining the base electrical parameters (R, L, C, G) of cable components within the LineCableModels.DataModel module. This includes implementations of standard engineering formulas for resistance, inductance, and geometric parameters of various conductor configurations.
Overview
- Implements basic electrical engineering formulas for calculating DC resistance and inductance of different conductor geometries (tubular, strip, wire arrays).
- Implements basic formulas for capacitance and dielectric losses in insulators and semiconductors.
- Provides functions for temperature correction of material properties.
- Calculates geometric mean radii for different conductor configurations.
- Includes functions for determining the effective length for helical wire arrangements.
- Calculates equivalent electrical parameters and correction factors for different geometries and configurations.
Dependencies
BaseLineCableModels.CommonsMeasurements
Exports
calc_circstrands_coordscalc_circstrands_gmrcalc_equivalent_alphacalc_equivalent_epscalc_equivalent_gmrcalc_equivalent_lossfactcalc_equivalent_mucalc_equivalent_rhocalc_gmdcalc_helical_paramscalc_inductance_trifoilcalc_parallel_equivalentcalc_shunt_capacitancecalc_shunt_conductancecalc_sigma_lossfactcalc_solenoid_correctioncalc_strip_resistancecalc_temperature_correctioncalc_tubular_gmrcalc_tubular_inductancecalc_tubular_resistance
LineCableModels.DataModel.BaseParams.calc_circstrands_coords — Method
calc_circstrands_coords(
num_wires::Int64,
radius_wire::Union{Float64, Measurements.Measurement{Float64}},
r_in::Union{Float64, Measurements.Measurement{Float64}},
C::Tuple{T<:Union{Float64, Measurements.Measurement{Float64}}, T<:Union{Float64, Measurements.Measurement{Float64}}}
) -> Vector{<:Tuple{Union{Float64, Measurements.Measurement{Float64}}, Union{Float64, Measurements.Measurement{Float64}}}}
Calculates the center coordinates of wires arranged in a circular pattern.
Arguments
num_wires: Number of wires in the circular arrangement [dimensionless].radius_wire: Radius of each individual wire [m].r_in: Inner radius of the wire array (to wire centers) [m].C: Optional tuple representing the center coordinates of the circular arrangement [m]. Default is (0.0, 0.0).
Returns
- Vector of tuples, where each tuple contains the
(x, y)coordinates [m] of the center of a wire.
Examples
# Create a 7-wire array with 2mm wire radius and 1cm inner radius
wire_coords = calc_circstrands_coords(7, 0.002, 0.01)
println(wire_coords[1]) # Output: First wire coordinates
# Create a wire array with custom center position
wire_coords = calc_circstrands_coords(7, 0.002, 0.01, C=(0.5, 0.3))See also
LineCableModels.DataModel.BaseParams.calc_circstrands_gmr — Method
calc_circstrands_gmr(
lay_rad::Union{Float64, Measurements.Measurement{Float64}},
N::Int64,
rad_wire::Union{Float64, Measurements.Measurement{Float64}},
mu_r::Union{Float64, Measurements.Measurement{Float64}}
) -> Union{Float64, Measurements.Measurement{Float64}}
Calculates the geometric mean radius (GMR) of a circular wire array, using formula (62), page 335, of the book by Edward Rosa [15]:
\[GMR = \sqrt[n] {r n a^{n-1}}\]
where $a$ is the layout radius, $n$ is the number of wires, and $r$ is the radius of each wire.
Arguments
lay_rad: Layout radius of the wire array [m].N: Number of wires in the array [dimensionless].rad_wire: Radius of an individual wire [m].mu_r: Relative permeability of the wire material [dimensionless].
Returns
- Geometric mean radius (GMR) of the wire array [m].
Examples
lay_rad = 0.05
N = 7
rad_wire = 0.002
mu_r = 1.0
gmr = calc_circstrands_gmr(lay_rad, N, rad_wire, mu_r)
println(gmr) # Expected output: 0.01187... [m]LineCableModels.DataModel.BaseParams.calc_equivalent_alpha — Method
calc_equivalent_alpha(
alpha1::Union{Float64, Measurements.Measurement{Float64}},
R1::Union{Float64, Measurements.Measurement{Float64}},
alpha2::Union{Float64, Measurements.Measurement{Float64}},
R2::Union{Float64, Measurements.Measurement{Float64}}
) -> Union{Float64, Measurements.Measurement{Float64}}
Calculates the equivalent temperature coefficient of resistance (alpha) when two conductors are connected in parallel, by cross-weighted-resistance averaging:
\[\alpha_{eq} = \frac{\alpha_1 R_2 + \alpha_2 R1}{R_1 + R_2}\]
where $\alpha_1$, $\alpha_2$ are the temperature coefficients of the conductors, and $R_1$, $R_2$ are the respective resistances.
Arguments
alpha1: Temperature coefficient of resistance of the first conductor [1/°C].R1: Resistance of the first conductor [Ω].alpha2: Temperature coefficient of resistance of the second conductor [1/°C].R2: Resistance of the second conductor [Ω].
Returns
- The equivalent temperature coefficient [1/°C] for the parallel combination.
Examples
alpha_conductor = 0.00393 # Copper
alpha_new_part = 0.00403 # Aluminum
R_conductor = 0.5
R_new_part = 1.0
alpha_eq = calc_equivalent_alpha(alpha_conductor, R_conductor, alpha_new_part, R_new_part)
println(alpha_eq) # Output: 0.00396 (approximately)LineCableModels.DataModel.BaseParams.calc_equivalent_eps — Method
calc_equivalent_eps(
C_eq::Union{Float64, Measurements.Measurement{Float64}},
r_ex::Union{Float64, Measurements.Measurement{Float64}},
r_in::Union{Float64, Measurements.Measurement{Float64}}
) -> Union{Float64, Measurements.Measurement{Float64}}
Calculates the equivalent permittivity for a coaxial cable insulation, using the formula [4]:
\[\varepsilon_{eq} = \frac{C_{eq} \log(\frac{r_{ext}}{r_{in}})}{2\pi \varepsilon_0}\]
where $\varepsilon_0$ is the permittivity of free space.
Arguments
C_eq: Equivalent capacitance of the insulation [F/m].r_ex: External radius of the insulation [m].r_in: Internal radius of the insulation [m].
Returns
- Equivalent relative permittivity of the insulation [dimensionless].
Examples
eps_eq = calc_equivalent_eps(1e-10, 0.01, 0.005) # Expected output: ~2.26 [dimensionless]See also
LineCableModels.DataModel.BaseParams.calc_equivalent_gmr — Method
calc_equivalent_gmr(
existing::LineCableModels.DataModel.AbstractCablePart,
new_layer::LineCableModels.DataModel.AbstractCablePart
) -> Any
Calculates the equivalent geometric mean radius (GMR) of a conductor after adding a new layer, by recursive application of the multizone stranded conductor defined as [3]:
\[GMR_{eq} = {GMR_{i-1}}^{\beta^2} \cdot {GMR_{i}}^{(1-\beta)^2} \cdot {GMD}^{2\beta(1-\beta)}\]
\[\beta = \frac{S_{i-1}}{S_{i-1} + S_{i}}\]
where:
- $S_{i-1}$ is the cumulative cross-sectional area of the existing cable part, $S_{i}$ is the total cross-sectional area after inclusion of the conducting layer ${i}$.
- $GMR_{i-1}$ is the cumulative GMR of the existing cable part, $GMR_{i}$ is the GMR of the conducting layer ${i}$.
- $GMD$ is the geometric mean distance between the existing cable part and the new layer, calculated using
calc_gmd.
Arguments
existing: The existing cable part (AbstractCablePart).new_layer: The new layer being added (AbstractCablePart).
Returns
- Updated equivalent GMR of the combined conductor [m].
Examples
material_props = Material(1.7241e-8, 1.0, 0.999994, 20.0, 0.00393)
conductor = Conductor(Strip(0.01, 0.002, 0.05, 10, material_props))
new_layer = CircStrands(0.02, 0.002, 7, 15, material_props)
equivalent_gmr = calc_equivalent_gmr(conductor, new_layer) # Expected output: Updated GMR value [m]See also
LineCableModels.DataModel.BaseParams.calc_equivalent_lossfact — Method
calc_equivalent_lossfact(
G_eq::Union{Float64, Measurements.Measurement{Float64}},
C_eq::Union{Float64, Measurements.Measurement{Float64}},
ω::Union{Float64, Measurements.Measurement{Float64}}
) -> Union{Float64, Measurements.Measurement{Float64}}
Calculates the equivalent loss factor (tangent) of a dielectric material:
\[\tan \delta = \frac{G_{eq}}{\omega \cdot C_{eq}}\]
where $\tan \delta$ is the loss factor (tangent).
Arguments
G_eq: Equivalent conductance of the material [S·m].C_eq: Equivalent capacitance of the material [F/m].ω: Angular frequency [rad/s].
Returns
- Equivalent loss factor of the dielectric material [dimensionless].
Examples
loss_factor = calc_equivalent_lossfact(1e-8, 1e-10, 2π*50) # Expected output: ~0.0318 [dimensionless]LineCableModels.DataModel.BaseParams.calc_equivalent_mu — Method
calc_equivalent_mu(
gmr::Union{Float64, Measurements.Measurement{Float64}},
r_ex::Union{Float64, Measurements.Measurement{Float64}},
r_in::Union{Float64, Measurements.Measurement{Float64}}
) -> Union{Float64, Measurements.Measurement{Float64}}
Calculates the relative permeability (mu_r) based on the geometric mean radius (GMR) and conductor dimensions, by executing the inverse of calc_tubular_gmr, and solving for mu_r:
\[\log GMR = \log r_2 - \mu_r \left[ \frac{r_1^4}{\left(r_2^2 - r_1^2\right)^2} \log\left(\frac{r_2}{r_1}\right) - \frac{3r_1^2 - r_2^2}{4\left(r_2^2 - r_1^2\right)} \right]\]
\[\mu_r = -\frac{\left(\log GMR - \log r_2\right)}{\frac{r_1^4}{\left(r_2^2 - r_1^2\right)^2} \log\left(\frac{r_2}{r_1}\right) - \frac{3r_1^2 - r_2^2}{4\left(r_2^2 - r_1^2\right)}}\]
where $r_1$ is the inner radius and $r_2$ is the outer radius.
Arguments
gmr: Geometric mean radius of the conductor [m].r_ex: External radius of the conductor [m].r_in: Internal radius of the conductor [m].
Returns
- Relative permeability (
mu_r) of the conductor material [dimensionless].
Errors
- Throws
ArgumentErrorifr_exis less thanr_in.
Notes
Assumes a tubular geometry for the conductor, reducing to the solid case if r_in is zero.
Examples
gmr = 0.015
r_ex = 0.02
r_in = 0.01
mu_r = calc_equivalent_mu(gmr, r_ex, r_in)
println(mu_r) # Expected output: ~1.7 [dimensionless]See also
LineCableModels.DataModel.BaseParams.calc_equivalent_rho — Method
calc_equivalent_rho(
R::Union{Float64, Measurements.Measurement{Float64}},
radius_ext_con::Union{Float64, Measurements.Measurement{Float64}},
radius_in_con::Union{Float64, Measurements.Measurement{Float64}}
) -> Union{Float64, Measurements.Measurement{Float64}}
Calculates the equivalent resistivity of a solid tubular conductor, using the formula [4]:
\[\rho_{eq} = R_{eq} S_{eff} = R_{eq} \pi (r_{ext}^2 - r_{in}^2)\]
where $S_{eff}$ is the effective cross-sectional area of the tubular conductor.
Arguments
R: Resistance of the conductor [Ω].radius_ext_con: External radius of the tubular conductor [m].radius_in_con: Internal radius of the tubular conductor [m].
Returns
- Equivalent resistivity of the tubular conductor [Ω·m].
Examples
rho_eq = calc_equivalent_rho(0.01, 0.02, 0.01) # Expected output: ~9.42e-4 [Ω·m]LineCableModels.DataModel.BaseParams.calc_gmd — Method
calc_gmd(
co1::LineCableModels.DataModel.AbstractCablePart,
co2::LineCableModels.DataModel.AbstractCablePart
) -> Any
Calculates the geometric mean distance (GMD) between two cable parts, by using the definition described in Grover [16]:
\[\log GMD = \left(\frac{\sum_{i=1}^{n_1}\sum_{j=1}^{n_2} (s_1 \cdot s_2) \cdot \log(d_{ij})}{\sum_{i=1}^{n_1}\sum_{j=1}^{n_2} (s_1 \cdot s_2)}\right)\]
where:
- $d_{ij}$ is the Euclidean distance between elements $i$ and $j$.
- $s_1$ and $s_2$ are the cross-sectional areas of the respective elements.
- $n_1$ and $n_2$ are the number of sub-elements in each cable part.
Arguments
co1: First cable part (AbstractCablePart).co2: Second cable part (AbstractCablePart).
Returns
- Geometric mean distance between the cable parts [m].
Notes
For concentric structures, the GMD converges to the external radii of the outermost element.
This implementation uses a weighted sum of logarithms rather than the traditional product formula $\Pi(d_{ij})^{(1/n)}$ found in textbooks. The logarithmic approach prevents numerical underflow/overflow when dealing with many conductors or extreme distance ratios, making it significantly more stable for practical calculations.
Examples
material_props = Material(1.7241e-8, 1.0, 0.999994, 20.0, 0.00393)
circstrands1 = CircStrands(0.01, 0.002, 7, 10, material_props)
circstrands2 = CircStrands(0.02, 0.002, 7, 15, material_props)
gmd = calc_gmd(circstrands1, circstrands2) # Expected output: GMD value [m]
strip = Strip(0.01, 0.002, 0.05, 10, material_props)
tubular = Tubular(0.01, 0.02, material_props)
gmd = calc_gmd(strip, tubular) # Expected output: GMD value [m]See also
LineCableModels.DataModel.BaseParams.calc_helical_params — Method
calc_helical_params(
r_in::Union{Float64, Measurements.Measurement{Float64}},
r_ex::Union{Float64, Measurements.Measurement{Float64}},
lay_ratio::Union{Float64, Measurements.Measurement{Float64}}
) -> Tuple{Union{Float64, Measurements.Measurement{Float64}}, Union{Float64, Measurements.Measurement{Float64}}, Union{Float64, Int64, Measurements.Measurement{Float64}}}
Calculates the mean diameter, pitch length, and overlength based on cable geometry parameters. The lay ratio is defined as the ratio of the pitch length $L_p$ to the external diameter $D_e$:
\[\lambda = \frac{L_p}{D_e}\]
where $D_e$ and $L_p$ are the dimensions represented in the figure.
Arguments
r_in: Inner radius of the cable layer [m].r_ex: Outer radius of the cable layer [m].lay_ratio: Ratio of the pitch (lay) length to the external diameter of the corresponding layer of wires [dimensionless].
Returns
mean_diameter: Mean diameter of the cable layer [m].pitch_length: The length over which the strands complete one full twist [m].overlength: Effective length increase resulting from the helical path [1/m].
Notes
Reference values for lay_ratio are given under standard EN 50182 [12]:
| Conductor type | Steel wires | Aluminum wires | Lay ratio - Steel | Lay ratio - Aluminum |
|---|---|---|---|---|
| AAAC 4 layers | - | 61 (1/6/12/18/24) | - | 15/13.5/12.5/11 |
| ACSR 3 layers | 7 (1/6) | 54 (12/18/24) | 19 | 15/13/11.5 |
| ACSR 2 layers | 7 (1/6) | 26 (10/16) | 19 | 14/11.5 |
| ACSR 1 layer | 7 (1/6) | 10 | 19 | 14 |
| ACCC/TW | - | 36 (8/12/16) | - | 15/13.5/11.5 |
Examples
r_in = 0.01
r_ex = 0.015
lay_ratio = 12
mean_diam, pitch, overlength = calc_helical_params(r_in, r_ex, lay_ratio)
# mean_diam ≈ 0.025 [m]
# pitch ≈ 0.3 [m]
# overlength > 1.0 [1/m]LineCableModels.DataModel.BaseParams.calc_inductance_trifoil — Method
calc_inductance_trifoil(
r_in_co::Union{Float64, Measurements.Measurement{Float64}},
r_ext_co::Union{Float64, Measurements.Measurement{Float64}},
rho_co::Union{Float64, Measurements.Measurement{Float64}},
mu_r_co::Union{Float64, Measurements.Measurement{Float64}},
r_in_scr::Union{Float64, Measurements.Measurement{Float64}},
r_ext_scr::Union{Float64, Measurements.Measurement{Float64}},
rho_scr::Union{Float64, Measurements.Measurement{Float64}},
mu_r_scr::Union{Float64, Measurements.Measurement{Float64}},
S::Union{Float64, Measurements.Measurement{Float64}},
rho_e::Union{Float64, Measurements.Measurement{Float64}},
f::Union{Float64, Measurements.Measurement{Float64}}
) -> Union{Float64, Measurements.Measurement{Float64}}
Calculates the positive-sequence inductance of a trifoil-configured cable system composed of core/screen assuming solid bonding, using the formula given under section 4.2.4.3 of CIGRE TB-531:
\[Z_d = \left[Z_a - Z_x\right] - \frac{\left( Z_m - Z_x \right)^2}{Z_s - Z_x}\]
\[L = \mathfrak{Im}\left(\frac{Z_d}{\omega}\right)\]
where $Z_a$, $Z_s$ are the self impedances of the core conductor and the screen, and $Z_m$, and $Z_x$ are the mutual impedances between core/screen and between cables, respectively, as per sections 4.2.3.4, 4.2.3.5, 4.2.3.6 and 4.2.3.8 of the same document [8].
Arguments
r_in_co: Internal radius of the phase conductor [m].r_ext_co: External radius of the phase conductor [m].rho_co: Electrical resistivity of the phase conductor material [Ω·m].mu_r_co: Relative permeability of the phase conductor material [dimensionless].r_in_scr: Internal radius of the metallic screen [m].r_ext_scr: External radius of the metallic screen [m].rho_scr: Electrical resistivity of the metallic screen material [Ω·m].mu_r_scr: Relative permeability of the screen conductor material [dimensionless].S: Spacing between conductors in trifoil configuration [m].rho_e: Soil resistivity [Ω·m]. Default: 100 Ω·m.f: Frequency [Hz]. Default:f₀.
Returns
- Positive-sequence inductance per unit length of the cable system [H/m].
Examples
L = calc_inductance_trifoil(0.01, 0.015, 1.72e-8, 1.0, 0.02, 0.025, 2.83e-8, 1.0, S=0.1, rho_e=50, f=50)
println(L) # Output: Inductance value in H/mSee also
LineCableModels.DataModel.BaseParams.calc_parallel_equivalent — Method
calc_parallel_equivalent(
Z1::Union{Float64, Complex{Measurements.Measurement{Float64}}, Measurements.Measurement{Float64}, ComplexF64},
Z2::Union{Float64, Complex{Measurements.Measurement{Float64}}, Measurements.Measurement{Float64}, ComplexF64}
) -> Any
Calculates the parallel equivalent of two impedances (or series equivalent of two admittances):
\[Z_{eq} = \frac{Z_1 Z_2}{Z_1 + Z_2}\]
This expression, when applied recursively to LineCableModels.DataModel.CircStrands objects, implements the formula for the hexagonal wiring pattern described in CIGRE TB-345 [1] [17]:
\[\frac{1}{R_{\text{dc}}} = \frac{\pi d^2}{4 \rho} \left( 1 + \sum_{1}^{n} \frac{6n}{k_n} \right)\]
\[k_n = \left[ 1 + \left( \pi \frac{D_n}{\lambda_n} \right)^2 \right]^{1/2}\]
where $R_{\text{dc}}$ is the DC resistance, $d$ is the diameter of each wire, $ho$ is the resistivity, $n$ is the number of layers following the hexagonal pattern, $D_n$ is the diameter of the $n$-th layer, and $\lambda_n$ is the pitch length of the $n$-th layer, obtained using calc_helical_params.
Arguments
Z1: The total impedance of the existing system [Ω].Z2: The impedance of the new layer being added [Ω].
Returns
- The parallel equivalent impedance [Ω].
Examples
Z1 = 5.0
Z2 = 10.0
Req = calc_parallel_equivalent(Z1, Z2)
println(Req) # Outputs: 3.3333333333333335See also
LineCableModels.DataModel.BaseParams.calc_shunt_capacitance — Method
calc_shunt_capacitance(
r_in::Union{Float64, Measurements.Measurement{Float64}},
r_ex::Union{Float64, Measurements.Measurement{Float64}},
epsr::Union{Float64, Measurements.Measurement{Float64}}
) -> Union{Float64, Measurements.Measurement{Float64}}
Calculates the shunt capacitance per unit length of a coaxial structure, using the standard formula for the capacitance of a coaxial structure [8] [4] [14]:
\[C = \frac{2 \pi \varepsilon_0 \varepsilon_r}{\log \left(\frac{r_{ext}}{r_{in}}\right)}\]
where $\varepsilon_0$ is the vacuum permittivity, $\varepsilon_r$ is the relative permittivity of the dielectric material, and $r_{in}$ and $r_{ext}$ are the inner and outer radii of the coaxial structure, respectively.
Arguments
r_in: Internal radius of the coaxial structure [m].r_ex: External radius of the coaxial structure [m].epsr: Relative permittivity of the dielectric material [dimensionless].
Returns
- Shunt capacitance per unit length [F/m].
Examples
r_in = 0.01
r_ex = 0.02
epsr = 2.3
capacitance = calc_shunt_capacitance(r_in, r_ex, epsr)
println(capacitance) # Expected output: ~1.24e-10 [F/m]LineCableModels.DataModel.BaseParams.calc_shunt_conductance — Method
calc_shunt_conductance(
r_in::Union{Float64, Measurements.Measurement{Float64}},
r_ex::Union{Float64, Measurements.Measurement{Float64}},
rho::Union{Float64, Measurements.Measurement{Float64}}
) -> Union{Float64, Measurements.Measurement{Float64}}
Calculates the shunt conductance per unit length of a coaxial structure, using the improved model reported in [4] [13] [18]:
\[G = \frac{2\pi\sigma}{\log(\frac{r_{ext}}{r_{in}})}\]
where $\sigma = \frac{1}{\rho}$ is the conductivity of the dielectric/semiconducting material, $r_{in}$ is the internal radius, and $r_{ext}$ is the external radius of the coaxial structure.
Arguments
r_in: Internal radius of the coaxial structure [m].r_ex: External radius of the coaxial structure [m].rho: Resistivity of the dielectric/semiconducting material [Ω·m].
Returns
- Shunt conductance per unit length [S·m].
Examples
r_in = 0.01
r_ex = 0.02
rho = 1e9
g = calc_shunt_conductance(r_in, r_ex, rho)
println(g) # Expected output: 2.7169e-9 [S·m]LineCableModels.DataModel.BaseParams.calc_sigma_lossfact — Method
calc_sigma_lossfact(
G_eq::Union{Float64, Measurements.Measurement{Float64}},
r_in::Union{Float64, Measurements.Measurement{Float64}},
r_ex::Union{Float64, Measurements.Measurement{Float64}}
) -> Union{Float64, Measurements.Measurement{Float64}}
Calculates the effective conductivity of a dielectric material from the known conductance (related to the loss factor $\tan \delta$) via [4] [13] [18]:
\[\sigma_{eq} = \frac{G_{eq}}{2\pi} \log(\frac{r_{ext}}{r_{in}})\]
where $\sigma_{eq} = \frac{1}{\rho_{eq}}$ is the conductivity of the dielectric/semiconducting material, $G_{eq}$ is the shunt conductance per unit length, $r_{in}$ is the internal radius, and $r_{ext}$ is the external radius of the coaxial structure.
Arguments
G_eq: Equivalent conductance of the material [S·m].r_in: Internal radius of the coaxial structure [m].r_ex: External radius of the coaxial structure [m].
Returns
- Effective material conductivity per unit length [S·m].
Examples
Geq = 2.7169e-9
sigma_eq = calc_sigma_lossfact(G_eq, r_in, r_ex)LineCableModels.DataModel.BaseParams.calc_solenoid_correction — Method
calc_solenoid_correction(
num_turns::Union{Float64, Measurements.Measurement{Float64}},
radius_ext_con::Union{Float64, Measurements.Measurement{Float64}},
radius_ext_ins::Union{Float64, Measurements.Measurement{Float64}}
) -> Union{Float64, Measurements.Measurement{Float64}}
Calculates the solenoid correction factor for magnetic permeability in insulated cables with helical conductors (CircStrands), using the formula from Gudmundsdottir et al. [5]:
\[\mu_{r, sol} = 1 + \frac{2 \pi^2 N^2 (r_{ins, ext}^2 - r_{con, ext}^2)}{\log(r_{ins, ext}/r_{con, ext})}\]
where:
- $N$ is the number of turns per unit length.
- $r_{con, ext}$ is the conductor external radius.
- $r_{ins, ext}$ is the insulator external radius.
Arguments
num_turns: Number of turns per unit length [1/m].radius_ext_con: External radius of the conductor [m].radius_ext_ins: External radius of the insulator [m].
Returns
- Correction factor for the insulator magnetic permeability [dimensionless].
Examples
# Cable with 10 turns per meter, conductor radius 5 mm, insulator radius 10 mm
correction = calc_solenoid_correction(10, 0.005, 0.01) # Expected output: > 1.0 [dimensionless]
# Non-helical cable (straight conductor)
correction = calc_solenoid_correction(NaN, 0.005, 0.01) # Expected output: 1.0 [dimensionless]LineCableModels.DataModel.BaseParams.calc_strip_resistance — Method
calc_strip_resistance(
thickness::Union{Float64, Measurements.Measurement{Float64}},
width::Union{Float64, Measurements.Measurement{Float64}},
rho::Union{Float64, Measurements.Measurement{Float64}},
alpha::Union{Float64, Measurements.Measurement{Float64}},
T0::Union{Float64, Measurements.Measurement{Float64}},
Top::Union{Float64, Measurements.Measurement{Float64}}
) -> Union{Float64, Measurements.Measurement{Float64}}
Calculates the DC resistance of a strip conductor based on its geometric and material properties, using the basic resistance formula in terms of the resistivity and cross-sectional area:
\[R = \rho \frac{\ell}{W T}\]
where $\ell$ is the length of the strip, $W$ is the width, and $T$ is the thickness. The length is assumed to be infinite in the direction of current flow, so the resistance is calculated per unit length.
Arguments
thickness: Thickness of the strip [m].width: Width of the strip [m].rho: Electrical resistivity of the conductor material [Ω·m].alpha: Temperature coefficient of resistivity [1/°C].T0: Reference temperature for the material properties [°C].Top: Operating temperature of the conductor [°C].
Returns
- DC resistance of the strip conductor [Ω].
Examples
thickness = 0.002
width = 0.05
rho = 1.7241e-8
alpha = 0.00393
T0 = 20
T = 25
resistance = calc_strip_resistance(thickness, width, rho, alpha, T0, T)
# Output: ~0.0001758 ΩSee also
LineCableModels.DataModel.BaseParams.calc_temperature_correction — Method
calc_temperature_correction(
alpha::Union{Float64, Measurements.Measurement{Float64}},
Top::Union{Float64, Measurements.Measurement{Float64}}
) -> Union{Float64, Measurements.Measurement{Float64}}
calc_temperature_correction(
alpha::Union{Float64, Measurements.Measurement{Float64}},
Top::Union{Float64, Measurements.Measurement{Float64}},
T0::Union{Float64, Measurements.Measurement{Float64}}
) -> Union{Float64, Measurements.Measurement{Float64}}
Calculates the temperature correction factor for material properties based on the standard linear temperature model [17]:
\[k(T) = 1 + \alpha (T - T_0)\]
where $\alpha$ is the temperature coefficient of the material resistivity, $T$ is the operating temperature, and $T_0$ is the reference temperature.
Arguments
alpha: Temperature coefficient of the material property [1/°C].T: Current temperature [°C].T0: Reference temperature at which the base material property was measured [°C]. Defaults to T₀.
Returns
- Temperature correction factor to be applied to the material property [dimensionless].
Examples
# Copper resistivity correction (alpha = 0.00393 [1/°C])
k = calc_temperature_correction(0.00393, 75.0, 20.0) # Expected output: 1.2161LineCableModels.DataModel.BaseParams.calc_tubular_gmr — Method
calc_tubular_gmr(
r_ex::Union{Float64, Measurements.Measurement{Float64}},
r_in::Union{Float64, Measurements.Measurement{Float64}},
mu_r::Union{Float64, Measurements.Measurement{Float64}}
) -> Union{Float64, Measurements.Measurement{Float64}}
Calculates the geometric mean radius (GMR) of a tubular conductor, using [2]:
\[\log GMR = \log r_2 - \mu_r \left[ \frac{r_1^4}{\left(r_2^2 - r_1^2\right)^2} \log\left(\frac{r_2}{r_1}\right) - \frac{3r_1^2 - r_2^2}{4\left(r_2^2 - r_1^2\right)} \right]\]
where $\mu_r$ is the material magnetic permeability (relative to free space), $r_1$ and $r_2$ are the inner and outer radii of the tubular conductor, respectively. If $r_2$ is approximately equal to $r_1$ , the tube collapses into a thin shell, and the GMR is equal to $r_2$. If the tube becomes infinitely thick (e.g., $r_2 \gg r_1$), the GMR diverges to infinity.
Arguments
r_ex: External radius of the tubular conductor [m].r_in: Internal radius of the tubular conductor [m].mu_r: Relative permeability of the conductor material [dimensionless].
Returns
- Geometric mean radius (GMR) of the tubular conductor [m].
Errors
- Throws
ArgumentErrorifr_exis less thanr_in.
Examples
r_ex = 0.02
r_in = 0.01
mu_r = 1.0
gmr = calc_tubular_gmr(r_ex, r_in, mu_r)
println(gmr) # Expected output: ~0.0135 [m]LineCableModels.DataModel.BaseParams.calc_tubular_inductance — Method
calc_tubular_inductance(
r_in::Union{Float64, Measurements.Measurement{Float64}},
r_ex::Union{Float64, Measurements.Measurement{Float64}},
mu_r::Union{Float64, Measurements.Measurement{Float64}}
) -> Union{Float64, Measurements.Measurement{Float64}}
Calculates the inductance of a tubular conductor per unit length, disregarding skin-effects (DC approximation) [4] [17] [14]:
\[L = \frac{\mu_r \mu_0}{2 \pi} \log \left( \frac{r_{ext}}{r_{in}} \right)\]
where $\mu_r$ is the relative permeability of the conductor material, $\mu_0$ is the vacuum permeability, and $r_{in}$ and $r_{ext}$ are the inner and outer radii of the conductor, respectively.
Arguments
r_in: Internal radius of the tubular conductor [m].r_ex: External radius of the tubular conductor [m].mu_r: Relative permeability of the conductor material [dimensionless].
Returns
- Internal inductance of the tubular conductor per unit length [H/m].
Examples
r_in = 0.01
r_ex = 0.02
mu_r = 1.0
L = calc_tubular_inductance(r_in, r_ex, mu_r)
# Output: ~2.31e-7 H/mSee also
LineCableModels.DataModel.BaseParams.calc_tubular_resistance — Method
calc_tubular_resistance(
r_in::Union{Float64, Measurements.Measurement{Float64}},
r_ex::Union{Float64, Measurements.Measurement{Float64}},
rho::Union{Float64, Measurements.Measurement{Float64}},
alpha::Union{Float64, Measurements.Measurement{Float64}},
T0::Union{Float64, Measurements.Measurement{Float64}},
Top::Union{Float64, Measurements.Measurement{Float64}}
) -> Union{Float64, Measurements.Measurement{Float64}}
Calculates the DC resistance of a tubular conductor based on its geometric and material properties, using the resistivity and cross-sectional area of a hollow cylinder with radii $r_{in}$ and $r_{ext}$:
\[R = \rho \frac{\ell}{\pi (r_{ext}^2 - r_{in}^2)}\]
where $\ell$ is the length of the conductor, $r_{in}$ and $r_{ext}$ are the inner and outer radii, respectively. The length is assumed to be infinite in the direction of current flow, so the resistance is calculated per unit length.
Arguments
r_in: Internal radius of the tubular conductor [m].r_ex: External radius of the tubular conductor [m].rho: Electrical resistivity of the conductor material [Ω·m].alpha: Temperature coefficient of resistivity [1/°C].T0: Reference temperature for the material properties [°C].Top: Operating temperature of the conductor [°C].
Returns
- DC resistance of the tubular conductor [Ω].
Examples
r_in = 0.01
r_ex = 0.02
rho = 1.7241e-8
alpha = 0.00393
T0 = 20
T = 25
resistance = calc_tubular_resistance(r_in, r_ex, rho, alpha, T0, T)
# Output: ~9.10e-8 ΩSee also
Earth properties
LineCableModels.EarthProps — Module
LineCableModels.EarthPropsThe EarthProps module provides functionality for modeling and computing earth properties within the LineCableModels.jl package. This module includes definitions for homogeneous and layered earth models, and formulations for frequency-dependent earth properties, to be used in impedance/admittance calculations.
Overview
- Defines the
EarthModelobject for representing horizontally or vertically multi-layered earth models with frequency-dependent properties (ρ, ε, μ). - Provides the
EarthLayertype for representing individual soil layers with electromagnetic properties. - Implements a multi-dispatch framework to allow different formulations of frequency-dependent earth models with
AbstractFDEMFormulation. - Contains utility functions for building complex multi-layered earth models and generating data summaries.
Dependencies
BaseLineCableModels.CommonsMeasurements
Exports
LineCableModels.EarthProps.CPEarth — Type
struct CPEarth <: LineCableModels.EarthProps.AbstractFDEMFormulationRepresents an earth model with constant properties (CP), i.e. frequency-invariant electromagnetic properties.
LineCableModels.EarthProps.CPEarth — Method
Functor implementation for CPEarth.
Computes frequency-dependent earth properties using the CPEarth formulation, which assumes frequency-invariant values for resistivity, permittivity, and permeability.
Arguments
frequencies: Vector of frequency values [Hz].base_rho_g: Base (DC) electrical resistivity of the soil [Ω·m].base_epsr_g: Base (DC) relative permittivity of the soil [dimensionless].base_mur_g: Base (DC) relative permeability of the soil [dimensionless].formulation: Instance of a subtype ofAbstractFDEMFormulationdefining the computation method.
Returns
rho: Vector of resistivity values [Ω·m] at the given frequencies.epsilon: Vector of permittivity values [F/m] at the given frequencies.mu: Vector of permeability values [H/m] at the given frequencies.
Examples
frequencies = [1e3, 1e4, 1e5]
# Using the CP model
rho, epsilon, mu = CPEarth(frequencies, 100, 10, 1, CPEarth())
println(rho) # Output: [100, 100, 100]
println(epsilon) # Output: [8.854e-11, 8.854e-11, 8.854e-11]
println(mu) # Output: [1.2566e-6, 1.2566e-6, 1.2566e-6]See also
LineCableModels.EarthProps.EarthLayer — Type
struct EarthLayer{T<:Union{Float64, Measurements.Measurement{Float64}}}Represents one single earth layer in an EarthModel object, with base and frequency-dependent properties, and attributes:
base_rho_g::Union{Float64, Measurements.Measurement{Float64}}: Base (DC) electrical resistivity [Ω·m].base_epsr_g::Union{Float64, Measurements.Measurement{Float64}}: Base (DC) relative permittivity [dimensionless].base_mur_g::Union{Float64, Measurements.Measurement{Float64}}: Base (DC) relative permeability [dimensionless].t::Union{Float64, Measurements.Measurement{Float64}}: Thickness of the layer [m].rho_g::Vector{T} where T<:Union{Float64, Measurements.Measurement{Float64}}: Computed resistivity values [Ω·m] at given frequencies.eps_g::Vector{T} where T<:Union{Float64, Measurements.Measurement{Float64}}: Computed permittivity values [F/m] at given frequencies.mu_g::Vector{T} where T<:Union{Float64, Measurements.Measurement{Float64}}: Computed permeability values [H/m] at given frequencies.
LineCableModels.EarthProps.EarthLayer — Method
Constructs an EarthLayer instance with specified base and frequency-dependent properties.
LineCableModels.EarthProps.EarthLayer — Method
EarthLayer(
frequencies::Array{T<:Union{Float64, Measurements.Measurement{Float64}}, 1},
base_rho_g::Union{Float64, Measurements.Measurement{Float64}},
base_epsr_g::Union{Float64, Measurements.Measurement{Float64}},
base_mur_g::Union{Float64, Measurements.Measurement{Float64}},
t::Union{Float64, Measurements.Measurement{Float64}},
freq_dependence::LineCableModels.EarthProps.AbstractFDEMFormulation
) -> LineCableModels.EarthProps.EarthLayer
Constructs an EarthLayer instance with specified base properties and computes its frequency-dependent values.
Arguments
frequencies: Vector of frequency values [Hz].base_rho_g: Base (DC) electrical resistivity of the layer [Ω·m].base_epsr_g: Base (DC) relative permittivity of the layer [dimensionless].base_mur_g: Base (DC) relative permeability of the layer [dimensionless].t: Thickness of the layer [m].freq_dependence: Instance of a subtype ofAbstractFDEMFormulationdefining the computation method for frequency-dependent properties.
Returns
- An
EarthLayerinstance with computed frequency-dependent properties.
Examples
frequencies = [1e3, 1e4, 1e5]
layer = EarthLayer(frequencies, 100, 10, 1, 5, CPEarth())
println(layer.rho_g) # Output: [100, 100, 100]
println(layer.eps_g) # Output: [8.854e-11, 8.854e-11, 8.854e-11]
println(layer.mu_g) # Output: [1.2566e-6, 1.2566e-6, 1.2566e-6]See also
LineCableModels.EarthProps.EarthModel — Type
struct EarthModel{T<:Union{Float64, Measurements.Measurement{Float64}}}Represents a multi-layered earth model with frequency-dependent properties, and attributes:
freq_dependence::LineCableModels.EarthProps.AbstractFDEMFormulation: Selected frequency-dependent formulation for earth properties.vertical_layers::Bool: Boolean flag indicating whether the model is treated as vertically layered.layers::Array{LineCableModels.EarthProps.EarthLayer{T}, 1} where T<:Union{Float64, Measurements.Measurement{Float64}}: Vector ofEarthLayerobjects, starting with an air layer and the specified first earth layer.
LineCableModels.EarthProps.EarthModel — Method
Constructs an EarthModel instance with specified attributes.
LineCableModels.EarthProps.EarthModel — Method
EarthModel(
frequencies::Array{T<:Union{Float64, Measurements.Measurement{Float64}}, 1},
rho_g::Union{Float64, Measurements.Measurement{Float64}},
epsr_g::Union{Float64, Measurements.Measurement{Float64}},
mur_g::Union{Float64, Measurements.Measurement{Float64}};
t,
freq_dependence,
vertical_layers,
air_layer
) -> EarthModel
Constructs an EarthModel instance with a specified first earth layer. A semi-infinite air layer is always added before the first earth layer.
Arguments
frequencies: Vector of frequency values [Hz].rho_g: Base (DC) electrical resistivity of the first earth layer [Ω·m].epsr_g: Base (DC) relative permittivity of the first earth layer [dimensionless].mur_g: Base (DC) relative permeability of the first earth layer [dimensionless].t: Thickness of the first earth layer [m]. For homogeneous earth models (or the bottommost layer), sett = Inf.freq_dependence: Instance of a subtype ofAbstractFDEMFormulationdefining the computation method for frequency-dependent properties (default:CPEarth).vertical_layers: Boolean flag indicating whether the model should be treated as vertically-layered (default:false).air_layer: optionalEarthLayerobject representing the semi-infinite air layer (default:EarthLayer(frequencies, Inf, 1.0, 1.0, Inf, freq_dependence)).
Returns
- An
EarthModelinstance with the specified attributes and computed frequency-dependent properties.
Examples
frequencies = [1e3, 1e4, 1e5]
earth_model = EarthModel(frequencies, 100, 10, 1, t=Inf)
println(length(earth_model.layers)) # Output: 2 (air + top layer)
println(earth_model.rho_eff) # Output: missingSee also
Import & export
LineCableModels.ImportExport — Module
LineCableModels.ImportExportThe ImportExport module provides methods for serializing and deserializing data structures in LineCableModels.jl, and data exchange with external programs.
Overview
This module provides functionality for:
- Saving and loading cable designs and material libraries to/from JSON and other formats.
- Exporting cable system models to PSCAD and ATP formats.
- Serializing custom types with special handling for measurements and complex numbers.
The module implements a generic serialization framework with automatic type reconstruction and proper handling of Julia-specific types like Measurement objects and Inf/NaN values.
Dependencies
BaseDataFramesDatesEzXMLJSON3LineCableModels.CommonsLinearAlgebraMeasurementsPrintfSerializationTablesXLSX
Exports
LineCableModels.ImportExport.export_data — Method
export_data(
::Val{:atp},
cable_system::LineCableSystem,
earth_props::EarthModel;
base_freq,
file_name
) -> Union{Nothing, String}
Export a LineCableSystem to an ATPDraw‑compatible XML file (LCC component with input data).
This routine serializes the cable system geometry (positions and outer radii) and the already‑computed, frequency‑specific equivalent parameters of each cable component to the ATPDraw XML schema. The result is written to disk and the absolute file path is returned on success.
Arguments
::Val{:atp}: Backend selector for the ATP/ATPDraw exporter.cable_system::LineCableSystem: The system to export. Each entry incable_system.cablesprovides one phase position and its associatedCableDesign. The number of phases exported equalslength(cable_system.cables).earth_props::EarthModel: Ground model used to populate ATP soil parameters. The exporter
uses the last layer’s base resistivity as Grnd resis.
base_freq::Number = f₀[Hz]: System frequency written to ATP (SysFreq) and stored in component metadata. *This exporter does not recompute R/L/C/G; it writes the values as
present in the groups/components at the time of export.*
file_name::String = "*_export.xml": Output file name or path. If a relative path is given, it is resolved against the exporter’s source directory. The absolute path of the saved file is returned.
Behavior
Create the ATPDraw
<project>root and header and insert a single LCC component withNumPhases = length(cable_system.cables).For each
CablePositionincable_system.cables:- Write a
<cable>element with:
- Write a
* `NumCond` = number of [`CableComponent`](@ref)s in the design,
* `Rout` = outermost radius of the design (m),
* `PosX`, `PosY` = cable coordinates (m).For each
CableComponentinside a cable:- Write one
<conductor>element with fields (all per unit length):
- Write one
* `Rin`, `Rout` — from the component’s conductor group,
* `rho` — conductor equivalence via [`calc_equivalent_rho`](@ref),
* `muC` — conductor relative permeability via [`calc_equivalent_mu`](@ref),
* `muI` — insulator relative permeability (taken from the first insulating layer’s material),
* `epsI` — insulation relative permittivity via [`calc_equivalent_eps`](@ref),
* `Cext`, `Gext` — shunt capacitance and conductance from the component’s insulator group.- Soil resistivity is written as Grnd resis using
earth_props.layers[end].base_rho_g. - The XML is pretty‑printed and written to
file_name. On I/O error, the function logs an error and returnsnothing.
Units
Units are printed in the XML file according to the ATPDraw specifications:
- Radii (
Rin,Rout,Routof cable): [m] - Coordinates (
PosX,PosY): [m] - Length (
Lengthtag): [m] - Frequency (
SysFreq/Freq): [Hz] - Resistivity (
rho, Grnd resis): [Ω·m] - Relative permittivity (
epsI) / permeability (muC,muI): [dimensionless] - Shunt capacitance (
Cext): [F/m] - Shunt conductance (
Gext): [S/m]
Notes
- The exporter assumes each component’s equivalent parameters (R/G/C and derived ρ/ε/μ) were already computed by the design/group constructors at the operating conditions of interest.
- Mixed numeric types are supported; values are stringified for XML output. When using uncertainty types (e.g.,
Measurements.Measurement), the uncertainty is removed. - Overlap checks between cables are enforced when building the system, not during export.
Examples
# Build or load a system `sys` and an earth model `earth`
file = export_data(Val(:atp), sys, earth; base_freq = 50.0,
file_name = "system_id_export.xml")
println("Exported to: ", file)See also
LineCableModels.ImportExport.export_data — Method
export_data(
::Val{:atp},
line_params::LineParameters;
file_name,
cable_system
) -> Union{Nothing, String}
Export calculated LineParameters (series impedance Z and shunt admittance Y) to an compliant ZY XML file.
This routine writes the complex Z and Y matrices versus frequency into a compact XML structure understood by external tools. Rows are emitted as comma‑separated complex entries (R+Xi / G+Bi) with one <Z>/<Y> block per frequency sample.
Arguments
::Val{:atp}: Backend selector for the ATP/ATPDraw ZY exporter.line_params::LineParameters: Object holding the frequency‑dependent matricesZ[:,:,k],Y[:,:,k], andf[k]inline_params.f.file_name::String = "ZY_export.xml": Output file name or path. If relative, it is resolved against the exporter’s source directory. The absolute path of the saved file is returned.cable_system::Union{LineCableSystem,Nothing} = nothing: Optional system used only to derive a default name. When provided andfile_nameis not overridden, the exporter uses"$(cable_system.system_id)_ZY_export.xml".
Behavior
The root tag
<ZY>includesNumPhases,Length(fixed to1.0), and format attributesZFmt="R+Xi",YFmt="G+Bi".For each frequency
fᵏ = line_params.f[k]:- Emit a
<Z Freq=...>block withnum_phaseslines, each line thek‑th slice of rowiformatted asreal(Z[i,j,k]) + imag(Z[i,j,k])i. - Emit a
<Y Freq=...>block in the same fashion (defaultG+Bi).
- Emit a
Close the
</ZY>element and write to disk. On I/O error the function logs and returnsnothing.
Units
Units are printed in the XML file according to the ATPDraw specifications:
freq(XMLFreqattribute): [Hz]Zentries: [Ω/km] (per unit length)Yentries: [S/km] (per unit length) whenYFmt = "G+Bi"- XML
Lengthattribute: [m]
Notes
- The exporter assumes
size(line_params.Z, 1) == size(line_params.Z, 2) == size(line_params.Y, 1) == size(line_params.Y, 2)andlength(line_params.f) == size(Z,3) == size(Y,3). - Numeric types are stringified; mixed numeric backends (e.g., with uncertainties) are acceptable as long as they can be printed via
@sprintf. - This exporter does not modify or recompute matrices; it serializes exactly what is in
line_params.
Examples
# Z, Y, f have already been computed into `lp::LineParameters`
file = export_data(:atp, lp; file_name = "ZY_export.xml")
println("Exported ZY to: ", file)
# Naming based on a cable system
file2 = export_data(:atp, lp; cable_system = sys)
println("Exported ZY to: ", file2) # => "$(sys.system_id)_ZY_export.xml"See also
LineParametersLineCableSystemexport_data(::Val{:atp}, cable_system, ...)— exporter that writes full LCC input data
LineCableModels.ImportExport.export_data — Method
export_data(
::Val{:pscad},
cable_system::LineCableSystem,
earth_props::EarthModel;
base_freq,
file_name
) -> Union{Nothing, String}
Exports a LineCableSystem to a PSCAD-compatible file format.
Arguments
cable_system: ALineCableSystemobject representing the cable system to be exported.earth_props: AnEarthModelobject containing the earth properties.base_freq: The base frequency [Hz] used for the PSCAD export.file_name: The path to the output file (default: "*_export.pscx")
Returns
- The absolute path of the saved file, or
nothingon failure.
Examples
cable_system = LineCableSystem(...)
earth_model = EarthModel(...)
export_data(cable_system, earth_model, base_freq=50)See also
LineCableModels.ImportExport.load! — Method
load!(library::CablesLibrary; file_name) -> CablesLibrary
Loads cable designs from a file into an existing CablesLibrary object. Modifies the library in-place. The format is determined by the file extension:
.json: Loads using the custom JSON deserialization and reconstruction..jls: Loads using Julia's native binary deserialization.
Arguments
library: TheCablesLibraryinstance to populate (modified in-place).file_name: Path to the file to load (default: "cables_library.json").
Returns
- The modified
CablesLibraryinstance.
LineCableModels.ImportExport.load! — Method
load!(
library::MaterialsLibrary;
file_name
) -> MaterialsLibrary
Loads materials from a JSON file into an existing MaterialsLibrary object. Modifies the library in-place.
Arguments
library: TheMaterialsLibraryinstance to populate (modified in-place).file_name: Path to the JSON file to load (default: "materials_library.json").
Returns
- The modified
MaterialsLibraryinstance.
See also
LineCableModels.ImportExport.save — Method
save(
library::CablesLibrary;
file_name
) -> Union{Nothing, String}
Saves a CablesLibrary to a file. The format is determined by the file extension:
.json: Saves using the custom JSON serialization..jls: Saves using Julia native binary serialization.
Arguments
library: TheCablesLibraryinstance to save.file_name: The path to the output file (default: "cables_library.json").
Returns
- The absolute path of the saved file, or
nothingon failure.
LineCableModels.ImportExport.save — Method
save(
library::MaterialsLibrary;
file_name
) -> Union{Nothing, String}
Saves a MaterialsLibrary to a JSON file.
Arguments
library: TheMaterialsLibraryinstance to save.file_name: The path to the output JSON file (default: "materials_library.json").
Returns
- The absolute path of the saved file, or
nothingon failure.
Materials library
LineCableModels.Materials — Module
LineCableModels.MaterialsThe Materials module provides functionality for managing and utilizing material properties within the LineCableModels.jl package. This module includes definitions for material properties, a library for storing and retrieving materials, and functions for manipulating material data.
Overview
- Defines the
Materialstruct representing fundamental physical properties of materials. - Provides the
MaterialsLibrarymutable struct for storing a collection of materials. - Includes functions for adding, removing, and retrieving materials from the library.
- Supports loading and saving material data from/to JSON files.
- Contains utility functions for displaying material data.
Dependencies
BaseLineCableModels.CommonsMeasurements
Exports
LineCableModels.Materials.Material — Type
struct Material{T<:Union{Float64, Measurements.Measurement{Float64}}}Defines electromagnetic and thermal properties of a material used in cable modeling:
rho::Union{Float64, Measurements.Measurement{Float64}}: Electrical resistivity of the material [Ω·m].eps_r::Union{Float64, Measurements.Measurement{Float64}}: Relative permittivity [dimensionless].mu_r::Union{Float64, Measurements.Measurement{Float64}}: Relative permeability [dimensionless].T0::Union{Float64, Measurements.Measurement{Float64}}: Reference temperature for property evaluations [°C].alpha::Union{Float64, Measurements.Measurement{Float64}}: Temperature coefficient of resistivity [1/°C].
LineCableModels.Materials.Material — Method
Material(rho, eps_r, mu_r, T0, alpha) -> Material
Weakly-typed constructor that infers the target scalar type T from the arguments, coerces values to T, and calls the strict numeric kernel.
Arguments
rho: Resistivity [Ω·m].eps_r: Relative permittivity [1].mu_r: Relative permeability [1].T0: Reference temperature [°C].alpha: Temperature coefficient of resistivity [1/°C].
Returns
Material{T}whereT = resolve_T(rho, eps_r, mu_r, T0, alpha).
LineCableModels.Materials.MaterialsLibrary — Type
mutable struct MaterialsLibrary <: AbstractDict{String, Material}Stores a collection of predefined materials for cable modeling, indexed by material name:
data::Dict{String, Material}: Dictionary mapping material names toMaterialobjects.
LineCableModels.Materials.MaterialsLibrary — Method
MaterialsLibrary(; add_defaults) -> MaterialsLibrary
Constructs an empty MaterialsLibrary instance and initializes with default materials.
Arguments
- None.
Returns
- A
MaterialsLibraryobject populated with default materials.
Examples
# Create a new, empty library
library = MaterialsLibrary()See also
Utilities
LineCableModels.Utils — Module
LineCableModels.UtilsThe Utils module provides utility functions for the LineCableModels.jl package. This module includes functions for handling measurements, numerical comparisons, and other common tasks.
Overview
- Provides general constants used throughout the package.
- Includes utility functions for numerical comparisons and handling measurements.
- Contains functions to compute uncertainties and bounds for measurements.
Dependencies
BaseDatesLineCableModels.CommonsLinearAlgebraLoggingPlotsPrintfStatistics
Exports
LineCableModels.Utils.bias_to_uncertain — Method
bias_to_uncertain(
nominal::Float64,
measurements::Vector{<:Measurements.Measurement}
) -> Any
Computes the uncertainty of a measurement by incorporating systematic bias.
Arguments
nominal: The deterministic nominal value (Float64).measurements: A vector ofMeasurementvalues from theMeasurements.jlpackage.
Returns
- A new
Measurementobject representing the mean measurement value with an uncertainty that accounts for both statistical variation and systematic bias.
Notes
- Computes the mean value and its associated uncertainty from the given
measurements. - Determines the bias as the absolute difference between the deterministic
nominalvalue and the mean measurement. - The final uncertainty is the sum of the standard uncertainty (
sigma_mean) and the systematic bias.
Examples
using Measurements
nominal = 10.0
measurements = [10.2 ± 0.1, 9.8 ± 0.2, 10.1 ± 0.15]
result = bias_to_uncertain(nominal, measurements)
println(result) # Output: Measurement with adjusted uncertaintyLineCableModels.Utils.coerce_to_T — Function
Public coercion API. Converts scalars and containers to a target type T, applying element‑wise coercion recursively. Complex numbers are handled by splitting into real and imaginary parts and coercing each side independently.
Arguments
x: Input value (scalar or container) [dimensionless].::Type{T}: Target type [dimensionless].
Returns
- Value coerced to the target type:
- For
Real → Complex{P}: constructsComplex{P}(coerce_to_T(re, P), coerce_to_T(im, P))(imaginary part from0). - For
Complex → Real: discards the imaginary part and coerces the real part. - For
AbstractArray,Tuple,NamedTuple: coerces each element recursively. - For other types: defers to
_coerce_elt_to_T.
- For
Examples
using Measurements
# Scalar
coerce_to_T(1.2, Float32) # 1.2f0
coerce_to_T(1.2, Measurement{Float64}) # 1.2 ± 0.0
coerce_to_T(1 + 2im, Complex{Float32}) # 1.0f0 + 2.0f0im
coerce_to_T(1 + 2im, Float64) # 1.0
# Containers
coerce_to_T([1.0, 2.0], Measurement{Float64}) # measurement array
coerce_to_T((1.0, 2.0), Float32) # (1.0f0, 2.0f0)
coerce_to_T((; a=1.0, b=2.0), Float32) # (a = 1.0f0, b = 2.0f0)Methods
coerce_to_T(x, _)defined at typecoercion.jl:262.
coerce_to_T(x, _)defined at typecoercion.jl:266.
coerce_to_T(x, _)defined at typecoercion.jl:269.
coerce_to_T(x, _)defined at typecoercion.jl:272.
coerce_to_T(x, _)defined at typecoercion.jl:276.
coerce_to_T(x, _)defined at typecoercion.jl:279.
coerce_to_T(A, _)defined at typecoercion.jl:283.
coerce_to_T(A, _)defined at typecoercion.jl:284.
coerce_to_T(t, _)defined at typecoercion.jl:287.
coerce_to_T(t, _)defined at typecoercion.jl:288.
coerce_to_T(nt, _)defined at typecoercion.jl:292.
coerce_to_T(nt, _)defined at typecoercion.jl:294.
coerce_to_T(x, _)defined at typecoercion.jl:298.
See also
LineCableModels.Utils.is_headless — Method
is_headless() -> Bool
Determines if the current execution environment is headless (without display capability).
Returns
trueif running in a continuous integration environment or without display access.falseotherwise when a display is available.
Examples
if is_headless()
# Use non-graphical backend
gr()
else
# Use interactive backend
plotlyjs()
endLineCableModels.Utils.is_in_testset — Method
is_in_testset() -> Any
Checks if the code is running inside a @testset by checking if Test is loaded in the current session and then calling get_testset_depth().
LineCableModels.Utils.percent_error — Method
percent_error(m::Number) -> Any
Computes the percentage uncertainty of a measurement.
Arguments
m: A numerical value, expected to be of typeMeasurementfrom theMeasurements.jlpackage.
Returns
- The percentage uncertainty, computed as
100 * uncertainty(m) / value(m), ifmis aMeasurement. NaNifmis not aMeasurement.
Examples
using Measurements
m = 10.0 ± 2.0
percent_err = percent_error(m) # Output: 20.0
not_a_measurement = 5.0
percent_err_invalid = percent_error(not_a_measurement) # Output: NaNLineCableModels.Utils.percent_to_uncertain — Method
percent_to_uncertain(val, perc) -> Any
Converts a value to a measurement with uncertainty based on percentage.
Arguments
val: The nominal value.perc: The percentage uncertainty (0 to 100).
Returns
- A
Measurementtype with the given value and calculated uncertainty.
Examples
using Measurements
percent_to_uncertain(100.0, 5) # Output: 100.0 ± 5.0
percent_to_uncertain(10.0, 10) # Output: 10.0 ± 1.0LineCableModels.Utils.resolve_T — Method
resolve_T(args...) -> Type
Resolves the promotion target type to be used by constructors and coercion utilities based on the runtime arguments. The decision uses structure‑aware predicates for Measurement and Complex:
- If any argument contains
Measurementand any containsComplex, returnsComplex{Measurement{BASE_FLOAT}}. - Else if any contains
Measurement, returnsMeasurement{BASE_FLOAT}. - Else if any contains
Complex, returnsComplex{BASE_FLOAT}. - Otherwise returns
BASE_FLOAT.
Arguments
args...: Values whose types will drive the promotion decision [dimensionless].
Returns
- A
Typesuitable for numeric promotion in subsequent coercion.
Examples
using Measurements
T = resolve_T(1.0, 2.0) # BASE_FLOAT
T = resolve_T(1 + 0im, 2.0) # Complex{BASE_FLOAT}
T = resolve_T(measurement(1.0, 0.1), 2.0) # Measurement{BASE_FLOAT}
T = resolve_T(measurement(1.0, 0.1), 2 + 0im) # Complex{Measurement{BASE_FLOAT}}LineCableModels.Utils.to_certain — Method
to_certain(value) -> Any
Converts a measurement to a value with zero uncertainty, retaining the numeric type Measurement.
Arguments
value: Input value that may be aMeasurementtype or another type.
Returns
- If input is a
Measurement, returns the same value with zero uncertainty; otherwise returns the original value unchanged.
Examples
x = 5.0 ± 0.1
result = to_certain(x) # Output: 5.0 ± 0.0
y = 10.0
result = to_certain(y) # Output: 10.0LineCableModels.Utils.to_lower — Method
to_lower(m::Number) -> Any
Computes the lower bound of a measurement value.
Arguments
m: A numerical value, expected to be of typeMeasurementfrom theMeasurements.jlpackage.
Returns
- The lower bound, computed as
value(m) - uncertainty(m)ifmis aMeasurement. NaNifmis not aMeasurement.
Examples
using Measurements
m = 10.0 ± 2.0
lower = to_lower(m) # Output: 8.0
not_a_measurement = 5.0
lower_invalid = to_lower(not_a_measurement) # Output: NaNLineCableModels.Utils.to_nominal — Method
to_nominal(x::Measurements.Measurement) -> AbstractFloat
Returns the nominal (deterministic) value of inputs that may contain Measurements.Measurement numbers, recursively handling Complex and arrays.
Arguments
x: Input value which can be aMeasurementtype or any other type.
Returns
- Measurement → its
value - Complex →
complex(to_nominal(real(z)), to_nominal(imag(z))) - AbstractArray → broadcasts
to_nominalelementwise - Anything else → returned unchanged
Examples
using Measurements
to_nominal(1.0) # Output: 1.0
to_nominal(5.2 ± 0.3) # Output: 5.2LineCableModels.Utils.to_upper — Method
to_upper(m::Number) -> Any
Computes the upper bound of a measurement value.
Arguments
m: A numerical value, expected to be of typeMeasurementfrom theMeasurements.jlpackage.
Returns
- The upper bound of
m, computed asvalue(m) + uncertainty(m)ifmis aMeasurement. NaNifmis not aMeasurement.
Examples
using Measurements
m = 10.0 ± 2.0
upper = to_upper(m) # Output: 12.0
not_a_measurement = 5.0
upper_invalid = to_upper(not_a_measurement) # Output: NaNUncertain Bessels
LineCableModels.UncertainBessels — Module
LineCableModels.UncertainBesselsUncertainty-aware wrappers for Bessel functions.
UncertainBessels lifts selected functions from SpecialFunctions so they accept Measurement and Complex{Measurement} inputs. The wrapper evaluates the underlying function at the nominal complex argument and propagates uncertainty via first-order finite differences using the four partial derivatives $\frac{\partial \mathrm{Re} \, f}{\partial x}, \frac{\partial \mathrm{Re} \, f}{\partial y}, \frac{\partial \mathrm{Im} \, f}{\partial x}, \frac{\partial \mathrm{Im} \, f}{\partial y}$ with $x = \mathrm{Re}(z)$ and $y = \mathrm{Im}(z)$. No new Bessel algorithms are implemented: for plain numeric inputs, results and numerical behaviour are those of SpecialFunctions.
Numerical scaling (as defined by SpecialFunctions) is supported for the “x” variants (e.g. besselix, besselkx, besseljx, …) to improve stability for large or complex arguments. In particular, the modified functions use exponential factors to temper growth along $\mathrm{Re}(z)$ (e.g. $I_\nu$ and $K_\nu$); other scaled variants follow conventions in SpecialFunctions and DLMF guidance for complex arguments. See [19] and [20].
Overview
- Thin, uncertainty-aware wrappers around
SpecialFunctions(besselj,bessely,besseli,besselk,besselh) and their scaled counterparts (…x). - For
Complex{Measurement}inputs, uncertainty is propagated using the 4-component gradient with respect to $\mathrm{Re}(z)$ and $\mathrm{Im}(z)$. - For
Measurement(real) inputs, a 1-D finite-difference derivative is used. - No change in semantics for
Real/Complexinputs: calls delegate toSpecialFunctions.
Dependencies
BaseLineCableModels.Commons
Exports
Usage
# do not import SpecialFunctions directly
using LineCableModels.UncertainBessels
z = complex(1.0, 1.0 ± 0.5)
J0_cpl = besselj(0, z) # Complex{Measurement}
J0_nom = besselj(0, value(z)) # nominal comparison
I1 = besselix(1, z) # scaled I1 with uncertaintyNumerical notes
- Scaled modified Bessels remove large exponential factors along $\mathrm{Re}(z)$ (e.g., $I_\nu$ and $K_\nu$ are scaled by opposite signs of $|\mathrm{Re}(z)|$), improving conditioning. Scaled forms for the other families follow the definitions in
SpecialFunctionsand DLMF. - Uncertainty propagation is first order (linearization at the nominal point). Large uncertainties or strong nonlinearity may reduce accuracy.
See also
Private API
Data model
DataFrames.DataFrame — Type
DataFrame(design::CableDesign; ...) -> DataFrames.DataFrame
DataFrame(
design::CableDesign,
format::Symbol;
S,
rho_e
) -> DataFrames.DataFrame
Extracts and displays data from a CableDesign.
Arguments
design: ACableDesignobject to extract data from.format: Symbol indicating the level of detail::baseparams: Basic RLC parameters with nominal value comparison (default).:components: Component-level equivalent properties.:detailed: Individual cable part properties with layer-by-layer breakdown.
S: Separation distance between cables [m] (only used for:baseparamsformat). Default: outermost cable diameter.rho_e: Resistivity of the earth [Ω·m] (only used for:baseparamsformat). Default: 100.
Returns
- A
DataFramecontaining the requested cable data in the specified format.
Examples
# Get basic RLC parameters
data = DataFrame(design) # Default is :baseparams format
# Get component-level data
comp_data = DataFrame(design, :components)
# Get detailed part-by-part breakdown
detailed_data = DataFrame(design, :detailed)
# Specify earth parameters for core calculations
core_data = DataFrame(design, :baseparams, S=0.5, rho_e=150)See also
DataFrames.DataFrame — Method
DataFrame(library::CablesLibrary) -> DataFrames.DataFrame
Lists the cable designs in a CablesLibrary object as a DataFrame.
Arguments
library: An instance ofCablesLibrarywhose cable designs are to be displayed.
Returns
- A
DataFrameobject with the following columns:cable_id: The unique identifier for each cable design.nominal_data: A string representation of the nominal data for each cable design.components: A comma-separated string listing the components of each cable design.
Examples
library = CablesLibrary()
design1 = CableDesign("example1", nominal_data=NominalData(...), components=Dict("A"=>...))
design2 = CableDesign("example2", nominal_data=NominalData(...), components=Dict("C"=>...))
add!(library, design1)
add!(library, design2)
# Display the library as a DataFrame
df = DataFrame(library)
first(df, 5) # Show the first 5 rows of the DataFrameSee also
DataFrames.DataFrame — Method
DataFrame(system::LineCableSystem) -> DataFrames.DataFrame
Generates a summary DataFrame for cable positions and phase mappings within a LineCableSystem.
Arguments
system: ALineCableSystemobject containing the cable definitions and their configurations.
Returns
- A
DataFramecontaining:cable_id: Identifier of each cable design.horz: Horizontal coordinate of each cable [m].vert: Vertical coordinate of each cable [m].phase_mapping: Human-readable string representation mapping each cable component to its assigned phase.
Examples
df = DataFrame(cable_system)
println(df)
# Output:
# │ cable_id │ horz │ vert │ phase_mapping │
# │------------│------│-------│-------------------------│
# │ "Cable1" │ 0.0 │ -0.5 │ core: 1, sheath: 0 │
# │ "Cable2" │ 0.35 │ -1.25 │ core: 2, sheath: 0 │See also
LineCableModels.DataModel.AbstractCablePart — Type
abstract type AbstractCablePart{T}Abstract type representing a generic cable part.
LineCableModels.DataModel.AbstractConductorPart — Type
abstract type AbstractConductorPart{T} <: LineCableModels.DataModel.AbstractCablePart{T}Abstract type representing a conductive part of a cable.
Subtypes implement specific configurations:
LineCableModels.DataModel.AbstractInsulatorPart — Type
abstract type AbstractInsulatorPart{T} <: LineCableModels.DataModel.AbstractCablePart{T}Abstract type representing an insulating part of a cable.
Subtypes implement specific configurations:
LineCableModels.DataModel.AbstractStrandsLayer — Type
abstract type AbstractStrandsLayer{T} <: LineCableModels.DataModel.AbstractConductorPart{T}Abstract type representing all stranded configurations composed of grouped discrete geometric shapes.
Subtypes implement specific configurations:
LineCableModels.DataModel.RectStrandsShape — Type
struct RectStrandsShape{T<:Union{Float64, Measurements.Measurement{Float64}}, U<:Int64} <: LineCableModels.DataModel.AbstractShapeGeometryHolds the pure geometric layout for a concentric layer of rectangular/flat strands.
LineCableModels.DataModel.SectorGeometryValid — Type
struct SectorGeometryValid <: LineCableModels.Validation.RuleRule that enforces specific Sector geometry constraints (discriminant, angle limits, no overlap).
Base.delete! — Method
delete!(library::CablesLibrary, cable_id::String)
Removes a cable design from a CablesLibrary object by its ID.
Arguments
library: An instance ofCablesLibraryfrom which the cable design will be removed.cable_id: The ID of the cable design to remove.
Returns
- Nothing. Modifies the
datafield of theCablesLibraryobject in-place by removing the specified cable design if it exists.
Examples
library = CablesLibrary()
design = CableDesign("example", ...) # Initialize a CableDesign
add!(library, design)
# Remove the cable design
delete!(library, "example")
haskey(library, "example") # Returns falseSee also
Base.get — Function
get(
library::CablesLibrary,
cable_id::String
) -> Union{Nothing, CableDesign}
get(
library::CablesLibrary,
cable_id::String,
default
) -> Any
Retrieves a cable design from a CablesLibrary object by its ID.
Arguments
library: An instance ofCablesLibraryfrom which the cable design will be retrieved.cable_id: The ID of the cable design to retrieve.
Returns
- A
CableDesignobject corresponding to the givencable_idif found, otherwisenothing.
Examples
library = CablesLibrary()
design = CableDesign("example", ...) # Initialize a CableDesign
add!(library, design)
# Retrieve the cable design
retrieved_design = get(library, "cable1")
println(retrieved_design.id) # Prints "example"
# Attempt to retrieve a non-existent design
missing_design = get(library, "nonexistent_id")
println(missing_design === nothing) # Prints trueSee also
Base.show — Method
show(
io::IO,
_::MIME{Symbol("text/plain")},
component::CableComponent
)
Defines the display representation of a CableComponent object for REPL or text output.
Arguments
io: Output stream.::MIME"text/plain": MIME type for plain text output.component: TheCableComponentobject to be displayed.
Returns
- Nothing. Modifies
ioby writing text representation of the object.
Base.show — Method
show(
io::IO,
_::MIME{Symbol("text/plain")},
design::CableDesign
)
Defines the display representation of a CableDesign object for REPL or text output.
Arguments
io: Output stream.::MIME"text/plain": MIME type for plain text output.design: TheCableDesignobject to be displayed.
Returns
- Nothing. Modifies
ioby writing text representation of the object.
Base.show — Method
show(
io::IO,
_::MIME{Symbol("text/plain")},
group::Union{ConductorGroup, InsulatorGroup}
)
Defines the display representation of a ConductorGroup or InsulatorGroupobjects for REPL or text output.
Arguments
io: Output stream.::MIME"text/plain": MIME type for plain text output.group: TheConductorGrouporInsulatorGroupinstance to be displayed.
Returns
- Nothing. Modifies
ioby writing text representation of the object.
Base.show — Method
show(
io::IO,
_::MIME{Symbol("text/plain")},
part::LineCableModels.DataModel.AbstractCablePart
)
Defines the display representation of an AbstractCablePart object for REPL or text output.
Arguments
io: Output stream.::MIME"text/plain": MIME type for plain text output.part: TheAbstractCablePartinstance to be displayed.
Returns
- Nothing. Modifies
ioby writing text representation of the object.
LineCableModels.Commons.add! — Method
Stores a cable design in a CablesLibrary object.
Arguments
library: An instance ofCablesLibraryto which the cable design will be added.design: ACableDesignobject representing the cable design to be added. This object must have acable_idfield to uniquely identify it.
Returns
- None. Modifies the
datafield of theCablesLibraryobject in-place by adding the new cable design.
Examples
library = CablesLibrary()
design = CableDesign("example", ...) # Initialize CableDesign with required fields
add!(library, design)
println(library) # Prints the updated dictionary containing the new cable designSee also
LineCableModels.Commons.add! — Method
add!(
group::ConductorGroup{T},
part_type::Type{C<:LineCableModels.DataModel.AbstractConductorPart},
args...;
kwargs...
) -> ConductorGroup{Tg} where Tg
Add a new conductor part to a ConductorGroup, validating raw inputs, normalizing proxies, and promoting the group’s numeric type if required.
Behavior:
- Apply part-level keyword defaults.
- Default
r_intogroup.r_exif absent. - Compute
Tnew = resolve_T(group, r_in, args..., values(kwargs)...). - If
Tnew === T, mutate in place; elsecoerce_to_T(group, Tnew)then mutate and return the promoted group.
Arguments
group:ConductorGroupobject to which the new part will be added.part_type: Type of conductor part to add (AbstractConductorPart).args...: Positional arguments specific to the constructor of thepart_type(AbstractConductorPart) [various].kwargs...: Named arguments for the constructor including optional values specific to the constructor of thepart_type(AbstractConductorPart) [various].
Returns
- The function modifies the
ConductorGroupinstance in place and does not return a value.
Notes
- Updates
gmr,resistance,alpha,r_ex,cross_section, andnum_wiresto account for the new part. - The
temperatureof the new part defaults to the temperature of the first layer if not specified. - The
r_inof the new part defaults to the external radius of the existing conductor if not specified.
- When an
AbstractCablePartis provided asr_in, the constructor retrieves itsr_exvalue, allowing the new cable part to be placed directly over the existing part in a layered cable design. - In case of uncertain measurements, if the added cable part is of a different type than the existing one, the uncertainty is removed from the radius value before being passed to the new component. This ensures that measurement uncertainties do not inappropriately cascade across different cable parts.
Examples
material_props = Material(1.7241e-8, 1.0, 0.999994, 20.0, 0.00393)
conductor = ConductorGroup(Strip(0.01, 0.002, 0.05, 10, material_props))
add!(conductor, CircStrands, 0.02, 0.002, 7, 15, material_props, temperature = 25)See also
LineCableModels.Commons.add! — Method
add!(
group::InsulatorGroup{T},
part_type::Type{C<:LineCableModels.DataModel.AbstractInsulatorPart},
args...;
f,
kwargs...
) -> InsulatorGroup{Tg} where Tg
Adds a new part to an existing InsulatorGroup object and updates its equivalent electrical parameters.
Behavior:
- Apply part-level keyword defaults (from
Validation.keyword_defaults). - Default
r_intogroup.r_exif absent. - Compute
Tnew = resolve_T(group, r_in, args..., values(kwargs)..., f). - If
Tnew === T, mutate in place; elsecoerce_to_T(group, Tnew)then mutate and return the promoted group.
Arguments
group:InsulatorGroupobject to which the new part will be added.part_type: Type of insulator part to add (AbstractInsulatorPart).args...: Positional arguments specific to the constructor of thepart_type(AbstractInsulatorPart) [various].kwargs...: Named arguments for the constructor including optional values specific to the constructor of thepart_type(AbstractInsulatorPart) [various].
Returns
- The function modifies the
InsulatorGroupinstance in place and does not return a value.
Notes
- Updates
shunt_capacitance,shunt_conductance,r_ex, andcross_sectionto account for the new part. - The
r_inof the new part defaults to the external radius of the existing insulator group if not specified.
- When an
AbstractCablePartis provided asr_in, the constructor retrieves itsr_exvalue, allowing the new cable part to be placed directly over the existing part in a layered cable design. - In case of uncertain measurements, if the added cable part is of a different type than the existing one, the uncertainty is removed from the radius value before being passed to the new component. This ensures that measurement uncertainties do not inappropriately cascade across different cable parts.
Examples
material_props = Material(1e10, 3.0, 1.0, 20.0, 0.0)
insulator_group = InsulatorGroup(Insulator(0.01, 0.015, material_props))
add!(insulator_group, Semicon, 0.015, 0.018, material_props)See also
LineCableModels.Commons.add! — Method
add!(
system::LineCableSystem{T},
cable::CableDesign,
horz::Number,
vert::Number
) -> LineCableSystem{T} where T
add!(
system::LineCableSystem{T},
cable::CableDesign,
horz::Number,
vert::Number,
conn::Union{Nothing, Dict{String, Int64}}
) -> LineCableSystem{T} where T
Convenience add! that accepts a cable design and coordinates (and optional mapping). Builds a CablePosition and forwards to add!(system, pos).
LineCableModels.Commons.add! — Method
add!(
system::LineCableSystem{T},
pos::CablePosition
) -> LineCableSystem{T} where T
Adds a new cable position to an existing LineCableSystem, updating its phase mapping and cable count. If adding the position introduces a different numeric scalar type, the system is promoted and the promoted system is returned. Otherwise, mutation happens in place.
Arguments
system: Instance ofLineCableSystemto which the cable will be added.cable: ACableDesignobject defining the cable structure.horz: Horizontal coordinate [m].vert: Vertical coordinate [m].conn: Dictionary mapping component names to phase indices, ornothingfor automatic assignment.
Returns
- The modified
LineCableSystemobject with the new cable added.
Examples
cable_design = CableDesign("example", nominal_data, components_dict)
# Define coordinates for two cables
xa, ya = 0.0, -1.0
xb, yb = 1.0, -2.0
# Create initial system with one cable
cablepos1 = CablePosition(cable_design, xa, ya, Dict("core" => 1))
cable_system = LineCableSystem("test_case_1", 1000.0, cablepos1)
# Add second cable to system
add!(cable_system, cable_design, xb, yb, Dict("core" => 2))
println(cable_system.num_cables) # Prints: 2See also
LineCableModels.DataModel._coerced_args — Method
_coerced_args(_::Type{C}, ntv, Tp, order::Tuple) -> Tuple
Builds the positional argument tuple to feed the typed core constructor, coercing only the fields returned by coercive_fields(C) to type Tp. Non‑coercive fields (e.g., integer flags) are passed through unchanged. Field order is controlled by order (a tuple of symbols), typically (required_fields(C)..., keyword_fields(C)...).
Arguments
::Type{C}: Component type [dimensionless].ntv: NormalizedNamedTuplereturned byvalidate![dimensionless].Tp: Target element type for numeric coercion [dimensionless].order::Tuple: Field order used to assemble the positional tuple [dimensionless].
Returns
- A
Tupleof arguments in the requested order, with coercions applied where configured.
Examples
args = _coerced_args(Tubular, ntv, Float64, (:r_in, :r_ex, :material_props, :temperature))LineCableModels.DataModel._ctor_materialize — Method
_ctor_materialize(mod, x) -> Any
Utility for the constructor macro to materialize input tuples from either:
- A tuple literal expression (e.g.,
(:a, :b, :c)), or - A bound constant tuple name (e.g.,
_REQ_TUBULAR).
Used to keep macro call sites short while allowing both styles.
Arguments
mod: Module where constants are resolved [dimensionless].x: Expression or symbol representing a tuple [dimensionless].
Returns
- A standard Julia
Tuple(of symbols or defaults).
Errors
ErrorExceptionifxis neither a tuple literal nor a bound constant name.
Examples
syms = _ctor_materialize(@__MODULE__, :( :a, :b ))
syms = _ctor_materialize(@__MODULE__, :_REQ_TUBULAR)LineCableModels.DataModel._do_add! — Method
_do_add!(
group::ConductorGroup{Tg},
C::Type{<:LineCableModels.DataModel.AbstractConductorPart},
args...;
kwargs...
)
Internal, in-place insertion (no promotion logic). Assumes :r_in was materialized. Runs Validation → parsing, then coerces fields to the group’s T and updates equivalent properties and book-keeping.
LineCableModels.DataModel._do_add! — Method
_do_add!(
group::InsulatorGroup{Tg},
C::Type{<:LineCableModels.DataModel.AbstractInsulatorPart},
args...;
f,
kwargs...
)
Do the actual insertion for InsulatorGroup with the group already at the correct scalar type. Validates/parses the part, coerces to the group’s T, constructs the strict numeric core, and updates geometry and admittances at the provided frequency.
Returns the mutated group (same object).
LineCableModels.DataModel._extract_part_properties — Method
_extract_part_properties(part, properties) -> Any
Helper function to extract properties from a part for detailed format.
Arguments
part: An instance ofAbstractCablePartfrom which to extract properties.properties: A vector of symbols indicating which properties to extract (not used in the current implementation).
Returns
- A vector containing the extracted properties in the following order:
type: The lowercase string representation of the part's type.r_in: The inner radius of the part, if it exists, otherwisemissing.r_ex: The outer radius of the part, if it exists, otherwisemissing.diameter_in: The inner diameter of the part (2 * rin), if `rinexists, otherwisemissing`.diameter_ext: The outer diameter of the part (2 * rex), if `rexexists, otherwisemissing`.thickness: The difference betweenr_exandr_in, if both exist, otherwisemissing.cross_section: The cross-sectional area of the part, if it exists, otherwisemissing.num_wires: The number of wires in the part, if it exists, otherwisemissing.resistance: The resistance of the part, if it exists, otherwisemissing.alpha: The temperature coefficient of resistivity of the part or its material, if it exists, otherwisemissing.gmr: The geometric mean radius of the part, if it exists, otherwisemissing.gmr_ratio: The ratio ofgmrtor_ex, if both exist, otherwisemissing.shunt_capacitance: The shunt capacitance of the part, if it exists, otherwisemissing.shunt_conductance: The shunt conductance of the part, if it exists, otherwisemissing.
Notes
This function is used to create a standardized format for displaying detailed information about cable parts.
Examples
part = Conductor(...)
properties = [:r_in, :r_ex, :resistance] # Example of properties to extract
extracted_properties = _extract_part_properties(part, properties)
println(extracted_properties)LineCableModels.DataModel._normalize_radii — Method
_normalize_radii(_::Type{T}, rin, rex) -> Any
Resolves radius parameters for cable components, converting from various input formats to standardized inner radius, outer radius, and thickness values.
This function serves as a high-level interface to the radius resolution system. It processes inputs through a two-stage pipeline:
- First normalizes input parameters to consistent forms using
_parse_radius_operand. - Then delegates to specialized implementations via
_do_resolve_radiusbased on the component type.
Arguments
param_in: Inner boundary parameter (defaults to radius) [m]. Can be a number, aDiameter, aThickness, or anAbstractCablePart.param_ext: Outer boundary parameter (defaults to radius) [m]. Can be a number, aDiameter, aThickness, or anAbstractCablePart.object_type: Type associated to the constructor of the newAbstractCablePart.
Returns
r_in: Normalized inner radius [m].r_ex: Normalized outer radius [m].thickness: Computed thickness or specialized dimension depending on the method [m]. ForCircStrandscomponents, this value represents the wire radius instead of thickness.
See also
LineCableModels.DataModel._parse_radius_operand — Function
Parses input values into radius representation based on object type and input type.
Arguments
x: Input value that can be a raw number, aDiameter, aThickness, or other convertible type [m].object_type: Type parameter used for dispatch.
Returns
- Parsed radius value in appropriate units [m].
Examples
radius = _parse_radius_operand(10.0, ...) # Direct radius value
radius = _parse_radius_operand(Diameter(20.0), ...) # From diameter object
radius = _parse_radius_operand(Thickness(5.0), ...) # From thickness objectMethods
_parse_radius_operand(x, _)defined at radii.jl:68.
_parse_radius_operand(d, _)defined at radii.jl:69.
_parse_radius_operand(p, _)defined at radii.jl:70.
_parse_radius_operand(p, _)defined at radii.jl:71.
_parse_radius_operand(x, _)defined at radii.jl:75.
_parse_radius_operand(x, _)defined at radii.jl:81.
See also
LineCableModels.DataModel._print_fields — Method
_print_fields(
io::IO,
obj,
fields_to_show::Vector{Symbol};
sigdigits
) -> Int64
Print the specified fields of an object in a compact format.
Arguments
io: The output stream.obj: The object whose fields will be displayed.fields_to_show: Vector of field names (as Symbols) to display.sigdigits: Number of significant digits for rounding numeric values.
Returns
- Number of fields that were actually displayed.
LineCableModels.DataModel._promotion_T — Method
_promotion_T(_::Type{C}, ntv, _order::Tuple) -> Type
Determines the promoted numeric element type for convenience constructors of component C. The promotion is computed across the values of coercive_fields(C), extracted from the normalized NamedTuple ntv produced by validate!. This ensures all numeric fields that participate in calculations share a common element type (e.g., Float64, Measurement{Float64}).
Arguments
::Type{C}: Component type [dimensionless].ntv: NormalizedNamedTuplereturned byvalidate![dimensionless]._order::Tuple: Ignored by this method; present for arity symmetry with_coerced_args[dimensionless].
Returns
- The promoted numeric element type [dimensionless].
Examples
Tp = _promotion_T(Tubular, (r_in=0.01, r_ex=0.02, material_props=mat, temperature=20.0), ())LineCableModels.DataModel._with_kwdefaults — Method
_with_kwdefaults(
_::Type{C},
kwargs::NamedTuple
) -> NamedTuple
Merge per-part keyword defaults declared via Validation.keyword_defaults with user-provided kwargs and return a NamedTuple suitable for forwarding.
Defaults may be a NamedTuple or a Tuple zipped against Validation.keyword_fields(::Type{C}). User keys always win.
LineCableModels.DataModel.get_material_color_makie — Method
get_material_color_makie(material_props; mu_scale=1.0, eps_scale=1.0)Piecewise ρ→base color (metals→silver, semiconductors→amber, etc.) with blue/purple magnetic overlay (μr) and teal/cyan permittivity overlay (εr). mu_scale and eps_scale scale overlay strength (1.0 = default).
LineCableModels.DataModel.nonsensify — Method
nonsensify(originaldesign::CableDesign; newid::String="")::CableDesign
Recreates a cable design by bulldozing reality into a "simplified" shape with only the so-called "main" material properties.
Translation: if you wanted physics, you came to the wrong neighborhood.
For each component, this abomination does:
ConductorGroup(Tubular(...))with radii stolen from the first and last conductor layers, and material blindly copied from the first conductor layer. Because high-fidelity is for losers.InsulatorGroup(Insulator(...))spanning from the new conductor outer radius to the original insulator group's outer radius; material is taken from the firstInsulatorlayer available (or whatever warm body it can find).
⚠ WARNING: This is deliberately nonsensical. It laughs in the face of proper equivalent property corrections and just slaps the "main" props on like duct tape. Use only when you don’t give a damn about accuracy and just want something that looks cable-ish, e.g., never.
LineCableModels.DataModel.vdeparse — Method
vdeparse(code::AbstractString) -> Dict{Symbol,String}Parses VDE/DIN 0271/0276 cable codes:
- stub (first non-space token): designation → conductormaterial (default copper) → insulation (default paper) → screen → waterblocking → innersheath → armouring → sheath → grounding
- tail: cores × cross-section (optional screen csa), voltage, conductor type (R/S/O + E/M/H, optional
/V⇒ compact)
Only parsed keys are returned; defaults are materialized when omitted.
LineCableModels.Utils.coerce_to_T — Method
Identity: no allocation when already at T.
LineCableModels.Validation.is_radius_input — Method
is_radius_input(
_::Type{T},
_::Val{:r_ex},
_::LineCableModels.DataModel.AbstractCablePart
) -> Bool
Default policy for outer radius raw inputs (annular shells): reject AbstractCablePart proxies. Outer radius must be numeric or a Thickness wrapper to avoid creating zero‑thickness layers.
Arguments
::Type{T}: Component type [dimensionless].::Val{:r_ex}: Field tag for the outer radius [dimensionless].::AbstractCablePart: Proxy object [dimensionless].
Returns
falsealways.
Examples
Validation.is_radius_input(Tubular, Val(:r_ex), prev_layer) # falseLineCableModels.Validation.is_radius_input — Method
is_radius_input(
_::Type{T},
_::Val{:r_ex},
_::Thickness
) -> Bool
Default policy for outer radius raw inputs (annular shells): accept Thickness as a convenience wrapper. The thickness is expanded to an outer radius during parsing.
Arguments
::Type{T}: Component type [dimensionless].::Val{:r_ex}: Field tag for the outer radius [dimensionless].::Thickness: Thickness wrapper [dimensionless].
Returns
Boolindicating acceptance (true).
Examples
Validation.is_radius_input(Tubular, Val(:r_ex), Thickness(1e-3)) # trueLineCableModels.Validation.is_radius_input — Method
is_radius_input(
_::Type{T},
_::Val{:r_in},
p::LineCableModels.DataModel.AbstractCablePart
) -> Bool
Default policy for inner radius raw inputs: accept proxies that expose an outer radius. This permits stacking by hijacking p.r_ex during parsing.
Arguments
::Type{T}: Component type [dimensionless].::Val{:r_in}: Field tag for the inner radius [dimensionless].p::AbstractCablePart: Proxy object [dimensionless].
Returns
Boolindicating acceptance (trueifhasproperty(p, :r_ex)).
Examples
Validation.is_radius_input(Tubular, Val(:r_in), prev_layer) # true if prev_layer has :r_exLineCableModels.Validation.maxfill — Method
maxfill(_::Type{T}, args...) -> Any
Fallback method for the maxfill interface. Throws an explicit error indicating that the component type T has not implemented its physical capacity limit.
Arguments
::Type{T}: Component type [dimensionless].args...: Geometric parameters required for the calculation [dimensionless].
Returns
- Nothing. Always throws an
ArgumentError.
LineCableModels.DataModel.@construct — Macro
Generates a weakly‑typed convenience constructor for a component T. The generated method:
- Accepts exactly the positional fields listed in
REQ. - Accepts keyword arguments listed in
OPTwith defaultsDEFS. - Calls
validate!(T, ...)forwarding variables (not defaults), - Computes the promotion type via
_promotion_T(T, ntv, order), - Coerces only
coercive_fields(T)via_coerced_args(T, ntv, Tp, order), - Delegates to the numeric core
T(...)with the coerced positional tuple.
REQ, OPT, and DEFS can be provided as tuple literals or as names of bound constant tuples. order is implicitly (REQ..., OPT...).
Arguments
T: Component type (bare name) [dimensionless].REQ: Tuple of required positional field names [dimensionless].OPT: Tuple of optional keyword field names [dimensionless]. Defaults to().DEFS: Tuple of default values matchingOPT[dimensionless]. Defaults to().
Returns
- A method definition for the weakly‑typed constructor.
Examples
const _REQ_TUBULAR = (:r_in, :r_ex, :material_props)
const _OPT_TUBULAR = (:temperature,)
const _DEFS_TUBULAR = (T₀,)
@construct Tubular _REQ_TUBULAR _OPT_TUBULAR _DEFS_TUBULAR
# Expands roughly to:
# function Tubular(r_in, r_ex, material_props; temperature=T₀)
# ntv = validate!(Tubular, r_in, r_ex, material_props; temperature=temperature)
# Tp = _promotion_T(Tubular, ntv, (:r_in, :r_ex, :material_props, :temperature))
# args = _coerced_args(Tubular, ntv, Tp, (:r_in, :r_ex, :material_props, :temperature))
# return Tubular(args...)
# endNotes
- Defaults supplied in
DEFSare escaped into the method signature (evaluated at macro expansion time). - Forwarding into
validate!always uses variables (e.g.,temperature=temperature), never literal defaults. - The macro is hygiene‑aware; identifiers
validate!,_promotion_T,_coerced_args, and the type name are properly escaped.
Errors
ErrorExceptioniflength(OPT) != length(DEFS).
Earth properties
DataFrames.DataFrame — Method
DataFrame(earth_model::EarthModel) -> DataFrames.DataFrame
Generates a DataFrame summarizing basic properties of earth layers from an EarthModel.
Arguments
earth_model: Instance ofEarthModelcontaining earth layers.
Returns
- A
DataFramewith columns:rho_g: Base (DC) resistivity of each layer [Ω·m].epsr_g: Base (DC) relative permittivity of each layer [dimensionless].mur_g: Base (DC) relative permeability of each layer [dimensionless].thickness: Thickness of each layer [m].
Examples
df = DataFrame(earth_model)
println(df)LineCableModels.EarthProps.AbstractFDEMFormulation — Type
abstract type AbstractFDEMFormulationAbstract type representing different frequency-dependent earth models (FDEM). Used in the multi-dispatch implementations in modules LineCableModels.EarthProps and LineCableModels.Engine.
Currently available formulations
CPEarth: Constant properties (CP) model.
Base.show — Method
show(
io::IO,
_::MIME{Symbol("text/plain")},
model::EarthModel
)
Defines the display representation of a EarthModel object for REPL or text output.
Arguments
io: The output stream to write the representation to [IO].mime: The MIME type for plain text output [MIME"text/plain"].model: TheEarthModelinstance to be displayed.
Returns
- Nothing. Modifies
ioto format the output.
LineCableModels.Commons.add! — Method
add!(
model::EarthModel{T<:Union{Float64, Measurements.Measurement{Float64}}},
frequencies::Array{T<:Union{Float64, Measurements.Measurement{Float64}}, 1},
base_rho_g::Union{Float64, Measurements.Measurement{Float64}},
base_epsr_g::Union{Float64, Measurements.Measurement{Float64}},
base_mur_g::Union{Float64, Measurements.Measurement{Float64}};
t
) -> EarthModel
Adds a new earth layer to an existing EarthModel.
Arguments
model: Instance ofEarthModelto which the new layer will be added.frequencies: Vector of frequency values [Hz].base_rho_g: Base electrical resistivity of the new earth layer [Ω·m].base_epsr_g: Base relative permittivity of the new earth layer [dimensionless].base_mur_g: Base relative permeability of the new earth layer [dimensionless].t: Thickness of the new earth layer [m] (default:Inf).
Returns
- Modifies
modelin place by appending a newEarthLayer.
Notes
For horizontal layering (vertical_layers = false):
- Layer 1 (air) is always infinite (
t = Inf). - Layer 2 (first earth layer) can be infinite if modeling a homogeneous half-space.
- If adding a third layer (
length(EarthModel.layers) == 3), it can be infinite only if the previous layer is finite. - No two successive earth layers (
length(EarthModel.layers) > 2) can have infinite thickness.
For vertical layering (vertical_layers = true):
- Layer 1 (air) is always horizontal and infinite at
z > 0. - Layer 2 (first vertical layer) is always infinite in
z < 0andy < 0. The first vertical layer is assumed to always end aty = 0. - Layer 3 (second vertical layer) can be infinite (establishing a vertical interface at
y = 0). - Subsequent layers can be infinite only if the previous is finite.
- No two successive vertical layers (
length(EarthModel.layers) > 3) can both be infinite.
Examples
frequencies = [1e3, 1e4, 1e5]
# Define a horizontal model with finite thickness for the first earth layer
horz_earth_model = EarthModel(frequencies, 100, 10, 1, t=5)
# Add a second horizontal earth layer
add!(horz_earth_model, frequencies, 200, 15, 1, t=10)
println(length(horz_earth_model.layers)) # Output: 3
# The bottom layer should be set to infinite thickness
add!(horz_earth_model, frequencies, 300, 15, 1, t=Inf)
println(length(horz_earth_model.layers)) # Output: 4
# Initialize a vertical-layered model with first interface at y = 0.
vert_earth_model = EarthModel(frequencies, 100, 10, 1, t=Inf, vertical_layers=true)
# Add a second vertical layer at y = 0 (this can also be infinite)
add!(vert_earth_model, frequencies, 150, 12, 1, t=Inf)
println(length(vert_earth_model.layers)) # Output: 3
# Attempt to add a third infinite layer (invalid case)
try
add!(vert_earth_model, frequencies, 120, 12, 1, t=Inf)
catch e
println(e) # Error: Cannot add consecutive vertical layers with infinite thickness.
end
# Fix: Set a finite thickness to the currently rightmost layer
vert_earth_model.layers[end].t = 3
# Add the third layer with infinite thickness now
add!(vert_earth_model, frequencies, 120, 12, 1, t=Inf)
println(length(vert_earth_model.layers)) # Output: 4See also
LineCableModels.Utils.coerce_to_T — Method
coerce_to_T(model::EarthModel, _::Type{T}) -> EarthModel
Converts an EarthModel{S} to EarthModel{T} by reconstructing the model with all layers coerced to the new scalar type T. Layer conversion is delegated to coerce_to_T(::EarthLayer, ::Type), and non-numeric metadata are forwarded unchanged.
Arguments
model: Source Earth model [dimensionless].::Type{T}: Target element type for numeric fields [dimensionless].
Returns
EarthModel{T}rebuilt with each layer and numeric payload converted toT.
Examples
m64 = coerce_to_T(model, Float64)
mM = coerce_to_T(model, Measurement{Float64})See also
LineCableModels.Utils.coerce_to_T — Method
coerce_to_T(
layer::LineCableModels.EarthProps.EarthLayer,
_::Type{T}
) -> LineCableModels.EarthProps.EarthLayer
Converts an EarthLayer{S} to EarthLayer{T} by coercing each stored field to the target element type T and rebuilding the layer via its inner constructor. Scalar and array fields are converted using the generic coerce_to_T machinery.
Arguments
layer: Source Earth layer [dimensionless].::Type{T}: Target element type for numeric fields [dimensionless].
Returns
EarthLayer{T}with all numeric state converted toT.
Examples
ℓ64 = coerce_to_T(layer, Float64)
ℓM = coerce_to_T(layer, Measurement{Float64})See also
Materials library
DataFrames.DataFrame — Method
DataFrame(library::MaterialsLibrary) -> DataFrames.DataFrame
Lists the contents of a MaterialsLibrary as a DataFrame.
Arguments
library: Instance ofMaterialsLibraryto be displayed.
Returns
- A
DataFramecontaining the material properties.
Examples
library = MaterialsLibrary()
df = DataFrame(library)See also
Base.delete! — Method
delete!(library::MaterialsLibrary, name::String)
Removes a material from a MaterialsLibrary.
Arguments
library: Instance ofMaterialsLibraryfrom which the material will be removed.name: Name of the material to be removed.
Returns
- The modified instance of
MaterialsLibrarywithout the specified material.
Errors
Throws an error if the material does not exist in the library.
Examples
library = MaterialsLibrary()
delete!(library, "copper")See also
Base.get — Function
get(
library::MaterialsLibrary,
name::String
) -> Union{Nothing, Material}
get(library::MaterialsLibrary, name::String, default) -> Any
Retrieves a material from a MaterialsLibrary by name.
Arguments
library: Instance ofMaterialsLibrarycontaining the materials.name: Name of the material to retrieve.
Returns
- The requested
Materialif found, otherwisenothing.
Examples
library = MaterialsLibrary()
material = get(library, "copper")See also
Base.show — Method
show(
io::IO,
_::MIME{Symbol("text/plain")},
dict::Dict{String, Material}
)
Defines the display representation of a MaterialsLibrary object for REPL or text output.
Arguments
io: Output stream.::MIME"text/plain": MIME type for plain text output.dict: TheMaterialsLibrarycontents to be displayed.
Returns
- Nothing. Modifies
ioby writing text representation of the library.
Base.show — Method
show(
io::IO,
_::MIME{Symbol("text/plain")},
library::MaterialsLibrary
)
Defines the display representation of a MaterialsLibrary object for REPL or text output.
Arguments
io: Output stream.::MIME"text/plain": MIME type for plain text output.library: TheMaterialsLibraryinstance to be displayed.
Returns
- Nothing. Modifies
ioby writing text representation of the library.
Base.show — Method
show(
io::IO,
_::MIME{Symbol("text/plain")},
material::Material
)
Defines the display representation of a Material object for REPL or text output.
Arguments
io: Output stream.::MIME"text/plain": MIME type for plain text output.material: TheMaterialinstance to be displayed.
Returns
- Nothing. Modifies
ioby writing text representation of the material.
LineCableModels.Commons.add! — Method
add!(
library::MaterialsLibrary,
name::AbstractString,
material::Material
) -> MaterialsLibrary
Adds a new material to a MaterialsLibrary.
Arguments
library: Instance ofMaterialsLibrarywhere the material will be added.name: Name of the material.material: Instance ofMaterialcontaining its properties.
Returns
- The modified instance of
MaterialsLibrarywith the new material added.
Errors
Throws an error if a material with the same name already exists in the library.
Examples
library = MaterialsLibrary()
material = Material(1.7241e-8, 1.0, 0.999994, 20.0, 0.00393)
add!(library, "copper", material)LineCableModels.Materials._add_default_materials! — Method
_add_default_materials!(
library::MaterialsLibrary
) -> MaterialsLibrary
Populates a MaterialsLibrary with commonly used materials, assigning predefined electrical and thermal properties.
Arguments
library: Instance ofMaterialsLibraryto be populated.
Returns
- The modified instance of
MaterialsLibrarycontaining the predefined materials.
Examples
library = MaterialsLibrary()
_add_default_materials!(library)See also
Utilities
LineCableModels.Utils._coerce_elt_to_T — Function
Element‑wise coercion kernel. Converts a single leaf value to the target type T while preserving semantics for Measurement, numeric types, and sentinels.
Arguments
x: Input leaf value [dimensionless].::Type{T}: Target type [dimensionless].
Returns
- Value coerced to the target, according to the rules below.
Notes
Number → R<:AbstractFloat: usesconvert(R, x).Number → M<:Measurement: embeds the number as a zero‑uncertainty measurement (i.e.,zero(M) + x).Measurement → M<:Measurement: recreates with the target inner type (value and uncertainty cast to_meas_inner(M)).Measurement → R<:AbstractFloat: drops uncertainty and converts the nominal value.nothingandmissingpass through unchanged.Bool,Symbol,String,Function,DataType: passed through unchanged for measurement/real targets.- Fallback: returns
xunchanged.
Examples
using Measurements
_coerce_elt_to_T(1.2, Float32) # 1.2f0
_coerce_elt_to_T(1.2, Measurement{Float64}) # 1.2 ± 0.0
_coerce_elt_to_T(measurement(2.0, 0.1), Float32) # 2.0f0
_coerce_elt_to_T(measurement(2.0, 0.1), Measurement{Float32}) # 2.0 ± 0.1 (Float32 inner)
_coerce_elt_to_T(missing, Float64) # missingMethods
_coerce_elt_to_T(x, _)defined at typecoercion.jl:205.
_coerce_elt_to_T(x, _)defined at typecoercion.jl:206.
_coerce_elt_to_T(m, _)defined at typecoercion.jl:207.
_coerce_elt_to_T(m, _)defined at typecoercion.jl:209.
_coerce_elt_to_T(_, _)defined at typecoercion.jl:210.
_coerce_elt_to_T(_, _)defined at typecoercion.jl:211.
_coerce_elt_to_T(x, _)defined at typecoercion.jl:212.
_coerce_elt_to_T(x, _)defined at typecoercion.jl:213.
_coerce_elt_to_T(x, _)defined at typecoercion.jl:214.
_coerce_elt_to_T(x, _)defined at typecoercion.jl:215.
LineCableModels.Utils._hascomplex_type — Function
Determines whether a Type contains or is a Complex number type somewhere in its structure. The check is recursive over arrays, tuples (including variadic tuples), named tuples, and union types. For concrete struct types, the predicate descends into field types. Guards are present to avoid infinite recursion through known self‑referential types (e.g., Measurements.Measurement).
Arguments
::Type: Type to inspect [dimensionless].
Returns
Boolindicating whether aComplextype occurs anywhere within the type structure.
Notes
- For
AbstractArray{S}, only the element typeSis inspected. - For
TupleandNamedTuple{N,T}, the parameters are traversed. - For
Union, both branches are inspected. - Concrete
Measurementtypes are treated as terminal and are not descended.
Examples
_hascomplex_type(Float64) # false
_hascomplex_type(Complex{Float64}) # true
_hascomplex_type(Vector{ComplexF64}) # true
_hascomplex_type(Tuple{Int, ComplexF64}) # trueMethods
_hascomplex_type(_)defined at typecoercion.jl:81.
_hascomplex_type(_)defined at typecoercion.jl:82.
_hascomplex_type(_)defined at typecoercion.jl:83.
_hascomplex_type(_)defined at typecoercion.jl:84.
_hascomplex_type(_)defined at typecoercion.jl:86.
_hascomplex_type(T)defined at typecoercion.jl:87.
_hascomplex_type(T)defined at typecoercion.jl:88.
_hascomplex_type(_)defined at typecoercion.jl:94.
LineCableModels.Utils._hasmeas_type — Method
_hasmeas_type(_::Type{<:Measurements.Measurement}) -> Bool
Determines whether a Type contains or is a Measurements.Measurement somewhere in its structure. The check is recursive over arrays, tuples (including variadic tuples), named tuples, and union types. For concrete struct types, the predicate descends into field types. Guards are present to avoid infinite recursion through known self‑contained types (e.g., Complex).
Arguments
::Type: Type to inspect [dimensionless].
Returns
Boolindicating whether aMeasurementoccurs anywhere within the type structure.
Notes
- For
AbstractArray{S}, only the element typeSis inspected. - For
TupleandNamedTuple{N,T}, the parameters are traversed. - For
Union, both branches are inspected. - Concrete
Complextypes are treated as terminal and are not descended.
Examples
using Measurements
_hasmeas_type(Float64) # false
_hasmeas_type(Measurement{Float64}) # true
_hasmeas_type(Vector{Measurement{Float64}}) # true
_hasmeas_type(Tuple{Int, Float64}) # false
_hasmeas_type(Union{Int, Measurement{Float64}}) # trueLineCableModels.Utils._meas_inner — Method
_meas_inner(_::Type{Measurements.Measurement{S}}) -> Any
Extracts the real inner type S from Measurement{S}.
Arguments
::Type{Measurement{S}}: Measurement type wrapper [dimensionless].
Returns
- The inner floating‐point type
S[dimensionless].
Examples
using Measurements
S = _meas_inner(Measurement{Float64}) # Float64LineCableModels.Utils.block_transform! — Method
Apply f to every square block of M defined by map, in-place.
M: n×n or n×n×nf (any eltype).map: length-n; equal ids => same block (non-contiguous ok).f: function likef(B::AbstractMatrix, args...) -> k×k Matrix.args...: extra positional args passed tof.slice_positions: positions inargsthat should be indexed asargs[i][idx]per block (useful forphase_map).
Returns M.
Index
LineCableModels.DataModelLineCableModels.DataModel.BaseParamsLineCableModels.EarthPropsLineCableModels.ImportExportLineCableModels.MaterialsLineCableModels.UncertainBesselsLineCableModels.UtilsDataFrames.DataFrameDataFrames.DataFrameDataFrames.DataFrameDataFrames.DataFrameDataFrames.DataFrameLineCableModels.DataModel.AbstractCablePartLineCableModels.DataModel.AbstractConductorPartLineCableModels.DataModel.AbstractInsulatorPartLineCableModels.DataModel.AbstractStrandsLayerLineCableModels.DataModel.CableComponentLineCableModels.DataModel.CableComponentLineCableModels.DataModel.CableComponentLineCableModels.DataModel.CableDesignLineCableModels.DataModel.CableDesignLineCableModels.DataModel.CableDesignLineCableModels.DataModel.CableDesignLineCableModels.DataModel.CablePositionLineCableModels.DataModel.CablePositionLineCableModels.DataModel.CablePositionLineCableModels.DataModel.CablesLibraryLineCableModels.DataModel.CablesLibraryLineCableModels.DataModel.CircStrandsLineCableModels.DataModel.CircStrandsLineCableModels.DataModel.ConductorGroupLineCableModels.DataModel.ConductorGroupLineCableModels.DataModel.DiameterLineCableModels.DataModel.InsulatorLineCableModels.DataModel.InsulatorLineCableModels.DataModel.InsulatorLineCableModels.DataModel.InsulatorGroupLineCableModels.DataModel.InsulatorGroupLineCableModels.DataModel.InsulatorGroupLineCableModels.DataModel.LineCableSystemLineCableModels.DataModel.LineCableSystemLineCableModels.DataModel.LineCableSystemLineCableModels.DataModel.LineCableSystemLineCableModels.DataModel.LineCableSystemLineCableModels.DataModel.NominalDataLineCableModels.DataModel.NominalDataLineCableModels.DataModel.RectStrandsLineCableModels.DataModel.RectStrandsLineCableModels.DataModel.RectStrandsShapeLineCableModels.DataModel.SectorLineCableModels.DataModel.SectorGeometryValidLineCableModels.DataModel.SectorInsulatorLineCableModels.DataModel.SectorParamsLineCableModels.DataModel.SemiconLineCableModels.DataModel.SemiconLineCableModels.DataModel.StripLineCableModels.DataModel.StripLineCableModels.DataModel.ThicknessLineCableModels.DataModel.TubularLineCableModels.DataModel.TubularLineCableModels.DataModel.TubularLineCableModels.EarthProps.AbstractFDEMFormulationLineCableModels.EarthProps.CPEarthLineCableModels.EarthProps.CPEarthLineCableModels.EarthProps.EarthLayerLineCableModels.EarthProps.EarthLayerLineCableModels.EarthProps.EarthLayerLineCableModels.EarthProps.EarthModelLineCableModels.EarthProps.EarthModelLineCableModels.EarthProps.EarthModelLineCableModels.Materials.MaterialLineCableModels.Materials.MaterialLineCableModels.Materials.MaterialsLibraryLineCableModels.Materials.MaterialsLibraryBase.delete!Base.delete!Base.getBase.getBase.showBase.showBase.showBase.showBase.showBase.showBase.showBase.showLineCableModels.Commons.add!LineCableModels.Commons.add!LineCableModels.Commons.add!LineCableModels.Commons.add!LineCableModels.Commons.add!LineCableModels.Commons.add!LineCableModels.Commons.add!LineCableModels.DataModel.BaseParams.calc_circstrands_coordsLineCableModels.DataModel.BaseParams.calc_circstrands_gmrLineCableModels.DataModel.BaseParams.calc_equivalent_alphaLineCableModels.DataModel.BaseParams.calc_equivalent_epsLineCableModels.DataModel.BaseParams.calc_equivalent_gmrLineCableModels.DataModel.BaseParams.calc_equivalent_lossfactLineCableModels.DataModel.BaseParams.calc_equivalent_muLineCableModels.DataModel.BaseParams.calc_equivalent_rhoLineCableModels.DataModel.BaseParams.calc_gmdLineCableModels.DataModel.BaseParams.calc_helical_paramsLineCableModels.DataModel.BaseParams.calc_inductance_trifoilLineCableModels.DataModel.BaseParams.calc_parallel_equivalentLineCableModels.DataModel.BaseParams.calc_shunt_capacitanceLineCableModels.DataModel.BaseParams.calc_shunt_conductanceLineCableModels.DataModel.BaseParams.calc_sigma_lossfactLineCableModels.DataModel.BaseParams.calc_solenoid_correctionLineCableModels.DataModel.BaseParams.calc_strip_resistanceLineCableModels.DataModel.BaseParams.calc_temperature_correctionLineCableModels.DataModel.BaseParams.calc_tubular_gmrLineCableModels.DataModel.BaseParams.calc_tubular_inductanceLineCableModels.DataModel.BaseParams.calc_tubular_resistanceLineCableModels.DataModel._coerced_argsLineCableModels.DataModel._ctor_materializeLineCableModels.DataModel._do_add!LineCableModels.DataModel._do_add!LineCableModels.DataModel._extract_part_propertiesLineCableModels.DataModel._normalize_radiiLineCableModels.DataModel._parse_radius_operandLineCableModels.DataModel._print_fieldsLineCableModels.DataModel._promotion_TLineCableModels.DataModel._with_kwdefaultsLineCableModels.DataModel.equivalentLineCableModels.DataModel.flat_formationLineCableModels.DataModel.get_material_color_makieLineCableModels.DataModel.nonsensifyLineCableModels.DataModel.trifoil_formationLineCableModels.DataModel.vdeparseLineCableModels.ImportExport.export_dataLineCableModels.ImportExport.export_dataLineCableModels.ImportExport.export_dataLineCableModels.ImportExport.load!LineCableModels.ImportExport.load!LineCableModels.ImportExport.saveLineCableModels.ImportExport.saveLineCableModels.Materials._add_default_materials!LineCableModels.Utils._coerce_elt_to_TLineCableModels.Utils._hascomplex_typeLineCableModels.Utils._hasmeas_typeLineCableModels.Utils._meas_innerLineCableModels.Utils.bias_to_uncertainLineCableModels.Utils.block_transform!LineCableModels.Utils.coerce_to_TLineCableModels.Utils.coerce_to_TLineCableModels.Utils.coerce_to_TLineCableModels.Utils.coerce_to_TLineCableModels.Utils.is_headlessLineCableModels.Utils.is_in_testsetLineCableModels.Utils.percent_errorLineCableModels.Utils.percent_to_uncertainLineCableModels.Utils.resolve_TLineCableModels.Utils.to_certainLineCableModels.Utils.to_lowerLineCableModels.Utils.to_nominalLineCableModels.Utils.to_upperLineCableModels.Validation.is_radius_inputLineCableModels.Validation.is_radius_inputLineCableModels.Validation.is_radius_inputLineCableModels.Validation.maxfillLineCableModels.DataModel.@construct