Buildings.Electrical.AC.OnePhase.Conversion

Package with models for AC/AC and AC/DC conversion

Information

This package contains models for AC/AC and AC/DC conversion.

Extends from Modelica.Icons.Package (Icon for standard packages).

Package Content

Name Description
Buildings.Electrical.AC.OnePhase.Conversion.ACACConverter ACACConverter AC AC converter single phase systems
Buildings.Electrical.AC.OnePhase.Conversion.ACACTransformer ACACTransformer AC AC transformer simplified equivalent circuit
Buildings.Electrical.AC.OnePhase.Conversion.ACACTransformerFull ACACTransformerFull AC AC transformer with detailed equivalent circuit
Buildings.Electrical.AC.OnePhase.Conversion.ACDCConverter ACDCConverter AC DC converter
Buildings.Electrical.AC.OnePhase.Conversion.Examples Examples Package with example models

Buildings.Electrical.AC.OnePhase.Conversion.ACACConverter Buildings.Electrical.AC.OnePhase.Conversion.ACACConverter

AC AC converter single phase systems

Buildings.Electrical.AC.OnePhase.Conversion.ACACConverter

Information

This is an AC/AC converter, based on a power balance between both circuit sides. The parameter conversionFactor defines the ratio between the RMS voltages as

V2 = conversionFactor V1

where V1 and V2 are the RMS voltages at the primary and secondary sides of the transformer, i.e., the connector N and P, respectively.

The loss of the converter is proportional to the power transmitted. The parameter eta is the efficiency of the transfer. The loss is computed as

Ploss = (1-η) Ptr,

where Ptr is the power transmitted. The model is bi-directional and the power can flow from the primary to the secondary side and vice-versa. Furthermore, reactive power on both side are set to zero.

Note:

This model is derived from Modelica.Electrical.QuasiStationary.SinglePhase.Utilities.IdealACDCConverter.

Extends from Buildings.Electrical.Icons.RefAngleConversion (Icon that represents if the angle symble should be displayed or not), Buildings.Electrical.Interfaces.PartialConversion (Model representing a generic two port system for conversion).

Parameters

TypeNameDefaultDescription
replaceable package PhaseSystem_pPartialPhaseSystemPhase system of terminal p
replaceable package PhaseSystem_nPartialPhaseSystemPhase system of terminal n
RealconversionFactor Ratio of QS rms voltage on side 2 / QS rms voltage on side 1
Realeta Converter efficiency, pLoss = (1-eta) * Ptr
Ground
side 1
Booleanground_1falseIf true, connect side 1 of converter to ground
side 2
Booleanground_2trueIf true, connect side 2 of converter to ground

Connectors

TypeNameDescription
replaceable package PhaseSystem_pPhase system of terminal p
replaceable package PhaseSystem_nPhase system of terminal n

Modelica definition

model ACACConverter "AC AC converter single phase systems" extends Buildings.Electrical.Icons.RefAngleConversion; extends Buildings.Electrical.Interfaces.PartialConversion( redeclare package PhaseSystem_p = PhaseSystems.OnePhase, redeclare package PhaseSystem_n = PhaseSystems.OnePhase, redeclare replaceable Interfaces.Terminal_n terminal_n constrainedby Interfaces.Terminal_n( i(start = zeros(PhaseSystem_n.n), each stateSelect = StateSelect.prefer)), redeclare replaceable Interfaces.Terminal_p terminal_p constrainedby Interfaces.Terminal_p( i(start = zeros(PhaseSystem_p.n), each stateSelect = StateSelect.prefer))); parameter Real conversionFactor(min = Modelica.Constants.eps) "Ratio of QS rms voltage on side 2 / QS rms voltage on side 1"; parameter Real eta(min=0, max=1) "Converter efficiency, pLoss = (1-eta) * Ptr"; parameter Boolean ground_1 = false "If true, connect side 1 of converter to ground"; parameter Boolean ground_2 = true "If true, connect side 2 of converter to ground"; Modelica.SIunits.Power LossPower[2] "Loss power"; protected Modelica.SIunits.Power P_p[2] = PhaseSystem_p.phasePowers_vi(terminal_p.v, terminal_p.i) "Power transmitted at pin p"; Modelica.SIunits.Power P_n[2] = PhaseSystem_n.phasePowers_vi(terminal_n.v, terminal_n.i) "Power transmitted at pin n"; equation // Ideal transformation terminal_p.v = conversionFactor*terminal_n.v; // Power loss term terminal_p.i[1] = terminal_n.i[1]/conversionFactor* Buildings.Utilities.Math.Functions.spliceFunction(eta-2, 1/(eta-2), P_p[1], deltax=0.1); terminal_p.i[2] = terminal_n.i[2]/conversionFactor* Buildings.Utilities.Math.Functions.spliceFunction(eta-2, 1/(eta-2), P_p[1], deltax=0.1); LossPower = P_p + P_n; // The two sides have the same reference angle terminal_p.theta = terminal_n.theta; if ground_1 then Connections.potentialRoot(terminal_n.theta); end if; if ground_2 then Connections.root(terminal_p.theta); end if; end ACACConverter;

Buildings.Electrical.AC.OnePhase.Conversion.ACACTransformer Buildings.Electrical.AC.OnePhase.Conversion.ACACTransformer

AC AC transformer simplified equivalent circuit

Buildings.Electrical.AC.OnePhase.Conversion.ACACTransformer

Information

This is a simplified equivalent transformer model. The model accounts for winding Joule losses and leakage reactances that are represented by a series of a resistance R and an inductance L. The resistance and the inductance represent both the effects of the secondary and primary side of the transformer.

The model is parameterized using the following parameters

Given the nominal conditions,the model computes the values of the resistance and the inductance.

Extends from Buildings.Electrical.Icons.RefAngleConversion (Icon that represents if the angle symble should be displayed or not), Buildings.Electrical.Interfaces.PartialConversion (Model representing a generic two port system for conversion).

Parameters

TypeNameDefaultDescription
replaceable package PhaseSystem_pPartialPhaseSystemPhase system of terminal p
replaceable package PhaseSystem_nPartialPhaseSystemPhase system of terminal n
VoltageVHigh Rms voltage on side 1 of the transformer (primary side) [V]
VoltageVLow Rms voltage on side 2 of the transformer (secondary side) [V]
ApparentPowerVABase Nominal power of the transformer [V.A]
RealXoverR Ratio between the complex and real components of the impedance (XL/R)
RealZperc Short circuit impedance
Ground
side 1
Booleanground_1falseIf true, connect side 1 of converter to ground
side 2
Booleanground_2trueIf true, connect side 2 of converter to ground
Initialization
Anglephi_10Angle of the voltage side 1 at initialization [rad]
Anglephi_2phi_1Angle of the voltage side 2 at initialization [rad]

Connectors

TypeNameDescription
replaceable package PhaseSystem_pPhase system of terminal p
replaceable package PhaseSystem_nPhase system of terminal n

Modelica definition

model ACACTransformer "AC AC transformer simplified equivalent circuit" extends Buildings.Electrical.Icons.RefAngleConversion; extends Buildings.Electrical.Interfaces.PartialConversion( redeclare package PhaseSystem_p = PhaseSystems.OnePhase, redeclare package PhaseSystem_n = PhaseSystems.OnePhase, redeclare replaceable Interfaces.Terminal_n terminal_n constrainedby Interfaces.Terminal_n( i(start = zeros(PhaseSystem_n.n), each stateSelect = StateSelect.prefer)), redeclare replaceable Interfaces.Terminal_p terminal_p constrainedby Interfaces.Terminal_p( i(start = zeros(PhaseSystem_p.n), each stateSelect = StateSelect.prefer))); parameter Modelica.SIunits.Voltage VHigh "Rms voltage on side 1 of the transformer (primary side)"; parameter Modelica.SIunits.Voltage VLow "Rms voltage on side 2 of the transformer (secondary side)"; parameter Modelica.SIunits.ApparentPower VABase "Nominal power of the transformer"; parameter Real XoverR "Ratio between the complex and real components of the impedance (XL/R)"; parameter Real Zperc "Short circuit impedance"; parameter Boolean ground_1 = false "If true, connect side 1 of converter to ground"; parameter Boolean ground_2 = true "If true, connect side 2 of converter to ground"; parameter Modelica.SIunits.Angle phi_1 = 0 "Angle of the voltage side 1 at initialization"; parameter Modelica.SIunits.Angle phi_2 = phi_1 "Angle of the voltage side 2 at initialization"; Modelica.SIunits.Efficiency eta "Efficiency"; Modelica.SIunits.Power PLoss[2] "Loss power"; Modelica.SIunits.Voltage V1[2]( start = PhaseSystem_n.phaseVoltages(VHigh, phi_1)) "Voltage at the winding - primary side"; Modelica.SIunits.Voltage V2[2]( start = PhaseSystem_p.phaseVoltages(VLow, phi_2)) "Voltage at the winding - secondary side"; protected Real N = VHigh/VLow "Winding ratio"; Modelica.SIunits.Current IHigh = VABase/VHigh "Nominal current on primary side"; Modelica.SIunits.Current ILow = VABase/VLow "Nominal current on secondary side"; Modelica.SIunits.Current IscHigh = IHigh/Zperc "Short circuit current on primary side"; Modelica.SIunits.Current IscLow = ILow/Zperc "Short circuit current on secondary side"; Modelica.SIunits.Impedance Zp = VHigh/IscHigh "Impedance of the primary side (module)"; Modelica.SIunits.Impedance Z1[2]= {Zp*cos(atan(XoverR)), Zp*sin(atan(XoverR))} "Impedance of the primary side of the transformer"; Modelica.SIunits.Impedance Zs = VLow/IscLow "Impedance of the secondary side (module)"; Modelica.SIunits.Impedance Z2[2]= {Zs*cos(atan(XoverR)), Zs*sin(atan(XoverR))} "Impedance of the secondary side of the transformer"; Modelica.SIunits.Power P_p[2]= PhaseSystem_p.phasePowers_vi(terminal_p.v, terminal_p.i) "Power transmitted at pin p (secondary)"; Modelica.SIunits.Power P_n[2]= PhaseSystem_n.phasePowers_vi(terminal_n.v, terminal_n.i) "Power transmitted at pin n (primary)"; Modelica.SIunits.Power S_p = Modelica.Fluid.Utilities.regRoot(P_p[1]^2 + P_p[2]^2, delta=0.1) "Apparent power at terminal p"; Modelica.SIunits.Power S_n = Modelica.Fluid.Utilities.regRoot(P_n[1]^2 + P_n[2]^2, delta=0.1) "Apparent power at terminal n"; equation // Efficiency eta = Buildings.Utilities.Math.Functions.smoothMin( x1= Modelica.Fluid.Utilities.regRoot(P_p[1]^2 + P_p[2]^2, delta=0.01)/ Modelica.Fluid.Utilities.regRoot(P_n[1]^2 + P_n[2]^2 + 1e-6, delta=0.01), x2= Modelica.Fluid.Utilities.regRoot(P_n[1]^2 + P_n[2]^2, delta=0.01)/ Modelica.Fluid.Utilities.regRoot(P_p[1]^2 + P_p[2]^2 + 1e-6, delta=0.01), deltaX = 0.01); // Ideal transformation V2 = V1/N; terminal_p.i[1] + terminal_n.i[1]*N = 0; terminal_p.i[2] + terminal_n.i[2]*N = 0; // Losses due to the impedance terminal_n.v = V1 + Buildings.Electrical.PhaseSystems.OnePhase.product( terminal_n.i, Z1); V2 = terminal_p.v; // Loss of power PLoss = P_p + P_n; // The two sides have the same reference angle terminal_p.theta = terminal_n.theta; if ground_1 then Connections.potentialRoot(terminal_n.theta); end if; if ground_2 then Connections.root(terminal_p.theta); end if; end ACACTransformer;

Buildings.Electrical.AC.OnePhase.Conversion.ACACTransformerFull Buildings.Electrical.AC.OnePhase.Conversion.ACACTransformerFull

AC AC transformer with detailed equivalent circuit

Buildings.Electrical.AC.OnePhase.Conversion.ACACTransformerFull

Information

This is a detailed transformer model that takes into account the winding Joule losses and the leakage reactances on both primary and secondary side. The model also describes the core or iron losses and the losses due to magnetization effects.

The losses are represented by a series of resistances R1, R2, Rm and inductances L1, L2, and Lm.

The model is parameterized using the following parameters

Given the nominal conditions, the model computes the values of the nominal impedances at both primary and secondary side. Given these values, the per unit values are transformed into the actual values of the resistances and inductances.

The magnetization losses can be enabled or disabled using the boolean flag magEffects.

Extends from Buildings.Electrical.Icons.RefAngleConversion (Icon that represents if the angle symble should be displayed or not), Buildings.Electrical.Interfaces.PartialConversion (Model representing a generic two port system for conversion).

Parameters

TypeNameDefaultDescription
replaceable package PhaseSystem_pPartialPhaseSystemPhase system of terminal p
replaceable package PhaseSystem_nPartialPhaseSystemPhase system of terminal n
VoltageVHigh RMS voltage on side 1 of the transformer (primary side) [V]
VoltageVLow RMS voltage on side 2 of the transformer (secondary side) [V]
ApparentPowerVABase Nominal power of the transformer [V.A]
Frequencyf Nominal frequency [Hz]
PerUnitR1 Resistance on side 1 of the transformer (pu) [1]
PerUnitL1 Inductance on side 1 of the transformer (pu) [1]
PerUnitR2 Resistance on side 2 of the transformer (pu) [1]
PerUnitL2 Inductance on side 2 of the transformer (pu) [1]
Magnetization
BooleanmagEffectsfalseIf true, introduce magnetization effects
PerUnitRm Magnetization resistance (pu) [1]
PerUnitLm Magnetization inductance (pu) [1]
Ground
side 1
Booleanground_1falseConnect side 1 of converter to ground
side 2
Booleanground_2trueConnect side 2 of converter to ground
Initialization
Anglephi_10Angle of the voltage side 1 at initialization [rad]
Anglephi_2phi_1Angle of the voltage side 2 at initialization [rad]

Connectors

TypeNameDescription
replaceable package PhaseSystem_pPhase system of terminal p
replaceable package PhaseSystem_nPhase system of terminal n

Modelica definition

model ACACTransformerFull "AC AC transformer with detailed equivalent circuit" extends Buildings.Electrical.Icons.RefAngleConversion; extends Buildings.Electrical.Interfaces.PartialConversion( redeclare package PhaseSystem_p = PhaseSystems.OnePhase, redeclare package PhaseSystem_n = PhaseSystems.OnePhase, redeclare replaceable Interfaces.Terminal_n terminal_n constrainedby Interfaces.Terminal_n( i(start = zeros(PhaseSystem_n.n), each stateSelect = StateSelect.prefer)), redeclare replaceable Interfaces.Terminal_p terminal_p constrainedby Interfaces.Terminal_p( i(start = zeros(PhaseSystem_p.n), each stateSelect = StateSelect.prefer))); parameter Modelica.SIunits.Voltage VHigh "RMS voltage on side 1 of the transformer (primary side)"; parameter Modelica.SIunits.Voltage VLow "RMS voltage on side 2 of the transformer (secondary side)"; parameter Modelica.SIunits.ApparentPower VABase "Nominal power of the transformer"; parameter Modelica.SIunits.Frequency f(start=60) "Nominal frequency"; parameter Buildings.Electrical.Types.PerUnit R1(min=0) "Resistance on side 1 of the transformer (pu)"; parameter Buildings.Electrical.Types.PerUnit L1(min=0) "Inductance on side 1 of the transformer (pu)"; parameter Buildings.Electrical.Types.PerUnit R2(min=0) "Resistance on side 2 of the transformer (pu)"; parameter Buildings.Electrical.Types.PerUnit L2(min=0) "Inductance on side 2 of the transformer (pu)"; parameter Boolean magEffects = false "If true, introduce magnetization effects"; parameter Buildings.Electrical.Types.PerUnit Rm(min=0) "Magnetization resistance (pu)"; parameter Buildings.Electrical.Types.PerUnit Lm(min=0) "Magnetization inductance (pu)"; parameter Boolean ground_1 = false "Connect side 1 of converter to ground"; parameter Boolean ground_2 = true "Connect side 2 of converter to ground"; parameter Modelica.SIunits.Angle phi_1 = 0 "Angle of the voltage side 1 at initialization"; parameter Modelica.SIunits.Angle phi_2 = phi_1 "Angle of the voltage side 2 at initialization"; Modelica.SIunits.Efficiency eta "Efficiency"; Modelica.SIunits.Power PLoss[2] "Loss power"; Modelica.SIunits.Voltage V1[2](start = PhaseSystem_n.phaseVoltages(VHigh, phi_1)) "Voltage at the winding - primary side"; Modelica.SIunits.Voltage V2[2](start = PhaseSystem_n.phaseVoltages(VLow, phi_2)) "Voltage at the winding - secondary side"; protected parameter Modelica.SIunits.AngularVelocity omega_n = 2*Modelica.Constants.pi*f; parameter Real N = VHigh/VLow "Winding ratio"; parameter Modelica.SIunits.Resistance RBaseHigh = VHigh^2/VABase "Base impedance of the primary side"; parameter Modelica.SIunits.Resistance RBaseLow = VLow^2/VABase "Base impedance of the secondary side"; Modelica.SIunits.Impedance Z1[2] = {RBaseHigh*R1, omega*L1*RBaseHigh/omega_n} "Impedance of the primary side of the transformer"; Modelica.SIunits.Impedance Z2[2] = {RBaseLow*R2, omega*L2*RBaseLow/omega_n} "Impedance of the secondary side of the transformer"; Modelica.SIunits.Impedance Zrm[2] = {RBaseHigh*Rm, 0} "Magnetization impedance - resistance"; Modelica.SIunits.Impedance Zlm[2] = {0, omega*Lm*RBaseHigh/omega_n} "Magnetization impedance - impedence"; Modelica.SIunits.Power P_p[2] = PhaseSystem_p.phasePowers_vi(terminal_p.v, terminal_p.i) "Power transmitted at pin p (secondary)"; Modelica.SIunits.Power P_n[2] = PhaseSystem_n.phasePowers_vi(terminal_n.v, terminal_n.i) "Power transmitted at pin n (primary)"; Modelica.SIunits.Power S_p = Modelica.Fluid.Utilities.regRoot(P_p[1]^2 + P_p[2]^2, delta=0.1) "Apparent power at terminal p"; Modelica.SIunits.Power S_n = Modelica.Fluid.Utilities.regRoot(P_n[1]^2 + P_n[2]^2, delta=0.1) "Apparent power at terminal n"; Modelica.SIunits.AngularVelocity omega "Angular velocity"; Modelica.SIunits.Current Im[2] "Magnetization current"; Modelica.SIunits.Angle theRef "Absolute angle of rotating reference system"; equation assert(sqrt(P_p[1]^2 + P_p[2]^2) <= VABase*1.01, "The load power of the transformer is higher than VABase"); // Angular velocity theRef = PhaseSystem_p.thetaRef(terminal_p.theta); omega = der(theRef); // Efficiency eta = Buildings.Utilities.Math.Functions.smoothMin( x1= Modelica.Fluid.Utilities.regRoot(P_p[1]^2 + P_p[2]^2, delta=0.01)/ Modelica.Fluid.Utilities.regRoot(P_n[1]^2 + P_n[2]^2 + 1e-6, delta=0.01), x2= Modelica.Fluid.Utilities.regRoot(P_n[1]^2 + P_n[2]^2, delta=0.01)/ Modelica.Fluid.Utilities.regRoot(P_p[1]^2 + P_p[2]^2 + 1e-6, delta=0.01), deltaX = 0.01); // Ideal transformation V2 = V1/N; terminal_p.i[1] + (terminal_n.i[1] - Im[1])*N = 0; terminal_p.i[2] + (terminal_n.i[2] - Im[2])*N = 0; // Magnetization current if magEffects then Im = Buildings.Electrical.PhaseSystems.OnePhase.divide(V1, Zrm) + Buildings.Electrical.PhaseSystems.OnePhase.divide(V1, Zlm); else Im = zeros(2); end if; // Losses due to the impedance - primary side terminal_n.v = V1 + Buildings.Electrical.PhaseSystems.OnePhase.product( terminal_n.i, Z1); // Losses due to the impedance - secondary side terminal_p.v = V2 + Buildings.Electrical.PhaseSystems.OnePhase.product( terminal_p.i, Z2); // Loss of power PLoss = P_p + P_n; // The two sides have the same reference angle terminal_p.theta = terminal_n.theta; if ground_1 then Connections.potentialRoot(terminal_n.theta); end if; if ground_2 then Connections.root(terminal_p.theta); end if; end ACACTransformerFull;

Buildings.Electrical.AC.OnePhase.Conversion.ACDCConverter Buildings.Electrical.AC.OnePhase.Conversion.ACDCConverter

AC DC converter

Buildings.Electrical.AC.OnePhase.Conversion.ACDCConverter

Information

This is an AC/DC converter, based on a power balance between both circuit sides. The parameter conversionFactor defines the ratio between the RMS voltages as

VDC = conversionFactor VAC,

where VDC is the voltage of the DC circuit and VAC is the RMS voltage at the primary side of the transformer.

The loss of the converter is proportional to the power transmitted. The parameter eta is the efficiency of the transfer. The loss is computed as

Ploss = (1-η) Ptr

where Ptr is the power transmitted. The model is bi-directional and the power can flow from both the primary to the secondary side and vice-versa. Furthermore, reactive power on both side are set to 0.

Note:

This model is derived from Modelica.Electrical.QuasiStationary.SinglePhase.Utilities.IdealACDCConverter.

Extends from Buildings.Electrical.Interfaces.PartialConversion (Model representing a generic two port system for conversion).

Parameters

TypeNameDefaultDescription
replaceable package PhaseSystem_pPartialPhaseSystemPhase system of terminal p
replaceable package PhaseSystem_nPartialPhaseSystemPhase system of terminal n
RealconversionFactor Ratio of DC voltage / AC RMS voltage
Realeta Converter efficiency, pLoss = (1-eta) * Ptr
Ground
AC side
Booleanground_ACfalseConnect AC side of converter to ground
DC side
Booleanground_DCtrueConnect DC side of converter to ground

Connectors

TypeNameDescription
replaceable package PhaseSystem_pPhase system of terminal p
replaceable package PhaseSystem_nPhase system of terminal n

Modelica definition

model ACDCConverter "AC DC converter" extends Buildings.Electrical.Interfaces.PartialConversion( redeclare package PhaseSystem_p = PhaseSystems.TwoConductor, redeclare package PhaseSystem_n = PhaseSystems.OnePhase, redeclare replaceable Interfaces.Terminal_n terminal_n constrainedby Interfaces.Terminal_n( i(start = zeros(PhaseSystem_n.n), each stateSelect = StateSelect.prefer)), redeclare DC.Interfaces.Terminal_p terminal_p( i(start = zeros(PhaseSystem_p.n), each stateSelect = StateSelect.prefer))); parameter Real conversionFactor(min = Modelica.Constants.eps) "Ratio of DC voltage / AC RMS voltage"; parameter Real eta(min=0, max=1) "Converter efficiency, pLoss = (1-eta) * Ptr"; Modelica.SIunits.Power PLoss "Loss power"; parameter Boolean ground_AC = false "Connect AC side of converter to ground"; parameter Boolean ground_DC = true "Connect DC side of converter to ground"; protected PhaseSystem_p.Current i_dc "DC current"; PhaseSystem_p.Voltage v_dc "DC voltage"; Modelica.SIunits.Power P_p[2] = PhaseSystem_p.phasePowers_vi(terminal_p.v, terminal_p.i) "Power transmitted at pin p (secondary)"; Modelica.SIunits.Power P_n[2](each start=0)= PhaseSystem_n.phasePowers_vi(terminal_n.v, terminal_n.i) "Power transmitted at pin n (primary)"; equation //voltage relation v_p = v_n*conversionFactor; // Power losses PLoss = (1-eta)* Buildings.Utilities.Math.Functions.spliceFunction(P_p[1], P_n[1], i_p, deltax=0.1); P_n + P_p = {PLoss, 0}; if ground_AC then Connections.potentialRoot(terminal_n.theta); end if; if ground_DC then v_dc = 0; Connections.root(terminal_p.theta); else i_dc = 0; Connections.potentialRoot(terminal_p.theta); end if; v_dc = terminal_p.v[2]; sum(terminal_p.i) + i_dc = 0; end ACDCConverter;