xref: /petsc/src/snes/tutorials/network/power/PFReadData.c (revision 31d78bcd2b98084dc1368b20eb1129c8b9fb39fe)
1c4762a1bSJed Brown #include "petscmat.h"
2c4762a1bSJed Brown #include "power.h"
3c4762a1bSJed Brown #include <string.h>
4c4762a1bSJed Brown #include <ctype.h>
5c4762a1bSJed Brown 
PFReadMatPowerData(PFDATA * pf,char * filename)6d71ae5a4SJacob Faibussowitsch PetscErrorCode PFReadMatPowerData(PFDATA *pf, char *filename)
7d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown   FILE        *fp;
9c4762a1bSJed Brown   VERTEX_Power Bus;
10c4762a1bSJed Brown   LOAD         Load;
11c4762a1bSJed Brown   GEN          Gen;
12c4762a1bSJed Brown   EDGE_Power   Branch;
13c4762a1bSJed Brown   PetscInt     line_counter   = 0;
14c4762a1bSJed Brown   PetscInt     bus_start_line = -1, bus_end_line = -1; /* xx_end_line points to the next line after the record ends */
15c4762a1bSJed Brown   PetscInt     gen_start_line = -1, gen_end_line = -1;
16c4762a1bSJed Brown   PetscInt     br_start_line = -1, br_end_line = -1;
17c4762a1bSJed Brown   char         line[MAXLINE];
18c4762a1bSJed Brown   PetscInt     loadi = 0, geni = 0, bri = 0, busi = 0, i, j;
19c4762a1bSJed Brown   int          extbusnum, bustype_i;
20c4762a1bSJed Brown   double       Pd, Qd;
21c4762a1bSJed Brown   PetscInt     maxbusnum = -1, intbusnum, *busext2intmap, genj, loadj;
22c4762a1bSJed Brown   GEN          newgen;
23c4762a1bSJed Brown   LOAD         newload;
24c4762a1bSJed Brown 
25c4762a1bSJed Brown   PetscFunctionBegin;
26c4762a1bSJed Brown   fp = fopen(filename, "r");
27c4762a1bSJed Brown   /* Check for valid file */
2828b400f6SJacob Faibussowitsch   PetscCheck(fp, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Can't open Matpower data file %s", filename);
29c4762a1bSJed Brown   pf->nload = 0;
30c4762a1bSJed Brown   while (fgets(line, MAXLINE, fp)) {
31c4762a1bSJed Brown     if (strstr(line, "mpc.bus = [")) bus_start_line = line_counter + 1;                     /* Bus data starts from next line */
32c4762a1bSJed Brown     if (strstr(line, "mpc.gen") && gen_start_line == -1) gen_start_line = line_counter + 1; /* Generator data starts from next line */
33c4762a1bSJed Brown     if (strstr(line, "mpc.branch")) br_start_line = line_counter + 1;                       /* Branch data starts from next line */
34c4762a1bSJed Brown     if (strstr(line, "];")) {
35c4762a1bSJed Brown       if (bus_start_line != -1 && bus_end_line == -1) bus_end_line = line_counter;
36c4762a1bSJed Brown       if (gen_start_line != -1 && gen_end_line == -1) gen_end_line = line_counter;
37c4762a1bSJed Brown       if (br_start_line != -1 && br_end_line == -1) br_end_line = line_counter;
38c4762a1bSJed Brown     }
39c4762a1bSJed Brown 
40c4762a1bSJed Brown     /* Count the number of pq loads */
41c4762a1bSJed Brown     if (bus_start_line != -1 && line_counter >= bus_start_line && bus_end_line == -1) {
42c4762a1bSJed Brown       sscanf(line, "%d %d %lf %lf", &extbusnum, &bustype_i, &Pd, &Qd);
43c4762a1bSJed Brown       if (!((Pd == 0.0) && (Qd == 0.0))) pf->nload++;
44c4762a1bSJed Brown       if (extbusnum > maxbusnum) maxbusnum = extbusnum;
45c4762a1bSJed Brown     }
46c4762a1bSJed Brown     line_counter++;
47c4762a1bSJed Brown   }
48c4762a1bSJed Brown   fclose(fp);
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   pf->nbus    = bus_end_line - bus_start_line;
51c4762a1bSJed Brown   pf->ngen    = gen_end_line - gen_start_line;
52c4762a1bSJed Brown   pf->nbranch = br_end_line - br_start_line;
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(pf->nbus, &pf->bus));
559566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(pf->ngen, &pf->gen));
569566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(pf->nload, &pf->load));
579566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(pf->nbranch, &pf->branch));
589371c9d4SSatish Balay   Bus    = pf->bus;
599371c9d4SSatish Balay   Gen    = pf->gen;
609371c9d4SSatish Balay   Load   = pf->load;
619371c9d4SSatish Balay   Branch = pf->branch;
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   /* Setting pf->sbase to 100 */
64c4762a1bSJed Brown   pf->sbase = 100.0;
65c4762a1bSJed Brown 
669566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxbusnum + 1, &busext2intmap));
67c4762a1bSJed Brown   for (i = 0; i < maxbusnum + 1; i++) busext2intmap[i] = -1;
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   fp = fopen(filename, "r");
70c4762a1bSJed Brown   /* Reading data */
71c4762a1bSJed Brown   for (i = 0; i < line_counter; i++) {
7208401ef6SPierre Jolivet     PetscCheck(fgets(line, MAXLINE, fp), PETSC_COMM_SELF, PETSC_ERR_SUP, "File is incorrectly formatted");
73c4762a1bSJed Brown 
74c4762a1bSJed Brown     if ((i >= bus_start_line) && (i < bus_end_line)) {
75c4762a1bSJed Brown       double gl, bl, vm, va, basekV;
76c4762a1bSJed Brown       int    bus_i, ide, area;
77c4762a1bSJed Brown       /* Bus data */
78c4762a1bSJed Brown       sscanf(line, "%d %d %lf %lf %lf %lf %d %lf %lf %lf", &bus_i, &ide, &Pd, &Qd, &gl, &bl, &area, &vm, &va, &basekV);
799371c9d4SSatish Balay       Bus[busi].bus_i                = bus_i;
809371c9d4SSatish Balay       Bus[busi].ide                  = ide;
819371c9d4SSatish Balay       Bus[busi].area                 = area;
829371c9d4SSatish Balay       Bus[busi].gl                   = gl;
839371c9d4SSatish Balay       Bus[busi].bl                   = bl;
849371c9d4SSatish Balay       Bus[busi].vm                   = vm;
859371c9d4SSatish Balay       Bus[busi].va                   = va;
869371c9d4SSatish Balay       Bus[busi].basekV               = basekV;
87c4762a1bSJed Brown       Bus[busi].internal_i           = busi;
88c4762a1bSJed Brown       busext2intmap[Bus[busi].bus_i] = busi;
89c4762a1bSJed Brown 
90c4762a1bSJed Brown       if (!((Pd == 0.0) && (Qd == 0.0))) {
91c4762a1bSJed Brown         Load[loadi].bus_i                 = Bus[busi].bus_i;
92c4762a1bSJed Brown         Load[loadi].status                = 1;
93c4762a1bSJed Brown         Load[loadi].pl                    = Pd;
94c4762a1bSJed Brown         Load[loadi].ql                    = Qd;
95c4762a1bSJed Brown         Load[loadi].area                  = Bus[busi].area;
96c4762a1bSJed Brown         Load[loadi].internal_i            = busi;
97c4762a1bSJed Brown         Bus[busi].lidx[Bus[busi].nload++] = loadi;
9808401ef6SPierre Jolivet         PetscCheck(Bus[busi].nload <= NLOAD_AT_BUS_MAX, PETSC_COMM_SELF, PETSC_ERR_SUP, "Exceeded maximum number of loads allowed at bus");
99c4762a1bSJed Brown         loadi++;
100c4762a1bSJed Brown       }
101c4762a1bSJed Brown       busi++;
102c4762a1bSJed Brown     }
103c4762a1bSJed Brown 
104c4762a1bSJed Brown     /* Read generator data */
105c4762a1bSJed Brown     if (i >= gen_start_line && i < gen_end_line) {
106c4762a1bSJed Brown       double pg, qg, qt, qb, vs, mbase, pt, pb;
107c4762a1bSJed Brown       int    bus_i, status;
108c4762a1bSJed Brown       sscanf(line, "%d %lf %lf %lf %lf %lf %lf %d %lf %lf", &bus_i, &pg, &qg, &qt, &qb, &vs, &mbase, &status, &pt, &pb);
1099371c9d4SSatish Balay       Gen[geni].bus_i  = bus_i;
1109371c9d4SSatish Balay       Gen[geni].status = status;
1119371c9d4SSatish Balay       Gen[geni].pg     = pg;
1129371c9d4SSatish Balay       Gen[geni].qg     = qg;
1139371c9d4SSatish Balay       Gen[geni].qt     = qt;
1149371c9d4SSatish Balay       Gen[geni].qb     = qb;
1159371c9d4SSatish Balay       Gen[geni].vs     = vs;
1169371c9d4SSatish Balay       Gen[geni].mbase  = mbase;
1179371c9d4SSatish Balay       Gen[geni].pt     = pt;
1189371c9d4SSatish Balay       Gen[geni].pb     = pb;
119c4762a1bSJed Brown 
120c4762a1bSJed Brown       intbusnum                                  = busext2intmap[Gen[geni].bus_i];
121c4762a1bSJed Brown       Gen[geni].internal_i                       = intbusnum;
122c4762a1bSJed Brown       Bus[intbusnum].gidx[Bus[intbusnum].ngen++] = geni;
123c4762a1bSJed Brown 
124c4762a1bSJed Brown       Bus[intbusnum].vm = Gen[geni].vs;
125c4762a1bSJed Brown 
12608401ef6SPierre Jolivet       PetscCheck(Bus[intbusnum].ngen <= NGEN_AT_BUS_MAX, PETSC_COMM_SELF, PETSC_ERR_SUP, "Exceeded maximum number of generators allowed at bus");
127c4762a1bSJed Brown       geni++;
128c4762a1bSJed Brown     }
129c4762a1bSJed Brown 
130c4762a1bSJed Brown     if (i >= br_start_line && i < br_end_line) {
131c4762a1bSJed Brown       PetscScalar R, X, Bc, B, G, Zm, tap, shift, tap2, tapr, tapi;
132c4762a1bSJed Brown       double      r, x, b, rateA, rateB, rateC, tapratio, phaseshift;
133c4762a1bSJed Brown       int         fbus, tbus, status;
134c4762a1bSJed Brown       sscanf(line, "%d %d %lf %lf %lf %lf %lf %lf %lf %lf %d", &fbus, &tbus, &r, &x, &b, &rateA, &rateB, &rateC, &tapratio, &phaseshift, &status);
1359371c9d4SSatish Balay       Branch[bri].fbus       = fbus;
1369371c9d4SSatish Balay       Branch[bri].tbus       = tbus;
1379371c9d4SSatish Balay       Branch[bri].status     = status;
1389371c9d4SSatish Balay       Branch[bri].r          = r;
1399371c9d4SSatish Balay       Branch[bri].x          = x;
1409371c9d4SSatish Balay       Branch[bri].b          = b;
1419371c9d4SSatish Balay       Branch[bri].rateA      = rateA;
1429371c9d4SSatish Balay       Branch[bri].rateB      = rateB;
1439371c9d4SSatish Balay       Branch[bri].rateC      = rateC;
1449371c9d4SSatish Balay       Branch[bri].tapratio   = tapratio;
1459371c9d4SSatish Balay       Branch[bri].phaseshift = phaseshift;
146c4762a1bSJed Brown 
147c4762a1bSJed Brown       if (Branch[bri].tapratio == 0.0) Branch[bri].tapratio = 1.0;
148c4762a1bSJed Brown       Branch[bri].phaseshift *= PETSC_PI / 180.0;
149c4762a1bSJed Brown 
150c4762a1bSJed Brown       intbusnum              = busext2intmap[Branch[bri].fbus];
151c4762a1bSJed Brown       Branch[bri].internal_i = intbusnum;
152c4762a1bSJed Brown 
153c4762a1bSJed Brown       intbusnum              = busext2intmap[Branch[bri].tbus];
154c4762a1bSJed Brown       Branch[bri].internal_j = intbusnum;
155c4762a1bSJed Brown 
156c4762a1bSJed Brown       /* Compute self and transfer admittances */
157c4762a1bSJed Brown       R  = Branch[bri].r;
158c4762a1bSJed Brown       X  = Branch[bri].x;
159c4762a1bSJed Brown       Bc = Branch[bri].b;
160c4762a1bSJed Brown 
161c4762a1bSJed Brown       Zm = R * R + X * X;
162c4762a1bSJed Brown       G  = R / Zm;
163c4762a1bSJed Brown       B  = -X / Zm;
164c4762a1bSJed Brown 
165c4762a1bSJed Brown       tap   = Branch[bri].tapratio;
166c4762a1bSJed Brown       shift = Branch[bri].phaseshift;
167c4762a1bSJed Brown       tap2  = tap * tap;
168c4762a1bSJed Brown       tapr  = tap * PetscCosScalar(shift);
169c4762a1bSJed Brown       tapi  = tap * PetscSinScalar(shift);
170c4762a1bSJed Brown 
171c4762a1bSJed Brown       Branch[bri].yff[0] = G / tap2;
172c4762a1bSJed Brown       Branch[bri].yff[1] = (B + Bc / 2.0) / tap2;
173c4762a1bSJed Brown 
174c4762a1bSJed Brown       Branch[bri].yft[0] = -(G * tapr - B * tapi) / tap2;
175c4762a1bSJed Brown       Branch[bri].yft[1] = -(B * tapr + G * tapi) / tap2;
176c4762a1bSJed Brown 
177c4762a1bSJed Brown       Branch[bri].ytf[0] = -(G * tapr + B * tapi) / tap2;
178c4762a1bSJed Brown       Branch[bri].ytf[1] = -(B * tapr - G * tapi) / tap2;
179c4762a1bSJed Brown 
180c4762a1bSJed Brown       Branch[bri].ytt[0] = G;
181c4762a1bSJed Brown       Branch[bri].ytt[1] = B + Bc / 2.0;
182c4762a1bSJed Brown 
183c4762a1bSJed Brown       bri++;
184c4762a1bSJed Brown     }
185c4762a1bSJed Brown   }
186c4762a1bSJed Brown   fclose(fp);
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   /* Reorder the generator data structure according to bus numbers */
1899371c9d4SSatish Balay   genj  = 0;
1909371c9d4SSatish Balay   loadj = 0;
1919566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(pf->ngen, &newgen));
1929566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(pf->nload, &newload));
193c4762a1bSJed Brown   for (i = 0; i < pf->nbus; i++) {
19448a46eb9SPierre Jolivet     for (j = 0; j < pf->bus[i].ngen; j++) PetscCall(PetscMemcpy(&newgen[genj++], &pf->gen[pf->bus[i].gidx[j]], sizeof(struct _p_GEN)));
19548a46eb9SPierre Jolivet     for (j = 0; j < pf->bus[i].nload; j++) PetscCall(PetscMemcpy(&newload[loadj++], &pf->load[pf->bus[i].lidx[j]], sizeof(struct _p_LOAD)));
196c4762a1bSJed Brown   }
1979566063dSJacob Faibussowitsch   PetscCall(PetscFree(pf->gen));
1989566063dSJacob Faibussowitsch   PetscCall(PetscFree(pf->load));
199c4762a1bSJed Brown   pf->gen  = newgen;
200c4762a1bSJed Brown   pf->load = newload;
201c4762a1bSJed Brown 
2029566063dSJacob Faibussowitsch   PetscCall(PetscFree(busext2intmap));
203*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
204c4762a1bSJed Brown }
205