2.2.2. Running simulations
2.2.2.1. Initializing
Create an instance of pyramses.sim:
import pyramses
ram = pyramses.sim()
See Running a basic simulation for a complete example and Full list of functions for the full API reference.
2.2.2.2. Start, pause, and continue
A properly configured test case (see Building a scenario and loading data) is required before running a simulation.
Run the simulation to completion:
ram.execSim(case)
Start the simulation and pause at a specific time:
ram.execSim(case, 10.0) # start and pause at t = 10 s
Continue a paused simulation:
ram.contSim(20.0) # resume until t = 20 s
ram.contSim(ram.getSimTime() + 60.0) # advance by 60 s from current time
ram.contSim(ram.getInfTime()) # run to the end of the time horizon
2.2.2.3. End simulation early
To terminate the simulation before reaching the time horizon:
ram.endSim()
2.2.2.4. Querying simulation state
When the simulation is paused, the following methods query the current system state.
Time information:
t = ram.getSimTime() # current simulation time (seconds)
t_inf = ram.getInfTime() # use as "infinity" to run to end
Component names:
buses = ram.getAllCompNames('BUS') # list of all bus names
gens = ram.getAllCompNames('SYNC') # list of all generator names
injs = ram.getAllCompNames('INJ') # list of all injector names
Supported component types: BUS, SYNC, INJ, DCTL, BRANCH, TWOP, SHUNT, LOAD.
Bus voltages and branch flows:
ram.execSim(case, 10.0)
busNames = ['g1', 'g2', '4032']
voltages = ram.getBusVolt(busNames) # voltage magnitudes (pu)
phases = ram.getBusPha(busNames) # voltage phase angles (rad)
powers = ram.getBranchPow(['1041-01']) # [[P_from, Q_from, P_to, Q_to]] (MW, Mvar)
Component observables and parameters:
comp_type = ['INJ', 'EXC', 'TOR']
comp_name = ['L_11', 'g2', 'g3']
obs_name = ['P', 'vf', 'Pm']
obs = ram.getObs(comp_type, comp_name, obs_name)
comp_type = ['EXC', 'EXC']
comp_name = ['g1', 'g2']
prm_name = ['V0', 'KPSS']
prms = ram.getPrm(comp_type, comp_name, prm_name)
Supported model types for observables and parameters: EXC (exciter), TOR (governor), INJ (injector), TWOP (two-port), DCTL (discrete controller), SYN (synchronous generator).
2.2.2.5. Adding disturbances at runtime
Disturbances can be added dynamically while the simulation is paused, enabling interactive scenario analysis:
ram.execSim(case, 80.0)
ram.addDisturb(100.0, 'BREAKER SYNC_MACH g7 0') # trip generator g7 at t=100 s
ram.addDisturb(100.0, 'FAULT BUS 4032 0. 0.') # 3-phase fault at bus 4032
ram.addDisturb(100.1, 'CLEAR BUS 4032') # clear fault
ram.addDisturb(100.0, 'CHGPRM DCTL 1-1041 Vsetpt -0.05 0') # step change in setpoint
ram.contSim(ram.getInfTime())
For disturbance syntax, see Disturbance description.
2.2.2.6. Exporting the Jacobian matrix
The system Jacobian in descriptor form (E, A matrices) can be exported at any pause point:
ram.execSim(case, 10.0)
ram.getJac() # writes jac_val.dat, jac_eqs.dat, jac_var.dat, jac_struc.dat
These files can then be used for small-signal stability analysis. See the RAMSES Eigenanalysis tool.
Note
Set $OMEGA_REF SYN ; in the solver settings data file when exporting the Jacobian for eigenanalysis.
2.2.2.7. Full list of functions
- class pyramses.sim(custLibDir=None)[source]
Simulation class.
- addDisturb(t_dist, disturb)[source]
Add a new disturbance at a specific time. Follows the same structure as the disturbances in the dst files.
- Parameters
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> case = pyramses.cfg("cmd.txt") >>> ram.execSim(case, 80.0) # simulate until 80 seconds and pause >>> ram.addDisturb(100.000, 'CHGPRM DCTL 1-1041 Vsetpt -0.05 0') # Decrease the setpoint of the DCTL by 0.015 pu, at t=100 s >>> ram.addDisturb(100.000, 'CHGPRM DCTL 2-1042 Vsetpt -0.05 0') >>> ram.addDisturb(100.000, 'CHGPRM DCTL 3-1043 Vsetpt -0.05 0') >>> ram.addDisturb(100.000, 'CHGPRM DCTL 4-1044 Vsetpt -0.05 0') >>> ram.addDisturb(100.000, 'CHGPRM DCTL 5-1045 Vsetpt -0.05 0') >>> ram.contSim(ram.getInfTime()) # continue the simulation
- addObserv(string)[source]
Add an element to be observed
- Parameters
string (str) – the element with proper format
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> case = pyramses.cfg("cmd.txt") # command file without any observables >>> ram.execSim(case, 0.0) # start >>> traj_filenm = 'obs.trj' >>> ram.initObserv(traj_filenm) >>> string = 'BUS *' # monitor all buses >>> ram.addObserv(string)
- contSim(pause=None)[source]
Continue the simulation until the given time.
The pause argument is optional and it gives the time until the simulation will stop.
- defineSS(ssID, filter1, filter2, filter3)[source]
Define a subsytem using three filters. The resulting list is an intersection of the filters.
- Parameters
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> case = pyramses.cfg("cmd.txt") >>> ram.execSim(case, 0.0) >>> ram.defineSS(1, ['735'], [], []) # SS 1 with all buses at 735 kV, no zones, no list of buses
Note
An empty filter means it is deactivated and discarded.
- execSim(cmd, pause=None)[source]
Run a simulation.
- Parameters
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> case = pyramses.cfg("cmd.txt") >>> ram.execSim(case) # start simulation
Note
If you have an existing command file, you can pass it to the simulator as pyramses.cfg(“cmd.txt”)
- finalObserv()[source]
Finalize observable selection, allocate buffer, and write header of trajectory file
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> case = pyramses.cfg("cmd.txt") # command file without any observables >>> ram.execSim(case, 0.0) # start >>> traj_filenm = 'obs.trj' >>> ram.initObserv(traj_filenm) >>> string = 'BUS *' # monitor all buses >>> ram.addObserv(string) >>> ram.finalObserv()
- getAllCompNames(comp_type)[source]
Get list of all components of type comp_type
- Parameters
comp_type (str) – the component type (BUS, SYNC, INJ, DCTL, BRANCH, TWOP, SHUNT, LOAD)
- Returns
list of component names
- Return type
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> case = pyramses.cfg("cmd.txt") >>> ram.execSim(case, 0) # start simulation paused >>> ram.getAllCompNames("BUS") ['B1', 'B2', 'B3']
- getBranchCur(branchName)[source]
Return the currents of a list of branches
- Parameters
- Returns
list of branch currents. These are x-y components at the origin and extremity respectively (ix_orig, iy_orig, ix_extr, iy_extr)
- Return type
list of floats
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> case = pyramses.cfg("cmd.txt") >>> ram.execSim(case, 10.0) # simulate until 10 seconds and pause >>> branchName = ['1011-1013','1012-1014','1021-1022'] >>> ram.getBranchCur(branchName)
- getBranchPow(branchName)[source]
Return the active and reactive powers of a list of branches
- Parameters
- Returns
list of branch powers. These are active and reactive power at the origin and extremity respectively (p_orig, q_orig, p_extr, q_extr)
- Return type
list of floats
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> case = pyramses.cfg("cmd.txt") >>> ram.execSim(case, 10.0) # simulate until 10 seconds and pause >>> branchName = ['1011-1013','1012-1014','1021-1022'] >>> ram.getBranchPow(branchName)
- getBusPha(busNames)[source]
Return the voltage phase of a list of buses
- Parameters
- Returns
list of bus voltage phase
- Return type
list of floats
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> case = pyramses.cfg("cmd.txt") >>> ram.execSim(case, 10.0) # simulate until 10 seconds and pause >>> buses = ['g1','g2','g3'] >>> ram.getBusPha(buses) [0.0000000000000000, 10.0615442362180327, 11.064702686997689]
- getBusVolt(busNames)[source]
Return the voltage magnitude of a list of buses
- Parameters
- Returns
list of bus voltage magnitudes
- Return type
list of floats
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> case = pyramses.cfg("cmd.txt") >>> ram.execSim(case, 10.0) # simulate until 10 seconds and pause >>> buses = ['g1','g2','g3'] >>> ram.getBusVolt(buses) [1.0736851673414456, 1.0615442362180327, 1.064702686997689]
- getCompName(comp_type, num)[source]
Return the name of the num:superscript:th component of type comp_type.
- Parameters
- Returns
The component name
- Return type
srt
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> case = pyramses.cfg("cmd.txt") >>> ram.execSim(case, 0) # start simulation paused >>> ram.getCompName("BUS",1) 'B1'
- getEndSim()[source]
Check if the simulation has ended
- Returns
0 -> simulation is still working, 1 -> simulation has ended
- Return type
integer
- getInfTime()[source]
Get the maximum representable double from the simulator ()
- Returns
maximum representable double
- Return type
- getJac()[source]
Gets the Jacobian matrix written in files
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> case = pyramses.cfg("cmd.txt") >>> ram.execSim(case, 0) # start simulation paused >>> ram.getJac()
- getLastErr()[source]
Return the last error message issued by RAMSES.
- Return type
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> case = pyramses.cfg("cmd.txt") >>> ram.execSim(case, 0) # start simulation paused >>> ram.getLastErr() 'This is an error msg'
- getObs(comp_type, comp_name, obs_name)[source]
Get the value of a named observable.
- Parameters
- Returns
list of observable values
- Return type
list of floats
Note
For the synchronous generator (‘SYN’) the accepted obs_name values are - ‘P’: Active power (MW) - ‘Q’: Reactive power (MVAr) -‘Omega’: Machine speed (pu) - ‘S’: Apparent power (MVA) - ‘SNOM’: Nominal apparent power (MVA) - ‘PNOM’: Nominal active power (MW)
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> case = pyramses.cfg("cmd.txt") >>> ram.execSim(case, 10.0) # simulate until 10 seconds and pause >>> comp_type = ['INJ','EXC','TOR'] >>> comp_name = ['L_11','g2','g3'] >>> obs_name = ['P','vf','Pm'] >>> ram.getObs(comp_type,comp_name, obs_name) [199.73704732259995, 1.3372945282585218, 0.816671075203357]
- getPrm(comp_type, comp_name, prm_name)[source]
Get the value of a named parameter
- Parameters
- Returns
list of parameter values
- Return type
list of floats
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> case = pyramses.cfg("cmd.txt") >>> ram.execSim(case, 10.0) # simulate until 10 seconds and pause >>> comp_type = ['EXC','EXC','EXC'] >>> comp_name = ['g1','g2','g3'] >>> prm_name = ['V0','V0','V0'] >>> ram.getPrm(comp_type,comp_name, prm_name) [1.0736851673414456, 1.0615442362180327, 1.064702686997689]
- getPrmNames(comp_type, comp_name)[source]
Get the named parameters of a model
- Parameters
- Returns
list of parameter names
- Return type
list of lists of strings
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> case = pyramses.cfg("cmd.txt") >>> ram.execSim(case, 0.0) # initialize and wait >>> comp_type = ['EXC','EXC','EXC'] >>> comp_name = ['g1','g2','g3'] # name of synchronous machines >>> ram.getPrmNames(comp_type,comp_name)
- getSS(ssID)[source]
Retrieve the buses of a subsytem.
>>> import pyramses >>> ram = pyramses.sim() >>> case = pyramses.cfg("cmd.txt") >>> ram.execSim(case, 0.0) >>> ram.defineSS(1, ['735'], [], []) # SS 1 with all buses at 735 kV >>> ram.getSS(1) # get list of buses in SS 1
- getSimTime()[source]
Get the current simulation time
- Returns
current simulation time in RAMSES.
- Return type
- getTrfoSS(ssID, location, in_service, rettype)[source]
Retrieve transformer information of subsystem after applying some filters
- Parameters
- Returns
list of transformer names
- Return type
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> case = pyramses.cfg("cmd.txt") >>> ram.execSim(case, 0.0) >>> ram.defineSS(1, ['735'], [], []) # SS 1 with all buses at 735 kV >>> ram.getTrfoSS(1,3,2,'Status')
Note
Tap is not implemented yet.
Todo
Implement Taps after discussing with Lester.
- get_MDL_no()[source]
Unload external DLL file with user defined models. Should be in current directory or absolute path.
- Returns
list of observable values
- Return type
list of floats
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> ram = pyramses.load_MDL("MDLs.dll") >>> case = pyramses.cfg("cmd.txt") >>> ram.execSim(case) # simulate >>> ram = pyramses.get_MDL_no()
- initObserv(traj_filenm)[source]
Initialize observable selection (structure and trajectory file)
- Parameters
traj_filenm (str) – the file to save the trajectory at
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> case = pyramses.cfg("cmd.txt") # command file without any observables >>> ram.execSim(case, 0.0) # start >>> traj_filenm = 'obs.trj' >>> ram.initObserv(traj_filenm)
- load_MDL()[source]
Load external DLL file with user defined models. Should be in current directory or absolute path.
- Parameters
MDLName (str) – path to file
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> ram = pyramses.load_MDL("MDLs.dll") >>> case = pyramses.cfg("cmd.txt") >>> ram.execSim(case) # simulate
- pauseSim(t_pause)[source]
Pause the simulation at the t_pause time
- Parameters
t_pause (float) – pause time
- unload_MDL()[source]
Unload external DLL file with user defined models. Should be in current directory or absolute path.
- Parameters
MDLName (str) – path to file
- Example
>>> import pyramses >>> ram = pyramses.sim() >>> ram = pyramses.load_MDL("MDLs.dll") >>> case = pyramses.cfg("cmd.txt") >>> ram.execSim(case) # simulate >>> ram = pyramses.unload_MDL("MDLs.dll")