xref: /petsc/src/snes/tutorials/network/water/waterfunctions.c (revision ed4474eb4a370bd0b975935bc736ff701abad52a)
1c4762a1bSJed Brown /* function subroutines used by water.c */
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include "water.h"
4c4762a1bSJed Brown #include <petscdmnetwork.h>
5c4762a1bSJed Brown 
Flow_Pipe(Pipe * pipe,PetscScalar hf,PetscScalar ht)6d71ae5a4SJacob Faibussowitsch PetscScalar Flow_Pipe(Pipe *pipe, PetscScalar hf, PetscScalar ht)
7d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown   PetscScalar flow_pipe;
9c4762a1bSJed Brown 
1057508eceSPierre Jolivet   flow_pipe = PetscSign(hf - ht) * PetscPowScalar(PetscAbsScalar(hf - ht) / pipe->k, 1 / pipe->n);
11c4762a1bSJed Brown   return flow_pipe;
12c4762a1bSJed Brown }
13c4762a1bSJed Brown 
Flow_Pump(Pump * pump,PetscScalar hf,PetscScalar ht)14d71ae5a4SJacob Faibussowitsch PetscScalar Flow_Pump(Pump *pump, PetscScalar hf, PetscScalar ht)
15d71ae5a4SJacob Faibussowitsch {
16c4762a1bSJed Brown   PetscScalar flow_pump;
17*097b9c6cSBarry Smith   flow_pump = PetscSign(hf - ht + pump->h0) * PetscPowScalar(PetscAbsScalar(hf - ht + pump->h0) / pump->r, 1 / pump->n);
18c4762a1bSJed Brown   return flow_pump;
19c4762a1bSJed Brown }
20c4762a1bSJed Brown 
FormFunction_Water(DM networkdm,Vec localX,Vec localF,PetscInt nv,PetscInt ne,const PetscInt * vtx,const PetscInt * edges,void * appctx)21d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction_Water(DM networkdm, Vec localX, Vec localF, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx)
22d71ae5a4SJacob Faibussowitsch {
23c4762a1bSJed Brown   const PetscScalar *xarr;
24c4762a1bSJed Brown   const PetscInt    *cone;
25c4762a1bSJed Brown   PetscScalar       *farr, hf, ht, flow;
262bf73ac6SHong Zhang   PetscInt           i, key, vnode1, vnode2, offsetnode1, offsetnode2, offset, ncomp;
27c4762a1bSJed Brown   PetscBool          ghostvtex;
28c4762a1bSJed Brown   VERTEX_Water       vertex, vertexnode1, vertexnode2;
29c4762a1bSJed Brown   EDGE_Water         edge;
30c4762a1bSJed Brown   Pipe              *pipe;
31c4762a1bSJed Brown   Pump              *pump;
32c4762a1bSJed Brown   Reservoir         *reservoir;
33c4762a1bSJed Brown   Tank              *tank;
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   PetscFunctionBegin;
36c4762a1bSJed Brown   /* Get arrays for the vectors */
379566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(localX, &xarr));
389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localF, &farr));
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   for (i = 0; i < ne; i++) { /* for each edge */
41c4762a1bSJed Brown     /* Get the offset and the key for the component for edge number e[i] */
429566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(networkdm, edges[i], 0, &key, (void **)&edge, NULL));
43c4762a1bSJed Brown 
44c4762a1bSJed Brown     /* Get the numbers for the vertices covering this edge */
459566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetConnectedVertices(networkdm, edges[i], &cone));
46c4762a1bSJed Brown     vnode1 = cone[0];
47c4762a1bSJed Brown     vnode2 = cone[1];
48c4762a1bSJed Brown 
492bf73ac6SHong Zhang     /* Get the components at the two vertices, their variable offsets */
509566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetNumComponents(networkdm, vnode1, &ncomp));
519566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(networkdm, vnode1, ncomp - 1, &key, (void **)&vertexnode1, NULL));
529566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(networkdm, vnode1, ncomp - 1, &offsetnode1));
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetNumComponents(networkdm, vnode2, &ncomp));
559566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(networkdm, vnode2, ncomp - 1, &key, (void **)&vertexnode2, NULL));
569566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(networkdm, vnode2, ncomp - 1, &offsetnode2));
57c4762a1bSJed Brown 
58c4762a1bSJed Brown     /* Variables at node1 and node 2 */
59c4762a1bSJed Brown     hf = xarr[offsetnode1];
60c4762a1bSJed Brown     ht = xarr[offsetnode2];
61c4762a1bSJed Brown 
62c4762a1bSJed Brown     flow = 0.0;
63c4762a1bSJed Brown     if (edge->type == EDGE_TYPE_PIPE) {
64c4762a1bSJed Brown       pipe = &edge->pipe;
65c4762a1bSJed Brown       flow = Flow_Pipe(pipe, hf, ht);
66c4762a1bSJed Brown     } else if (edge->type == EDGE_TYPE_PUMP) {
67c4762a1bSJed Brown       pump = &edge->pump;
68c4762a1bSJed Brown       flow = Flow_Pump(pump, hf, ht);
69c4762a1bSJed Brown     }
70c4762a1bSJed Brown     /* Convention: Node 1 has outgoing flow and Node 2 has incoming flow */
71c4762a1bSJed Brown     if (vertexnode1->type == VERTEX_TYPE_JUNCTION) farr[offsetnode1] -= flow;
72c4762a1bSJed Brown     if (vertexnode2->type == VERTEX_TYPE_JUNCTION) farr[offsetnode2] += flow;
73c4762a1bSJed Brown   }
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   /* Subtract Demand flows at the vertices */
76c4762a1bSJed Brown   for (i = 0; i < nv; i++) {
779566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[i], &ghostvtex));
78c4762a1bSJed Brown     if (ghostvtex) continue;
79c4762a1bSJed Brown 
809566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(networkdm, vtx[i], ALL_COMPONENTS, &offset));
819566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetNumComponents(networkdm, vtx[i], &ncomp));
829566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(networkdm, vtx[i], ncomp - 1, &key, (void **)&vertex, NULL));
83c4762a1bSJed Brown 
84c4762a1bSJed Brown     if (vertex->type == VERTEX_TYPE_JUNCTION) {
85c4762a1bSJed Brown       farr[offset] -= vertex->junc.demand;
86c4762a1bSJed Brown     } else if (vertex->type == VERTEX_TYPE_RESERVOIR) {
87c4762a1bSJed Brown       reservoir    = &vertex->res;
88c4762a1bSJed Brown       farr[offset] = xarr[offset] - reservoir->head;
89c4762a1bSJed Brown     } else if (vertex->type == VERTEX_TYPE_TANK) {
90c4762a1bSJed Brown       tank         = &vertex->tank;
91c4762a1bSJed Brown       farr[offset] = xarr[offset] - (tank->elev + tank->initlvl);
92c4762a1bSJed Brown     }
93c4762a1bSJed Brown   }
94c4762a1bSJed Brown 
959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(localX, &xarr));
969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localF, &farr));
973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
98c4762a1bSJed Brown }
99c4762a1bSJed Brown 
WaterFormFunction(SNES snes,Vec X,Vec F,void * user)100d71ae5a4SJacob Faibussowitsch PetscErrorCode WaterFormFunction(SNES snes, Vec X, Vec F, void *user)
101d71ae5a4SJacob Faibussowitsch {
102c4762a1bSJed Brown   DM              networkdm;
103c4762a1bSJed Brown   Vec             localX, localF;
104c4762a1bSJed Brown   const PetscInt *v, *e;
105c4762a1bSJed Brown   PetscInt        nv, ne;
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   PetscFunctionBegin;
108c4762a1bSJed Brown   /* Get the DM attached with the SNES */
1099566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &networkdm));
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   /* Get local vertices and edges */
1129566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetSubnetwork(networkdm, 0, &nv, &ne, &v, &e));
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   /* Get local vectors */
1159566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(networkdm, &localX));
1169566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(networkdm, &localF));
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   /* Scatter values from global to local vector */
1199566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX));
1209566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX));
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   /* Initialize residual */
1239566063dSJacob Faibussowitsch   PetscCall(VecSet(F, 0.0));
1249566063dSJacob Faibussowitsch   PetscCall(VecSet(localF, 0.0));
125c4762a1bSJed Brown 
1269566063dSJacob Faibussowitsch   PetscCall(FormFunction_Water(networkdm, localX, localF, nv, ne, v, e, NULL));
127c4762a1bSJed Brown 
1289566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(networkdm, &localX));
1299566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(networkdm, localF, ADD_VALUES, F));
1309566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(networkdm, localF, ADD_VALUES, F));
131c4762a1bSJed Brown 
1329566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(networkdm, &localF));
1333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
134c4762a1bSJed Brown }
135c4762a1bSJed Brown 
WaterSetInitialGuess(DM networkdm,Vec X)136d71ae5a4SJacob Faibussowitsch PetscErrorCode WaterSetInitialGuess(DM networkdm, Vec X)
137d71ae5a4SJacob Faibussowitsch {
138c4762a1bSJed Brown   PetscInt        nv, ne;
139c4762a1bSJed Brown   const PetscInt *vtx, *edges;
140c4762a1bSJed Brown   Vec             localX;
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   PetscFunctionBegin;
1439566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(networkdm, &localX));
144c4762a1bSJed Brown 
1459566063dSJacob Faibussowitsch   PetscCall(VecSet(X, 0.0));
1469566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX));
1479566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX));
148c4762a1bSJed Brown 
1492bf73ac6SHong Zhang   /* Get subnetwork */
1509566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetSubnetwork(networkdm, 0, &nv, &ne, &vtx, &edges));
1519566063dSJacob Faibussowitsch   PetscCall(SetInitialGuess_Water(networkdm, localX, nv, ne, vtx, edges, NULL));
152c4762a1bSJed Brown 
1539566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(networkdm, localX, ADD_VALUES, X));
1549566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(networkdm, localX, ADD_VALUES, X));
1559566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(networkdm, &localX));
1563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
157c4762a1bSJed Brown }
158c4762a1bSJed Brown 
GetListofEdges_Water(WATERDATA * water,PetscInt * edgelist)159d71ae5a4SJacob Faibussowitsch PetscErrorCode GetListofEdges_Water(WATERDATA *water, PetscInt *edgelist)
160d71ae5a4SJacob Faibussowitsch {
161c4762a1bSJed Brown   PetscInt  i, j, node1, node2;
162c4762a1bSJed Brown   Pipe     *pipe;
163c4762a1bSJed Brown   Pump     *pump;
164c4762a1bSJed Brown   PetscBool netview = PETSC_FALSE;
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   PetscFunctionBegin;
1679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-water_view", &netview));
168c4762a1bSJed Brown   for (i = 0; i < water->nedge; i++) {
169c4762a1bSJed Brown     if (water->edge[i].type == EDGE_TYPE_PIPE) {
170c4762a1bSJed Brown       pipe  = &water->edge[i].pipe;
171c4762a1bSJed Brown       node1 = pipe->node1;
172c4762a1bSJed Brown       node2 = pipe->node2;
17348a46eb9SPierre Jolivet       if (netview) PetscCall(PetscPrintf(PETSC_COMM_SELF, "edge %" PetscInt_FMT ", pipe v[%" PetscInt_FMT "] -> v[%" PetscInt_FMT "]\n", i, node1, node2));
174c4762a1bSJed Brown     } else {
175c4762a1bSJed Brown       pump  = &water->edge[i].pump;
176c4762a1bSJed Brown       node1 = pump->node1;
177c4762a1bSJed Brown       node2 = pump->node2;
17848a46eb9SPierre Jolivet       if (netview) PetscCall(PetscPrintf(PETSC_COMM_SELF, "edge %" PetscInt_FMT ", pump v[%" PetscInt_FMT "] -> v[%" PetscInt_FMT "]\n", i, node1, node2));
179c4762a1bSJed Brown     }
180c4762a1bSJed Brown 
181c4762a1bSJed Brown     for (j = 0; j < water->nvertex; j++) {
182c4762a1bSJed Brown       if (water->vertex[j].id == node1) {
183c4762a1bSJed Brown         edgelist[2 * i] = j;
184c4762a1bSJed Brown         break;
185c4762a1bSJed Brown       }
186c4762a1bSJed Brown     }
187c4762a1bSJed Brown 
188c4762a1bSJed Brown     for (j = 0; j < water->nvertex; j++) {
189c4762a1bSJed Brown       if (water->vertex[j].id == node2) {
190c4762a1bSJed Brown         edgelist[2 * i + 1] = j;
191c4762a1bSJed Brown         break;
192c4762a1bSJed Brown       }
193c4762a1bSJed Brown     }
194c4762a1bSJed Brown   }
1953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
196c4762a1bSJed Brown }
197c4762a1bSJed Brown 
SetInitialGuess_Water(DM networkdm,Vec localX,PetscInt nv,PetscInt ne,const PetscInt * vtx,const PetscInt * edges,void * appctx)198d71ae5a4SJacob Faibussowitsch PetscErrorCode SetInitialGuess_Water(DM networkdm, Vec localX, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx)
199d71ae5a4SJacob Faibussowitsch {
200c4762a1bSJed Brown   PetscInt     i, offset, key;
2012bf73ac6SHong Zhang   PetscBool    ghostvtex, sharedv;
202c4762a1bSJed Brown   VERTEX_Water vertex;
203c4762a1bSJed Brown   PetscScalar *xarr;
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   PetscFunctionBegin;
2069566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX, &xarr));
207c4762a1bSJed Brown   for (i = 0; i < nv; i++) {
2089566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[i], &ghostvtex));
2099566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsSharedVertex(networkdm, vtx[i], &sharedv));
2102bf73ac6SHong Zhang     if (ghostvtex || sharedv) continue;
211c4762a1bSJed Brown 
2129566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(networkdm, vtx[i], 0, &key, (void **)&vertex, NULL));
2139566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(networkdm, vtx[i], 0, &offset));
214c4762a1bSJed Brown     if (vertex->type == VERTEX_TYPE_JUNCTION) {
215c4762a1bSJed Brown       xarr[offset] = 100;
216c4762a1bSJed Brown     } else if (vertex->type == VERTEX_TYPE_RESERVOIR) {
217c4762a1bSJed Brown       xarr[offset] = vertex->res.head;
218c4762a1bSJed Brown     } else {
219c4762a1bSJed Brown       xarr[offset] = vertex->tank.initlvl + vertex->tank.elev;
220c4762a1bSJed Brown     }
221c4762a1bSJed Brown   }
2229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX, &xarr));
2233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
224c4762a1bSJed Brown }
225