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