xref: /petsc/src/snes/tutorials/network/water/waterreaddata.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
1c4762a1bSJed Brown #include "water.h"
2c4762a1bSJed Brown #include <string.h>
3c4762a1bSJed Brown #include <ctype.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown PetscErrorCode PumpHeadCurveResidual(SNES snes,Vec X, Vec F,void *ctx)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   const PetscScalar *x;
8c4762a1bSJed Brown   PetscScalar *f;
9c4762a1bSJed Brown   Pump        *pump=(Pump*)ctx;
10c4762a1bSJed Brown   PetscScalar *head=pump->headcurve.head,*flow=pump->headcurve.flow;
11c4762a1bSJed Brown   PetscInt i;
12c4762a1bSJed Brown 
13c4762a1bSJed Brown   PetscFunctionBegin;
145f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
155f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   f[0] = f[1] = f[2] = 0;
18c4762a1bSJed Brown   for (i=0; i < pump->headcurve.npt;i++) {
19c4762a1bSJed Brown     f[0] +=   x[0] - x[1]*PetscPowScalar(flow[i],x[2]) - head[i]; /* Partial w.r.t x[0] */
20c4762a1bSJed Brown     f[1] +=  (x[0] - x[1]*PetscPowScalar(flow[i],x[2]) - head[i])*-1*PetscPowScalar(flow[i],x[2]); /*Partial w.r.t x[1] */
21c4762a1bSJed Brown     f[2] +=  (x[0] - x[1]*PetscPowScalar(flow[i],x[2]) - head[i])*-1*x[1]*x[2]*PetscPowScalar(flow[i],x[2]-1); /*Partial w.r.t x[2] */
22c4762a1bSJed Brown   }
23c4762a1bSJed Brown 
245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   PetscFunctionReturn(0);
28c4762a1bSJed Brown }
29c4762a1bSJed Brown 
30c4762a1bSJed Brown PetscErrorCode SetPumpHeadCurveParams(Pump *pump)
31c4762a1bSJed Brown {
32c4762a1bSJed Brown   SNES           snes;
33c4762a1bSJed Brown   Vec            X,F;
34c4762a1bSJed Brown   PetscScalar   *head,*flow,*x;
35c4762a1bSJed Brown   SNESConvergedReason reason;
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   PetscFunctionBegin;
38c4762a1bSJed Brown   head = pump->headcurve.head;
39c4762a1bSJed Brown   flow = pump->headcurve.flow;
40c4762a1bSJed Brown   if (pump->headcurve.npt == 1) {
41c4762a1bSJed Brown     /* Single point head curve, set the other two data points */
42c4762a1bSJed Brown     flow[1] = 0;
43c4762a1bSJed Brown     head[1] = 1.33*head[0]; /* 133% of design head -- From EPANET manual */
44c4762a1bSJed Brown     flow[2] = 2*flow[0];    /* 200% of design flow -- From EPANET manual */
45c4762a1bSJed Brown     head[2] = 0;
46c4762a1bSJed Brown     pump->headcurve.npt += 2;
47c4762a1bSJed Brown   }
48c4762a1bSJed Brown 
495f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(PETSC_COMM_SELF,&snes));
50c4762a1bSJed Brown 
515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_SELF,&X));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(X,3,3));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(X));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(X,&F));
55c4762a1bSJed Brown 
565f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFunction(snes,F,PumpHeadCurveResidual,(void*)pump));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetJacobian(snes,NULL,NULL,SNESComputeJacobianDefault,NULL));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
59c4762a1bSJed Brown 
605f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(X,&x));
61c4762a1bSJed Brown   x[0] = head[1]; x[1] = 10; x[2] = 3;
625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(X,&x));
63c4762a1bSJed Brown 
645f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes,NULL,X));
65c4762a1bSJed Brown 
665f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetConvergedReason(snes,&reason));
67c4762a1bSJed Brown   if (reason < 0) {
68546078acSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"Pump head curve did not converge");
69c4762a1bSJed Brown   }
70c4762a1bSJed Brown 
715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(X,&x));
72c4762a1bSJed Brown   pump->h0 = x[0];
73c4762a1bSJed Brown   pump->r  = x[1];
74c4762a1bSJed Brown   pump->n  = x[2];
755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(X,&x));
76c4762a1bSJed Brown 
775f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&F));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESDestroy(&snes));
80c4762a1bSJed Brown   PetscFunctionReturn(0);
81c4762a1bSJed Brown }
82c4762a1bSJed Brown 
83c4762a1bSJed Brown int LineStartsWith(const char *a, const char *b)
84c4762a1bSJed Brown {
85c4762a1bSJed Brown   if (strncmp(a, b, strlen(b)) == 0) return 1;
86c4762a1bSJed Brown   return 0;
87c4762a1bSJed Brown }
88c4762a1bSJed Brown 
89c4762a1bSJed Brown int CheckDataSegmentEnd(const char *line)
90c4762a1bSJed Brown {
91c4762a1bSJed Brown   if (LineStartsWith(line,"[JUNCTIONS]") || \
92c4762a1bSJed Brown      LineStartsWith(line,"[RESERVOIRS]") || \
93c4762a1bSJed Brown      LineStartsWith(line,"[TANKS]") || \
94c4762a1bSJed Brown      LineStartsWith(line,"[PIPES]") || \
95c4762a1bSJed Brown      LineStartsWith(line,"[PUMPS]") || \
96c4762a1bSJed Brown      LineStartsWith(line,"[CURVES]") || \
97c4762a1bSJed Brown      LineStartsWith(line,"[VALVES]") || \
98c4762a1bSJed Brown      LineStartsWith(line,"[PATTERNS]") || \
99c4762a1bSJed Brown      LineStartsWith(line,"[VALVES]") || \
100c4762a1bSJed Brown      LineStartsWith(line,"[QUALITY]") || \
101c4762a1bSJed Brown      LineStartsWith(line,"\n") || LineStartsWith(line,"\r\n")) {
102c4762a1bSJed Brown     return 1;
103c4762a1bSJed Brown   }
104c4762a1bSJed Brown   return 0;
105c4762a1bSJed Brown }
106c4762a1bSJed Brown 
107c4762a1bSJed Brown /* Gets the file pointer positiion for the start of the data segment and the
108c4762a1bSJed Brown    number of data segments (lines) read
109c4762a1bSJed Brown */
110c4762a1bSJed Brown PetscErrorCode GetDataSegment(FILE *fp,char *line,fpos_t *data_segment_start_pos,PetscInt *ndatalines)
111c4762a1bSJed Brown {
112c4762a1bSJed Brown   PetscInt data_segment_end;
113c4762a1bSJed Brown   PetscInt nlines=0;
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   PetscFunctionBegin;
116c4762a1bSJed Brown   data_segment_end = 0;
117c4762a1bSJed Brown   fgetpos(fp,data_segment_start_pos);
1182c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!fgets(line,MAXLINE,fp),PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Cannot read data segment from file");
119c4762a1bSJed Brown   while (LineStartsWith(line,";")) {
120c4762a1bSJed Brown     fgetpos(fp,data_segment_start_pos);
1212c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!fgets(line,MAXLINE,fp),PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Cannot read data segment from file");
122c4762a1bSJed Brown   }
123c4762a1bSJed Brown   while (!data_segment_end) {
1242c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!fgets(line,MAXLINE,fp),PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Cannot read data segment from file");
125c4762a1bSJed Brown     nlines++;
126c4762a1bSJed Brown     data_segment_end = CheckDataSegmentEnd(line);
127c4762a1bSJed Brown   }
128c4762a1bSJed Brown   *ndatalines = nlines;
129c4762a1bSJed Brown   PetscFunctionReturn(0);
130c4762a1bSJed Brown }
131c4762a1bSJed Brown 
132c4762a1bSJed Brown PetscErrorCode WaterReadData(WATERDATA *water,char *filename)
133c4762a1bSJed Brown {
134c4762a1bSJed Brown   FILE           *fp=NULL;
135c4762a1bSJed Brown   VERTEX_Water   vert;
136c4762a1bSJed Brown   EDGE_Water     edge;
137c4762a1bSJed Brown   fpos_t         junc_start_pos,res_start_pos,tank_start_pos,pipe_start_pos,pump_start_pos;
138c4762a1bSJed Brown   fpos_t         curve_start_pos,title_start_pos;
139c4762a1bSJed Brown   char           line[MAXLINE];
140c4762a1bSJed Brown   PetscInt       i,j,nv=0,ne=0,ncurve=0,ntitle=0,nlines,ndata,curve_id;
141c4762a1bSJed Brown   Junction       *junction=NULL;
142c4762a1bSJed Brown   Reservoir      *reservoir=NULL;
143c4762a1bSJed Brown   Tank           *tank=NULL;
144c4762a1bSJed Brown   Pipe           *pipe=NULL;
145c4762a1bSJed Brown   Pump           *pump=NULL;
146c4762a1bSJed Brown   PetscScalar    curve_x,curve_y;
147c4762a1bSJed Brown   double         v1,v2,v3,v4,v5,v6;
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   PetscFunctionBegin;
150c4762a1bSJed Brown   water->nvertex = water->nedge = 0;
151c4762a1bSJed Brown   fp = fopen(filename,"rb");
152c4762a1bSJed Brown   /* Check for valid file */
153*28b400f6SJacob Faibussowitsch   PetscCheck(fp,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Can't open EPANET data file %s",filename);
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   /* Read file and get line numbers for different data segments */
156c4762a1bSJed Brown   while (fgets(line,MAXLINE,fp)) {
157c4762a1bSJed Brown 
158c4762a1bSJed Brown     if (strstr(line,"[TITLE]")) {
159c4762a1bSJed Brown       GetDataSegment(fp,line,&title_start_pos,&ntitle);
160c4762a1bSJed Brown     }
161c4762a1bSJed Brown 
162c4762a1bSJed Brown     if (strstr(line,"[JUNCTIONS]")) {
163c4762a1bSJed Brown       GetDataSegment(fp,line,&junc_start_pos,&nlines);
164c4762a1bSJed Brown       water->nvertex += nlines;
165c4762a1bSJed Brown       water->njunction = nlines;
166c4762a1bSJed Brown     }
167c4762a1bSJed Brown 
168c4762a1bSJed Brown     if (strstr(line,"[RESERVOIRS]")) {
169c4762a1bSJed Brown       GetDataSegment(fp,line,&res_start_pos,&nlines);
170c4762a1bSJed Brown       water->nvertex += nlines;
171c4762a1bSJed Brown       water->nreservoir = nlines;
172c4762a1bSJed Brown     }
173c4762a1bSJed Brown 
174c4762a1bSJed Brown     if (strstr(line,"[TANKS]")) {
175c4762a1bSJed Brown       GetDataSegment(fp,line,&tank_start_pos,&nlines);
176c4762a1bSJed Brown       water->nvertex += nlines;
177c4762a1bSJed Brown       water->ntank = nlines;
178c4762a1bSJed Brown     }
179c4762a1bSJed Brown 
180c4762a1bSJed Brown     if (strstr(line,"[PIPES]")) {
181c4762a1bSJed Brown       GetDataSegment(fp,line,&pipe_start_pos,&nlines);
182c4762a1bSJed Brown       water->nedge += nlines;
183c4762a1bSJed Brown       water->npipe = nlines;
184c4762a1bSJed Brown     }
185c4762a1bSJed Brown 
186c4762a1bSJed Brown     if (strstr(line,"[PUMPS]")) {
187c4762a1bSJed Brown       GetDataSegment(fp,line,&pump_start_pos,&nlines);
188c4762a1bSJed Brown       water->nedge += nlines;
189c4762a1bSJed Brown       water->npump = nlines;
190c4762a1bSJed Brown     }
191c4762a1bSJed Brown 
192c4762a1bSJed Brown     if (strstr(line,"[CURVES]")) {
193c4762a1bSJed Brown       GetDataSegment(fp,line,&curve_start_pos,&ncurve);
194c4762a1bSJed Brown     }
195c4762a1bSJed Brown   }
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   /* Allocate vertex and edge data structs */
1985f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(water->nvertex,&water->vertex));
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(water->nedge,&water->edge));
200c4762a1bSJed Brown   vert = water->vertex;
201c4762a1bSJed Brown   edge = water->edge;
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   /* Junctions */
204c4762a1bSJed Brown   fsetpos(fp,&junc_start_pos);
205c4762a1bSJed Brown   for (i=0; i < water->njunction; i++) {
206c4762a1bSJed Brown     int id=0,pattern=0;
2072c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!fgets(line,MAXLINE,fp),PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Cannot read junction from file");
208c4762a1bSJed Brown     vert[nv].type = VERTEX_TYPE_JUNCTION;
209c4762a1bSJed Brown     junction = &vert[nv].junc;
2102c71b3e2SJacob Faibussowitsch     ndata = sscanf(line,"%d %lf %lf %d",&id,&v1,&v2,&pattern);PetscCheckFalse(ndata < 3,PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Unable to read junction data");
211c4762a1bSJed Brown     vert[nv].id          = id;
212c4762a1bSJed Brown     junction->dempattern = pattern;
213c4762a1bSJed Brown     junction->elev   = (PetscScalar)v1;
214c4762a1bSJed Brown     junction->demand = (PetscScalar)v2;
215c4762a1bSJed Brown     junction->demand *= GPM_CFS;
216c4762a1bSJed Brown     junction->id = vert[nv].id;
217c4762a1bSJed Brown     nv++;
218c4762a1bSJed Brown   }
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   /* Reservoirs */
221c4762a1bSJed Brown   fsetpos(fp,&res_start_pos);
222c4762a1bSJed Brown   for (i=0; i < water->nreservoir; i++) {
223c4762a1bSJed Brown     int id=0,pattern=0;
2242c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!fgets(line,MAXLINE,fp),PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Cannot read reservoir from file");
225c4762a1bSJed Brown     vert[nv].type = VERTEX_TYPE_RESERVOIR;
226c4762a1bSJed Brown     reservoir = &vert[nv].res;
2272c71b3e2SJacob Faibussowitsch     ndata = sscanf(line,"%d %lf %d",&id,&v1,&pattern);PetscCheckFalse(ndata < 2,PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Unable to read reservoir data");
228c4762a1bSJed Brown     vert[nv].id            = id;
229c4762a1bSJed Brown     reservoir->headpattern = pattern;
230c4762a1bSJed Brown     reservoir->head = (PetscScalar)v1;
231c4762a1bSJed Brown     reservoir->id   = vert[nv].id;
232c4762a1bSJed Brown     nv++;
233c4762a1bSJed Brown   }
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   /* Tanks */
236c4762a1bSJed Brown   fsetpos(fp,&tank_start_pos);
237c4762a1bSJed Brown   for (i=0; i < water->ntank; i++) {
238c4762a1bSJed Brown     int id=0,curve=0;
2392c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!fgets(line,MAXLINE,fp),PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Cannot read data tank from file");
240c4762a1bSJed Brown     vert[nv].type = VERTEX_TYPE_TANK;
241c4762a1bSJed Brown     tank = &vert[nv].tank;
2422c71b3e2SJacob Faibussowitsch     ndata = sscanf(line,"%d %lf %lf %lf %lf %lf %lf %d",&id,&v1,&v2,&v3,&v4,&v5,&v6,&curve);PetscCheckFalse(ndata < 7,PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Unable to read tank data");
243c4762a1bSJed Brown     vert[nv].id       = id;
244c4762a1bSJed Brown     tank->volumecurve = curve;
245c4762a1bSJed Brown     tank->elev      = (PetscScalar)v1;
246c4762a1bSJed Brown     tank->initlvl   = (PetscScalar)v2;
247c4762a1bSJed Brown     tank->minlvl    = (PetscScalar)v3;
248c4762a1bSJed Brown     tank->maxlvl    = (PetscScalar)v4;
249c4762a1bSJed Brown     tank->diam      = (PetscScalar)v5;
250c4762a1bSJed Brown     tank->minvolume = (PetscScalar)v6;
251c4762a1bSJed Brown     tank->id        = vert[nv].id;
252c4762a1bSJed Brown     nv++;
253c4762a1bSJed Brown   }
254c4762a1bSJed Brown 
255c4762a1bSJed Brown   /* Pipes */
256c4762a1bSJed Brown   fsetpos(fp,&pipe_start_pos);
257c4762a1bSJed Brown   for (i=0; i < water->npipe; i++) {
258c4762a1bSJed Brown     int id=0,node1=0,node2=0;
2592c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!fgets(line,MAXLINE,fp),PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Cannot read data pipe from file");
260c4762a1bSJed Brown     edge[ne].type = EDGE_TYPE_PIPE;
261c4762a1bSJed Brown     pipe = &edge[ne].pipe;
262c4762a1bSJed Brown     ndata = sscanf(line,"%d %d %d %lf %lf %lf %lf %s",&id,&node1,&node2,&v1,&v2,&v3,&v4,pipe->stat);
263c4762a1bSJed Brown     pipe->id        = id;
264c4762a1bSJed Brown     pipe->node1     = node1;
265c4762a1bSJed Brown     pipe->node2     = node2;
266c4762a1bSJed Brown     pipe->length    = (PetscScalar)v1;
267c4762a1bSJed Brown     pipe->diam      = (PetscScalar)v2;
268c4762a1bSJed Brown     pipe->roughness = (PetscScalar)v3;
269c4762a1bSJed Brown     pipe->minorloss = (PetscScalar)v4;
270c4762a1bSJed Brown     edge[ne].id     = pipe->id;
271c4762a1bSJed Brown     if (strcmp(pipe->stat,"OPEN") == 0) pipe->status = PIPE_STATUS_OPEN;
272c4762a1bSJed Brown     if (ndata < 8) {
273c4762a1bSJed Brown       strcpy(pipe->stat,"OPEN"); /* default OPEN */
274c4762a1bSJed Brown       pipe->status = PIPE_STATUS_OPEN;
275c4762a1bSJed Brown     }
276c4762a1bSJed Brown     if (ndata < 7) pipe->minorloss = 0.;
277c4762a1bSJed Brown     pipe->n = 1.85;
278c4762a1bSJed Brown     pipe->k = 4.72*pipe->length/(PetscPowScalar(pipe->roughness,pipe->n)*PetscPowScalar(0.0833333*pipe->diam,4.87));
279c4762a1bSJed Brown     ne++;
280c4762a1bSJed Brown   }
281c4762a1bSJed Brown 
282c4762a1bSJed Brown   /* Pumps */
283c4762a1bSJed Brown   fsetpos(fp,&pump_start_pos);
284c4762a1bSJed Brown   for (i=0; i < water->npump; i++) {
285c4762a1bSJed Brown     int id=0,node1=0,node2=0,paramid=0;
2862c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!fgets(line,MAXLINE,fp),PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Cannot read data pump from file");
287c4762a1bSJed Brown     edge[ne].type = EDGE_TYPE_PUMP;
288c4762a1bSJed Brown     pump = &edge[ne].pump;
2892c71b3e2SJacob Faibussowitsch     ndata = sscanf(line,"%d %d %d %s %d",&id,&node1,&node2,pump->param,&paramid);PetscCheckFalse(ndata != 5,PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Unable to read pump data");
290c4762a1bSJed Brown     pump->id      = id;
291c4762a1bSJed Brown     pump->node1   = node1;
292c4762a1bSJed Brown     pump->node2   = node2;
293c4762a1bSJed Brown     pump->paramid = paramid;
294c4762a1bSJed Brown     edge[ne].id   = pump->id;
295c4762a1bSJed Brown     ne++;
296c4762a1bSJed Brown   }
297c4762a1bSJed Brown 
298c4762a1bSJed Brown   /* Curves */
299c4762a1bSJed Brown   fsetpos(fp,&curve_start_pos);
300c4762a1bSJed Brown   for (i=0; i < ncurve; i++) {
301c4762a1bSJed Brown     int icurve_id=0;
3022c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!fgets(line,MAXLINE,fp),PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Cannot read data curve from file");
3032c71b3e2SJacob Faibussowitsch     ndata = sscanf(line,"%d %lf %lf",&icurve_id,&v1,&v2);PetscCheckFalse(ndata != 3,PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Unable to read curve data");
304c4762a1bSJed Brown     curve_id = icurve_id;
305c4762a1bSJed Brown     curve_x  = (PetscScalar)v1;
306c4762a1bSJed Brown     curve_y  = (PetscScalar)v2;
307c4762a1bSJed Brown     /* Check for pump with the curve_id */
308c4762a1bSJed Brown     for (j=water->npipe;j < water->npipe+water->npump;j++) {
309c4762a1bSJed Brown       if (water->edge[j].pump.paramid == curve_id) {
310c4762a1bSJed Brown         if (pump->headcurve.npt == 3) {
31198921bdaSJacob Faibussowitsch           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Pump %d [%d --> %d]: No support for more than 3-pt head-flow curve",pump->id,pump->node1,pump->node2);
312c4762a1bSJed Brown         }
313c4762a1bSJed Brown         pump = &water->edge[j].pump;
314c4762a1bSJed Brown         pump->headcurve.flow[pump->headcurve.npt] = curve_x*GPM_CFS;
315c4762a1bSJed Brown         pump->headcurve.head[pump->headcurve.npt] = curve_y;
316c4762a1bSJed Brown         pump->headcurve.npt++;
317c4762a1bSJed Brown         break;
318c4762a1bSJed Brown       }
319c4762a1bSJed Brown     }
320c4762a1bSJed Brown   }
321c4762a1bSJed Brown 
322c4762a1bSJed Brown   fclose(fp);
323c4762a1bSJed Brown 
324c4762a1bSJed Brown   /* Get pump curve parameters */
325c4762a1bSJed Brown   for (j=water->npipe;j < water->npipe+water->npump;j++) {
326c4762a1bSJed Brown     pump = &water->edge[j].pump;
327c4762a1bSJed Brown     if (strcmp(pump->param,"HEAD") == 0) {
328c4762a1bSJed Brown       /* Head-flow curve */
3295f80ce2aSJacob Faibussowitsch       CHKERRQ(SetPumpHeadCurveParams(pump));
330c4762a1bSJed Brown     }
331c4762a1bSJed Brown   }
332c4762a1bSJed Brown   PetscFunctionReturn(0);
333c4762a1bSJed Brown }
334