PowerModelsDistribution.jl Library

PowerModelsDistribution.constraint_SWL_psdMethod

Take a multi-conductor voltage variable V and a current variable I. The associated power is then defined as S = VI^H Define the lifted variables as W and L as W = VV^H, L = I*I^H Then, it is equally valid that [W S; S^H L] ∈ PSDCone, rank([W S; S^H L])=1 This function adds this PSD constraint for the rectangular coordinates of S, W and L.

source
PowerModelsDistribution.constraint_mc_loadMethod

CONSTANT POWER Fixes the load power sd. sd = [sd1, sd2, sd3] What is actually fixed, depends on whether the load is connected in delta or wye. When connected in wye, the load power equals the per-phase power sn drawn at the bus to which the load is connected. sd1 = va.conj(ia) = sn_a

CONSTANT CURRENT Sets the active and reactive load power sd to be proportional to the the voltage magnitude. pd = cp.|vm| qd = cq.|vm| sd = cp.|vm| + j.cq.|vm|

CONSTANT IMPEDANCE Sets the active and reactive power drawn by the load to be proportional to the square of the voltage magnitude. pd = cp.|vm|^2 qd = cq.|vm|^2 sd = cp.|vm|^2 + j.cq.|vm|^2

DELTA When connected in delta, the load power gives the reference in the delta reference frame. This means sd1 = vab.conj(iab) = (va-vb).conj(iab) We can relate this to the per-phase power by sna = va.conj(ia) = va.conj(iab-ica) = va.conj(conj(sab/vab) - conj(sca/vca)) = va.(sab/(va-vb) - sca/(vc-va)) So for delta, sn is constrained indirectly.

source
PowerModelsDistribution.constraint_mc_load_power_deltaMethod

For a delta load, sd = (sab, sbc, sca), but we want to fix s = (sa, sb, sc) s is a non-linear transform of v and sd, s=f(v,sd) sa = vaconj(sab/(va-vb) - sca/(vc-va)) sb = vbconj(sab/(va-vb) - sca/(vc-va)) sc = vc*conj(sab/(va-vb) - sca/(vc-va))

source
PowerModelsDistribution.constraint_mc_load_power_deltaMethod

For a delta load, sd = (sab, sbc, sca), but we want to fix s = (sa, sb, sc) s is a non-linear transform of v and sd, s=f(v,sd) sa = vaconj(sab/(va-vb) - sca/(vc-va)) sb = vbconj(sab/(va-vb) - sca/(vc-va)) sc = vc*conj(sab/(va-vb) - sca/(vc-va))

source
PowerModelsDistribution.constraint_mc_ohms_yt_fromMethod

Creates Ohms constraints (yt post fix indicates that Y and T values are in rectangular form)

p_fr ==     g[c,c] * vm_fr[c]^2 +
            sum( g[c,d]*vm_fr[c]*vm_fr[d]*cos(va_fr[c]-va_fr[d]) +
                 b[c,d]*vm_fr[c]*vm_fr[d]*sin(va_fr[c]-va_fr[d]) for d in _PMs.conductor_ids(pm) if d != c) +
            sum(-g[c,d]*vm_fr[c]*vm_to[d]*cos(va_fr[c]-va_to[d]) +
                -b[c,d]*vm_fr[c]*vm_to[d]*sin(va_fr[c]-va_to[d]) for d in _PMs.conductor_ids(pm))
            + g_fr[c,c] * vm_fr[c]^2 +
            sum( g_fr[c,d]*vm_fr[c]*vm_fr[d]*cos(va_fr[c]-va_fr[d]) +
                 b_fr[c,d]*vm_fr[c]*vm_fr[d]*sin(va_fr[c]-va_fr[d]) for d in _PMs.conductor_ids(pm) if d != c)
            )
q_fr == -b[c,c] *vm_fr[c]^2 -
            sum( b[c,d]*vm_fr[c]*vm_fr[d]*cos(va_fr[c]-va_fr[d]) -
                 g[c,d]*vm_fr[c]*vm_fr[d]*sin(va_fr[c]-va_fr[d]) for d in _PMs.conductor_ids(pm) if d != c) -
            sum(-b[c,d]*vm_fr[c]*vm_to[d]*cos(va_fr[c]-va_to[d]) +
                 g[c,d]*vm_fr[c]*vm_to[d]*sin(va_fr[c]-va_to[d]) for d in _PMs.conductor_ids(pm))
            -b_fr[c,c] *vm_fr[c]^2 -
            sum( b_fr[c,d]*vm_fr[c]*vm_fr[d]*cos(va_fr[c]-va_fr[d]) -
                 g_fr[c,d]*vm_fr[c]*vm_fr[d]*sin(va_fr[c]-va_fr[d]) for d in _PMs.conductor_ids(pm) if d != c)
            )
source
PowerModelsDistribution.constraint_mc_ohms_yt_fromMethod

Creates Ohms constraints (yt post fix indicates that Y and T values are in rectangular form)

p_fr ==  sum(
             vr_fr[c]*(g[c,d]*(vr_fr[d]-vr_to[d])-b[c,d]*(vi_fr[d]-vi_to[d]))
            -vi_fr[c]*(-b[c,d]*(vr_fr[d]-vr_to[d])-g[c,d]*(vi_fr[d]-vi_to[d]))
            for d in _PMs.conductor_ids(pm))
         + sum(
             vr_fr[c]*(g_fr[c,d]*vr_fr[d]-b_fr[c,d]*vi_fr[d])
            -vi_fr[c]*(-b_fr[c,d]*vr_fr[d]-g_fr[c,d]*vi_fr[d])
            for d in _PMs.conductor_ids(pm))
q_fr ==  sum(
            -vr_fr[c]*(b[c,d]*(vr_fr[d]-vr_to[d])+g[c,d]*(vi_fr[d]-vi_to[d]))
            +vi_fr[c]*(g[c,d]*(vr_fr[d]-vr_to[d])-b[c,d]*(vi_fr[d]-vi_to[d]))
            for d in _PMs.conductor_ids(pm))
          + sum(
            -vr_fr[c]*(b_fr[c,d]*vr_fr[d]+g_fr[c,d]*vi_fr[d])
            +vi_fr[c]*(g_fr[c,d]*vr_fr[d]-b_fr[c,d]*vi_fr[d])
            for d in _PMs.conductor_ids(pm))
source
PowerModelsDistribution.constraint_mc_ohms_yt_toMethod

Creates Ohms constraints (yt post fix indicates that Y and T values are in rectangular form)

p[t_idx] ==  (g+g_to)*v[t_bus]^2 + (-g*tr-b*ti)/tm*(v[t_bus]*v[f_bus]*cos(t[t_bus]-t[f_bus])) + (-b*tr+g*ti)/tm*(v[t_bus]*v[f_bus]*sin(t[t_bus]-t[f_bus]))
q[t_idx] == -(b+b_to)*v[t_bus]^2 - (-b*tr+g*ti)/tm*(v[t_bus]*v[f_bus]*cos(t[f_bus]-t[t_bus])) + (-g*tr-b*ti)/tm*(v[t_bus]*v[f_bus]*sin(t[t_bus]-t[f_bus]))
source
PowerModelsDistribution.constraint_mc_ohms_yt_toMethod

Creates Ohms constraints (yt post fix indicates that Y and T values are in rectangular form)

p[t_idx] ==  (g+g_to)*v[t_bus]^2 + (-g*tr-b*ti)/tm*(v[t_bus]*v[f_bus]*cos(t[t_bus]-t[f_bus])) + (-b*tr+g*ti)/tm*(v[t_bus]*v[f_bus]*sin(t[t_bus]-t[f_bus]))
q[t_idx] == -(b+b_to)*v[t_bus]^2 - (-b*tr+g*ti)/tm*(v[t_bus]*v[f_bus]*cos(t[f_bus]-t[t_bus])) + (-g*tr-b*ti)/tm*(v[t_bus]*v[f_bus]*sin(t[t_bus]-t[f_bus]))
source
PowerModelsDistribution.constraint_mc_voltage_balanceMethod

Impose all balance related constraints for which key present in data model of bus. For a discussion of sequence components and voltage unbalance factor (VUF), see @INPROCEEDINGS{girigoudarmolzahnroald-2019, author={K. Girigoudar and D. K. Molzahn and L. A. Roald}, booktitle={submitted}, title={{Analytical and Empirical Comparisons of Voltage Unbalance Definitions}}, year={2019}, month={}, url={https://molzahn.github.io/pubs/girigoudarmolzahnroald-2019.pdf} }

source
PowerModelsDistribution.constraint_pqwMethod

Creates the constraints modelling the (relaxed) voltage-dependency of the power consumed in each phase, s=p+jq. This is completely symmetrical for p and q, with appropriate substitutions of the variables and parameters: p->q, a->b, alpha->beta, pmin->qmin, pmax->qmax

source
PowerModelsDistribution.parse_dssMethod
parse_dss(filename)

Parses a OpenDSS file given by filename into a Dict{Array{Dict}}. Only supports components and options, but not commands, e.g. "plot" or "solve". Will also parse files defined inside of the originating DSS file via the "compile", "redirect" or "buscoords" commands.

source
PowerModelsDistribution.relaxation_psd_to_socMethod

See section 4.3 for complex to real PSD constraint transformation: @article{Fazel2001, author = {Fazel, M. and Hindi, H. and Boyd, S.P.}, title = {{A rank minimization heuristic with application to minimum order system approximation}}, doi = {10.1109/ACC.2001.945730}, journal = {Proc. American Control Conf.}, number = {2}, pages = {4734–4739}, url = {http://ieeexplore.ieee.org/lpdocs/epic03/wrapper.htm?arnumber=945730}, volume = {6}, year = {2001} }

source
PowerModelsDistribution.relaxation_psd_to_soc_complexMethod

SDP to SOC relaxation of type 2, applied to complex-value matrix, as described in:

@article{Kim2003,
author = {Kim, S and Kojima, M and Yamashita, M},
title = {{Second order cone programming relaxation of a positive semidefinite constraint}},
doi = {10.1080/1055678031000148696},
journal = {Optimization Methods and Software},
number = {5},
pages = {535--541},
volume = {18},
year = {2003}
}
source
PowerModelsDistribution.relaxation_psd_to_soc_complex_conicMethod

SDP to SOC relaxation of type 2, applied to complex-value matrix, as described in:

@article{Kim2003,
author = {Kim, S and Kojima, M and Yamashita, M},
title = {{Second order cone programming relaxation of a positive semidefinite constraint}},
doi = {10.1080/1055678031000148696},
journal = {Optimization Methods and Software},
number = {5},
pages = {535--541},
volume = {18},
year = {2003}
}
source
PowerModelsDistribution.relaxation_psd_to_soc_conicMethod

See section 4.3 for complex to real PSD constraint transformation: @article{Fazel2001, author = {Fazel, M. and Hindi, H. and Boyd, S.P.}, title = {{A rank minimization heuristic with application to minimum order system approximation}}, doi = {10.1109/ACC.2001.945730}, journal = {Proc. American Control Conf.}, number = {2}, pages = {4734–4739}, url = {http://ieeexplore.ieee.org/lpdocs/epic03/wrapper.htm?arnumber=945730}, volume = {6}, year = {2001} }

source
PowerModelsDistribution.relaxation_psd_to_soc_realMethod

SDP to SOC relaxation of type 2, applied to real-value matrix, as described in:

@article{Kim2003,
author = {Kim, S and Kojima, M and Yamashita, M},
title = {{Second order cone programming relaxation of a positive semidefinite constraint}},
doi = {10.1080/1055678031000148696},
journal = {Optimization Methods and Software},
number = {5},
pages = {535--541},
volume = {18},
year = {2003}
}
source
PowerModelsDistribution.relaxation_psd_to_soc_real_conicMethod

SDP to SOC relaxation of type 2, applied to real-value matrix, as described in:

@article{Kim2003,
author = {Kim, S and Kojima, M and Yamashita, M},
title = {{Second order cone programming relaxation of a positive semidefinite constraint}},
doi = {10.1080/1055678031000148696},
journal = {Optimization Methods and Software},
number = {5},
pages = {535--541},
volume = {18},
year = {2003}
}
source
PowerModelsDistribution.variable_mc_loadMethod

The variable creation for the loads is rather complicated because Expressions are used wherever possible instead of explicit variables. All loads need a current variable; for wye loads, this variable will be in the wye reference frame whilst for delta currents it will be in the delta reference frame. All loads need variables for the off-diagonals of the nodal power variables. In some cases, the diagonals elements can be created as Expressions. Delta loads only need a current variable and auxilary power variable (X), and all other load model variables are then linear transformations of these (linear Expressions).

source
PowerModelsDistribution.variable_mc_loadMethod

The variable creation for the loads is rather complicated because Expressions are used wherever possible instead of explicit variables. Delta loads always need a current variable and auxilary power variable (X), and all other load model variables are then linear transformations of these (linear Expressions). Wye loads however, don't need any variables when the load is modelled as constant power or constant impedance. In all other cases (e.g. when a cone is used to constrain the power), variables need to be created.

source
PowerModelsDistribution.variable_mc_load_delta_auxMethod

Creates power matrix variable X for delta windings; this defines both the wye-side power Sy and the delta-side power Sd through the lin. transformations Sy = X.Td, Sd = Td.X with Td=[1 -1 0; 0 1 -1; -1 0 1]

See the paper by Zhao et al. for the first convex relaxation of delta transformations. @INPROCEEDINGS{zhaooptimal2017, author={C. Zhao, E. Dall'Anese and S. Low}, booktitle={IREP 2017 Bulk Power Systems Dynamics and Control Symposium}, title={{Optimal Power Flow in Multiphase Radial Networks with Delta Connections}}, year={2017}, month={}, url={https://www.nrel.gov/docs/fy18osti/67852.pdf} }

See upcoming paper for discussion of bounds. [reference added when accepted]

source
PowerModelsDistribution.variable_mc_load_power_busMethod

The bus qualifier denotes that this is the power withdrawn at the bus; Only for grounded wye-connected loads, this is the same as the power consumed by the multi-phase load. The off-diagonals only need to be created for the matrix KCL formulation.

source
PowerModelsDistribution.variable_mx_complexMethod

Shorthand to create two real matrix variables, where the first is the real part and the second the imaginary part. If the name argument is a String, it will be suffixed with 're' and 'im'. It is possible to specify the names of the real and imaginary part directly as a Tuple as well (to achieve P and Q instead of Sre and Sim for example).

source
PowerModelsDistribution.variable_mx_complex_with_diagMethod

Same as variablemxcomplex, but square and the diagonal of the matrix variables consists of the constants passed as the diagre and diagim argument. The diag argument is a dictionary of (index, 1d-array) pairs. Useful for power matrices with specified diagonals (constant power wye loads).

source
PowerModelsDistribution.variable_mx_realMethod

This function creates a set of real matrix variables of size NxM, indexed over the elements of the indices argument. The upper and lower bounds have to be specified, and are dictionaries with the indices as keys and the matrix bounds as values. The name and prefix arguments will be combined into the base_name argument for JuMP; the prefix will typically be the network number nw. Instead of sequentially creating the matrix variables, the elements of the matrices are created sequentially for all matrices at once. I.e., we loop over the elements, and not over the indices. This is needed so that the variable names printed by JuMP are in line with the current design.

Returns a dictionary of (index, matrix variable) pairs

source
PowerModelsDistribution.variable_mx_real_with_diagMethod

Same as variablemxreal, but has to be square and the diagonal of the matrix variables consists of the elements passed as the diag argument. The diag argument is a dictionary of (index, 1d-array) pairs. Useful for power matrices with specified diagonals (constant power wye loads). If not specified, the diagonal elements are set to zero.

source
PowerModelsDistribution._add_component!Method
_add_component!(dss_data, ctype_name, compDict)

Adds a component of type ctype_name with properties given by compDict to the existing dss_data structure. If a component of the same type has already been added to dss_data, the new component is appeneded to the existing array of components of that type, otherwise a new array is created.

source
PowerModelsDistribution._add_propertyMethod
_add_property(compDict, key, value)

Adds a property to an existing component properties dictionary compDict given the key and value of the property. If a property of the same name already exists inside compDict, the original value is converted to an array, and the new value is appended to the end.

source
PowerModelsDistribution._add_shunt!Method

Helper function to add a new shunt. The shunt element is always inserted at the internal bus of the second winding in OpenDSS. If one of the branches of the loss model connected to this bus, has zero impedance (for example, if XHL==0 or XLT==0 or R[3]==0), then this bus might be removed by rmredundantpdelements!, in which case a new shunt should be inserted at the remaining bus of the removed branch.

source
PowerModelsDistribution._adjust_base!Method
function _adjust_base!(pmd_data)

Updates the voltage base at each bus, so that the ratios of the voltage bases across a transformer are consistent with the ratios of voltage ratings of the windings. Default behaviour is to start at the primary winding of the first transformer, and to propagate from there. Branches are updated; the impedances and addmittances are rescaled to be consistent with the new voltage bases.

source
PowerModelsDistribution._adjust_base_rec!Method

This is the recursive code that goes with adjustbase!; adjustbase! initializes arrays and other data that is passed along in the calls to this recursive function. For very large networks, this might have to be rewritten to not rely on recursion.

source
PowerModelsDistribution._adjust_sourcegen_bounds!Method
_adjust_sourcegen_bounds!(pmd_data)

Changes the bounds for the sourcebus generator by checking the emergamps of all of the branches attached to the sourcebus and taking the sum of non-infinite values. Defaults to Inf if all emergamps connected to sourcebus are also Inf. This method was updated to include connected transformers as well. It know has to occur after the call to InfrastructureModels.arraystodicts, so the code was adjusted to accomodate that.

source
PowerModelsDistribution._calc_branch_power_ub_frtoMethod

Returns a total (shunt+series) power magnitude bound for the from and to side of a branch. The total current rating also implies a current bound through the upper bound on the voltage magnitude of the connected buses.

source
PowerModelsDistribution._calc_bus_vm_ll_boundsMethod

Returns bounds in line-to-line bounds on the voltage magnitude. If these are not part of the problem specification, then a valid upper bound is implied by the line-to-neutral bounds, but a lower bound (greater than zero) is not. Therefore, a default lower bound is then used, specified by the keyword argument vdmin_eps. The returned bounds are for the pairs 1->2, 2->3, 3->1

source
PowerModelsDistribution._calc_load_vboundsMethod

Returns the voltage magnitude bounds for the individual load elements in a multiphase load. These are inferred from vmin/vmax for wye loads and from calcbusvmll_bounds for delta loads.

source
PowerModelsDistribution._create_capacitorFunction
_create_capacitor(bus1, name, bus2=0; kwargs)

Creates a Dict{String,Any} containing all of the expected properties for a Capacitor. If bus2 is not specified, the capacitor will be treated as a shunt. See OpenDSS documentation for valid fields and ways to specify the different properties.

source
PowerModelsDistribution._create_generatorFunction
_create_generator(bus1, name; kwargs...)

Creates a Dict{String,Any} containing all of the expected properties for a Generator. See OpenDSS documentation for valid fields and ways to specify the different properties.

source
PowerModelsDistribution._create_lineFunction
_create_line(bus1, bus2, name; kwargs...)

Creates a Dict{String,Any} containing all of the properties for a Line. See OpenDSS documentation for valid fields and ways to specify the different properties.

source
PowerModelsDistribution._create_linecodeFunction
_create_linecode(name; kwargs...)

Creates a Dict{String,Any} containing all of the properties of a Linecode. See OpenDSS documentation for valid fields and ways to specify the different properties. DEPRECIATED: Calculation all done inside of createline() due to Rg, Xg. Merge linecode values into line kwargs values BEFORE calling createline(). This is now mainly used for parsing linecode dicts into correct data types.

source
PowerModelsDistribution._create_loadFunction
_create_load(bus1, name; kwargs...)

Creates a Dict{String,Any} containing all of the expected properties for a Load. See OpenDSS documentation for valid fields and ways to specify the different properties.

source
PowerModelsDistribution._create_pvsystemFunction
_create_pvsystem(bus1, name; kwargs...)

Creates a Dict{String,Any} containing all of the expected properties for a PVSystem. See OpenDSS document https://github.com/tshort/OpenDSS/blob/master/Doc/OpenDSS%20PVSystem%20Model.doc for valid fields and ways to specify the different properties.

source
PowerModelsDistribution._create_reactorFunction
_create_reactor(bus1, name, bus2=0; kwargs...)

Creates a Dict{String,Any} containing all of the expected properties for a Reactor. If bus2 is not specified Reactor is treated like a shunt. See OpenDSS documentation for valid fields and ways to specify the different properties.

source
PowerModelsDistribution._create_storageFunction
_create_storage(bus1, name; kwargs...)

Creates a Dict{String,Any} containing all expected properties for a storage element. See OpenDSS documentation for valid fields and ways to specify the different properties.

source
PowerModelsDistribution._create_transformerFunction
_create_transformer(name; kwargs...)

Creates a Dict{String,Any} containing all of the expected properties for a Transformer. See OpenDSS documentation for valid fields and ways to specify the different properties.

source
PowerModelsDistribution._create_vbranch!Method

This function adds a new branch to the data model and returns its dictionary. It is virtual in the sense that it does not correspond to a branch in the network, but is part of the decomposition of the transformer.

source
PowerModelsDistribution._create_vbus!Method

This function adds a new bus to the data model and returns its dictionary. It is virtual in the sense that it does not correspond to a bus in the network, but is part of the decomposition of the transformer.

source
PowerModelsDistribution._create_vsourceFunction
_create_vsource(bus1, name, bus2=0; kwargs...)

Creates a Dict{String,Any} containing all of the expected properties for a Voltage Source. If bus2 is not specified, VSource will be treated like a generator. Mostly used as sourcebus which represents the circuit. See OpenDSS documentation for valid fields and ways to specify the different properties.

source
PowerModelsDistribution._load_expmodel_paramsMethod

Returns the exponential load model parameters for a load. For an exponential load it simply returns certain data model properties, whilst for constantpower, constantcurrent and constant_impedance it returns the equivalent exponential model parameters.

source
PowerModelsDistribution._make_matrix_variable_elementMethod

Sometimes we want to bound only a subset of the elements of a matrix variable. For example, an unbounded Hermitian variable usually still has a lower bound of zero on the real diagonal elements. When there is a mix of bounded and unbounded elements, the unboundedness is encoded as 'Inf' and '-Inf' in the bound parameters. This cannot be passed directlty to JuMP, because it would lead to an error in Mosek for example. Instead, this method checks whether all bounds for an element (n,m) are Inf, and if so, does not pass a bound to JuMP.

source
PowerModelsDistribution._merge_dss!Method
_merge_dss!(dss_prime, dss_to_add)

Merges two (partially) parsed OpenDSS files to the same dictionary dss_prime. Used in cases where files are referenced via the "compile" or "redirect" OpenDSS commands inside the originating file.

source
PowerModelsDistribution._parse_arrayMethod
_parse_array(dtype, data)

Parses a OpenDSS style array string data into a one dimensional array of type dtype. Array strings are capped by either brackets, single quotes, or double quotes, and elements are separated by spaces.

source
PowerModelsDistribution._parse_buscoordsMethod
_parse_buscoords(file)

Parses a Bus Coordinate file, in either "dat" or "csv" formats, where in "dat", columns are separated by spaces, and in "csv" by commas. File expected to contain "bus,x,y" on each line.

source
PowerModelsDistribution._parse_componentFunction
_parse_component(component, properies, compDict=Dict{String,Any}())

Parses a component with properties into a compDict. If compDict is not defined, an empty dictionary will be used. Assumes that unnamed properties are given in order, but named properties can be given anywhere.

source
PowerModelsDistribution._parse_lineFunction
_parse_line(elements, curCompDict=Dict{String,Any}())

Parses an already separated line given by elements (an array) of an OpenDSS file into curCompDict. If not defined, curCompDict is an empty dictionary.

source
PowerModelsDistribution._parse_matrixMethod
_parse_matrix(dtype, data)

Parses a OpenDSS style triangular matrix string data into a two dimensional array of type dtype. Matrix strings are capped by either parenthesis or brackets, rows are separated by "|", and columns are separated by spaces.

source
PowerModelsDistribution._parse_propertiesMethod
_parse_properties(properties)

Parses a string of properties of a component type, character by character into an array with each element containing (if present) the property name, "=", and the property value.

source
PowerModelsDistribution._rm_redundant_pd_elements!Method

This function removes zero impedance branches. Only for transformer loss model! Branches with zero impedances are deleted, and one of the buses it connects. For now, the implementation should only be used on the loss model of transformers. When deleting buses, references at shunts, loads... should be updated accordingly. In the current implementation, that is only done for shunts. The other elements, such as loads, do not appear in the transformer loss model.

source
PowerModelsDistribution._sc2br_impedanceMethod

Converts a set of short-circuit tests to an equivalent reactance network. Reference: R. C. Dugan, “A perspective on transformer modeling for distribution system analysis,” in 2003 IEEE Power Engineering Society General Meeting (IEEE Cat. No.03CH37491), 2003, vol. 1, pp. 114-119 Vol. 1.

source