1c4762a1bSJed Brown #include "water.h"
2c4762a1bSJed Brown #include <string.h>
3c4762a1bSJed Brown #include <ctype.h>
4c4762a1bSJed Brown
PumpHeadCurveResidual(SNES snes,Vec X,Vec F,PetscCtx ctx)5*2a8381b2SBarry Smith PetscErrorCode PumpHeadCurveResidual(SNES snes, Vec X, Vec F, PetscCtx ctx)
6d71ae5a4SJacob Faibussowitsch {
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;
149566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
159566063dSJacob Faibussowitsch PetscCall(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
249566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
259566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
27c4762a1bSJed Brown }
28c4762a1bSJed Brown
SetPumpHeadCurveParams(Pump * pump)29d71ae5a4SJacob Faibussowitsch PetscErrorCode SetPumpHeadCurveParams(Pump *pump)
30d71ae5a4SJacob Faibussowitsch {
31c4762a1bSJed Brown SNES snes;
32c4762a1bSJed Brown Vec X, F;
33c4762a1bSJed Brown PetscScalar *head, *flow, *x;
34c4762a1bSJed Brown SNESConvergedReason reason;
35c4762a1bSJed Brown
36c4762a1bSJed Brown PetscFunctionBegin;
37c4762a1bSJed Brown head = pump->headcurve.head;
38c4762a1bSJed Brown flow = pump->headcurve.flow;
39c4762a1bSJed Brown if (pump->headcurve.npt == 1) {
40c4762a1bSJed Brown /* Single point head curve, set the other two data points */
41c4762a1bSJed Brown flow[1] = 0;
42c4762a1bSJed Brown head[1] = 1.33 * head[0]; /* 133% of design head -- From EPANET manual */
43c4762a1bSJed Brown flow[2] = 2 * flow[0]; /* 200% of design flow -- From EPANET manual */
44c4762a1bSJed Brown head[2] = 0;
45c4762a1bSJed Brown pump->headcurve.npt += 2;
46c4762a1bSJed Brown }
47c4762a1bSJed Brown
489566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_SELF, &snes));
49c4762a1bSJed Brown
509566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &X));
519566063dSJacob Faibussowitsch PetscCall(VecSetSizes(X, 3, 3));
529566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(X));
539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &F));
54c4762a1bSJed Brown
559566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, F, PumpHeadCurveResidual, (void *)pump));
569566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, NULL, NULL, SNESComputeJacobianDefault, NULL));
579566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes));
58c4762a1bSJed Brown
599566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x));
609371c9d4SSatish Balay x[0] = head[1];
619371c9d4SSatish Balay x[1] = 10;
629371c9d4SSatish Balay x[2] = 3;
639566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x));
64c4762a1bSJed Brown
659566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, X));
66c4762a1bSJed Brown
679566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes, &reason));
687a46b595SBarry Smith PetscCheck(reason >= 0, PETSC_COMM_SELF, PETSC_ERR_CONV_FAILED, "Pump head curve did not converge");
69c4762a1bSJed Brown
709566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x));
71c4762a1bSJed Brown pump->h0 = x[0];
72c4762a1bSJed Brown pump->r = x[1];
73c4762a1bSJed Brown pump->n = x[2];
749566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x));
75c4762a1bSJed Brown
769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X));
779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&F));
789566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes));
793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
80c4762a1bSJed Brown }
81c4762a1bSJed Brown
LineStartsWith(const char * a,const char * b)82d71ae5a4SJacob Faibussowitsch int LineStartsWith(const char *a, const char *b)
83d71ae5a4SJacob Faibussowitsch {
84c4762a1bSJed Brown if (strncmp(a, b, strlen(b)) == 0) return 1;
85c4762a1bSJed Brown return 0;
86c4762a1bSJed Brown }
87c4762a1bSJed Brown
CheckDataSegmentEnd(const char * line)88d71ae5a4SJacob Faibussowitsch int CheckDataSegmentEnd(const char *line)
89d71ae5a4SJacob Faibussowitsch {
909371c9d4SSatish 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")) {
91c4762a1bSJed Brown return 1;
92c4762a1bSJed Brown }
93c4762a1bSJed Brown return 0;
94c4762a1bSJed Brown }
95c4762a1bSJed Brown
96c4762a1bSJed Brown /* Gets the file pointer positiion for the start of the data segment and the
97c4762a1bSJed Brown number of data segments (lines) read
98c4762a1bSJed Brown */
GetDataSegment(FILE * fp,char * line,fpos_t * data_segment_start_pos,PetscInt * ndatalines)99d71ae5a4SJacob Faibussowitsch PetscErrorCode GetDataSegment(FILE *fp, char *line, fpos_t *data_segment_start_pos, PetscInt *ndatalines)
100d71ae5a4SJacob Faibussowitsch {
101c4762a1bSJed Brown PetscInt data_segment_end;
102c4762a1bSJed Brown PetscInt nlines = 0;
103c4762a1bSJed Brown
104c4762a1bSJed Brown PetscFunctionBegin;
105c4762a1bSJed Brown data_segment_end = 0;
106c4762a1bSJed Brown fgetpos(fp, data_segment_start_pos);
10708401ef6SPierre Jolivet PetscCheck(fgets(line, MAXLINE, fp), PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "Cannot read data segment from file");
108c4762a1bSJed Brown while (LineStartsWith(line, ";")) {
109c4762a1bSJed Brown fgetpos(fp, data_segment_start_pos);
11008401ef6SPierre Jolivet PetscCheck(fgets(line, MAXLINE, fp), PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "Cannot read data segment from file");
111c4762a1bSJed Brown }
112c4762a1bSJed Brown while (!data_segment_end) {
11308401ef6SPierre Jolivet PetscCheck(fgets(line, MAXLINE, fp), PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "Cannot read data segment from file");
114c4762a1bSJed Brown nlines++;
115c4762a1bSJed Brown data_segment_end = CheckDataSegmentEnd(line);
116c4762a1bSJed Brown }
117c4762a1bSJed Brown *ndatalines = nlines;
1183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
119c4762a1bSJed Brown }
120c4762a1bSJed Brown
WaterReadData(WATERDATA * water,char * filename)121d71ae5a4SJacob Faibussowitsch PetscErrorCode WaterReadData(WATERDATA *water, char *filename)
122d71ae5a4SJacob Faibussowitsch {
123c4762a1bSJed Brown FILE *fp = NULL;
124c4762a1bSJed Brown VERTEX_Water vert;
125c4762a1bSJed Brown EDGE_Water edge;
126c4762a1bSJed Brown fpos_t junc_start_pos, res_start_pos, tank_start_pos, pipe_start_pos, pump_start_pos;
127c4762a1bSJed Brown fpos_t curve_start_pos, title_start_pos;
128c4762a1bSJed Brown char line[MAXLINE];
129c4762a1bSJed Brown PetscInt i, j, nv = 0, ne = 0, ncurve = 0, ntitle = 0, nlines, ndata, curve_id;
130c4762a1bSJed Brown Junction *junction = NULL;
131c4762a1bSJed Brown Reservoir *reservoir = NULL;
132c4762a1bSJed Brown Tank *tank = NULL;
133c4762a1bSJed Brown Pipe *pipe = NULL;
134c4762a1bSJed Brown Pump *pump = NULL;
135c4762a1bSJed Brown PetscScalar curve_x, curve_y;
136c4762a1bSJed Brown double v1, v2, v3, v4, v5, v6;
137c4762a1bSJed Brown
138c4762a1bSJed Brown PetscFunctionBegin;
139c4762a1bSJed Brown water->nvertex = water->nedge = 0;
140c4762a1bSJed Brown fp = fopen(filename, "rb");
141c4762a1bSJed Brown /* Check for valid file */
14228b400f6SJacob Faibussowitsch PetscCheck(fp, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Can't open EPANET data file %s", filename);
143c4762a1bSJed Brown
144c4762a1bSJed Brown /* Read file and get line numbers for different data segments */
145c4762a1bSJed Brown while (fgets(line, MAXLINE, fp)) {
1463ba16761SJacob Faibussowitsch if (strstr(line, "[TITLE]")) PetscCall(GetDataSegment(fp, line, &title_start_pos, &ntitle));
147c4762a1bSJed Brown
148c4762a1bSJed Brown if (strstr(line, "[JUNCTIONS]")) {
1493ba16761SJacob Faibussowitsch PetscCall(GetDataSegment(fp, line, &junc_start_pos, &nlines));
150c4762a1bSJed Brown water->nvertex += nlines;
151c4762a1bSJed Brown water->njunction = nlines;
152c4762a1bSJed Brown }
153c4762a1bSJed Brown
154c4762a1bSJed Brown if (strstr(line, "[RESERVOIRS]")) {
1553ba16761SJacob Faibussowitsch PetscCall(GetDataSegment(fp, line, &res_start_pos, &nlines));
156c4762a1bSJed Brown water->nvertex += nlines;
157c4762a1bSJed Brown water->nreservoir = nlines;
158c4762a1bSJed Brown }
159c4762a1bSJed Brown
160c4762a1bSJed Brown if (strstr(line, "[TANKS]")) {
1613ba16761SJacob Faibussowitsch PetscCall(GetDataSegment(fp, line, &tank_start_pos, &nlines));
162c4762a1bSJed Brown water->nvertex += nlines;
163c4762a1bSJed Brown water->ntank = nlines;
164c4762a1bSJed Brown }
165c4762a1bSJed Brown
166c4762a1bSJed Brown if (strstr(line, "[PIPES]")) {
1673ba16761SJacob Faibussowitsch PetscCall(GetDataSegment(fp, line, &pipe_start_pos, &nlines));
168c4762a1bSJed Brown water->nedge += nlines;
169c4762a1bSJed Brown water->npipe = nlines;
170c4762a1bSJed Brown }
171c4762a1bSJed Brown
172c4762a1bSJed Brown if (strstr(line, "[PUMPS]")) {
1733ba16761SJacob Faibussowitsch PetscCall(GetDataSegment(fp, line, &pump_start_pos, &nlines));
174c4762a1bSJed Brown water->nedge += nlines;
175c4762a1bSJed Brown water->npump = nlines;
176c4762a1bSJed Brown }
177c4762a1bSJed Brown
1783ba16761SJacob Faibussowitsch if (strstr(line, "[CURVES]")) PetscCall(GetDataSegment(fp, line, &curve_start_pos, &ncurve));
179c4762a1bSJed Brown }
180c4762a1bSJed Brown
181c4762a1bSJed Brown /* Allocate vertex and edge data structs */
1829566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(water->nvertex, &water->vertex));
1839566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(water->nedge, &water->edge));
184c4762a1bSJed Brown vert = water->vertex;
185c4762a1bSJed Brown edge = water->edge;
186c4762a1bSJed Brown
187c4762a1bSJed Brown /* Junctions */
188c4762a1bSJed Brown fsetpos(fp, &junc_start_pos);
189c4762a1bSJed Brown for (i = 0; i < water->njunction; i++) {
190c4762a1bSJed Brown int id = 0, pattern = 0;
19108401ef6SPierre Jolivet PetscCheck(fgets(line, MAXLINE, fp), PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "Cannot read junction from file");
192c4762a1bSJed Brown vert[nv].type = VERTEX_TYPE_JUNCTION;
193c4762a1bSJed Brown junction = &vert[nv].junc;
1949371c9d4SSatish Balay ndata = sscanf(line, "%d %lf %lf %d", &id, &v1, &v2, &pattern);
1959371c9d4SSatish Balay PetscCheck(ndata >= 3, PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "Unable to read junction data");
196c4762a1bSJed Brown vert[nv].id = id;
197c4762a1bSJed Brown junction->dempattern = pattern;
198c4762a1bSJed Brown junction->elev = (PetscScalar)v1;
199c4762a1bSJed Brown junction->demand = (PetscScalar)v2;
200c4762a1bSJed Brown junction->demand *= GPM_CFS;
201c4762a1bSJed Brown junction->id = vert[nv].id;
202c4762a1bSJed Brown nv++;
203c4762a1bSJed Brown }
204c4762a1bSJed Brown
205c4762a1bSJed Brown /* Reservoirs */
206c4762a1bSJed Brown fsetpos(fp, &res_start_pos);
207c4762a1bSJed Brown for (i = 0; i < water->nreservoir; i++) {
208c4762a1bSJed Brown int id = 0, pattern = 0;
20908401ef6SPierre Jolivet PetscCheck(fgets(line, MAXLINE, fp), PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "Cannot read reservoir from file");
210c4762a1bSJed Brown vert[nv].type = VERTEX_TYPE_RESERVOIR;
211c4762a1bSJed Brown reservoir = &vert[nv].res;
2129371c9d4SSatish Balay ndata = sscanf(line, "%d %lf %d", &id, &v1, &pattern);
2139371c9d4SSatish Balay PetscCheck(ndata >= 2, PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "Unable to read reservoir data");
214c4762a1bSJed Brown vert[nv].id = id;
215c4762a1bSJed Brown reservoir->headpattern = pattern;
216c4762a1bSJed Brown reservoir->head = (PetscScalar)v1;
217c4762a1bSJed Brown reservoir->id = vert[nv].id;
218c4762a1bSJed Brown nv++;
219c4762a1bSJed Brown }
220c4762a1bSJed Brown
221c4762a1bSJed Brown /* Tanks */
222c4762a1bSJed Brown fsetpos(fp, &tank_start_pos);
223c4762a1bSJed Brown for (i = 0; i < water->ntank; i++) {
224c4762a1bSJed Brown int id = 0, curve = 0;
22508401ef6SPierre Jolivet PetscCheck(fgets(line, MAXLINE, fp), PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "Cannot read data tank from file");
226c4762a1bSJed Brown vert[nv].type = VERTEX_TYPE_TANK;
227c4762a1bSJed Brown tank = &vert[nv].tank;
2289371c9d4SSatish Balay ndata = sscanf(line, "%d %lf %lf %lf %lf %lf %lf %d", &id, &v1, &v2, &v3, &v4, &v5, &v6, &curve);
2299371c9d4SSatish Balay PetscCheck(ndata >= 7, PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "Unable to read tank data");
230c4762a1bSJed Brown vert[nv].id = id;
231c4762a1bSJed Brown tank->volumecurve = curve;
232c4762a1bSJed Brown tank->elev = (PetscScalar)v1;
233c4762a1bSJed Brown tank->initlvl = (PetscScalar)v2;
234c4762a1bSJed Brown tank->minlvl = (PetscScalar)v3;
235c4762a1bSJed Brown tank->maxlvl = (PetscScalar)v4;
236c4762a1bSJed Brown tank->diam = (PetscScalar)v5;
237c4762a1bSJed Brown tank->minvolume = (PetscScalar)v6;
238c4762a1bSJed Brown tank->id = vert[nv].id;
239c4762a1bSJed Brown nv++;
240c4762a1bSJed Brown }
241c4762a1bSJed Brown
242c4762a1bSJed Brown /* Pipes */
243c4762a1bSJed Brown fsetpos(fp, &pipe_start_pos);
244c4762a1bSJed Brown for (i = 0; i < water->npipe; i++) {
245c4762a1bSJed Brown int id = 0, node1 = 0, node2 = 0;
24608401ef6SPierre Jolivet PetscCheck(fgets(line, MAXLINE, fp), PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "Cannot read data pipe from file");
247c4762a1bSJed Brown edge[ne].type = EDGE_TYPE_PIPE;
248c4762a1bSJed Brown pipe = &edge[ne].pipe;
249c4762a1bSJed Brown ndata = sscanf(line, "%d %d %d %lf %lf %lf %lf %s", &id, &node1, &node2, &v1, &v2, &v3, &v4, pipe->stat);
250c4762a1bSJed Brown pipe->id = id;
251c4762a1bSJed Brown pipe->node1 = node1;
252c4762a1bSJed Brown pipe->node2 = node2;
253c4762a1bSJed Brown pipe->length = (PetscScalar)v1;
254c4762a1bSJed Brown pipe->diam = (PetscScalar)v2;
255c4762a1bSJed Brown pipe->roughness = (PetscScalar)v3;
256c4762a1bSJed Brown pipe->minorloss = (PetscScalar)v4;
257c4762a1bSJed Brown edge[ne].id = pipe->id;
258c4762a1bSJed Brown if (strcmp(pipe->stat, "OPEN") == 0) pipe->status = PIPE_STATUS_OPEN;
259c4762a1bSJed Brown if (ndata < 8) {
260c4762a1bSJed Brown strcpy(pipe->stat, "OPEN"); /* default OPEN */
261c4762a1bSJed Brown pipe->status = PIPE_STATUS_OPEN;
262c4762a1bSJed Brown }
263c4762a1bSJed Brown if (ndata < 7) pipe->minorloss = 0.;
264c4762a1bSJed Brown pipe->n = 1.85;
265c4762a1bSJed Brown pipe->k = 4.72 * pipe->length / (PetscPowScalar(pipe->roughness, pipe->n) * PetscPowScalar(0.0833333 * pipe->diam, 4.87));
266c4762a1bSJed Brown ne++;
267c4762a1bSJed Brown }
268c4762a1bSJed Brown
269c4762a1bSJed Brown /* Pumps */
270c4762a1bSJed Brown fsetpos(fp, &pump_start_pos);
271c4762a1bSJed Brown for (i = 0; i < water->npump; i++) {
272c4762a1bSJed Brown int id = 0, node1 = 0, node2 = 0, paramid = 0;
27308401ef6SPierre Jolivet PetscCheck(fgets(line, MAXLINE, fp), PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "Cannot read data pump from file");
274c4762a1bSJed Brown edge[ne].type = EDGE_TYPE_PUMP;
275c4762a1bSJed Brown pump = &edge[ne].pump;
2769371c9d4SSatish Balay ndata = sscanf(line, "%d %d %d %s %d", &id, &node1, &node2, pump->param, ¶mid);
2779371c9d4SSatish Balay PetscCheck(ndata == 5, PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "Unable to read pump data");
278c4762a1bSJed Brown pump->id = id;
279c4762a1bSJed Brown pump->node1 = node1;
280c4762a1bSJed Brown pump->node2 = node2;
281c4762a1bSJed Brown pump->paramid = paramid;
282c4762a1bSJed Brown edge[ne].id = pump->id;
283c4762a1bSJed Brown ne++;
284c4762a1bSJed Brown }
285c4762a1bSJed Brown
286c4762a1bSJed Brown /* Curves */
287c4762a1bSJed Brown fsetpos(fp, &curve_start_pos);
288c4762a1bSJed Brown for (i = 0; i < ncurve; i++) {
289c4762a1bSJed Brown int icurve_id = 0;
29008401ef6SPierre Jolivet PetscCheck(fgets(line, MAXLINE, fp), PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "Cannot read data curve from file");
2919371c9d4SSatish Balay ndata = sscanf(line, "%d %lf %lf", &icurve_id, &v1, &v2);
2929371c9d4SSatish Balay PetscCheck(ndata == 3, PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "Unable to read curve data");
293c4762a1bSJed Brown curve_id = icurve_id;
294c4762a1bSJed Brown curve_x = (PetscScalar)v1;
295c4762a1bSJed Brown curve_y = (PetscScalar)v2;
296c4762a1bSJed Brown /* Check for pump with the curve_id */
297c4762a1bSJed Brown for (j = water->npipe; j < water->npipe + water->npump; j++) {
298c4762a1bSJed Brown if (water->edge[j].pump.paramid == curve_id) {
2997a46b595SBarry 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);
300c4762a1bSJed Brown pump = &water->edge[j].pump;
301c4762a1bSJed Brown pump->headcurve.flow[pump->headcurve.npt] = curve_x * GPM_CFS;
302c4762a1bSJed Brown pump->headcurve.head[pump->headcurve.npt] = curve_y;
303c4762a1bSJed Brown pump->headcurve.npt++;
304c4762a1bSJed Brown break;
305c4762a1bSJed Brown }
306c4762a1bSJed Brown }
307c4762a1bSJed Brown }
308c4762a1bSJed Brown
309c4762a1bSJed Brown fclose(fp);
310c4762a1bSJed Brown
311c4762a1bSJed Brown /* Get pump curve parameters */
312c4762a1bSJed Brown for (j = water->npipe; j < water->npipe + water->npump; j++) {
313c4762a1bSJed Brown pump = &water->edge[j].pump;
314c4762a1bSJed Brown if (strcmp(pump->param, "HEAD") == 0) {
315c4762a1bSJed Brown /* Head-flow curve */
3169566063dSJacob Faibussowitsch PetscCall(SetPumpHeadCurveParams(pump));
317c4762a1bSJed Brown }
318c4762a1bSJed Brown }
3193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
320c4762a1bSJed Brown }
321