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