xref: /petsc/src/snes/tutorials/network/power/PFReadData.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
1c4762a1bSJed Brown #include "petscmat.h"
2c4762a1bSJed Brown #include "power.h"
3c4762a1bSJed Brown #include <string.h>
4c4762a1bSJed Brown #include <ctype.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown PetscErrorCode PFReadMatPowerData(PFDATA *pf,char *filename)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   FILE           *fp;
9c4762a1bSJed Brown   PetscErrorCode ierr;
10c4762a1bSJed Brown   VERTEX_Power   Bus;
11c4762a1bSJed Brown   LOAD           Load;
12c4762a1bSJed Brown   GEN            Gen;
13c4762a1bSJed Brown   EDGE_Power     Branch;
14c4762a1bSJed Brown   PetscInt       line_counter=0;
15c4762a1bSJed Brown   PetscInt       bus_start_line=-1,bus_end_line=-1; /* xx_end_line points to the next line after the record ends */
16c4762a1bSJed Brown   PetscInt       gen_start_line=-1,gen_end_line=-1;
17c4762a1bSJed Brown   PetscInt       br_start_line=-1,br_end_line=-1;
18c4762a1bSJed Brown   char           line[MAXLINE];
19c4762a1bSJed Brown   PetscInt       loadi=0,geni=0,bri=0,busi=0,i,j;
20c4762a1bSJed Brown   int            extbusnum,bustype_i;
21c4762a1bSJed Brown   double         Pd,Qd;
22c4762a1bSJed Brown   PetscInt       maxbusnum=-1,intbusnum,*busext2intmap,genj,loadj;
23c4762a1bSJed Brown   GEN            newgen;
24c4762a1bSJed Brown   LOAD           newload;
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   PetscFunctionBegin;
27c4762a1bSJed Brown   fp = fopen(filename,"r");
28c4762a1bSJed Brown   /* Check for valid file */
29*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!fp,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Can't open Matpower data file %s",filename);
30c4762a1bSJed Brown   pf->nload=0;
31c4762a1bSJed Brown   while (fgets(line,MAXLINE,fp)) {
32c4762a1bSJed Brown     if (strstr(line,"mpc.bus = ["))    bus_start_line = line_counter+1; /* Bus data starts from next line */
33c4762a1bSJed Brown     if (strstr(line,"mpc.gen") && gen_start_line == -1)    gen_start_line = line_counter+1; /* Generator data starts from next line */
34c4762a1bSJed Brown     if (strstr(line,"mpc.branch")) br_start_line = line_counter+1; /* Branch data starts from next line */
35c4762a1bSJed Brown     if (strstr(line,"];")) {
36c4762a1bSJed Brown       if (bus_start_line != -1 && bus_end_line == -1) bus_end_line = line_counter;
37c4762a1bSJed Brown       if (gen_start_line != -1 && gen_end_line == -1) gen_end_line = line_counter;
38c4762a1bSJed Brown       if (br_start_line  != -1 && br_end_line == -1) br_end_line = line_counter;
39c4762a1bSJed Brown     }
40c4762a1bSJed Brown 
41c4762a1bSJed Brown     /* Count the number of pq loads */
42c4762a1bSJed Brown     if (bus_start_line != -1 && line_counter >= bus_start_line && bus_end_line == -1) {
43c4762a1bSJed Brown       sscanf(line,"%d %d %lf %lf",&extbusnum,&bustype_i,&Pd,&Qd);
44c4762a1bSJed Brown       if (!((Pd == 0.0) && (Qd == 0.0))) pf->nload++;
45c4762a1bSJed Brown       if (extbusnum > maxbusnum) maxbusnum = extbusnum;
46c4762a1bSJed Brown     }
47c4762a1bSJed Brown     line_counter++;
48c4762a1bSJed Brown   }
49c4762a1bSJed Brown   fclose(fp);
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   pf->nbus    = bus_end_line - bus_start_line;
52c4762a1bSJed Brown   pf->ngen    = gen_end_line - gen_start_line;
53c4762a1bSJed Brown   pf->nbranch = br_end_line  - br_start_line;
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   ierr = PetscCalloc1(pf->nbus,&pf->bus);CHKERRQ(ierr);
56c4762a1bSJed Brown   ierr = PetscCalloc1(pf->ngen,&pf->gen);CHKERRQ(ierr);
57c4762a1bSJed Brown   ierr = PetscCalloc1(pf->nload,&pf->load);CHKERRQ(ierr);
58c4762a1bSJed Brown   ierr = PetscCalloc1(pf->nbranch,&pf->branch);CHKERRQ(ierr);
59c4762a1bSJed Brown   Bus = pf->bus; Gen = pf->gen; Load = pf->load; Branch = pf->branch;
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* Setting pf->sbase to 100 */
62c4762a1bSJed Brown   pf->sbase = 100.0;
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   ierr = PetscMalloc1(maxbusnum+1,&busext2intmap);CHKERRQ(ierr);
65c4762a1bSJed Brown   for (i=0; i < maxbusnum+1; i++) busext2intmap[i] = -1;
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   fp = fopen(filename,"r");
68c4762a1bSJed Brown   /* Reading data */
69c4762a1bSJed Brown   for (i=0;i<line_counter;i++) {
70*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!fgets(line,MAXLINE,fp),PETSC_COMM_SELF,PETSC_ERR_SUP,"File is incorrectly formatted");
71c4762a1bSJed Brown 
72c4762a1bSJed Brown     if ((i >= bus_start_line) && (i < bus_end_line)) {
73c4762a1bSJed Brown       double gl,bl,vm,va,basekV;
74c4762a1bSJed Brown       int    bus_i,ide,area;
75c4762a1bSJed Brown       /* Bus data */
76c4762a1bSJed Brown       sscanf(line,"%d %d %lf %lf %lf %lf %d %lf %lf %lf",&bus_i,&ide,&Pd,&Qd,&gl,&bl,&area,&vm,&va,&basekV);
77c4762a1bSJed Brown       Bus[busi].bus_i = bus_i; Bus[busi].ide = ide; Bus[busi].area = area;
78c4762a1bSJed Brown       Bus[busi].gl = gl; Bus[busi].bl = bl;
79c4762a1bSJed Brown       Bus[busi].vm = vm; Bus[busi].va = va; Bus[busi].basekV = basekV;
80c4762a1bSJed Brown       Bus[busi].internal_i = busi;
81c4762a1bSJed Brown       busext2intmap[Bus[busi].bus_i] = busi;
82c4762a1bSJed Brown 
83c4762a1bSJed Brown       if (!((Pd == 0.0) && (Qd == 0.0))) {
84c4762a1bSJed Brown         Load[loadi].bus_i = Bus[busi].bus_i;
85c4762a1bSJed Brown         Load[loadi].status = 1;
86c4762a1bSJed Brown         Load[loadi].pl = Pd;
87c4762a1bSJed Brown         Load[loadi].ql = Qd;
88c4762a1bSJed Brown         Load[loadi].area = Bus[busi].area;
89c4762a1bSJed Brown         Load[loadi].internal_i = busi;
90c4762a1bSJed Brown         Bus[busi].lidx[Bus[busi].nload++] = loadi;
91*2c71b3e2SJacob Faibussowitsch         PetscCheckFalse(Bus[busi].nload > NLOAD_AT_BUS_MAX,PETSC_COMM_SELF,PETSC_ERR_SUP,"Exceeded maximum number of loads allowed at bus");
92c4762a1bSJed Brown         loadi++;
93c4762a1bSJed Brown       }
94c4762a1bSJed Brown       busi++;
95c4762a1bSJed Brown     }
96c4762a1bSJed Brown 
97c4762a1bSJed Brown     /* Read generator data */
98c4762a1bSJed Brown     if (i >= gen_start_line && i < gen_end_line) {
99c4762a1bSJed Brown       double pg,qg,qt,qb,vs,mbase,pt,pb;
100c4762a1bSJed Brown       int    bus_i,status;
101c4762a1bSJed Brown       sscanf(line,"%d %lf %lf %lf %lf %lf %lf %d %lf %lf",&bus_i,&pg,&qg,&qt,&qb,&vs,&mbase,&status,&pt,&pb);
102c4762a1bSJed Brown       Gen[geni].bus_i = bus_i; Gen[geni].status = status;
103c4762a1bSJed Brown       Gen[geni].pg = pg; Gen[geni].qg = qg; Gen[geni].qt = qt; Gen[geni].qb = qb;
104c4762a1bSJed Brown       Gen[geni].vs = vs; Gen[geni].mbase = mbase; Gen[geni].pt = pt; Gen[geni].pb = pb;
105c4762a1bSJed Brown 
106c4762a1bSJed Brown       intbusnum = busext2intmap[Gen[geni].bus_i];
107c4762a1bSJed Brown       Gen[geni].internal_i = intbusnum;
108c4762a1bSJed Brown       Bus[intbusnum].gidx[Bus[intbusnum].ngen++] = geni;
109c4762a1bSJed Brown 
110c4762a1bSJed Brown       Bus[intbusnum].vm = Gen[geni].vs;
111c4762a1bSJed Brown 
112*2c71b3e2SJacob Faibussowitsch       PetscCheckFalse(Bus[intbusnum].ngen > NGEN_AT_BUS_MAX,PETSC_COMM_SELF,PETSC_ERR_SUP,"Exceeded maximum number of generators allowed at bus");
113c4762a1bSJed Brown       geni++;
114c4762a1bSJed Brown     }
115c4762a1bSJed Brown 
116c4762a1bSJed Brown     if (i >= br_start_line && i < br_end_line) {
117c4762a1bSJed Brown       PetscScalar R,X,Bc,B,G,Zm,tap,shift,tap2,tapr,tapi;
118c4762a1bSJed Brown       double      r,x,b,rateA,rateB,rateC,tapratio,phaseshift;
119c4762a1bSJed Brown       int         fbus,tbus,status;
120c4762a1bSJed 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);
121c4762a1bSJed Brown       Branch[bri].fbus = fbus; Branch[bri].tbus = tbus; Branch[bri].status = status;
122c4762a1bSJed Brown       Branch[bri].r = r; Branch[bri].x = x; Branch[bri].b = b;
123c4762a1bSJed Brown       Branch[bri].rateA = rateA; Branch[bri].rateB = rateB; Branch[bri].rateC = rateC;
124c4762a1bSJed Brown       Branch[bri].tapratio = tapratio; Branch[bri].phaseshift = phaseshift;
125c4762a1bSJed Brown 
126c4762a1bSJed Brown       if (Branch[bri].tapratio == 0.0) Branch[bri].tapratio = 1.0;
127c4762a1bSJed Brown       Branch[bri].phaseshift *= PETSC_PI/180.0;
128c4762a1bSJed Brown 
129c4762a1bSJed Brown       intbusnum = busext2intmap[Branch[bri].fbus];
130c4762a1bSJed Brown       Branch[bri].internal_i = intbusnum;
131c4762a1bSJed Brown 
132c4762a1bSJed Brown       intbusnum = busext2intmap[Branch[bri].tbus];
133c4762a1bSJed Brown       Branch[bri].internal_j = intbusnum;
134c4762a1bSJed Brown 
135c4762a1bSJed Brown       /* Compute self and transfer admittances */
136c4762a1bSJed Brown       R = Branch[bri].r;
137c4762a1bSJed Brown       X = Branch[bri].x;
138c4762a1bSJed Brown       Bc = Branch[bri].b;
139c4762a1bSJed Brown 
140c4762a1bSJed Brown       Zm = R*R + X*X;
141c4762a1bSJed Brown       G  = R/Zm;
142c4762a1bSJed Brown       B  = -X/Zm;
143c4762a1bSJed Brown 
144c4762a1bSJed Brown       tap = Branch[bri].tapratio;
145c4762a1bSJed Brown       shift = Branch[bri].phaseshift;
146c4762a1bSJed Brown       tap2 = tap*tap;
147c4762a1bSJed Brown       tapr = tap*PetscCosScalar(shift);
148c4762a1bSJed Brown       tapi = tap*PetscSinScalar(shift);
149c4762a1bSJed Brown 
150c4762a1bSJed Brown       Branch[bri].yff[0] = G/tap2;
151c4762a1bSJed Brown       Branch[bri].yff[1] = (B+Bc/2.0)/tap2;
152c4762a1bSJed Brown 
153c4762a1bSJed Brown       Branch[bri].yft[0] = -(G*tapr - B*tapi)/tap2;
154c4762a1bSJed Brown       Branch[bri].yft[1] = -(B*tapr + G*tapi)/tap2;
155c4762a1bSJed Brown 
156c4762a1bSJed Brown       Branch[bri].ytf[0] = -(G*tapr + B*tapi)/tap2;
157c4762a1bSJed Brown       Branch[bri].ytf[1] = -(B*tapr - G*tapi)/tap2;
158c4762a1bSJed Brown 
159c4762a1bSJed Brown       Branch[bri].ytt[0] = G;
160c4762a1bSJed Brown       Branch[bri].ytt[1] = B+Bc/2.0;
161c4762a1bSJed Brown 
162c4762a1bSJed Brown       bri++;
163c4762a1bSJed Brown     }
164c4762a1bSJed Brown   }
165c4762a1bSJed Brown   fclose(fp);
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /* Reorder the generator data structure according to bus numbers */
168c4762a1bSJed Brown   genj=0; loadj=0;
169c4762a1bSJed Brown   ierr = PetscMalloc1(pf->ngen,&newgen);CHKERRQ(ierr);
170c4762a1bSJed Brown   ierr = PetscMalloc1(pf->nload,&newload);CHKERRQ(ierr);
171c4762a1bSJed Brown   for (i = 0; i < pf->nbus; i++) {
172c4762a1bSJed Brown     for (j = 0; j < pf->bus[i].ngen; j++) {
173c4762a1bSJed Brown       ierr = PetscMemcpy(&newgen[genj++],&pf->gen[pf->bus[i].gidx[j]],sizeof(struct _p_GEN));CHKERRQ(ierr);
174c4762a1bSJed Brown     }
175c4762a1bSJed Brown     for (j = 0; j < pf->bus[i].nload; j++) {
176c4762a1bSJed Brown       ierr = PetscMemcpy(&newload[loadj++],&pf->load[pf->bus[i].lidx[j]],sizeof(struct _p_LOAD));CHKERRQ(ierr);
177c4762a1bSJed Brown     }
178c4762a1bSJed Brown   }
179c4762a1bSJed Brown   ierr = PetscFree(pf->gen);CHKERRQ(ierr);
180c4762a1bSJed Brown   ierr = PetscFree(pf->load);CHKERRQ(ierr);
181c4762a1bSJed Brown   pf->gen = newgen;
182c4762a1bSJed Brown   pf->load = newload;
183c4762a1bSJed Brown 
184c4762a1bSJed Brown   ierr = PetscFree(busext2intmap);CHKERRQ(ierr);
185c4762a1bSJed Brown   PetscFunctionReturn(0);
186c4762a1bSJed Brown }
187