Documentation Index
Fetch the complete documentation index at: https://mintlify.com/GridOPTICS/GridPACK/llms.txt
Use this file to discover all available pages before exploring further.
GridPACK applications are built around three core class hierarchies — Bus, Branch, and Factory — that together describe network physics, drive matrix and vector assembly, and control the parallel execution loop. This page walks through each piece in detail, using the power flow application in src/applications/examples/powerflow as a concrete reference.
Overview of the three main tasks
When you build a new GridPACK application you must provide:
- A Bus class — inherits from
BaseBusComponent, holds per-bus state, and implements MatVecInterface methods so the mapper can assemble matrices and vectors.
- A Branch class — inherits from
BaseBranchComponent, holds per-branch state, and implements the off-diagonal matrix block methods.
- A Factory class — inherits from
BaseFactory, calls load() on every component to populate fields from parsed network data, and exposes any extra setup routines your application needs (e.g. setYBus).
The network type itself requires no hand-written code — it is a typedef over the templated BaseNetwork:
typedef gridpack::network::BaseNetwork<MyBus, MyBranch> MyNetwork;
Writing the Bus class
#ifndef MY_BUS_HPP
#define MY_BUS_HPP
#include "gridpack/include/gridpack.hpp"
namespace gridpack {
namespace myapp {
class MyBus : public gridpack::component::BaseBusComponent {
public:
MyBus();
~MyBus();
// Populate fields from the parsed DataCollection
void load(const boost::shared_ptr<
gridpack::component::DataCollection> &data) override;
// MatVecInterface — diagonal blocks
bool matrixDiagSize(int *isize, int *jsize) const override;
bool matrixDiagValues(gridpack::ComplexType *values) override;
// MatVecInterface — vector contribution
bool vectorSize(int *isize) const override;
bool vectorValues(gridpack::ComplexType *values) override;
void setValues(gridpack::ComplexType *values) override;
// Exchange buffer for ghost bus updates
int getXCBufSize() override;
void setXCBuf(void *buf) override;
private:
double p_voltage, p_angle;
double p_ybusr, p_ybusi;
int p_mode;
// Pointers into the exchange buffer
double *p_vMag_ptr;
double *p_vAng_ptr;
friend class boost::serialization::access;
template<class Archive>
void serialize(Archive &ar, const unsigned int version)
{
ar & boost::serialization::base_object<
gridpack::component::BaseBusComponent>(*this)
& p_voltage & p_angle;
}
};
} // namespace myapp
} // namespace gridpack
#endif
Implementing load()
load() is called by the base factory on every bus and branch after the network is partitioned. It receives the DataCollection object that was populated by the parser and must copy the relevant fields into the component’s private members.
void gridpack::myapp::MyBus::load(
const boost::shared_ptr<gridpack::component::DataCollection> &data)
{
data->getValue(CASE_SBASE, &p_sbase);
data->getValue(BUS_VOLTAGE_ANG, &p_angle);
data->getValue(BUS_VOLTAGE_MAG, &p_voltage);
// Convert angle from degrees to radians
double pi = 4.0 * atan(1.0);
p_angle = p_angle * pi / 180.0;
int itype;
data->getValue(BUS_TYPE, &itype);
if (itype == 3) setReferenceBus(true);
// Multi-valued fields: retrieve by index
int ngen;
if (data->getValue(GENERATOR_NUMBER, &ngen)) {
for (int i = 0; i < ngen; i++) {
double pg, qg;
int gstatus;
bool ok = data->getValue(GENERATOR_PG, &pg, i);
ok = ok && data->getValue(GENERATOR_QG, &qg, i);
ok = ok && data->getValue(GENERATOR_STAT, &gstatus, i);
if (ok) {
p_pg.push_back(pg);
p_qg.push_back(qg);
p_gstatus.push_back(gstatus);
}
}
}
}
The keys (CASE_SBASE, BUS_VOLTAGE_ANG, etc.) are preprocessor symbols defined in src/parser/dictionary.hpp. Include that header in your component file.
Ghost buses created during partitioning receive a copy of the original DataCollection, so load() runs correctly on them without any extra data exchange.
The MatVecInterface — matrix and vector methods
BaseBusComponent inherits MatVecInterface, which is the contract the mapper classes use to build distributed PETSc matrices and vectors. You must implement the following groups of methods.
Diagonal block methods (Bus)
matrixDiagSize() returns the number of rows and columns that this bus contributes to the matrix diagonal block. Return false if the bus makes no contribution (e.g. a reference bus in certain formulations).
bool gridpack::myapp::MyBus::matrixDiagSize(int *isize, int *jsize) const
{
if (p_mode == Jacobian) {
*isize = 2; *jsize = 2;
return true;
} else if (p_mode == YBus) {
*isize = 1; *jsize = 1;
return true;
}
return false;
}
matrixDiagValues() fills the values array in row-major order. For a 2×2 block, values[0] is position (0,0), values[1] is (0,1), values[2] is (1,0), and values[3] is (1,1). Values are always ComplexType; set the imaginary part to zero for real matrices.
bool gridpack::myapp::MyBus::matrixDiagValues(
gridpack::ComplexType *values)
{
if (p_mode == YBus) {
values[0] = gridpack::ComplexType(p_ybusr, p_ybusi);
return true;
} else if (p_mode == Jacobian) {
if (!getReferenceBus()) {
values[0] = -p_Qinj - p_ybusi * p_v * p_v;
values[1] = p_Pinj - p_ybusr * p_v * p_v;
values[2] = p_Pinj / p_v + p_ybusr * p_v;
values[3] = p_Qinj / p_v - p_ybusi * p_v;
return true;
} else {
// Identity block for reference bus
values[0] = 1.0; values[1] = 0.0;
values[2] = 0.0; values[3] = 1.0;
return true;
}
}
return false;
}
Off-diagonal block methods (Branch)
These are implemented in the Branch class. matrixForwardSize() and matrixForwardValues() describe the block at row = bus1’s index, column = bus2’s index. The corresponding Reverse functions describe the transposed position.
bool gridpack::myapp::MyBranch::matrixForwardSize(
int *isize, int *jsize) const
{
if (p_mode == Jacobian) {
auto *bus1 = dynamic_cast<MyBus*>(getBus1().get());
auto *bus2 = dynamic_cast<MyBus*>(getBus2().get());
if (!bus1->getReferenceBus() && !bus2->getReferenceBus()) {
*isize = 2; *jsize = 2;
return true;
}
return false;
} else if (p_mode == YBus) {
*isize = 1; *jsize = 1;
return true;
}
return false;
}
matrixForwardValues() fills the 2×2 off-diagonal Jacobian block. The branch must read current voltages from its terminal buses using accessor methods defined in the Bus class.
Vector methods
vectorSize() and vectorValues() follow the same convention as the diagonal matrix methods but for right-hand side or solution vectors. setValues() is the reverse path — it is called by the mapper after a linear solve to push the solution increments back into the bus state:
void gridpack::myapp::MyBus::setValues(gridpack::ComplexType *values)
{
// Decrement angle and voltage magnitude by the Newton step
p_a -= real(values[0]);
p_v -= real(values[1]);
// Write back into the exchange buffer so ghost buses get updated values
*p_vAng_ptr = p_a;
*p_vMag_ptr = p_v;
}
Exchange buffers
To propagate bus state to ghost copies on neighbouring processors, implement getXCBufSize() and setXCBuf():
int gridpack::myapp::MyBus::getXCBufSize() { return 2 * sizeof(double); }
void gridpack::myapp::MyBus::setXCBuf(void *buf)
{
p_vAng_ptr = static_cast<double*>(buf);
p_vMag_ptr = p_vAng_ptr + 1;
*p_vAng_ptr = p_a;
*p_vMag_ptr = p_v;
}
Writing the Factory class
The Factory inherits BaseFactory<MyNetwork>. The base class provides load(), setComponents(), and setExchange(). Add application-specific methods that loop over buses and branches to invoke component-level routines:
void gridpack::myapp::MyFactory::setYBus()
{
int numBus = p_network->numBuses();
int numBranch = p_network->numBranches();
for (int i = 0; i < numBus; i++) p_network->getBus(i).get()->setYBus();
for (int i = 0; i < numBranch; i++) p_network->getBranch(i).get()->setYBus();
}
The getBus() and getBranch() methods return typed pointers to your concrete MyBus / MyBranch objects, so you do not need a dynamic_cast to call application-specific methods from the factory.
The application loop
The top-level driver follows a fixed sequence: initialise the environment, read the network, partition, load, map to matrices/vectors, solve, and write results.
Initialise the environment
Create a gridpack::Environment as the first object in main(). Its constructor initialises MPI, Global Arrays, and the math library. Place all GridPACK objects inside a nested scope so their destructors run before env goes out of scope.int main(int argc, char **argv)
{
gridpack::Environment env(argc, argv);
{
gridpack::myapp::MyApp app;
app.execute();
}
// Destructors for all GridPACK objects called here,
// before MPI_Finalize is invoked by ~Environment()
}
Create the network and parse the input
gridpack::parallel::Communicator world;
boost::shared_ptr<MyNetwork> network(new MyNetwork(world));
gridpack::utility::Configuration *config =
gridpack::utility::Configuration::configuration();
config->open("input.xml", world);
auto *cursor = config->getCursor("Configuration.MyApp");
std::string filename;
cursor->get("networkConfiguration", &filename);
gridpack::parser::PTI23_parser<MyNetwork> parser(network);
parser.parse(filename.c_str());
Use PTI33_parser for PSS/E version 33 files. The parser key in input.xml controls which format version is active.
Partition and initialise components
network->partition();
MyFactory factory(network);
factory.load(); // calls MyBus::load() and MyBranch::load() on all components
factory.setComponents(); // establishes neighbour pointers and mapper indices
factory.setExchange(); // allocates exchange buffers (calls getXCBufSize / setXCBuf)
network->initBusUpdate();
factory.setYBus(); // application-specific pre-computation
Build matrices and vectors
factory.setMode(RHS);
gridpack::mapper::BusVectorMap<MyNetwork> vMap(network);
boost::shared_ptr<gridpack::math::Vector> PQ = vMap.mapToVector();
factory.setMode(Jacobian);
gridpack::mapper::FullMatrixMap<MyNetwork> jMap(network);
boost::shared_ptr<gridpack::math::Matrix> J = jMap.mapToMatrix();
Solve
gridpack::math::LinearSolver solver(*J);
solver.configure(cursor);
boost::shared_ptr<gridpack::math::Vector> X(PQ->clone());
double tol = 1.0;
int iter = 0, max_iter = 100;
double tolerance = 1.0e-6;
while (real(tol) > tolerance && iter < max_iter) {
factory.setMode(RHS);
vMap.mapToBus(X);
network->updateBuses();
vMap.mapToVector(PQ);
factory.setMode(Jacobian);
jMap.mapToMatrix(J);
X->zero();
solver.solve(*PQ, *X);
tol = PQ->normInfinity();
iter++;
}
Write results
gridpack::serial_io::SerialBusIO<MyNetwork> busIO(128, network);
busIO.header("\n Bus Voltages and Phase Angles\n");
busIO.write();
Each bus writes its own line via the serialWrite() method you define in MyBus.
CMakeLists.txt for a custom application
Use the installed GridPACK.cmake module to pick up all library and include paths automatically. A minimal CMakeLists.txt for a standalone application looks like this:
cmake_minimum_required(VERSION 3.10)
if (NOT GRIDPACK_DIR)
set(GRIDPACK_DIR "$ENV{GRIDPACK_DIR}"
CACHE PATH "GridPACK installation directory")
endif()
include("${GRIDPACK_DIR}/lib/GridPACK.cmake")
project(MyGridPACKApp)
enable_language(CXX)
gridpack_setup()
add_definitions(${GRIDPACK_DEFINITIONS})
include_directories(BEFORE ${CMAKE_CURRENT_SOURCE_DIR})
include_directories(BEFORE ${GRIDPACK_INCLUDE_DIRS})
add_executable(myapp.x
myapp_main.cpp
myapp_app.cpp
myapp_factory.cpp
)
target_link_libraries(myapp.x ${GRIDPACK_LIBS})
GridPACK must be installed (not just built in-tree) before this pattern works. Run cmake --install with the same CMAKE_INSTALL_PREFIX that you set as GRIDPACK_DIR.
If you are building inside the GridPACK source tree, follow the in-tree pattern from src/applications/examples/powerflow/CMakeLists.txt, which lists individual component libraries (gridpack_pfmatrix_components, gridpack_math, gridpack_parallel, etc.) rather than using the installed module.
Where to find working examples
| Example | Location |
|---|
| Hello World (minimal bus/branch) | src/applications/examples/hello_world/ |
| Resistor grid (linear solve) | src/applications/examples/resistor_grid/ |
| Power flow (full Newton-Raphson) | src/applications/examples/powerflow/ |
| Contingency analysis (multi-task) | src/applications/examples/contingency_analysis/ |
| Power flow module (reusable module) | src/applications/modules/powerflow/ |
| PF bus/branch components | src/applications/components/pf_matrix/ |