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