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