1c4762a1bSJed Brown /* function subroutines used by water.c */ 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include "water.h" 4c4762a1bSJed Brown #include <petscdmnetwork.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown PetscScalar Flow_Pipe(Pipe *pipe,PetscScalar hf,PetscScalar ht) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown PetscScalar flow_pipe; 9c4762a1bSJed Brown 10c4762a1bSJed Brown flow_pipe = PetscSign(hf-ht)*PetscPowScalar(PetscAbsScalar(hf - ht)/pipe->k,(1/pipe->n)); 11c4762a1bSJed Brown return flow_pipe; 12c4762a1bSJed Brown } 13c4762a1bSJed Brown 14c4762a1bSJed Brown PetscScalar Flow_Pump(Pump *pump,PetscScalar hf, PetscScalar ht) 15c4762a1bSJed Brown { 16c4762a1bSJed Brown PetscScalar flow_pump; 17c4762a1bSJed Brown flow_pump = PetscPowScalar((hf - ht + pump->h0)/pump->r,(1/pump->n)); 18c4762a1bSJed Brown return flow_pump; 19c4762a1bSJed Brown } 20c4762a1bSJed Brown 21c4762a1bSJed Brown PetscErrorCode FormFunction_Water(DM networkdm,Vec localX,Vec localF,PetscInt nv,PetscInt ne,const PetscInt* vtx,const PetscInt* edges,void* appctx) 22c4762a1bSJed Brown { 23c4762a1bSJed Brown PetscErrorCode ierr; 24c4762a1bSJed Brown const PetscScalar *xarr; 25c4762a1bSJed Brown const PetscInt *cone; 26c4762a1bSJed Brown PetscScalar *farr,hf,ht,flow; 27*2bf73ac6SHong Zhang PetscInt i,key,vnode1,vnode2,offsetnode1,offsetnode2,offset,ncomp; 28c4762a1bSJed Brown PetscBool ghostvtex; 29c4762a1bSJed Brown VERTEX_Water vertex,vertexnode1,vertexnode2; 30c4762a1bSJed Brown EDGE_Water edge; 31c4762a1bSJed Brown Pipe *pipe; 32c4762a1bSJed Brown Pump *pump; 33c4762a1bSJed Brown Reservoir *reservoir; 34c4762a1bSJed Brown Tank *tank; 35c4762a1bSJed Brown 36c4762a1bSJed Brown PetscFunctionBegin; 37c4762a1bSJed Brown /* Get arrays for the vectors */ 38c4762a1bSJed Brown ierr = VecGetArrayRead(localX,&xarr);CHKERRQ(ierr); 39c4762a1bSJed Brown ierr = VecGetArray(localF,&farr);CHKERRQ(ierr); 40c4762a1bSJed Brown 41c4762a1bSJed Brown for (i=0; i<ne; i++) { /* for each edge */ 42c4762a1bSJed Brown /* Get the offset and the key for the component for edge number e[i] */ 43*2bf73ac6SHong Zhang ierr = DMNetworkGetComponent(networkdm,edges[i],0,&key,(void**)&edge,NULL);CHKERRQ(ierr); 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* Get the numbers for the vertices covering this edge */ 46c4762a1bSJed Brown ierr = DMNetworkGetConnectedVertices(networkdm,edges[i],&cone);CHKERRQ(ierr); 47c4762a1bSJed Brown vnode1 = cone[0]; 48c4762a1bSJed Brown vnode2 = cone[1]; 49c4762a1bSJed Brown 50*2bf73ac6SHong Zhang /* Get the components at the two vertices, their variable offsets */ 51*2bf73ac6SHong Zhang ierr = DMNetworkGetNumComponents(networkdm,vnode1,&ncomp);CHKERRQ(ierr); 52*2bf73ac6SHong Zhang ierr = DMNetworkGetComponent(networkdm,vnode1,ncomp-1,&key,(void**)&vertexnode1,NULL);CHKERRQ(ierr); 53*2bf73ac6SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vnode1,ncomp-1,&offsetnode1);CHKERRQ(ierr); 54c4762a1bSJed Brown 55*2bf73ac6SHong Zhang ierr = DMNetworkGetNumComponents(networkdm,vnode2,&ncomp);CHKERRQ(ierr); 56*2bf73ac6SHong Zhang ierr = DMNetworkGetComponent(networkdm,vnode2,ncomp-1,&key,(void**)&vertexnode2,NULL);CHKERRQ(ierr); 57*2bf73ac6SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vnode2,ncomp-1,&offsetnode2);CHKERRQ(ierr); 58c4762a1bSJed Brown 59c4762a1bSJed Brown /* Variables at node1 and node 2 */ 60c4762a1bSJed Brown hf = xarr[offsetnode1]; 61c4762a1bSJed Brown ht = xarr[offsetnode2]; 62c4762a1bSJed Brown 63c4762a1bSJed Brown flow = 0.0; 64c4762a1bSJed Brown if (edge->type == EDGE_TYPE_PIPE) { 65c4762a1bSJed Brown pipe = &edge->pipe; 66c4762a1bSJed Brown flow = Flow_Pipe(pipe,hf,ht); 67c4762a1bSJed Brown } else if (edge->type == EDGE_TYPE_PUMP) { 68c4762a1bSJed Brown pump = &edge->pump; 69c4762a1bSJed Brown flow = Flow_Pump(pump,hf,ht); 70c4762a1bSJed Brown } 71c4762a1bSJed Brown /* Convention: Node 1 has outgoing flow and Node 2 has incoming flow */ 72c4762a1bSJed Brown if (vertexnode1->type == VERTEX_TYPE_JUNCTION) farr[offsetnode1] -= flow; 73c4762a1bSJed Brown if (vertexnode2->type == VERTEX_TYPE_JUNCTION) farr[offsetnode2] += flow; 74c4762a1bSJed Brown } 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* Subtract Demand flows at the vertices */ 77c4762a1bSJed Brown for (i=0; i<nv; i++) { 78c4762a1bSJed Brown ierr = DMNetworkIsGhostVertex(networkdm,vtx[i],&ghostvtex);CHKERRQ(ierr); 79c4762a1bSJed Brown if (ghostvtex) continue; 80c4762a1bSJed Brown 81*2bf73ac6SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vtx[i],ALL_COMPONENTS,&offset);CHKERRQ(ierr); 82*2bf73ac6SHong Zhang ierr = DMNetworkGetNumComponents(networkdm,vtx[i],&ncomp);CHKERRQ(ierr); 83*2bf73ac6SHong Zhang ierr = DMNetworkGetComponent(networkdm,vtx[i],ncomp-1,&key,(void**)&vertex,NULL);CHKERRQ(ierr); 84c4762a1bSJed Brown 85c4762a1bSJed Brown if (vertex->type == VERTEX_TYPE_JUNCTION) { 86c4762a1bSJed Brown farr[offset] -= vertex->junc.demand; 87c4762a1bSJed Brown } else if (vertex->type == VERTEX_TYPE_RESERVOIR) { 88c4762a1bSJed Brown reservoir = &vertex->res; 89c4762a1bSJed Brown farr[offset] = xarr[offset] - reservoir->head; 90c4762a1bSJed Brown } else if (vertex->type == VERTEX_TYPE_TANK) { 91c4762a1bSJed Brown tank = &vertex->tank; 92c4762a1bSJed Brown farr[offset] = xarr[offset] - (tank->elev + tank->initlvl); 93c4762a1bSJed Brown } 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 96c4762a1bSJed Brown ierr = VecRestoreArrayRead(localX,&xarr);CHKERRQ(ierr); 97c4762a1bSJed Brown ierr = VecRestoreArray(localF,&farr);CHKERRQ(ierr); 98c4762a1bSJed Brown PetscFunctionReturn(0); 99c4762a1bSJed Brown } 100c4762a1bSJed Brown 101c4762a1bSJed Brown PetscErrorCode WaterFormFunction(SNES snes,Vec X, Vec F, void *user) 102c4762a1bSJed Brown { 103c4762a1bSJed Brown PetscErrorCode ierr; 104c4762a1bSJed Brown DM networkdm; 105c4762a1bSJed Brown Vec localX,localF; 106c4762a1bSJed Brown const PetscInt *v,*e; 107c4762a1bSJed Brown PetscInt nv,ne; 108c4762a1bSJed Brown 109c4762a1bSJed Brown PetscFunctionBegin; 110c4762a1bSJed Brown /* Get the DM attached with the SNES */ 111c4762a1bSJed Brown ierr = SNESGetDM(snes,&networkdm);CHKERRQ(ierr); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* Get local vertices and edges */ 114*2bf73ac6SHong Zhang ierr = DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&v,&e);CHKERRQ(ierr); 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* Get local vectors */ 117c4762a1bSJed Brown ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr); 118c4762a1bSJed Brown ierr = DMGetLocalVector(networkdm,&localF);CHKERRQ(ierr); 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* Scatter values from global to local vector */ 121c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 123c4762a1bSJed Brown 124c4762a1bSJed Brown /* Initialize residual */ 125c4762a1bSJed Brown ierr = VecSet(F,0.0);CHKERRQ(ierr); 126c4762a1bSJed Brown ierr = VecSet(localF,0.0);CHKERRQ(ierr); 127c4762a1bSJed Brown 128c4762a1bSJed Brown ierr = FormFunction_Water(networkdm,localX,localF,nv,ne,v,e,NULL);CHKERRQ(ierr); 129c4762a1bSJed Brown 130c4762a1bSJed Brown ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr); 131c4762a1bSJed Brown ierr = DMLocalToGlobalBegin(networkdm,localF,ADD_VALUES,F);CHKERRQ(ierr); 132c4762a1bSJed Brown ierr = DMLocalToGlobalEnd(networkdm,localF,ADD_VALUES,F);CHKERRQ(ierr); 133c4762a1bSJed Brown 134c4762a1bSJed Brown ierr = DMRestoreLocalVector(networkdm,&localF);CHKERRQ(ierr); 135c4762a1bSJed Brown PetscFunctionReturn(0); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown 138c4762a1bSJed Brown PetscErrorCode WaterSetInitialGuess(DM networkdm,Vec X) 139c4762a1bSJed Brown { 140c4762a1bSJed Brown PetscErrorCode ierr; 141c4762a1bSJed Brown PetscInt nv,ne; 142c4762a1bSJed Brown const PetscInt *vtx,*edges; 143c4762a1bSJed Brown Vec localX; 144c4762a1bSJed Brown 145c4762a1bSJed Brown PetscFunctionBegin; 146c4762a1bSJed Brown ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr); 147c4762a1bSJed Brown 148c4762a1bSJed Brown ierr = VecSet(X,0.0);CHKERRQ(ierr); 149c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 150c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 151c4762a1bSJed Brown 152*2bf73ac6SHong Zhang /* Get subnetwork */ 153*2bf73ac6SHong Zhang ierr = DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 154c4762a1bSJed Brown ierr = SetInitialGuess_Water(networkdm,localX,nv,ne,vtx,edges,NULL);CHKERRQ(ierr); 155c4762a1bSJed Brown 156c4762a1bSJed Brown ierr = DMLocalToGlobalBegin(networkdm,localX,ADD_VALUES,X);CHKERRQ(ierr); 157c4762a1bSJed Brown ierr = DMLocalToGlobalEnd(networkdm,localX,ADD_VALUES,X);CHKERRQ(ierr); 158c4762a1bSJed Brown ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr); 159c4762a1bSJed Brown PetscFunctionReturn(0); 160c4762a1bSJed Brown } 161c4762a1bSJed Brown 162c4762a1bSJed Brown PetscErrorCode GetListofEdges_Water(WATERDATA *water,PetscInt *edgelist) 163c4762a1bSJed Brown { 164c4762a1bSJed Brown PetscErrorCode ierr; 165c4762a1bSJed Brown PetscInt i,j,node1,node2; 166c4762a1bSJed Brown Pipe *pipe; 167c4762a1bSJed Brown Pump *pump; 168c4762a1bSJed Brown PetscBool netview=PETSC_FALSE; 169c4762a1bSJed Brown 170c4762a1bSJed Brown PetscFunctionBegin; 171c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL, "-water_view",&netview);CHKERRQ(ierr); 172c4762a1bSJed Brown for (i=0; i < water->nedge; i++) { 173c4762a1bSJed Brown if (water->edge[i].type == EDGE_TYPE_PIPE) { 174c4762a1bSJed Brown pipe = &water->edge[i].pipe; 175c4762a1bSJed Brown node1 = pipe->node1; 176c4762a1bSJed Brown node2 = pipe->node2; 177c4762a1bSJed Brown if (netview) { 178c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"edge %d, pipe v[%d] -> v[%d]\n",i,node1,node2);CHKERRQ(ierr); 179c4762a1bSJed Brown } 180c4762a1bSJed Brown } else { 181c4762a1bSJed Brown pump = &water->edge[i].pump; 182c4762a1bSJed Brown node1 = pump->node1; 183c4762a1bSJed Brown node2 = pump->node2; 184c4762a1bSJed Brown if (netview) { 185c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"edge %d, pump v[%d] -> v[%d]\n",i,node1,node2);CHKERRQ(ierr); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 189c4762a1bSJed Brown for (j=0; j < water->nvertex; j++) { 190c4762a1bSJed Brown if (water->vertex[j].id == node1) { 191c4762a1bSJed Brown edgelist[2*i] = j; 192c4762a1bSJed Brown break; 193c4762a1bSJed Brown } 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196c4762a1bSJed Brown for (j=0; j < water->nvertex; j++) { 197c4762a1bSJed Brown if (water->vertex[j].id == node2) { 198c4762a1bSJed Brown edgelist[2*i+1] = j; 199c4762a1bSJed Brown break; 200c4762a1bSJed Brown } 201c4762a1bSJed Brown } 202c4762a1bSJed Brown } 203c4762a1bSJed Brown PetscFunctionReturn(0); 204c4762a1bSJed Brown } 205c4762a1bSJed Brown 206c4762a1bSJed Brown PetscErrorCode SetInitialGuess_Water(DM networkdm,Vec localX,PetscInt nv,PetscInt ne, const PetscInt *vtx, const PetscInt *edges,void* appctx) 207c4762a1bSJed Brown { 208c4762a1bSJed Brown PetscErrorCode ierr; 209c4762a1bSJed Brown PetscInt i,offset,key; 210*2bf73ac6SHong Zhang PetscBool ghostvtex,sharedv; 211c4762a1bSJed Brown VERTEX_Water vertex; 212c4762a1bSJed Brown PetscScalar *xarr; 213c4762a1bSJed Brown 214c4762a1bSJed Brown PetscFunctionBegin; 215c4762a1bSJed Brown ierr = VecGetArray(localX,&xarr);CHKERRQ(ierr); 216c4762a1bSJed Brown for (i=0; i < nv; i++) { 217c4762a1bSJed Brown ierr = DMNetworkIsGhostVertex(networkdm,vtx[i],&ghostvtex);CHKERRQ(ierr); 218*2bf73ac6SHong Zhang ierr = DMNetworkIsSharedVertex(networkdm,vtx[i],&sharedv);CHKERRQ(ierr); 219*2bf73ac6SHong Zhang if (ghostvtex || sharedv) continue; 220c4762a1bSJed Brown 221*2bf73ac6SHong Zhang ierr = DMNetworkGetComponent(networkdm,vtx[i],0,&key,(void**)&vertex,NULL);CHKERRQ(ierr); 222*2bf73ac6SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vtx[i],0,&offset);CHKERRQ(ierr); 223c4762a1bSJed Brown if (vertex->type == VERTEX_TYPE_JUNCTION) { 224c4762a1bSJed Brown xarr[offset] = 100; 225c4762a1bSJed Brown } else if (vertex->type == VERTEX_TYPE_RESERVOIR) { 226c4762a1bSJed Brown xarr[offset] = vertex->res.head; 227c4762a1bSJed Brown } else { 228c4762a1bSJed Brown xarr[offset] = vertex->tank.initlvl + vertex->tank.elev; 229c4762a1bSJed Brown } 230c4762a1bSJed Brown } 231c4762a1bSJed Brown ierr = VecRestoreArray(localX,&xarr);CHKERRQ(ierr); 232c4762a1bSJed Brown PetscFunctionReturn(0); 233c4762a1bSJed Brown } 234