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; 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)); 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 499566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_SELF,&snes)); 50c4762a1bSJed Brown 519566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF,&X)); 529566063dSJacob Faibussowitsch PetscCall(VecSetSizes(X,3,3)); 539566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(X)); 549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&F)); 55c4762a1bSJed Brown 569566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,F,PumpHeadCurveResidual,(void*)pump)); 579566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,NULL,NULL,SNESComputeJacobianDefault,NULL)); 589566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 59c4762a1bSJed Brown 609566063dSJacob Faibussowitsch PetscCall(VecGetArray(X,&x)); 61c4762a1bSJed Brown x[0] = head[1]; x[1] = 10; 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)); 67*7a46b595SBarry 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 81c4762a1bSJed Brown int LineStartsWith(const char *a, const char *b) 82c4762a1bSJed Brown { 83c4762a1bSJed Brown if (strncmp(a, b, strlen(b)) == 0) return 1; 84c4762a1bSJed Brown return 0; 85c4762a1bSJed Brown } 86c4762a1bSJed Brown 87c4762a1bSJed Brown int CheckDataSegmentEnd(const char *line) 88c4762a1bSJed Brown { 89c4762a1bSJed Brown if (LineStartsWith(line,"[JUNCTIONS]") || \ 90c4762a1bSJed Brown LineStartsWith(line,"[RESERVOIRS]") || \ 91c4762a1bSJed Brown LineStartsWith(line,"[TANKS]") || \ 92c4762a1bSJed Brown LineStartsWith(line,"[PIPES]") || \ 93c4762a1bSJed Brown LineStartsWith(line,"[PUMPS]") || \ 94c4762a1bSJed Brown LineStartsWith(line,"[CURVES]") || \ 95c4762a1bSJed Brown LineStartsWith(line,"[VALVES]") || \ 96c4762a1bSJed Brown LineStartsWith(line,"[PATTERNS]") || \ 97c4762a1bSJed Brown LineStartsWith(line,"[VALVES]") || \ 98c4762a1bSJed Brown LineStartsWith(line,"[QUALITY]") || \ 99c4762a1bSJed Brown LineStartsWith(line,"\n") || LineStartsWith(line,"\r\n")) { 100c4762a1bSJed Brown return 1; 101c4762a1bSJed Brown } 102c4762a1bSJed Brown return 0; 103c4762a1bSJed Brown } 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* Gets the file pointer positiion for the start of the data segment and the 106c4762a1bSJed Brown number of data segments (lines) read 107c4762a1bSJed Brown */ 108c4762a1bSJed Brown PetscErrorCode GetDataSegment(FILE *fp,char *line,fpos_t *data_segment_start_pos,PetscInt *ndatalines) 109c4762a1bSJed Brown { 110c4762a1bSJed Brown PetscInt data_segment_end; 111c4762a1bSJed Brown PetscInt nlines=0; 112c4762a1bSJed Brown 113c4762a1bSJed Brown PetscFunctionBegin; 114c4762a1bSJed Brown data_segment_end = 0; 115c4762a1bSJed Brown fgetpos(fp,data_segment_start_pos); 11608401ef6SPierre Jolivet PetscCheck(fgets(line,MAXLINE,fp),PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Cannot read data segment from file"); 117c4762a1bSJed Brown while (LineStartsWith(line,";")) { 118c4762a1bSJed Brown fgetpos(fp,data_segment_start_pos); 11908401ef6SPierre Jolivet PetscCheck(fgets(line,MAXLINE,fp),PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Cannot read data segment from file"); 120c4762a1bSJed Brown } 121c4762a1bSJed Brown while (!data_segment_end) { 12208401ef6SPierre Jolivet PetscCheck(fgets(line,MAXLINE,fp),PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Cannot read data segment from file"); 123c4762a1bSJed Brown nlines++; 124c4762a1bSJed Brown data_segment_end = CheckDataSegmentEnd(line); 125c4762a1bSJed Brown } 126c4762a1bSJed Brown *ndatalines = nlines; 127c4762a1bSJed Brown PetscFunctionReturn(0); 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 130c4762a1bSJed Brown PetscErrorCode WaterReadData(WATERDATA *water,char *filename) 131c4762a1bSJed Brown { 132c4762a1bSJed Brown FILE *fp=NULL; 133c4762a1bSJed Brown VERTEX_Water vert; 134c4762a1bSJed Brown EDGE_Water edge; 135c4762a1bSJed Brown fpos_t junc_start_pos,res_start_pos,tank_start_pos,pipe_start_pos,pump_start_pos; 136c4762a1bSJed Brown fpos_t curve_start_pos,title_start_pos; 137c4762a1bSJed Brown char line[MAXLINE]; 138c4762a1bSJed Brown PetscInt i,j,nv=0,ne=0,ncurve=0,ntitle=0,nlines,ndata,curve_id; 139c4762a1bSJed Brown Junction *junction=NULL; 140c4762a1bSJed Brown Reservoir *reservoir=NULL; 141c4762a1bSJed Brown Tank *tank=NULL; 142c4762a1bSJed Brown Pipe *pipe=NULL; 143c4762a1bSJed Brown Pump *pump=NULL; 144c4762a1bSJed Brown PetscScalar curve_x,curve_y; 145c4762a1bSJed Brown double v1,v2,v3,v4,v5,v6; 146c4762a1bSJed Brown 147c4762a1bSJed Brown PetscFunctionBegin; 148c4762a1bSJed Brown water->nvertex = water->nedge = 0; 149c4762a1bSJed Brown fp = fopen(filename,"rb"); 150c4762a1bSJed Brown /* Check for valid file */ 15128b400f6SJacob Faibussowitsch PetscCheck(fp,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Can't open EPANET data file %s",filename); 152c4762a1bSJed Brown 153c4762a1bSJed Brown /* Read file and get line numbers for different data segments */ 154c4762a1bSJed Brown while (fgets(line,MAXLINE,fp)) { 155c4762a1bSJed Brown 156c4762a1bSJed Brown if (strstr(line,"[TITLE]")) { 157c4762a1bSJed Brown GetDataSegment(fp,line,&title_start_pos,&ntitle); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown if (strstr(line,"[JUNCTIONS]")) { 161c4762a1bSJed Brown GetDataSegment(fp,line,&junc_start_pos,&nlines); 162c4762a1bSJed Brown water->nvertex += nlines; 163c4762a1bSJed Brown water->njunction = nlines; 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 166c4762a1bSJed Brown if (strstr(line,"[RESERVOIRS]")) { 167c4762a1bSJed Brown GetDataSegment(fp,line,&res_start_pos,&nlines); 168c4762a1bSJed Brown water->nvertex += nlines; 169c4762a1bSJed Brown water->nreservoir = nlines; 170c4762a1bSJed Brown } 171c4762a1bSJed Brown 172c4762a1bSJed Brown if (strstr(line,"[TANKS]")) { 173c4762a1bSJed Brown GetDataSegment(fp,line,&tank_start_pos,&nlines); 174c4762a1bSJed Brown water->nvertex += nlines; 175c4762a1bSJed Brown water->ntank = nlines; 176c4762a1bSJed Brown } 177c4762a1bSJed Brown 178c4762a1bSJed Brown if (strstr(line,"[PIPES]")) { 179c4762a1bSJed Brown GetDataSegment(fp,line,&pipe_start_pos,&nlines); 180c4762a1bSJed Brown water->nedge += nlines; 181c4762a1bSJed Brown water->npipe = nlines; 182c4762a1bSJed Brown } 183c4762a1bSJed Brown 184c4762a1bSJed Brown if (strstr(line,"[PUMPS]")) { 185c4762a1bSJed Brown GetDataSegment(fp,line,&pump_start_pos,&nlines); 186c4762a1bSJed Brown water->nedge += nlines; 187c4762a1bSJed Brown water->npump = nlines; 188c4762a1bSJed Brown } 189c4762a1bSJed Brown 190c4762a1bSJed Brown if (strstr(line,"[CURVES]")) { 191c4762a1bSJed Brown GetDataSegment(fp,line,&curve_start_pos,&ncurve); 192c4762a1bSJed Brown } 193c4762a1bSJed Brown } 194c4762a1bSJed Brown 195c4762a1bSJed Brown /* Allocate vertex and edge data structs */ 1969566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(water->nvertex,&water->vertex)); 1979566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(water->nedge,&water->edge)); 198c4762a1bSJed Brown vert = water->vertex; 199c4762a1bSJed Brown edge = water->edge; 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* Junctions */ 202c4762a1bSJed Brown fsetpos(fp,&junc_start_pos); 203c4762a1bSJed Brown for (i=0; i < water->njunction; i++) { 204c4762a1bSJed Brown int id=0,pattern=0; 20508401ef6SPierre Jolivet PetscCheck(fgets(line,MAXLINE,fp),PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Cannot read junction from file"); 206c4762a1bSJed Brown vert[nv].type = VERTEX_TYPE_JUNCTION; 207c4762a1bSJed Brown junction = &vert[nv].junc; 20808401ef6SPierre Jolivet ndata = sscanf(line,"%d %lf %lf %d",&id,&v1,&v2,&pattern);PetscCheck(ndata >= 3,PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Unable to read junction data"); 209c4762a1bSJed Brown vert[nv].id = id; 210c4762a1bSJed Brown junction->dempattern = pattern; 211c4762a1bSJed Brown junction->elev = (PetscScalar)v1; 212c4762a1bSJed Brown junction->demand = (PetscScalar)v2; 213c4762a1bSJed Brown junction->demand *= GPM_CFS; 214c4762a1bSJed Brown junction->id = vert[nv].id; 215c4762a1bSJed Brown nv++; 216c4762a1bSJed Brown } 217c4762a1bSJed Brown 218c4762a1bSJed Brown /* Reservoirs */ 219c4762a1bSJed Brown fsetpos(fp,&res_start_pos); 220c4762a1bSJed Brown for (i=0; i < water->nreservoir; i++) { 221c4762a1bSJed Brown int id=0,pattern=0; 22208401ef6SPierre Jolivet PetscCheck(fgets(line,MAXLINE,fp),PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Cannot read reservoir from file"); 223c4762a1bSJed Brown vert[nv].type = VERTEX_TYPE_RESERVOIR; 224c4762a1bSJed Brown reservoir = &vert[nv].res; 22508401ef6SPierre Jolivet ndata = sscanf(line,"%d %lf %d",&id,&v1,&pattern);PetscCheck(ndata >= 2,PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Unable to read reservoir data"); 226c4762a1bSJed Brown vert[nv].id = id; 227c4762a1bSJed Brown reservoir->headpattern = pattern; 228c4762a1bSJed Brown reservoir->head = (PetscScalar)v1; 229c4762a1bSJed Brown reservoir->id = vert[nv].id; 230c4762a1bSJed Brown nv++; 231c4762a1bSJed Brown } 232c4762a1bSJed Brown 233c4762a1bSJed Brown /* Tanks */ 234c4762a1bSJed Brown fsetpos(fp,&tank_start_pos); 235c4762a1bSJed Brown for (i=0; i < water->ntank; i++) { 236c4762a1bSJed Brown int id=0,curve=0; 23708401ef6SPierre Jolivet PetscCheck(fgets(line,MAXLINE,fp),PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Cannot read data tank from file"); 238c4762a1bSJed Brown vert[nv].type = VERTEX_TYPE_TANK; 239c4762a1bSJed Brown tank = &vert[nv].tank; 24008401ef6SPierre Jolivet ndata = sscanf(line,"%d %lf %lf %lf %lf %lf %lf %d",&id,&v1,&v2,&v3,&v4,&v5,&v6,&curve);PetscCheck(ndata >= 7,PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Unable to read tank data"); 241c4762a1bSJed Brown vert[nv].id = id; 242c4762a1bSJed Brown tank->volumecurve = curve; 243c4762a1bSJed Brown tank->elev = (PetscScalar)v1; 244c4762a1bSJed Brown tank->initlvl = (PetscScalar)v2; 245c4762a1bSJed Brown tank->minlvl = (PetscScalar)v3; 246c4762a1bSJed Brown tank->maxlvl = (PetscScalar)v4; 247c4762a1bSJed Brown tank->diam = (PetscScalar)v5; 248c4762a1bSJed Brown tank->minvolume = (PetscScalar)v6; 249c4762a1bSJed Brown tank->id = vert[nv].id; 250c4762a1bSJed Brown nv++; 251c4762a1bSJed Brown } 252c4762a1bSJed Brown 253c4762a1bSJed Brown /* Pipes */ 254c4762a1bSJed Brown fsetpos(fp,&pipe_start_pos); 255c4762a1bSJed Brown for (i=0; i < water->npipe; i++) { 256c4762a1bSJed Brown int id=0,node1=0,node2=0; 25708401ef6SPierre Jolivet PetscCheck(fgets(line,MAXLINE,fp),PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Cannot read data pipe from file"); 258c4762a1bSJed Brown edge[ne].type = EDGE_TYPE_PIPE; 259c4762a1bSJed Brown pipe = &edge[ne].pipe; 260c4762a1bSJed Brown ndata = sscanf(line,"%d %d %d %lf %lf %lf %lf %s",&id,&node1,&node2,&v1,&v2,&v3,&v4,pipe->stat); 261c4762a1bSJed Brown pipe->id = id; 262c4762a1bSJed Brown pipe->node1 = node1; 263c4762a1bSJed Brown pipe->node2 = node2; 264c4762a1bSJed Brown pipe->length = (PetscScalar)v1; 265c4762a1bSJed Brown pipe->diam = (PetscScalar)v2; 266c4762a1bSJed Brown pipe->roughness = (PetscScalar)v3; 267c4762a1bSJed Brown pipe->minorloss = (PetscScalar)v4; 268c4762a1bSJed Brown edge[ne].id = pipe->id; 269c4762a1bSJed Brown if (strcmp(pipe->stat,"OPEN") == 0) pipe->status = PIPE_STATUS_OPEN; 270c4762a1bSJed Brown if (ndata < 8) { 271c4762a1bSJed Brown strcpy(pipe->stat,"OPEN"); /* default OPEN */ 272c4762a1bSJed Brown pipe->status = PIPE_STATUS_OPEN; 273c4762a1bSJed Brown } 274c4762a1bSJed Brown if (ndata < 7) pipe->minorloss = 0.; 275c4762a1bSJed Brown pipe->n = 1.85; 276c4762a1bSJed Brown pipe->k = 4.72*pipe->length/(PetscPowScalar(pipe->roughness,pipe->n)*PetscPowScalar(0.0833333*pipe->diam,4.87)); 277c4762a1bSJed Brown ne++; 278c4762a1bSJed Brown } 279c4762a1bSJed Brown 280c4762a1bSJed Brown /* Pumps */ 281c4762a1bSJed Brown fsetpos(fp,&pump_start_pos); 282c4762a1bSJed Brown for (i=0; i < water->npump; i++) { 283c4762a1bSJed Brown int id=0,node1=0,node2=0,paramid=0; 28408401ef6SPierre Jolivet PetscCheck(fgets(line,MAXLINE,fp),PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Cannot read data pump from file"); 285c4762a1bSJed Brown edge[ne].type = EDGE_TYPE_PUMP; 286c4762a1bSJed Brown pump = &edge[ne].pump; 28708401ef6SPierre Jolivet ndata = sscanf(line,"%d %d %d %s %d",&id,&node1,&node2,pump->param,¶mid);PetscCheck(ndata == 5,PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Unable to read pump data"); 288c4762a1bSJed Brown pump->id = id; 289c4762a1bSJed Brown pump->node1 = node1; 290c4762a1bSJed Brown pump->node2 = node2; 291c4762a1bSJed Brown pump->paramid = paramid; 292c4762a1bSJed Brown edge[ne].id = pump->id; 293c4762a1bSJed Brown ne++; 294c4762a1bSJed Brown } 295c4762a1bSJed Brown 296c4762a1bSJed Brown /* Curves */ 297c4762a1bSJed Brown fsetpos(fp,&curve_start_pos); 298c4762a1bSJed Brown for (i=0; i < ncurve; i++) { 299c4762a1bSJed Brown int icurve_id=0; 30008401ef6SPierre Jolivet PetscCheck(fgets(line,MAXLINE,fp),PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Cannot read data curve from file"); 30108401ef6SPierre Jolivet ndata = sscanf(line,"%d %lf %lf",&icurve_id,&v1,&v2);PetscCheck(ndata == 3,PETSC_COMM_SELF,PETSC_ERR_FILE_READ,"Unable to read curve data"); 302c4762a1bSJed Brown curve_id = icurve_id; 303c4762a1bSJed Brown curve_x = (PetscScalar)v1; 304c4762a1bSJed Brown curve_y = (PetscScalar)v2; 305c4762a1bSJed Brown /* Check for pump with the curve_id */ 306c4762a1bSJed Brown for (j=water->npipe;j < water->npipe+water->npump;j++) { 307c4762a1bSJed Brown if (water->edge[j].pump.paramid == curve_id) { 308*7a46b595SBarry 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); 309c4762a1bSJed Brown pump = &water->edge[j].pump; 310c4762a1bSJed Brown pump->headcurve.flow[pump->headcurve.npt] = curve_x*GPM_CFS; 311c4762a1bSJed Brown pump->headcurve.head[pump->headcurve.npt] = curve_y; 312c4762a1bSJed Brown pump->headcurve.npt++; 313c4762a1bSJed Brown break; 314c4762a1bSJed Brown } 315c4762a1bSJed Brown } 316c4762a1bSJed Brown } 317c4762a1bSJed Brown 318c4762a1bSJed Brown fclose(fp); 319c4762a1bSJed Brown 320c4762a1bSJed Brown /* Get pump curve parameters */ 321c4762a1bSJed Brown for (j=water->npipe;j < water->npipe+water->npump;j++) { 322c4762a1bSJed Brown pump = &water->edge[j].pump; 323c4762a1bSJed Brown if (strcmp(pump->param,"HEAD") == 0) { 324c4762a1bSJed Brown /* Head-flow curve */ 3259566063dSJacob Faibussowitsch PetscCall(SetPumpHeadCurveParams(pump)); 326c4762a1bSJed Brown } 327c4762a1bSJed Brown } 328c4762a1bSJed Brown PetscFunctionReturn(0); 329c4762a1bSJed Brown } 330