xref: /petsc/src/ts/tutorials/network/pipes.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
142dc13f1SHong Zhang static char help[] = "This example demonstrates DMNetwork. It is used for testing parallel generation of dmnetwork, then redistribute. \n\\n";
242dc13f1SHong Zhang /*
342dc13f1SHong Zhang   Example: mpiexec -n <np> ./pipes -ts_max_steps 10
442dc13f1SHong Zhang */
542dc13f1SHong Zhang 
642dc13f1SHong Zhang #include "wash.h"
742dc13f1SHong Zhang 
842dc13f1SHong Zhang /*
942dc13f1SHong Zhang   WashNetworkDistribute - proc[0] distributes sequential wash object
1042dc13f1SHong Zhang    Input Parameters:
1142dc13f1SHong Zhang .  comm - MPI communicator
1242dc13f1SHong Zhang .  wash - wash context with all network data in proc[0]
1342dc13f1SHong Zhang 
1442dc13f1SHong Zhang    Output Parameter:
1542dc13f1SHong Zhang .  wash - wash context with nedge, nvertex and edgelist distributed
1642dc13f1SHong Zhang 
1742dc13f1SHong Zhang    Note: The routine is used for testing parallel generation of dmnetwork, then redistribute.
1842dc13f1SHong Zhang */
1942dc13f1SHong Zhang PetscErrorCode WashNetworkDistribute(MPI_Comm comm,Wash wash)
2042dc13f1SHong Zhang {
2142dc13f1SHong Zhang   PetscMPIInt    rank,size,tag=0;
2242dc13f1SHong Zhang   PetscInt       i,e,v,numEdges,numVertices,nedges,*eowners=NULL,estart,eend,*vtype=NULL,nvertices;
2342dc13f1SHong Zhang   PetscInt       *edgelist = wash->edgelist,*nvtx=NULL,*vtxDone=NULL;
2442dc13f1SHong Zhang 
2542dc13f1SHong Zhang   PetscFunctionBegin;
26*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm,&size));
2742dc13f1SHong Zhang   if (size == 1) PetscFunctionReturn(0);
2842dc13f1SHong Zhang 
29*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm,&rank));
3042dc13f1SHong Zhang   numEdges    = wash->nedge;
3142dc13f1SHong Zhang   numVertices = wash->nvertex;
3242dc13f1SHong Zhang 
3342dc13f1SHong Zhang   /* (1) all processes get global and local number of edges */
34*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Bcast(&numEdges,1,MPIU_INT,0,comm));
3542dc13f1SHong Zhang   nedges = numEdges/size; /* local nedges */
36dd400576SPatrick Sanan   if (rank == 0) {
3742dc13f1SHong Zhang     nedges += numEdges - size*(numEdges/size);
3842dc13f1SHong Zhang   }
3942dc13f1SHong Zhang   wash->Nedge = numEdges;
4042dc13f1SHong Zhang   wash->nedge = nedges;
41*5f80ce2aSJacob Faibussowitsch   /* CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[%d] nedges %d, numEdges %d\n",rank,nedges,numEdges)); */
4242dc13f1SHong Zhang 
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc3(size+1,&eowners,size,&nvtx,numVertices,&vtxDone));
44*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allgather(&nedges,1,MPIU_INT,eowners+1,1,MPIU_INT,PETSC_COMM_WORLD));
4542dc13f1SHong Zhang   eowners[0] = 0;
4642dc13f1SHong Zhang   for (i=2; i<=size; i++) {
4742dc13f1SHong Zhang     eowners[i] += eowners[i-1];
4842dc13f1SHong Zhang   }
4942dc13f1SHong Zhang 
5042dc13f1SHong Zhang   estart = eowners[rank];
5142dc13f1SHong Zhang   eend   = eowners[rank+1];
52*5f80ce2aSJacob Faibussowitsch   /* CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[%d] own lists row %d - %d\n",rank,estart,eend)); */
5342dc13f1SHong Zhang 
5442dc13f1SHong Zhang   /* (2) distribute row block edgelist to all processors */
55dd400576SPatrick Sanan   if (rank == 0) {
5642dc13f1SHong Zhang     vtype = wash->vtype;
5742dc13f1SHong Zhang     for (i=1; i<size; i++) {
5842dc13f1SHong Zhang       /* proc[0] sends edgelist to proc[i] */
59*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Send(edgelist+2*eowners[i],2*(eowners[i+1]-eowners[i]),MPIU_INT,i,tag,comm));
6042dc13f1SHong Zhang 
6142dc13f1SHong Zhang       /* proc[0] sends vtype to proc[i] */
62*5f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Send(vtype+2*eowners[i],2*(eowners[i+1]-eowners[i]),MPIU_INT,i,tag,comm));
6342dc13f1SHong Zhang     }
6442dc13f1SHong Zhang   } else {
6542dc13f1SHong Zhang     MPI_Status      status;
66*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(2*(eend-estart),&vtype));
67*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(2*(eend-estart),&edgelist));
6842dc13f1SHong Zhang 
69*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Recv(edgelist,2*(eend-estart),MPIU_INT,0,tag,comm,&status));
70*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Recv(vtype,2*(eend-estart),MPIU_INT,0,tag,comm,&status));
7142dc13f1SHong Zhang   }
7242dc13f1SHong Zhang 
7342dc13f1SHong Zhang   wash->edgelist = edgelist;
7442dc13f1SHong Zhang 
7542dc13f1SHong Zhang   /* (3) all processes get global and local number of vertices, without ghost vertices */
76dd400576SPatrick Sanan   if (rank == 0) {
7742dc13f1SHong Zhang     for (i=0; i<size; i++) {
7842dc13f1SHong Zhang       for (e=eowners[i]; e<eowners[i+1]; e++) {
7942dc13f1SHong Zhang         v = edgelist[2*e];
8042dc13f1SHong Zhang         if (!vtxDone[v]) {
8142dc13f1SHong Zhang           nvtx[i]++; vtxDone[v] = 1;
8242dc13f1SHong Zhang         }
8342dc13f1SHong Zhang         v = edgelist[2*e+1];
8442dc13f1SHong Zhang         if (!vtxDone[v]) {
8542dc13f1SHong Zhang           nvtx[i]++; vtxDone[v] = 1;
8642dc13f1SHong Zhang         }
8742dc13f1SHong Zhang       }
8842dc13f1SHong Zhang     }
8942dc13f1SHong Zhang   }
90*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Bcast(&numVertices,1,MPIU_INT,0,PETSC_COMM_WORLD));
91*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Scatter(nvtx,1,MPIU_INT,&nvertices,1,MPIU_INT,0,PETSC_COMM_WORLD));
92*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(eowners,nvtx,vtxDone));
9342dc13f1SHong Zhang 
9442dc13f1SHong Zhang   wash->Nvertex = numVertices;
9542dc13f1SHong Zhang   wash->nvertex = nvertices;
9642dc13f1SHong Zhang   wash->vtype   = vtype;
9742dc13f1SHong Zhang   PetscFunctionReturn(0);
9842dc13f1SHong Zhang }
9942dc13f1SHong Zhang 
10042dc13f1SHong Zhang PetscErrorCode WASHIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void* ctx)
10142dc13f1SHong Zhang {
10242dc13f1SHong Zhang   Wash           wash=(Wash)ctx;
10342dc13f1SHong Zhang   DM             networkdm;
10442dc13f1SHong Zhang   Vec            localX,localXdot,localF, localXold;
10542dc13f1SHong Zhang   const PetscInt *cone;
10642dc13f1SHong Zhang   PetscInt       vfrom,vto,offsetfrom,offsetto,varoffset;
10742dc13f1SHong Zhang   PetscInt       v,vStart,vEnd,e,eStart,eEnd;
10842dc13f1SHong Zhang   PetscInt       nend,type;
10942dc13f1SHong Zhang   PetscBool      ghost;
11042dc13f1SHong Zhang   PetscScalar    *farr,*juncf, *pipef;
11142dc13f1SHong Zhang   PetscReal      dt;
11242dc13f1SHong Zhang   Pipe           pipe;
11342dc13f1SHong Zhang   PipeField      *pipex,*pipexdot,*juncx;
11442dc13f1SHong Zhang   Junction       junction;
11542dc13f1SHong Zhang   DMDALocalInfo  info;
11642dc13f1SHong Zhang   const PetscScalar *xarr,*xdotarr, *xoldarr;
11742dc13f1SHong Zhang 
11842dc13f1SHong Zhang   PetscFunctionBegin;
11942dc13f1SHong Zhang   localX    = wash->localX;
12042dc13f1SHong Zhang   localXdot = wash->localXdot;
12142dc13f1SHong Zhang 
122*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolution(ts,&localXold));
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&networkdm));
124*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTimeStep(ts,&dt));
125*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(networkdm,&localF));
12642dc13f1SHong Zhang 
12742dc13f1SHong Zhang   /* Set F and localF as zero */
128*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(F,0.0));
129*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(localF,0.0));
13042dc13f1SHong Zhang 
13142dc13f1SHong Zhang   /* Update ghost values of locaX and locaXdot */
132*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX));
133*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX));
13442dc13f1SHong Zhang 
135*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(networkdm,Xdot,INSERT_VALUES,localXdot));
136*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(networkdm,Xdot,INSERT_VALUES,localXdot));
13742dc13f1SHong Zhang 
138*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(localX,&xarr));
139*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(localXdot,&xdotarr));
140*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(localXold,&xoldarr));
141*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(localF,&farr));
14242dc13f1SHong Zhang 
14342dc13f1SHong Zhang    /* junction->type == JUNCTION:
14442dc13f1SHong Zhang            juncf[0] = -qJ + sum(qin); juncf[1] = qJ - sum(qout)
14542dc13f1SHong Zhang        junction->type == RESERVOIR (upper stream):
14642dc13f1SHong Zhang            juncf[0] = -hJ + H0; juncf[1] = qJ - sum(qout)
14742dc13f1SHong Zhang        junction->type == VALVE (down stream):
14842dc13f1SHong Zhang            juncf[0] =  -qJ + sum(qin); juncf[1] = qJ
14942dc13f1SHong Zhang   */
15042dc13f1SHong Zhang   /* Vertex/junction initialization */
151*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkGetVertexRange(networkdm,&vStart,&vEnd));
15242dc13f1SHong Zhang   for (v=vStart; v<vEnd; v++) {
153*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkIsGhostVertex(networkdm,v,&ghost));
15442dc13f1SHong Zhang     if (ghost) continue;
15542dc13f1SHong Zhang 
156*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetComponent(networkdm,v,0,&type,(void**)&junction,NULL));
157*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetLocalVecOffset(networkdm,v,ALL_COMPONENTS,&varoffset));
15842dc13f1SHong Zhang     juncx      = (PipeField*)(xarr + varoffset);
15942dc13f1SHong Zhang     juncf      = (PetscScalar*)(farr + varoffset);
16042dc13f1SHong Zhang 
16142dc13f1SHong Zhang     juncf[0] = -juncx[0].q;
16242dc13f1SHong Zhang     juncf[1] =  juncx[0].q;
16342dc13f1SHong Zhang 
16442dc13f1SHong Zhang     if (junction->type == RESERVOIR) { /* upstream reservoir */
16542dc13f1SHong Zhang       juncf[0] = juncx[0].h - wash->H0;
16642dc13f1SHong Zhang     }
16742dc13f1SHong Zhang   }
16842dc13f1SHong Zhang 
16942dc13f1SHong Zhang   /* Edge/pipe */
170*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd));
17142dc13f1SHong Zhang   for (e=eStart; e<eEnd; e++) {
172*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetComponent(networkdm,e,0,&type,(void**)&pipe,NULL));
173*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetLocalVecOffset(networkdm,e,ALL_COMPONENTS,&varoffset));
17442dc13f1SHong Zhang     pipex    = (PipeField*)(xarr + varoffset);
17542dc13f1SHong Zhang     pipexdot = (PipeField*)(xdotarr + varoffset);
17642dc13f1SHong Zhang     pipef    = (PetscScalar*)(farr + varoffset);
17742dc13f1SHong Zhang 
17842dc13f1SHong Zhang     /* Get some data into the pipe structure: note, some of these operations
17942dc13f1SHong Zhang      * might be redundant. Will it consume too much time? */
18042dc13f1SHong Zhang     pipe->dt   = dt;
18142dc13f1SHong Zhang     pipe->xold = (PipeField*)(xoldarr + varoffset);
18242dc13f1SHong Zhang 
18342dc13f1SHong Zhang     /* Evaluate F over this edge/pipe: pipef[1], ...,pipef[2*nend] */
184*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetLocalInfo(pipe->da,&info));
185*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PipeIFunctionLocal_Lax(&info,t,pipex,pipexdot,pipef,pipe));
18642dc13f1SHong Zhang 
18742dc13f1SHong Zhang     /* Get boundary values from connected vertices */
188*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetConnectedVertices(networkdm,e,&cone));
18942dc13f1SHong Zhang     vfrom = cone[0]; /* local ordering */
19042dc13f1SHong Zhang     vto   = cone[1];
191*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetLocalVecOffset(networkdm,vfrom,ALL_COMPONENTS,&offsetfrom));
192*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetLocalVecOffset(networkdm,vto,ALL_COMPONENTS,&offsetto));
19342dc13f1SHong Zhang 
19442dc13f1SHong Zhang     /* Evaluate upstream boundary */
195*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetComponent(networkdm,vfrom,0,&type,(void**)&junction,NULL));
1962c71b3e2SJacob Faibussowitsch     PetscCheckFalse(junction->type != JUNCTION && junction->type != RESERVOIR,PETSC_COMM_SELF,PETSC_ERR_SUP,"junction type is not supported");
19742dc13f1SHong Zhang     juncx = (PipeField*)(xarr + offsetfrom);
19842dc13f1SHong Zhang     juncf = (PetscScalar*)(farr + offsetfrom);
19942dc13f1SHong Zhang 
20042dc13f1SHong Zhang     pipef[0] = pipex[0].h - juncx[0].h;
20142dc13f1SHong Zhang     juncf[1] -= pipex[0].q;
20242dc13f1SHong Zhang 
20342dc13f1SHong Zhang     /* Evaluate downstream boundary */
204*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetComponent(networkdm,vto,0,&type,(void**)&junction,NULL));
2052c71b3e2SJacob Faibussowitsch     PetscCheckFalse(junction->type != JUNCTION && junction->type != VALVE,PETSC_COMM_SELF,PETSC_ERR_SUP,"junction type is not supported");
20642dc13f1SHong Zhang     juncx = (PipeField*)(xarr + offsetto);
20742dc13f1SHong Zhang     juncf = (PetscScalar*)(farr + offsetto);
20842dc13f1SHong Zhang     nend  = pipe->nnodes - 1;
20942dc13f1SHong Zhang 
21042dc13f1SHong Zhang     pipef[2*nend + 1] = pipex[nend].h - juncx[0].h;
21142dc13f1SHong Zhang     juncf[0] += pipex[nend].q;
21242dc13f1SHong Zhang   }
21342dc13f1SHong Zhang 
214*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(localX,&xarr));
215*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(localXdot,&xdotarr));
216*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(localF,&farr));
21742dc13f1SHong Zhang 
218*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalBegin(networkdm,localF,ADD_VALUES,F));
219*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalEnd(networkdm,localF,ADD_VALUES,F));
220*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(networkdm,&localF));
22142dc13f1SHong Zhang   /*
222*5f80ce2aSJacob Faibussowitsch    CHKERRQ(PetscPrintf(PETSC_COMM_WORLD("F:\n"));
223*5f80ce2aSJacob Faibussowitsch    CHKERRQ(VecView(F,PETSC_VIEWER_STDOUT_WORLD));
22442dc13f1SHong Zhang    */
22542dc13f1SHong Zhang   PetscFunctionReturn(0);
22642dc13f1SHong Zhang }
22742dc13f1SHong Zhang 
22842dc13f1SHong Zhang PetscErrorCode WASHSetInitialSolution(DM networkdm,Vec X,Wash wash)
22942dc13f1SHong Zhang {
23042dc13f1SHong Zhang   PetscInt       k,nx,vkey,vfrom,vto,offsetfrom,offsetto;
23142dc13f1SHong Zhang   PetscInt       type,varoffset;
23242dc13f1SHong Zhang   PetscInt       e,eStart,eEnd;
23342dc13f1SHong Zhang   Vec            localX;
23442dc13f1SHong Zhang   PetscScalar    *xarr;
23542dc13f1SHong Zhang   Pipe           pipe;
23642dc13f1SHong Zhang   Junction       junction;
23742dc13f1SHong Zhang   const PetscInt *cone;
23842dc13f1SHong Zhang   const PetscScalar *xarray;
23942dc13f1SHong Zhang 
24042dc13f1SHong Zhang   PetscFunctionBegin;
241*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(X,0.0));
242*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(networkdm,&localX));
243*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(localX,&xarr));
24442dc13f1SHong Zhang 
24542dc13f1SHong Zhang   /* Edge */
246*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd));
24742dc13f1SHong Zhang   for (e=eStart; e<eEnd; e++) {
248*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetLocalVecOffset(networkdm,e,ALL_COMPONENTS,&varoffset));
249*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetComponent(networkdm,e,0,&type,(void**)&pipe,NULL));
25042dc13f1SHong Zhang 
25142dc13f1SHong Zhang     /* set initial values for this pipe */
252*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PipeComputeSteadyState(pipe,wash->Q0,wash->H0));
253*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetSize(pipe->x,&nx));
25442dc13f1SHong Zhang 
255*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArrayRead(pipe->x,&xarray));
25642dc13f1SHong Zhang     /* copy pipe->x to xarray */
25742dc13f1SHong Zhang     for (k=0; k<nx; k++) {
25842dc13f1SHong Zhang       (xarr+varoffset)[k] = xarray[k];
25942dc13f1SHong Zhang     }
26042dc13f1SHong Zhang 
26142dc13f1SHong Zhang     /* set boundary values into vfrom and vto */
262*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetConnectedVertices(networkdm,e,&cone));
26342dc13f1SHong Zhang     vfrom = cone[0]; /* local ordering */
26442dc13f1SHong Zhang     vto   = cone[1];
265*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetLocalVecOffset(networkdm,vfrom,ALL_COMPONENTS,&offsetfrom));
266*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetLocalVecOffset(networkdm,vto,ALL_COMPONENTS,&offsetto));
26742dc13f1SHong Zhang 
26842dc13f1SHong Zhang     /* if vform is a head vertex: */
269*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetComponent(networkdm,vfrom,0,&vkey,(void**)&junction,NULL));
27042dc13f1SHong Zhang     if (junction->type == RESERVOIR) {
27142dc13f1SHong Zhang       (xarr+offsetfrom)[1] = wash->H0; /* 1st H */
27242dc13f1SHong Zhang     }
27342dc13f1SHong Zhang 
27442dc13f1SHong Zhang     /* if vto is an end vertex: */
275*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetComponent(networkdm,vto,0,&vkey,(void**)&junction,NULL));
27642dc13f1SHong Zhang     if (junction->type == VALVE) {
27742dc13f1SHong Zhang       (xarr+offsetto)[0] = wash->QL; /* last Q */
27842dc13f1SHong Zhang     }
279*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(pipe->x,&xarray));
28042dc13f1SHong Zhang   }
28142dc13f1SHong Zhang 
282*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(localX,&xarr));
283*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalBegin(networkdm,localX,ADD_VALUES,X));
284*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalEnd(networkdm,localX,ADD_VALUES,X));
285*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(networkdm,&localX));
28642dc13f1SHong Zhang 
28742dc13f1SHong Zhang #if 0
28842dc13f1SHong Zhang   PetscInt N;
289*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(X,&N));
290*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"initial solution %d:\n",N));
291*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(X,PETSC_VIEWER_STDOUT_WORLD));
29242dc13f1SHong Zhang #endif
29342dc13f1SHong Zhang   PetscFunctionReturn(0);
29442dc13f1SHong Zhang }
29542dc13f1SHong Zhang 
29642dc13f1SHong Zhang PetscErrorCode TSDMNetworkMonitor(TS ts, PetscInt step, PetscReal t, Vec x, void *context)
29742dc13f1SHong Zhang {
29842dc13f1SHong Zhang   DMNetworkMonitor   monitor;
29942dc13f1SHong Zhang 
30042dc13f1SHong Zhang   PetscFunctionBegin;
30142dc13f1SHong Zhang   monitor = (DMNetworkMonitor)context;
302*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkMonitorView(monitor,x));
30342dc13f1SHong Zhang   PetscFunctionReturn(0);
30442dc13f1SHong Zhang }
30542dc13f1SHong Zhang 
30642dc13f1SHong Zhang PetscErrorCode PipesView(DM networkdm,PetscInt KeyPipe,Vec X)
30742dc13f1SHong Zhang {
30842dc13f1SHong Zhang   PetscInt       i,numkeys=1,*blocksize,*numselectedvariable,**selectedvariables,n;
30942dc13f1SHong Zhang   IS             isfrom_q,isfrom_h,isfrom;
31042dc13f1SHong Zhang   Vec            Xto;
31142dc13f1SHong Zhang   VecScatter     ctx;
31242dc13f1SHong Zhang   MPI_Comm       comm;
31342dc13f1SHong Zhang 
31442dc13f1SHong Zhang   PetscFunctionBegin;
315*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject)networkdm,&comm));
31642dc13f1SHong Zhang 
31742dc13f1SHong Zhang   /* 1. Create isfrom_q for q-variable of pipes */
318*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc3(numkeys,&blocksize,numkeys,&numselectedvariable,numkeys,&selectedvariables));
31942dc13f1SHong Zhang   for (i=0; i<numkeys; i++) {
32042dc13f1SHong Zhang     blocksize[i]           = 2;
32142dc13f1SHong Zhang     numselectedvariable[i] = 1;
322*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(numselectedvariable[i],&selectedvariables[i]));
32342dc13f1SHong Zhang     selectedvariables[i][0] = 0; /* q-variable */
32442dc13f1SHong Zhang   }
325*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkCreateIS(networkdm,numkeys,&KeyPipe,blocksize,numselectedvariable,selectedvariables,&isfrom_q));
32642dc13f1SHong Zhang 
32742dc13f1SHong Zhang   /* 2. Create Xto and isto */
328*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(isfrom_q, &n));
329*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(comm,&Xto));
330*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(Xto,n,PETSC_DECIDE));
331*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(Xto));
332*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(Xto,0.0));
33342dc13f1SHong Zhang 
33442dc13f1SHong Zhang   /* 3. Create scatter */
335*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(X,isfrom_q,Xto,NULL,&ctx));
33642dc13f1SHong Zhang 
33742dc13f1SHong Zhang   /* 4. Scatter to Xq */
338*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(ctx,X,Xto,INSERT_VALUES,SCATTER_FORWARD));
339*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(ctx,X,Xto,INSERT_VALUES,SCATTER_FORWARD));
340*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&ctx));
341*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isfrom_q));
34242dc13f1SHong Zhang 
343*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Xq:\n"));
344*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(Xto,PETSC_VIEWER_STDOUT_WORLD));
34542dc13f1SHong Zhang 
34642dc13f1SHong Zhang   /* 5. Create isfrom_h for h-variable of pipes; Create scatter; Scatter to Xh */
34742dc13f1SHong Zhang   for (i=0; i<numkeys; i++) {
34842dc13f1SHong Zhang     selectedvariables[i][0] = 1; /* h-variable */
34942dc13f1SHong Zhang   }
350*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkCreateIS(networkdm,numkeys,&KeyPipe,blocksize,numselectedvariable,selectedvariables,&isfrom_h));
35142dc13f1SHong Zhang 
352*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(X,isfrom_h,Xto,NULL,&ctx));
353*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(ctx,X,Xto,INSERT_VALUES,SCATTER_FORWARD));
354*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(ctx,X,Xto,INSERT_VALUES,SCATTER_FORWARD));
355*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&ctx));
356*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isfrom_h));
35742dc13f1SHong Zhang 
358*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Xh:\n"));
359*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(Xto,PETSC_VIEWER_STDOUT_WORLD));
360*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Xto));
36142dc13f1SHong Zhang 
36242dc13f1SHong Zhang   /* 6. Create isfrom for all pipe variables; Create scatter; Scatter to Xpipes */
36342dc13f1SHong Zhang   for (i=0; i<numkeys; i++) {
36442dc13f1SHong Zhang     blocksize[i] = -1; /* select all the variables of the i-th component */
36542dc13f1SHong Zhang   }
366*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkCreateIS(networkdm,numkeys,&KeyPipe,blocksize,NULL,NULL,&isfrom));
367*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isfrom));
368*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkCreateIS(networkdm,numkeys,&KeyPipe,NULL,NULL,NULL,&isfrom));
36942dc13f1SHong Zhang 
370*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(isfrom, &n));
371*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(comm,&Xto));
372*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(Xto,n,PETSC_DECIDE));
373*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(Xto));
374*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(Xto,0.0));
37542dc13f1SHong Zhang 
376*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(X,isfrom,Xto,NULL,&ctx));
377*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(ctx,X,Xto,INSERT_VALUES,SCATTER_FORWARD));
378*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(ctx,X,Xto,INSERT_VALUES,SCATTER_FORWARD));
379*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&ctx));
380*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isfrom));
38142dc13f1SHong Zhang 
382*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Xpipes:\n"));
383*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(Xto,PETSC_VIEWER_STDOUT_WORLD));
38442dc13f1SHong Zhang 
38542dc13f1SHong Zhang   /* 7. Free spaces */
38642dc13f1SHong Zhang   for (i=0; i<numkeys; i++) {
387*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(selectedvariables[i]));
38842dc13f1SHong Zhang   }
389*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree3(blocksize,numselectedvariable,selectedvariables));
390*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Xto));
39142dc13f1SHong Zhang   PetscFunctionReturn(0);
39242dc13f1SHong Zhang }
39342dc13f1SHong Zhang 
39442dc13f1SHong Zhang PetscErrorCode ISJunctionsView(DM networkdm,PetscInt KeyJunc)
39542dc13f1SHong Zhang {
39642dc13f1SHong Zhang   PetscInt       numkeys=1;
39742dc13f1SHong Zhang   IS             isfrom;
39842dc13f1SHong Zhang   MPI_Comm       comm;
39942dc13f1SHong Zhang   PetscMPIInt    rank;
40042dc13f1SHong Zhang 
40142dc13f1SHong Zhang   PetscFunctionBegin;
402*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject)networkdm,&comm));
403*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm,&rank));
40442dc13f1SHong Zhang 
40542dc13f1SHong Zhang   /* Create a global isfrom for all junction variables */
406*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkCreateIS(networkdm,numkeys,&KeyJunc,NULL,NULL,NULL,&isfrom));
407*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"ISJunctions:\n"));
408*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(isfrom,PETSC_VIEWER_STDOUT_WORLD));
409*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isfrom));
41042dc13f1SHong Zhang 
41142dc13f1SHong Zhang   /* Create a local isfrom for all junction variables */
412*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkCreateLocalIS(networkdm,numkeys,&KeyJunc,NULL,NULL,NULL,&isfrom));
41342dc13f1SHong Zhang   if (!rank) {
414*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[%d] ISLocalJunctions:\n",rank));
415*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISView(isfrom,PETSC_VIEWER_STDOUT_SELF));
41642dc13f1SHong Zhang   }
417*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isfrom));
41842dc13f1SHong Zhang   PetscFunctionReturn(0);
41942dc13f1SHong Zhang }
42042dc13f1SHong Zhang 
42142dc13f1SHong Zhang PetscErrorCode WashNetworkCleanUp(Wash wash)
42242dc13f1SHong Zhang {
42342dc13f1SHong Zhang   PetscMPIInt    rank;
42442dc13f1SHong Zhang 
42542dc13f1SHong Zhang   PetscFunctionBegin;
426*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(wash->comm,&rank));
427*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(wash->edgelist));
428*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(wash->vtype));
429dd400576SPatrick Sanan   if (rank == 0) {
430*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree2(wash->junction,wash->pipe));
43142dc13f1SHong Zhang   }
43242dc13f1SHong Zhang   PetscFunctionReturn(0);
43342dc13f1SHong Zhang }
43442dc13f1SHong Zhang 
43542dc13f1SHong Zhang PetscErrorCode WashNetworkCreate(MPI_Comm comm,PetscInt pipesCase,Wash *wash_ptr)
43642dc13f1SHong Zhang {
43742dc13f1SHong Zhang   PetscInt       npipes;
43842dc13f1SHong Zhang   PetscMPIInt    rank;
43942dc13f1SHong Zhang   Wash           wash=NULL;
44042dc13f1SHong Zhang   PetscInt       i,numVertices,numEdges,*vtype;
44142dc13f1SHong Zhang   PetscInt       *edgelist;
44242dc13f1SHong Zhang   Junction       junctions=NULL;
44342dc13f1SHong Zhang   Pipe           pipes=NULL;
44442dc13f1SHong Zhang   PetscBool      washdist=PETSC_TRUE;
44542dc13f1SHong Zhang 
44642dc13f1SHong Zhang   PetscFunctionBegin;
447*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm,&rank));
44842dc13f1SHong Zhang 
449*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(1,&wash));
45042dc13f1SHong Zhang   wash->comm = comm;
45142dc13f1SHong Zhang   *wash_ptr  = wash;
45242dc13f1SHong Zhang   wash->Q0   = 0.477432; /* RESERVOIR */
45342dc13f1SHong Zhang   wash->H0   = 150.0;
45442dc13f1SHong Zhang   wash->HL   = 143.488;  /* VALVE */
45542dc13f1SHong Zhang   wash->QL   = 0.0;
45642dc13f1SHong Zhang   wash->nnodes_loc = 0;
45742dc13f1SHong Zhang 
45842dc13f1SHong Zhang   numVertices = 0;
45942dc13f1SHong Zhang   numEdges    = 0;
46042dc13f1SHong Zhang   edgelist    = NULL;
46142dc13f1SHong Zhang 
46242dc13f1SHong Zhang   /* proc[0] creates a sequential wash and edgelist */
463*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Setup pipesCase %D\n",pipesCase));
46442dc13f1SHong Zhang 
46542dc13f1SHong Zhang   /* Set global number of pipes, edges, and junctions */
46642dc13f1SHong Zhang   /*-------------------------------------------------*/
46742dc13f1SHong Zhang   switch (pipesCase) {
46842dc13f1SHong Zhang   case 0:
46942dc13f1SHong Zhang     /* pipeCase 0: */
47042dc13f1SHong Zhang     /* =================================================
47142dc13f1SHong Zhang     (RESERVOIR) v0 --E0--> v1--E1--> v2 --E2-->v3 (VALVE)
47242dc13f1SHong Zhang     ====================================================  */
47342dc13f1SHong Zhang     npipes = 3;
474*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetInt(NULL,NULL, "-npipes", &npipes, NULL));
47542dc13f1SHong Zhang     wash->nedge   = npipes;
47642dc13f1SHong Zhang     wash->nvertex = npipes + 1;
47742dc13f1SHong Zhang 
47842dc13f1SHong Zhang     /* Set local edges and vertices -- proc[0] sets entire network, then distributes */
47942dc13f1SHong Zhang     numVertices = 0;
48042dc13f1SHong Zhang     numEdges    = 0;
48142dc13f1SHong Zhang     edgelist    = NULL;
482dd400576SPatrick Sanan     if (rank == 0) {
48342dc13f1SHong Zhang       numVertices = wash->nvertex;
48442dc13f1SHong Zhang       numEdges    = wash->nedge;
48542dc13f1SHong Zhang 
486*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscCalloc1(2*numEdges,&edgelist));
48742dc13f1SHong Zhang       for (i=0; i<numEdges; i++) {
48842dc13f1SHong Zhang         edgelist[2*i] = i; edgelist[2*i+1] = i+1;
48942dc13f1SHong Zhang       }
49042dc13f1SHong Zhang 
49142dc13f1SHong Zhang       /* Add network components */
49242dc13f1SHong Zhang       /*------------------------*/
493*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscCalloc2(numVertices,&junctions,numEdges,&pipes));
49442dc13f1SHong Zhang 
49542dc13f1SHong Zhang       /* vertex */
49642dc13f1SHong Zhang       for (i=0; i<numVertices; i++) {
49742dc13f1SHong Zhang         junctions[i].id = i;
49842dc13f1SHong Zhang         junctions[i].type = JUNCTION;
49942dc13f1SHong Zhang       }
50042dc13f1SHong Zhang 
50142dc13f1SHong Zhang       junctions[0].type             = RESERVOIR;
50242dc13f1SHong Zhang       junctions[numVertices-1].type = VALVE;
50342dc13f1SHong Zhang     }
50442dc13f1SHong Zhang     break;
50542dc13f1SHong Zhang   case 1:
50642dc13f1SHong Zhang     /* pipeCase 1: */
50742dc13f1SHong Zhang     /* ==========================
50842dc13f1SHong Zhang                 v2 (VALVE)
50942dc13f1SHong Zhang                 ^
51042dc13f1SHong Zhang                 |
51142dc13f1SHong Zhang                E2
51242dc13f1SHong Zhang                 |
51342dc13f1SHong Zhang     v0 --E0--> v3--E1--> v1
51442dc13f1SHong Zhang   (RESERVOIR)            (RESERVOIR)
51542dc13f1SHong Zhang     =============================  */
51642dc13f1SHong Zhang     npipes = 3;
51742dc13f1SHong Zhang     wash->nedge   = npipes;
51842dc13f1SHong Zhang     wash->nvertex = npipes + 1;
51942dc13f1SHong Zhang 
52042dc13f1SHong Zhang     /* Set local edges and vertices -- proc[0] sets entire network, then distributes */
521dd400576SPatrick Sanan     if (rank == 0) {
52242dc13f1SHong Zhang       numVertices = wash->nvertex;
52342dc13f1SHong Zhang       numEdges    = wash->nedge;
52442dc13f1SHong Zhang 
525*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscCalloc1(2*numEdges,&edgelist));
52642dc13f1SHong Zhang       edgelist[0] = 0; edgelist[1] = 3;  /* edge[0] */
52742dc13f1SHong Zhang       edgelist[2] = 3; edgelist[3] = 1;  /* edge[1] */
52842dc13f1SHong Zhang       edgelist[4] = 3; edgelist[5] = 2;  /* edge[2] */
52942dc13f1SHong Zhang 
53042dc13f1SHong Zhang       /* Add network components */
53142dc13f1SHong Zhang       /*------------------------*/
532*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscCalloc2(numVertices,&junctions,numEdges,&pipes));
53342dc13f1SHong Zhang       /* vertex */
53442dc13f1SHong Zhang       for (i=0; i<numVertices; i++) {
53542dc13f1SHong Zhang         junctions[i].id   = i;
53642dc13f1SHong Zhang         junctions[i].type = JUNCTION;
53742dc13f1SHong Zhang       }
53842dc13f1SHong Zhang 
53942dc13f1SHong Zhang       junctions[0].type = RESERVOIR;
54042dc13f1SHong Zhang       junctions[1].type = VALVE;
54142dc13f1SHong Zhang       junctions[2].type = VALVE;
54242dc13f1SHong Zhang     }
54342dc13f1SHong Zhang     break;
54442dc13f1SHong Zhang   case 2:
54542dc13f1SHong Zhang     /* pipeCase 2: */
54642dc13f1SHong Zhang     /* ==========================
54742dc13f1SHong Zhang     (RESERVOIR)  v2--> E2
54842dc13f1SHong Zhang                        |
54942dc13f1SHong Zhang             v0 --E0--> v3--E1--> v1
55042dc13f1SHong Zhang     (RESERVOIR)               (VALVE)
55142dc13f1SHong Zhang     =============================  */
55242dc13f1SHong Zhang 
55342dc13f1SHong Zhang     /* Set application parameters -- to be used in function evalutions */
55442dc13f1SHong Zhang     npipes = 3;
55542dc13f1SHong Zhang     wash->nedge   = npipes;
55642dc13f1SHong Zhang     wash->nvertex = npipes + 1;
55742dc13f1SHong Zhang 
55842dc13f1SHong Zhang     /* Set local edges and vertices -- proc[0] sets entire network, then distributes */
559dd400576SPatrick Sanan     if (rank == 0) {
56042dc13f1SHong Zhang       numVertices = wash->nvertex;
56142dc13f1SHong Zhang       numEdges    = wash->nedge;
56242dc13f1SHong Zhang 
563*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscCalloc1(2*numEdges,&edgelist));
56442dc13f1SHong Zhang       edgelist[0] = 0; edgelist[1] = 3;  /* edge[0] */
56542dc13f1SHong Zhang       edgelist[2] = 3; edgelist[3] = 1;  /* edge[1] */
56642dc13f1SHong Zhang       edgelist[4] = 2; edgelist[5] = 3;  /* edge[2] */
56742dc13f1SHong Zhang 
56842dc13f1SHong Zhang       /* Add network components */
56942dc13f1SHong Zhang       /*------------------------*/
570*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscCalloc2(numVertices,&junctions,numEdges,&pipes));
57142dc13f1SHong Zhang       /* vertex */
57242dc13f1SHong Zhang       for (i=0; i<numVertices; i++) {
57342dc13f1SHong Zhang         junctions[i].id = i;
57442dc13f1SHong Zhang         junctions[i].type = JUNCTION;
57542dc13f1SHong Zhang       }
57642dc13f1SHong Zhang 
57742dc13f1SHong Zhang       junctions[0].type = RESERVOIR;
57842dc13f1SHong Zhang       junctions[1].type = VALVE;
57942dc13f1SHong Zhang       junctions[2].type = RESERVOIR;
58042dc13f1SHong Zhang     }
58142dc13f1SHong Zhang     break;
58242dc13f1SHong Zhang   default:
58342dc13f1SHong Zhang     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"not done yet");
58442dc13f1SHong Zhang   }
58542dc13f1SHong Zhang 
58642dc13f1SHong Zhang   /* set edge global id */
58742dc13f1SHong Zhang   for (i=0; i<numEdges; i++) pipes[i].id = i;
58842dc13f1SHong Zhang 
589dd400576SPatrick Sanan   if (rank == 0) { /* set vtype for proc[0] */
59042dc13f1SHong Zhang     PetscInt v;
591*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(2*numEdges,&vtype));
59242dc13f1SHong Zhang     for (i=0; i<2*numEdges; i++) {
59342dc13f1SHong Zhang       v        = edgelist[i];
59442dc13f1SHong Zhang       vtype[i] = junctions[v].type;
59542dc13f1SHong Zhang     }
59642dc13f1SHong Zhang     wash->vtype = vtype;
59742dc13f1SHong Zhang   }
59842dc13f1SHong Zhang 
59942dc13f1SHong Zhang   *wash_ptr      = wash;
60042dc13f1SHong Zhang   wash->nedge    = numEdges;
60142dc13f1SHong Zhang   wash->nvertex  = numVertices;
60242dc13f1SHong Zhang   wash->edgelist = edgelist;
60342dc13f1SHong Zhang   wash->junction = junctions;
60442dc13f1SHong Zhang   wash->pipe     = pipes;
60542dc13f1SHong Zhang 
60642dc13f1SHong Zhang   /* Distribute edgelist to other processors */
607*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-wash_distribute",&washdist,NULL));
60842dc13f1SHong Zhang   if (washdist) {
60942dc13f1SHong Zhang     /*
610*5f80ce2aSJacob Faibussowitsch      CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," Distribute sequential wash ...\n"));
61142dc13f1SHong Zhang      */
612*5f80ce2aSJacob Faibussowitsch     CHKERRQ(WashNetworkDistribute(comm,wash));
61342dc13f1SHong Zhang   }
61442dc13f1SHong Zhang   PetscFunctionReturn(0);
61542dc13f1SHong Zhang }
61642dc13f1SHong Zhang 
61742dc13f1SHong Zhang /* ------------------------------------------------------- */
61842dc13f1SHong Zhang int main(int argc,char ** argv)
61942dc13f1SHong Zhang {
62042dc13f1SHong Zhang   PetscErrorCode    ierr;
62142dc13f1SHong Zhang   Wash              wash;
62242dc13f1SHong Zhang   Junction          junctions,junction;
62342dc13f1SHong Zhang   Pipe              pipe,pipes;
624f11a936eSBarry Smith   PetscInt          KeyPipe,KeyJunction,*edgelist = NULL,*vtype = NULL;
625f11a936eSBarry Smith   PetscInt          i,e,v,eStart,eEnd,vStart,vEnd,key,vkey,type;
626f11a936eSBarry Smith   PetscInt          steps=1,nedges,nnodes=6;
62742dc13f1SHong Zhang   const PetscInt    *cone;
62842dc13f1SHong Zhang   DM                networkdm;
62942dc13f1SHong Zhang   PetscMPIInt       size,rank;
63042dc13f1SHong Zhang   PetscReal         ftime;
63142dc13f1SHong Zhang   Vec               X;
63242dc13f1SHong Zhang   TS                ts;
63342dc13f1SHong Zhang   TSConvergedReason reason;
63442dc13f1SHong Zhang   PetscBool         viewpipes,viewjuncs,monipipes=PETSC_FALSE,userJac=PETSC_TRUE,viewdm=PETSC_FALSE,viewX=PETSC_FALSE;
63542dc13f1SHong Zhang   PetscInt          pipesCase=0;
63642dc13f1SHong Zhang   DMNetworkMonitor  monitor;
63742dc13f1SHong Zhang   MPI_Comm          comm;
63842dc13f1SHong Zhang 
63942dc13f1SHong Zhang   ierr = PetscInitialize(&argc,&argv,"pOption",help);if (ierr) return ierr;
64042dc13f1SHong Zhang 
64142dc13f1SHong Zhang   /* Read runtime options */
642*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL, "-case", &pipesCase, NULL));
643*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-user_Jac",&userJac,NULL));
644*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-pipe_monitor",&monipipes,NULL));
645*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-viewdm",&viewdm,NULL));
646*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-viewX",&viewX,NULL));
647*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL, "-npipenodes", &nnodes, NULL));
64842dc13f1SHong Zhang 
64942dc13f1SHong Zhang   /* Create networkdm */
65042dc13f1SHong Zhang   /*------------------*/
651*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkCreate(PETSC_COMM_WORLD,&networkdm));
652*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject)networkdm,&comm));
653*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm,&rank));
654*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm,&size));
65542dc13f1SHong Zhang 
65642dc13f1SHong Zhang   if (size == 1 && monipipes) {
657*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkMonitorCreate(networkdm,&monitor));
65842dc13f1SHong Zhang   }
65942dc13f1SHong Zhang 
66042dc13f1SHong Zhang   /* Register the components in the network */
661*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkRegisterComponent(networkdm,"junctionstruct",sizeof(struct _p_Junction),&KeyJunction));
662*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkRegisterComponent(networkdm,"pipestruct",sizeof(struct _p_Pipe),&KeyPipe));
66342dc13f1SHong Zhang 
66442dc13f1SHong Zhang   /* Create a distributed wash network (user-specific) */
665*5f80ce2aSJacob Faibussowitsch   CHKERRQ(WashNetworkCreate(comm,pipesCase,&wash));
66642dc13f1SHong Zhang   nedges      = wash->nedge;
66742dc13f1SHong Zhang   edgelist    = wash->edgelist;
66842dc13f1SHong Zhang   vtype       = wash->vtype;
66942dc13f1SHong Zhang   junctions   = wash->junction;
67042dc13f1SHong Zhang   pipes       = wash->pipe;
67142dc13f1SHong Zhang 
67242dc13f1SHong Zhang   /* Set up the network layout */
673*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkSetNumSubNetworks(networkdm,PETSC_DECIDE,1));
674*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkAddSubnetwork(networkdm,NULL,nedges,edgelist,NULL));
67542dc13f1SHong Zhang 
676*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkLayoutSetUp(networkdm));
67742dc13f1SHong Zhang 
678*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd));
679*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkGetVertexRange(networkdm,&vStart,&vEnd));
680*5f80ce2aSJacob Faibussowitsch   /* CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[%d] eStart/End: %d - %d; vStart/End: %d - %d\n",rank,eStart,eEnd,vStart,vEnd)); */
68142dc13f1SHong Zhang 
68242dc13f1SHong Zhang   if (rank) { /* junctions[] and pipes[] for proc[0] are allocated in WashNetworkCreate() */
68342dc13f1SHong Zhang     /* vEnd - vStart = nvertices + number of ghost vertices! */
684*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc2(vEnd - vStart,&junctions,nedges,&pipes));
68542dc13f1SHong Zhang   }
68642dc13f1SHong Zhang 
68742dc13f1SHong Zhang   /* Add Pipe component and number of variables to all local edges */
68842dc13f1SHong Zhang   for (e = eStart; e < eEnd; e++) {
68942dc13f1SHong Zhang     pipes[e-eStart].nnodes = nnodes;
690*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkAddComponent(networkdm,e,KeyPipe,&pipes[e-eStart],2*pipes[e-eStart].nnodes));
69142dc13f1SHong Zhang 
69242dc13f1SHong Zhang     if (size == 1 && monipipes) { /* Add monitor -- show Q_{pipes[e-eStart].id}? */
69342dc13f1SHong Zhang       pipes[e-eStart].length = 600.0;
694*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMNetworkMonitorAdd(monitor, "Pipe Q", e, pipes[e-eStart].nnodes, 0, 2, 0.0,pipes[e-eStart].length, -0.8, 0.8, PETSC_TRUE));
695*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMNetworkMonitorAdd(monitor, "Pipe H", e, pipes[e-eStart].nnodes, 1, 2, 0.0,pipes[e-eStart].length, -400.0, 800.0, PETSC_TRUE));
69642dc13f1SHong Zhang     }
69742dc13f1SHong Zhang   }
69842dc13f1SHong Zhang 
699eac198afSGetnet   /* Add Junction component and number of variables to all local vertices */
70042dc13f1SHong Zhang   for (v = vStart; v < vEnd; v++) {
701*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkAddComponent(networkdm,v,KeyJunction,&junctions[v-vStart],2));
70242dc13f1SHong Zhang   }
70342dc13f1SHong Zhang 
70442dc13f1SHong Zhang   if (size > 1) {  /* must be called before DMSetUp()???. Other partitioners do not work yet??? -- cause crash in proc[0]! */
70542dc13f1SHong Zhang     DM               plexdm;
70642dc13f1SHong Zhang     PetscPartitioner part;
707*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetPlex(networkdm,&plexdm));
708*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetPartitioner(plexdm, &part));
709*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPartitionerSetType(part,PETSCPARTITIONERSIMPLE));
710*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsSetValue(NULL,"-dm_plex_csr_alg","mat")); /* for parmetis */
71142dc13f1SHong Zhang   }
71242dc13f1SHong Zhang 
71342dc13f1SHong Zhang   /* Set up DM for use */
714*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(networkdm));
71542dc13f1SHong Zhang   if (viewdm) {
716*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nOriginal networkdm, DMView:\n"));
717*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMView(networkdm,PETSC_VIEWER_STDOUT_WORLD));
71842dc13f1SHong Zhang   }
71942dc13f1SHong Zhang 
72042dc13f1SHong Zhang   /* Set user physical parameters to the components */
72142dc13f1SHong Zhang   for (e = eStart; e < eEnd; e++) {
722*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetConnectedVertices(networkdm,e,&cone));
72342dc13f1SHong Zhang     /* vfrom */
724*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetComponent(networkdm,cone[0],0,&vkey,(void**)&junction,NULL));
72542dc13f1SHong Zhang     junction->type = (VertexType)vtype[2*e];
72642dc13f1SHong Zhang 
72742dc13f1SHong Zhang     /* vto */
728*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetComponent(networkdm,cone[1],0,&vkey,(void**)&junction,NULL));
72942dc13f1SHong Zhang     junction->type = (VertexType)vtype[2*e+1];
73042dc13f1SHong Zhang   }
73142dc13f1SHong Zhang 
732*5f80ce2aSJacob Faibussowitsch   CHKERRQ(WashNetworkCleanUp(wash));
73342dc13f1SHong Zhang 
73442dc13f1SHong Zhang   /* Network partitioning and distribution of data */
735*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkDistribute(&networkdm,0));
73642dc13f1SHong Zhang   if (viewdm) {
73742dc13f1SHong Zhang     PetscPrintf(PETSC_COMM_WORLD,"\nAfter DMNetworkDistribute, DMView:\n");CHKERRQ(ierr);
738*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMView(networkdm,PETSC_VIEWER_STDOUT_WORLD));
73942dc13f1SHong Zhang   }
74042dc13f1SHong Zhang 
74142dc13f1SHong Zhang   /* create vectors */
742*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(networkdm,&X));
743*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateLocalVector(networkdm,&wash->localX));
744*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateLocalVector(networkdm,&wash->localXdot));
74542dc13f1SHong Zhang 
74642dc13f1SHong Zhang   /* PipeSetUp -- each process only sets its own pipes */
74742dc13f1SHong Zhang   /*---------------------------------------------------*/
748*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkGetVertexRange(networkdm,&vStart,&vEnd));
74942dc13f1SHong Zhang 
75042dc13f1SHong Zhang   userJac = PETSC_TRUE;
751*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkHasJacobian(networkdm,userJac,userJac));
752*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd));
75342dc13f1SHong Zhang   for (e=eStart; e<eEnd; e++) { /* each edge has only one component, pipe */
754*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetComponent(networkdm,e,0,&type,(void**)&pipe,NULL));
75542dc13f1SHong Zhang 
75642dc13f1SHong Zhang     wash->nnodes_loc += pipe->nnodes; /* local total number of nodes, will be used by PipesView() */
75742dc13f1SHong Zhang     ierr = PipeSetParameters(pipe,
75842dc13f1SHong Zhang                              600.0,          /* length */
75942dc13f1SHong Zhang                              0.5,            /* diameter */
76042dc13f1SHong Zhang                              1200.0,         /* a */
76142dc13f1SHong Zhang                              0.018);CHKERRQ(ierr);    /* friction */
762*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PipeSetUp(pipe));
76342dc13f1SHong Zhang 
76442dc13f1SHong Zhang     if (userJac) {
76542dc13f1SHong Zhang       /* Create Jacobian matrix structures for a Pipe */
76642dc13f1SHong Zhang       Mat            *J;
767*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PipeCreateJacobian(pipe,NULL,&J));
768*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMNetworkEdgeSetMatrix(networkdm,e,J));
76942dc13f1SHong Zhang     }
77042dc13f1SHong Zhang   }
77142dc13f1SHong Zhang 
77242dc13f1SHong Zhang   if (userJac) {
773*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetVertexRange(networkdm,&vStart,&vEnd));
77442dc13f1SHong Zhang     for (v=vStart; v<vEnd; v++) {
77542dc13f1SHong Zhang       Mat            *J;
776*5f80ce2aSJacob Faibussowitsch       CHKERRQ(JunctionCreateJacobian(networkdm,v,NULL,&J));
777*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMNetworkVertexSetMatrix(networkdm,v,J));
77842dc13f1SHong Zhang 
779*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMNetworkGetComponent(networkdm,v,0,&vkey,(void**)&junction,NULL));
78042dc13f1SHong Zhang       junction->jacobian = J;
78142dc13f1SHong Zhang     }
78242dc13f1SHong Zhang   }
78342dc13f1SHong Zhang 
78442dc13f1SHong Zhang   /* Setup solver                                           */
78542dc13f1SHong Zhang   /*--------------------------------------------------------*/
786*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
78742dc13f1SHong Zhang 
788*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetDM(ts,(DM)networkdm));
789*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIFunction(ts,NULL,WASHIFunction,wash));
79042dc13f1SHong Zhang 
791*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxSteps(ts,steps));
792*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
793*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,0.1));
794*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSBEULER));
79542dc13f1SHong Zhang   if (size == 1 && monipipes) {
796*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSMonitorSet(ts, TSDMNetworkMonitor, monitor, NULL));
79742dc13f1SHong Zhang   }
798*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
79942dc13f1SHong Zhang 
800*5f80ce2aSJacob Faibussowitsch   CHKERRQ(WASHSetInitialSolution(networkdm,X,wash));
80142dc13f1SHong Zhang 
802*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,X));
80342dc13f1SHong Zhang 
804*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolveTime(ts,&ftime));
805*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&steps));
806*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetConvergedReason(ts,&reason));
807*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps));
80842dc13f1SHong Zhang   if (viewX) {
809*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(X,PETSC_VIEWER_STDOUT_WORLD));
81042dc13f1SHong Zhang   }
81142dc13f1SHong Zhang 
81242dc13f1SHong Zhang   viewpipes = PETSC_FALSE;
813*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL, "-Jac_view", &viewpipes,NULL));
81442dc13f1SHong Zhang   if (viewpipes) {
81542dc13f1SHong Zhang     SNES snes;
81642dc13f1SHong Zhang     Mat  Jac;
817*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetSNES(ts,&snes));
818*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESGetJacobian(snes,&Jac,NULL,NULL,NULL));
819*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(Jac,PETSC_VIEWER_DRAW_WORLD));
82042dc13f1SHong Zhang   }
82142dc13f1SHong Zhang 
82242dc13f1SHong Zhang   /* View solutions */
82342dc13f1SHong Zhang   /* -------------- */
82442dc13f1SHong Zhang   viewpipes = PETSC_FALSE;
825*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL, "-pipe_view", &viewpipes,NULL));
82642dc13f1SHong Zhang   if (viewpipes) {
827*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PipesView(networkdm,KeyPipe,X));
82842dc13f1SHong Zhang   }
82942dc13f1SHong Zhang 
83042dc13f1SHong Zhang   /* Test IS */
83142dc13f1SHong Zhang   viewjuncs = PETSC_FALSE;
832*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL, "-isjunc_view", &viewjuncs,NULL));
83342dc13f1SHong Zhang   if (viewjuncs) {
834*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISJunctionsView(networkdm,KeyJunction));
83542dc13f1SHong Zhang   }
83642dc13f1SHong Zhang 
83742dc13f1SHong Zhang   /* Free spaces */
83842dc13f1SHong Zhang   /* ----------- */
839*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
840*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X));
841*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&wash->localX));
842*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&wash->localXdot));
84342dc13f1SHong Zhang 
84442dc13f1SHong Zhang   /* Destroy objects from each pipe that are created in PipeSetUp() */
845*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMNetworkGetEdgeRange(networkdm,&eStart, &eEnd));
84642dc13f1SHong Zhang   for (i = eStart; i < eEnd; i++) {
847*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkGetComponent(networkdm,i,0,&key,(void**)&pipe,NULL));
848*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PipeDestroy(&pipe));
84942dc13f1SHong Zhang   }
85042dc13f1SHong Zhang   if (userJac) {
85142dc13f1SHong Zhang     for (v=vStart; v<vEnd; v++) {
852*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMNetworkGetComponent(networkdm,v,0,&vkey,(void**)&junction,NULL));
853*5f80ce2aSJacob Faibussowitsch       CHKERRQ(JunctionDestroyJacobian(networkdm,v,junction));
85442dc13f1SHong Zhang     }
85542dc13f1SHong Zhang   }
85642dc13f1SHong Zhang 
85742dc13f1SHong Zhang   if (size == 1 && monipipes) {
858*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMNetworkMonitorDestroy(&monitor));
85942dc13f1SHong Zhang   }
860*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&networkdm));
861*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(wash));
86242dc13f1SHong Zhang 
86342dc13f1SHong Zhang   if (rank) {
864*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree2(junctions,pipes));
86542dc13f1SHong Zhang   }
86642dc13f1SHong Zhang   ierr = PetscFinalize();
86742dc13f1SHong Zhang   return ierr;
86842dc13f1SHong Zhang }
86942dc13f1SHong Zhang 
87042dc13f1SHong Zhang /*TEST
87142dc13f1SHong Zhang 
87242dc13f1SHong Zhang    build:
87342dc13f1SHong Zhang      depends: pipeInterface.c pipeImpls.c
87442dc13f1SHong Zhang      requires: mumps
87542dc13f1SHong Zhang 
87642dc13f1SHong Zhang    test:
87742dc13f1SHong Zhang       args: -ts_monitor -case 1 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -options_left no -viewX
87842dc13f1SHong Zhang       localrunfiles: pOption
87942dc13f1SHong Zhang       output_file: output/pipes_1.out
88042dc13f1SHong Zhang 
88142dc13f1SHong Zhang    test:
88242dc13f1SHong Zhang       suffix: 2
88342dc13f1SHong Zhang       nsize: 2
88442dc13f1SHong Zhang       args: -ts_monitor -case 1 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -petscpartitioner_type simple -options_left no -viewX
88542dc13f1SHong Zhang       localrunfiles: pOption
88642dc13f1SHong Zhang       output_file: output/pipes_2.out
88742dc13f1SHong Zhang 
88842dc13f1SHong Zhang    test:
88942dc13f1SHong Zhang       suffix: 3
89042dc13f1SHong Zhang       nsize: 2
89142dc13f1SHong Zhang       args: -ts_monitor -case 0 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -petscpartitioner_type simple -options_left no -viewX
89242dc13f1SHong Zhang       localrunfiles: pOption
89342dc13f1SHong Zhang       output_file: output/pipes_3.out
89442dc13f1SHong Zhang 
89542dc13f1SHong Zhang    test:
89642dc13f1SHong Zhang       suffix: 4
89742dc13f1SHong Zhang       args: -ts_monitor -case 2 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -options_left no -viewX
89842dc13f1SHong Zhang       localrunfiles: pOption
89942dc13f1SHong Zhang       output_file: output/pipes_4.out
90042dc13f1SHong Zhang 
90142dc13f1SHong Zhang    test:
90242dc13f1SHong Zhang       suffix: 5
90342dc13f1SHong Zhang       nsize: 3
90442dc13f1SHong Zhang       args: -ts_monitor -case 2 -ts_max_steps 10 -pc_factor_mat_solver_type mumps -petscpartitioner_type simple -options_left no -viewX
90542dc13f1SHong Zhang       localrunfiles: pOption
90642dc13f1SHong Zhang       output_file: output/pipes_5.out
90742dc13f1SHong Zhang 
90842dc13f1SHong Zhang    test:
90942dc13f1SHong Zhang       suffix: 6
91042dc13f1SHong Zhang       nsize: 2
91142dc13f1SHong Zhang       args: -ts_monitor -case 1 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -petscpartitioner_type simple -options_left no -wash_distribute 0 -viewX
91242dc13f1SHong Zhang       localrunfiles: pOption
91342dc13f1SHong Zhang       output_file: output/pipes_6.out
91442dc13f1SHong Zhang 
91542dc13f1SHong Zhang    test:
91642dc13f1SHong Zhang       suffix: 7
91742dc13f1SHong Zhang       nsize: 2
91842dc13f1SHong Zhang       args: -ts_monitor -case 2 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -petscpartitioner_type simple -options_left no -wash_distribute 0 -viewX
91942dc13f1SHong Zhang       localrunfiles: pOption
92042dc13f1SHong Zhang       output_file: output/pipes_7.out
92142dc13f1SHong Zhang 
92242dc13f1SHong Zhang    test:
92342dc13f1SHong Zhang       suffix: 8
92442dc13f1SHong Zhang       nsize: 2
92542dc13f1SHong Zhang       requires: parmetis
92642dc13f1SHong Zhang       args: -ts_monitor -case 2 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -petscpartitioner_type parmetis -options_left no -wash_distribute 1
92742dc13f1SHong Zhang       localrunfiles: pOption
92842dc13f1SHong Zhang       output_file: output/pipes_8.out
92942dc13f1SHong Zhang 
93042dc13f1SHong Zhang    test:
93142dc13f1SHong Zhang       suffix: 9
93242dc13f1SHong Zhang       nsize: 2
93342dc13f1SHong Zhang       args: -case 0 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -petscpartitioner_type simple -options_left no -wash_distribute 0 -pipe_view
93442dc13f1SHong Zhang       localrunfiles: pOption
93542dc13f1SHong Zhang       output_file: output/pipes_9.out
93642dc13f1SHong Zhang 
93742dc13f1SHong Zhang    test:
93842dc13f1SHong Zhang       suffix: 10
93942dc13f1SHong Zhang       nsize: 2
94042dc13f1SHong Zhang       args: -case 0 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -petscpartitioner_type simple -options_left no -wash_distribute 0 -isjunc_view
94142dc13f1SHong Zhang       localrunfiles: pOption
94242dc13f1SHong Zhang       output_file: output/pipes_10.out
94342dc13f1SHong Zhang 
94442dc13f1SHong Zhang TEST*/
945