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