Power Flow

For further information on this topic, please see the AC Power Flow or DC Power Flow sections of the Manual. Below, we have provided a list of functions that can be utilized for power flow analysis.


AC Power Flow
DC Power Flow

To load power flow API functionalities into the current scope, utilize the following command:

using JuliaGrid

AC Power Flow

JuliaGrid.newtonRaphsonFunction
newtonRaphson(system::PowerSystem, factorization::Factorization)

The function sets up the Newton-Raphson method to solve the AC power flow.

Arguments

The function requires the PowerSystem composite type to establish the framework. Next, the Factorization argument, while optional, determines the method used to solve the linear system of equations within each iteration. It can take one of the following values:

  • LU: utilizes LU factorization (default);
  • QR: utilizes QR factorization.

Updates

If the AC model has not been created, the function automatically initiates an update within the ac field of the PowerSystem type. It also performs a check on bus types and rectifies any mistakes present.

Returns

The function returns an instance of the ACPowerFlow type, which includes the following fields:

  • voltage: the bus voltage magnitudes and angles;
  • power: the variable allocated to store the active and reactive powers;
  • current: the variable allocated to store the currents;
  • method: the Jacobian matrix, its factorization, mismatches, increments, and indices.

Examples

Set up the Newton-Raphson method utilizing LU factorization:

system = powerSystem("case14.h5")
acModel!(system)

analysis = newtonRaphson(system)

Set up the Newton-Raphson method utilizing QR factorization:

system = powerSystem("case14.h5")
acModel!(system)

analysis = newtonRaphson(system, QR)
source
JuliaGrid.fastNewtonRaphsonBXFunction
fastNewtonRaphsonBX(system::PowerSystem, factorization::Factorization)

The function sets up the fast Newton-Raphson method of version BX to solve the AC power flow.

Arguments

The function requires the PowerSystem composite type to establish the framework. Next, the Factorization argument, while optional, determines the method used to solve the linear system of equations within each iteration. It can take one of the following values:

  • LU: utilizes LU factorization (default);
  • QR: utilizes QR factorization.

Updates

If the AC model has not been created, the function automatically initiates an update within the ac field of the PowerSystem type. It also performs a check on bus types and rectifies any mistakes present.

Returns

The function returns an instance of the ACPowerFlow type, which includes the following fields:

  • voltage: the bus voltage magnitudes and angles;
  • power: the variable allocated to store the active and reactive powers;
  • current: the variable allocated to store the currents;
  • method: the Jacobian matrices, their factorizations, mismatches, increments, and indices.

Examples

Set up the fast Newton-Raphson method utilizing LU factorization:

system = powerSystem("case14.h5")
acModel!(system)

analysis = fastNewtonRaphsonBX(system)

Set up the fast Newton-Raphson method utilizing QR factorization:

system = powerSystem("case14.h5")
acModel!(system)

analysis = fastNewtonRaphsonBX(system, QR)
source
JuliaGrid.fastNewtonRaphsonXBFunction
fastNewtonRaphsonXB(system::PowerSystem, factorization::Factorization)

The function sets up the fast Newton-Raphson method of version XB to solve the AC power flow.

Arguments

The function requires the PowerSystem composite type to establish the framework. Next, the Factorization argument, while optional, determines the method used to solve the linear system of equations within each iteration. It can take one of the following values:

  • LU: utilizes LU factorization (default);
  • QR: utilizes QR factorization.

Updates

If the AC model has not been created, the function automatically initiates an update within the ac field of the PowerSystem type. It also performs a check on bus types and rectifies any mistakes present.

Returns

The function returns an instance of the ACPowerFlow type, which includes the following fields:

  • voltage: the bus voltage magnitudes and angles;
  • power: the variable allocated to store the active and reactive powers;
  • current: the variable allocated to store the currents;
  • method: the Jacobian matrices, their factorizations, mismatches, increments, and indices.

Examples

Set up the fast Newton-Raphson method utilizing LU factorization:

system = powerSystem("case14.h5")
acModel!(system)

analysis = fastNewtonRaphsonXB(system)

Set up the fast Newton-Raphson method utilizing QR factorization:

system = powerSystem("case14.h5")
acModel!(system)

analysis = fastNewtonRaphsonXB(system, QR)
source
JuliaGrid.gaussSeidelFunction
gaussSeidel(system::PowerSystem)

The function sets up the Gauss-Seidel to solve the AC power flow.

Arguments

The function requires the PowerSystem composite type to establish the framework.

Updates

If the AC model has not been created, the function automatically initiates an update within the ac field of the PowerSystem type. It also performs a check on bus types and rectifies any mistakes present.

Returns

The function returns an instance of the ACPowerFlow type, which includes the following fields:

  • voltage: the bus voltage magnitudes and angles;
  • power: the variable allocated to store the active and reactive powers;
  • current: the variable allocated to store the currents;
  • method: contains the bus complex voltages and indices.

Example

system = powerSystem("case14.h5")
acModel!(system)

analysis = gaussSeidel(system)
source
JuliaGrid.mismatch!Method
mismatch!(system::PowerSystem, analysis::ACPowerFlow)

The function calculates both active and reactive power injection mismatches.

Updates

This function updates the mismatch variables in the Newton-Raphson and fast Newton-Raphson methods. It should be employed during the iteration loop before invoking the solve! function.

Returns

The function returns maximum absolute values of the sctive and reactive power injection mismatches, which can be utilized to terminate the iteration loop of the Newton-Raphson, fast Newton-Raphson, or Gauss-Seidel methods employed to solve the AC power flow problem.

Example

system = powerSystem("case14.h5")
acModel!(system)

analysis = newtonRaphson(system)
mismatch!(system, analysis)
source
JuliaGrid.solve!Method
solve!(system::PowerSystem, analysis::ACPowerFlow)

The function employs the Newton-Raphson, fast Newton-Raphson, or Gauss-Seidel method to solve the AC power flow model and calculate bus voltage magnitudes and angles.

After the mismatch! function is called, this function should be executed to perform a single iteration of the method.

Updates

The calculated voltages are stored in the voltage field of the ACPowerFlow type.

Example

system = powerSystem("case14.h5")
acModel!(system)

analysis = newtonRaphson(system)
for i = 1:10
    stopping = mismatch!(system, analysis)
    if all(stopping .< 1e-8)
        break
    end
    solve!(system, analysis)
end
source
JuliaGrid.startingVoltage!Function
startingVoltage!(system::PowerSystem, analysis::ACPowerFlow)

The function extracts bus voltage magnitudes and angles from the PowerSystem composite type and assigns them to the ACPowerFlow type, enabling users to initialize voltage values as required.

Updates

This function only updates the voltage field of the ACPowerFlow type.

Example

system = powerSystem("case14.h5")
acModel!(system)

analysis = newtonRaphson(system)
for i = 1:10
    stopping = mismatch!(system, analysis)
    if all(stopping .< 1e-8)
        break
    end
    solve!(system, analysis)
end

updateBus!(system, analysis; label = 14, reactive = 0.13, magnitude = 1.2, angle = -0.17)

startingVoltage!(system, analysis)
for i = 1:10
    stopping = mismatch!(system, analysis)
    if all(stopping .< 1e-8)
        break
    end
    solve!(system, analysis)
end
source
JuliaGrid.reactiveLimit!Function
reactiveLimit!(system::PowerSystem, analysis::ACPowerFlow)

The function verifies whether the generators in a power system exceed their reactive power limits. This is done by setting the reactive power of the generators to within the limits if they are violated, after determining the bus voltage magnitudes and angles. If the limits are violated, the corresponding generator buses or the slack bus are converted to demand buses.

Updates

The function assigns values to the generator.output.active and bus.supply.active variables of the PowerSystem type. Additionally, it examines the reactive powers of the generator and adjusts them to their maximum or minimum values if they exceed the specified threshold. Subsequently, the generator.output.reactive variable of the PowerSystem type is modified accordingly. As a result of this adjustment, the bus.supply.reactive variable of the PowerSystem type is also updated, and the bus types specified in bus.layout.type are modified. If the slack bus is converted, the bus.layout.slack field is correspondingly adjusted.

Returns

The function returns the variable to indicate which buses violate the limits, with -1 indicating a violation of the minimum limits and 1 indicating a violation of the maximum limits.

Example

system = powerSystem("case14.h5")
acModel!(system)

analysis = newtonRaphson(system)
for i = 1:10
    stopping = mismatch!(system, analysis)
    if all(stopping .< 1e-8)
        break
    end
    solve!(system, analysis)
end

violate = reactiveLimit!(system, analysis)

analysis = newtonRaphson(system)
for i = 1:10
    stopping = mismatch!(system, analysis)
    if all(stopping .< 1e-8)
        break
    end
    solve!(system, analysis)
end
source
JuliaGrid.adjustAngle!Function
adjustAngle!(system::PowerSystem, analysis::ACPowerFlow; slack)

The function modifies the bus voltage angles based on a different slack bus than the one identified by the bus.layout.slack field.

For instance, if the reactive power of the generator exceeds the limit on the slack bus, the reactiveLimit! function will change that bus to the demand bus and designate the first generator bus in the sequence as the new slack bus. After obtaining the updated AC power flow solution based on the new slack bus, it is possible to adjust the voltage angles to align with the angle of the original slack bus. The slack keyword specifies the bus label of the original slack bus.

Updates

This function only updates the voltage.angle variable of the ACPowerFlow type.

Example

system = powerSystem("case14.h5")
acModel!(system)

analysis = newtonRaphson(system)
for i = 1:10
    stopping = mismatch!(system, analysis)
    if all(stopping .< 1e-8)
        break
    end
    solve!(system, analysis)
end

reactiveLimit!(system, analysis)

analysis = newtonRaphson(system)
for i = 1:10
    stopping = mismatch!(system, analysis)
    if all(stopping .< 1e-8)
        break
    end
    solve!(system, analysis)
end

adjustAngle!(system, analysis; slack = 1)
source

DC Power Flow

JuliaGrid.dcPowerFlowFunction
dcPowerFlow(system::PowerSystem, factorization::Factorization)

The function sets up the framework to solve the DC power flow.

Arguments

The function requires the PowerSystem composite type to establish the framework. Next, the Factorization argument, while optional, determines the method used to solve the linear system of equations. It can take one of the following values:

  • LU: utilizes LU factorization (default);
  • LDLt: utilizes LDLt factorization;
  • QR: utilizes QR factorization.

Updates

If the DC model was not created, the function will automatically initiate an update of the dc field within the PowerSystem composite type. Additionally, if the slack bus lacks an in-service generator, JuliaGrid considers it a mistake and defines a new slack bus as the first generator bus with an in-service generator in the bus type list.

Returns

The function returns an instance of the DCPowerFlow type, which includes the following fields:

  • voltage: the variable allocated to store the bus voltage angles;
  • power: the variable allocated to store the active powers;
  • method: the factorized nodal matrix.

Examples

Set up the DC power flow utilizing LU factorization:

system = powerSystem("case14.h5")
dcModel!(system)

analysis = dcPowerFlow(system)

Set up the DC power flow utilizing QR factorization:

system = powerSystem("case14.h5")
dcModel!(system)

analysis = dcPowerFlow(system, QR)
source
JuliaGrid.solve!Method
solve!(system::PowerSystem, analysis::DCPowerFlow)

The function solves the DC power flow model and calculates bus voltage angles.

Updates

The calculated voltage angles are stored in the voltage field of the DCPowerFlow type.

Example

system = powerSystem("case14.h5")
dcModel!(system)

analysis = dcPowerFlow(system)
solve!(system, analysis)
source