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