1 static char help[] = "This example demonstrates DMNetwork. It is used for testing parallel generation of dmnetwork, then redistribute. \n\\n";
2 /*
3   Example: mpiexec -n <np> ./pipes1 -ts_max_steps 10
4 */
5 
6 #include "wash.h"
7 #include <petscdmplex.h>
8 
9 /*
10   WashNetworkDistribute - proc[0] distributes sequential wash object
11    Input Parameters:
12 .  comm - MPI communicator
13 .  wash - wash context with all network data in proc[0]
14 
15    Output Parameter:
16 .  wash - wash context with nedge, nvertex and edgelist distributed
17 
18    Note: The routine is used for testing parallel generation of dmnetwork, then redistribute.
19 */
WashNetworkDistribute(MPI_Comm comm,Wash wash)20 PetscErrorCode WashNetworkDistribute(MPI_Comm comm,Wash wash)
21 {
22   PetscErrorCode ierr;
23   PetscMPIInt    rank,size,tag=0;
24   PetscInt       i,e,v,numEdges,numVertices,nedges,*eowners=NULL,estart,eend,*vtype=NULL,nvertices;
25   PetscInt       *edgelist = wash->edgelist,*nvtx=NULL,*vtxDone=NULL;
26 
27   PetscFunctionBegin;
28   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
29   if (size == 1) PetscFunctionReturn(0);
30 
31   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
32   numEdges    = wash->nedge;
33   numVertices = wash->nvertex;
34 
35   /* (1) all processes get global and local number of edges */
36   ierr = MPI_Bcast(&numEdges,1,MPIU_INT,0,comm);CHKERRQ(ierr);
37   nedges = numEdges/size; /* local nedges */
38   if (!rank) {
39     nedges += numEdges - size*(numEdges/size);
40   }
41   wash->Nedge = numEdges;
42   wash->nedge = nedges;
43   /* ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] nedges %d, numEdges %d\n",rank,nedges,numEdges);CHKERRQ(ierr); */
44 
45   ierr = PetscCalloc3(size+1,&eowners,size,&nvtx,numVertices,&vtxDone);CHKERRQ(ierr);
46   ierr = MPI_Allgather(&nedges,1,MPIU_INT,eowners+1,1,MPIU_INT,PETSC_COMM_WORLD);CHKERRQ(ierr);
47   eowners[0] = 0;
48   for (i=2; i<=size; i++) {
49     eowners[i] += eowners[i-1];
50   }
51 
52   estart = eowners[rank];
53   eend   = eowners[rank+1];
54   /* ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] own lists row %d - %d\n",rank,estart,eend);CHKERRQ(ierr); */
55 
56   /* (2) distribute row block edgelist to all processors */
57   if (!rank) {
58     vtype = wash->vtype;
59     for (i=1; i<size; i++) {
60       /* proc[0] sends edgelist to proc[i] */
61       ierr = MPI_Send(edgelist+2*eowners[i],2*(eowners[i+1]-eowners[i]),MPIU_INT,i,tag,comm);CHKERRQ(ierr);
62 
63       /* proc[0] sends vtype to proc[i] */
64       ierr = MPI_Send(vtype+2*eowners[i],2*(eowners[i+1]-eowners[i]),MPIU_INT,i,tag,comm);CHKERRQ(ierr);
65     }
66   } else {
67     MPI_Status      status;
68     ierr = PetscMalloc1(2*(eend-estart),&vtype);CHKERRQ(ierr);
69     ierr = PetscMalloc1(2*(eend-estart),&edgelist);CHKERRQ(ierr);
70 
71     ierr = MPI_Recv(edgelist,2*(eend-estart),MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr);
72     ierr = MPI_Recv(vtype,2*(eend-estart),MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr);
73   }
74 
75   wash->edgelist = edgelist;
76 
77   /* (3) all processes get global and local number of vertices, without ghost vertices */
78   if (!rank) {
79     for (i=0; i<size; i++) {
80       for (e=eowners[i]; e<eowners[i+1]; e++) {
81         v = edgelist[2*e];
82         if (!vtxDone[v]) {
83           nvtx[i]++; vtxDone[v] = 1;
84         }
85         v = edgelist[2*e+1];
86         if (!vtxDone[v]) {
87           nvtx[i]++; vtxDone[v] = 1;
88         }
89       }
90     }
91   }
92   ierr = MPI_Bcast(&numVertices,1,MPIU_INT,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
93   ierr = MPI_Scatter(nvtx,1,MPIU_INT,&nvertices,1,MPIU_INT,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
94   ierr = PetscFree3(eowners,nvtx,vtxDone);CHKERRQ(ierr);
95 
96   wash->Nvertex = numVertices;
97   wash->nvertex = nvertices;
98   wash->vtype   = vtype;
99   PetscFunctionReturn(0);
100 }
101 
WASHIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void * ctx)102 PetscErrorCode WASHIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void* ctx)
103 {
104   PetscErrorCode ierr;
105   Wash           wash=(Wash)ctx;
106   DM             networkdm;
107   Vec            localX,localXdot,localF, localXold;
108   const PetscInt *cone;
109   PetscInt       vfrom,vto,offsetfrom,offsetto,varoffset;
110   PetscInt       v,vStart,vEnd,e,eStart,eEnd;
111   PetscInt       nend,type;
112   PetscBool      ghost;
113   PetscScalar    *farr,*juncf, *pipef;
114   PetscReal      dt;
115   Pipe           pipe;
116   PipeField      *pipex,*pipexdot,*juncx;
117   Junction       junction;
118   DMDALocalInfo  info;
119   const PetscScalar *xarr,*xdotarr, *xoldarr;
120 
121   PetscFunctionBegin;
122   localX    = wash->localX;
123   localXdot = wash->localXdot;
124 
125   ierr = TSGetSolution(ts,&localXold);CHKERRQ(ierr);
126   ierr = TSGetDM(ts,&networkdm);CHKERRQ(ierr);
127   ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
128   ierr = DMGetLocalVector(networkdm,&localF);CHKERRQ(ierr);
129 
130   /* Set F and localF as zero */
131   ierr = VecSet(F,0.0);CHKERRQ(ierr);
132   ierr = VecSet(localF,0.0);CHKERRQ(ierr);
133 
134   /* Update ghost values of locaX and locaXdot */
135   ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
136   ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
137 
138   ierr = DMGlobalToLocalBegin(networkdm,Xdot,INSERT_VALUES,localXdot);CHKERRQ(ierr);
139   ierr = DMGlobalToLocalEnd(networkdm,Xdot,INSERT_VALUES,localXdot);CHKERRQ(ierr);
140 
141   ierr = VecGetArrayRead(localX,&xarr);CHKERRQ(ierr);
142   ierr = VecGetArrayRead(localXdot,&xdotarr);CHKERRQ(ierr);
143   ierr = VecGetArrayRead(localXold,&xoldarr);CHKERRQ(ierr);
144   ierr = VecGetArray(localF,&farr);CHKERRQ(ierr);
145 
146    /* junction->type == JUNCTION:
147            juncf[0] = -qJ + sum(qin); juncf[1] = qJ - sum(qout)
148        junction->type == RESERVOIR (upper stream):
149            juncf[0] = -hJ + H0; juncf[1] = qJ - sum(qout)
150        junction->type == VALVE (down stream):
151            juncf[0] =  -qJ + sum(qin); juncf[1] = qJ
152   */
153   /* Vertex/junction initialization */
154   ierr = DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);CHKERRQ(ierr);
155   for (v=vStart; v<vEnd; v++) {
156     ierr = DMNetworkIsGhostVertex(networkdm,v,&ghost);CHKERRQ(ierr);
157     if (ghost) continue;
158 
159     ierr = DMNetworkGetComponent(networkdm,v,0,&type,(void**)&junction);CHKERRQ(ierr);
160     ierr = DMNetworkGetVariableOffset(networkdm,v,&varoffset);CHKERRQ(ierr);
161     juncx      = (PipeField*)(xarr + varoffset);
162     juncf      = (PetscScalar*)(farr + varoffset);
163 
164     juncf[0] = -juncx[0].q;
165     juncf[1] =  juncx[0].q;
166 
167     if (junction->type == RESERVOIR) { /* upstream reservoir */
168       juncf[0] = juncx[0].h - wash->H0;
169     }
170   }
171 
172   /* Edge/pipe */
173   ierr = DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);CHKERRQ(ierr);
174   for (e=eStart; e<eEnd; e++) {
175     ierr = DMNetworkGetComponent(networkdm,e,0,&type,(void**)&pipe);CHKERRQ(ierr);
176     ierr = DMNetworkGetVariableOffset(networkdm,e,&varoffset);CHKERRQ(ierr);
177     pipex    = (PipeField*)(xarr + varoffset);
178     pipexdot = (PipeField*)(xdotarr + varoffset);
179     pipef    = (PetscScalar*)(farr + varoffset);
180 
181     /* Get some data into the pipe structure: note, some of these operations
182      * might be redundant. Will it consume too much time? */
183     pipe->dt   = dt;
184     pipe->xold = (PipeField*)(xoldarr + varoffset);
185 
186     /* Evaluate F over this edge/pipe: pipef[1], ...,pipef[2*nend] */
187     ierr = DMDAGetLocalInfo(pipe->da,&info);CHKERRQ(ierr);
188     ierr = PipeIFunctionLocal_Lax(&info,t,pipex,pipexdot,pipef,pipe);CHKERRQ(ierr);
189 
190     /* Get boundary values from connected vertices */
191     ierr = DMNetworkGetConnectedVertices(networkdm,e,&cone);CHKERRQ(ierr);
192     vfrom = cone[0]; /* local ordering */
193     vto   = cone[1];
194     ierr = DMNetworkGetVariableOffset(networkdm,vfrom,&offsetfrom);CHKERRQ(ierr);
195     ierr = DMNetworkGetVariableOffset(networkdm,vto,&offsetto);CHKERRQ(ierr);
196 
197     /* Evaluate upstream boundary */
198     ierr = DMNetworkGetComponent(networkdm,vfrom,0,&type,(void**)&junction);CHKERRQ(ierr);
199     if (junction->type != JUNCTION && junction->type != RESERVOIR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"junction type is not supported");
200     juncx = (PipeField*)(xarr + offsetfrom);
201     juncf = (PetscScalar*)(farr + offsetfrom);
202 
203     pipef[0] = pipex[0].h - juncx[0].h;
204     juncf[1] -= pipex[0].q;
205 
206     /* Evaluate downstream boundary */
207     ierr = DMNetworkGetComponent(networkdm,vto,0,&type,(void**)&junction);CHKERRQ(ierr);
208     if (junction->type != JUNCTION && junction->type != VALVE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"junction type is not supported");
209     juncx = (PipeField*)(xarr + offsetto);
210     juncf = (PetscScalar*)(farr + offsetto);
211     nend  = pipe->nnodes - 1;
212 
213     pipef[2*nend + 1] = pipex[nend].h - juncx[0].h;
214     juncf[0] += pipex[nend].q;
215   }
216 
217   ierr = VecRestoreArrayRead(localX,&xarr);CHKERRQ(ierr);
218   ierr = VecRestoreArrayRead(localXdot,&xdotarr);CHKERRQ(ierr);
219   ierr = VecRestoreArray(localF,&farr);CHKERRQ(ierr);
220 
221   ierr = DMLocalToGlobalBegin(networkdm,localF,ADD_VALUES,F);CHKERRQ(ierr);
222   ierr = DMLocalToGlobalEnd(networkdm,localF,ADD_VALUES,F);CHKERRQ(ierr);
223   ierr = DMRestoreLocalVector(networkdm,&localF);CHKERRQ(ierr);
224   /*
225    ierr = PetscPrintf(PETSC_COMM_WORLD("F:\n");CHKERRQ(ierr);
226    ierr = VecView(F,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
227    */
228   PetscFunctionReturn(0);
229 }
230 
WASHSetInitialSolution(DM networkdm,Vec X,Wash wash)231 PetscErrorCode WASHSetInitialSolution(DM networkdm,Vec X,Wash wash)
232 {
233   PetscErrorCode ierr;
234   PetscInt       k,nx,vkey,vfrom,vto,offsetfrom,offsetto;
235   PetscInt       type,varoffset;
236   PetscInt       e,eStart,eEnd;
237   Vec            localX;
238   PetscScalar    *xarr;
239   Pipe           pipe;
240   Junction       junction;
241   const PetscInt *cone;
242   const PetscScalar *xarray;
243 
244   PetscFunctionBegin;
245   ierr = VecSet(X,0.0);CHKERRQ(ierr);
246   ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr);
247   ierr = VecGetArray(localX,&xarr);CHKERRQ(ierr);
248 
249   /* Edge */
250   ierr = DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);CHKERRQ(ierr);
251   for (e=eStart; e<eEnd; e++) {
252     ierr = DMNetworkGetVariableOffset(networkdm,e,&varoffset);CHKERRQ(ierr);
253     ierr = DMNetworkGetComponent(networkdm,e,0,&type,(void**)&pipe);CHKERRQ(ierr);
254 
255     /* set initial values for this pipe */
256     ierr = PipeComputeSteadyState(pipe,wash->Q0,wash->H0);CHKERRQ(ierr);
257     ierr = VecGetSize(pipe->x,&nx);CHKERRQ(ierr);
258 
259     ierr = VecGetArrayRead(pipe->x,&xarray);CHKERRQ(ierr);
260     /* copy pipe->x to xarray */
261     for (k=0; k<nx; k++) {
262       (xarr+varoffset)[k] = xarray[k];
263     }
264 
265     /* set boundary values into vfrom and vto */
266     ierr = DMNetworkGetConnectedVertices(networkdm,e,&cone);CHKERRQ(ierr);
267     vfrom = cone[0]; /* local ordering */
268     vto   = cone[1];
269     ierr = DMNetworkGetVariableOffset(networkdm,vfrom,&offsetfrom);CHKERRQ(ierr);
270     ierr = DMNetworkGetVariableOffset(networkdm,vto,&offsetto);CHKERRQ(ierr);
271 
272     /* if vform is a head vertex: */
273     ierr = DMNetworkGetComponent(networkdm,vfrom,0,&vkey,(void**)&junction);CHKERRQ(ierr);
274     if (junction->type == RESERVOIR) {
275       (xarr+offsetfrom)[1] = wash->H0; /* 1st H */
276     }
277 
278     /* if vto is an end vertex: */
279     ierr = DMNetworkGetComponent(networkdm,vto,0,&vkey,(void**)&junction);CHKERRQ(ierr);
280     if (junction->type == VALVE) {
281       (xarr+offsetto)[0] = wash->QL; /* last Q */
282     }
283     ierr = VecRestoreArrayRead(pipe->x,&xarray);CHKERRQ(ierr);
284   }
285 
286   ierr = VecRestoreArray(localX,&xarr);CHKERRQ(ierr);
287   ierr = DMLocalToGlobalBegin(networkdm,localX,ADD_VALUES,X);CHKERRQ(ierr);
288   ierr = DMLocalToGlobalEnd(networkdm,localX,ADD_VALUES,X);CHKERRQ(ierr);
289   ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr);
290 
291 #if 0
292   PetscInt N;
293   ierr = VecGetSize(X,&N);CHKERRQ(ierr);
294   ierr = PetscPrintf(PETSC_COMM_WORLD,"initial solution %d:\n",N);CHKERRQ(ierr);
295   ierr = VecView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
296 #endif
297   PetscFunctionReturn(0);
298 }
299 
TSDMNetworkMonitor(TS ts,PetscInt step,PetscReal t,Vec x,void * context)300 PetscErrorCode TSDMNetworkMonitor(TS ts, PetscInt step, PetscReal t, Vec x, void *context)
301 {
302   PetscErrorCode     ierr;
303   DMNetworkMonitor   monitor;
304 
305   PetscFunctionBegin;
306   monitor = (DMNetworkMonitor)context;
307   ierr = DMNetworkMonitorView(monitor,x);CHKERRQ(ierr);
308   PetscFunctionReturn(0);
309 }
310 
PipesView(Vec X,DM networkdm,Wash wash)311 PetscErrorCode PipesView(Vec X,DM networkdm,Wash wash)
312 {
313   PetscErrorCode       ierr;
314   Pipe                 pipe;
315   PetscInt             key,Start,End;
316   PetscMPIInt          rank;
317   PetscInt             nx,nnodes,nidx,*idx1,*idx2,*idx1_h,*idx2_h,idx_start,i,k,k1,xstart,j1;
318   Vec                  Xq,Xh,localX;
319   IS                   is1_q,is2_q,is1_h,is2_h;
320   VecScatter           ctx_q,ctx_h;
321 
322   PetscFunctionBegin;
323   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
324 
325   /* get num of local and global total nnodes */
326   nidx = wash->nnodes_loc;
327   ierr = MPIU_Allreduce(&nidx,&nx,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD);CHKERRQ(ierr);
328 
329   ierr = VecCreate(PETSC_COMM_WORLD,&Xq);CHKERRQ(ierr);
330   if (rank == 0) { /* all entries of Xq are in proc[0] */
331     ierr = VecSetSizes(Xq,nx,PETSC_DECIDE);CHKERRQ(ierr);
332   } else {
333     ierr = VecSetSizes(Xq,0,PETSC_DECIDE);CHKERRQ(ierr);
334   }
335   ierr = VecSetFromOptions(Xq);CHKERRQ(ierr);
336   ierr = VecSet(Xq,0.0);CHKERRQ(ierr);
337   ierr = VecDuplicate(Xq,&Xh);CHKERRQ(ierr);
338 
339   ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr);
340 
341   /* set idx1 and idx2 */
342   ierr = PetscCalloc4(nidx,&idx1,nidx,&idx2,nidx,&idx1_h,nidx,&idx2_h);CHKERRQ(ierr);
343 
344   ierr = DMNetworkGetEdgeRange(networkdm,&Start, &End);CHKERRQ(ierr);
345 
346   ierr = VecGetOwnershipRange(X,&xstart,NULL);CHKERRQ(ierr);
347   k1 = 0;
348   j1 = 0;
349   for (i = Start; i < End; i++) {
350     ierr = DMNetworkGetComponent(networkdm,i,0,&key,(void**)&pipe);CHKERRQ(ierr);
351     nnodes = pipe->nnodes;
352     idx_start = pipe->id*nnodes;
353     for (k=0; k<nnodes; k++) {
354       idx1[k1] = xstart + j1*2*nnodes + 2*k;
355       idx2[k1] = idx_start + k;
356 
357       idx1_h[k1] = xstart + j1*2*nnodes + 2*k + 1;
358       idx2_h[k1] = idx_start + k;
359       k1++;
360     }
361     j1++;
362   }
363 
364   ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx1,PETSC_COPY_VALUES,&is1_q);CHKERRQ(ierr);
365   ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx2,PETSC_COPY_VALUES,&is2_q);CHKERRQ(ierr);
366   ierr = VecScatterCreate(X,is1_q,Xq,is2_q,&ctx_q);CHKERRQ(ierr);
367   ierr = VecScatterBegin(ctx_q,X,Xq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
368   ierr = VecScatterEnd(ctx_q,X,Xq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
369 
370   ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx1_h,PETSC_COPY_VALUES,&is1_h);CHKERRQ(ierr);
371   ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx2_h,PETSC_COPY_VALUES,&is2_h);CHKERRQ(ierr);
372   ierr = VecScatterCreate(X,is1_h,Xh,is2_h,&ctx_h);CHKERRQ(ierr);
373   ierr = VecScatterBegin(ctx_h,X,Xh,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
374   ierr = VecScatterEnd(ctx_h,X,Xh,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
375 
376   ierr = PetscPrintf(PETSC_COMM_WORLD,"Xq: \n");CHKERRQ(ierr);
377   ierr = VecView(Xq,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
378   ierr = PetscPrintf(PETSC_COMM_WORLD,"Xh: \n");CHKERRQ(ierr);
379   ierr = VecView(Xh,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
380 
381   ierr = VecScatterDestroy(&ctx_q);CHKERRQ(ierr);
382   ierr = PetscFree4(idx1,idx2,idx1_h,idx2_h);CHKERRQ(ierr);
383   ierr = ISDestroy(&is1_q);CHKERRQ(ierr);
384   ierr = ISDestroy(&is2_q);CHKERRQ(ierr);
385 
386   ierr = VecScatterDestroy(&ctx_h);CHKERRQ(ierr);
387   ierr = ISDestroy(&is1_h);CHKERRQ(ierr);
388   ierr = ISDestroy(&is2_h);CHKERRQ(ierr);
389 
390   ierr = VecDestroy(&Xq);CHKERRQ(ierr);
391   ierr = VecDestroy(&Xh);CHKERRQ(ierr);
392   ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr);
393   PetscFunctionReturn(0);
394 }
395 
WashNetworkCleanUp(Wash wash)396 PetscErrorCode WashNetworkCleanUp(Wash wash)
397 {
398   PetscErrorCode ierr;
399   PetscMPIInt    rank;
400 
401   PetscFunctionBegin;
402   ierr = MPI_Comm_rank(wash->comm,&rank);CHKERRQ(ierr);
403   ierr = PetscFree(wash->edgelist);CHKERRQ(ierr);
404   ierr = PetscFree(wash->vtype);CHKERRQ(ierr);
405   if (!rank) {
406     ierr = PetscFree2(wash->junction,wash->pipe);CHKERRQ(ierr);
407   }
408   PetscFunctionReturn(0);
409 }
410 
WashNetworkCreate(MPI_Comm comm,PetscInt pipesCase,Wash * wash_ptr)411 PetscErrorCode WashNetworkCreate(MPI_Comm comm,PetscInt pipesCase,Wash *wash_ptr)
412 {
413   PetscErrorCode ierr;
414   PetscInt       npipes;
415   PetscMPIInt    rank;
416   Wash           wash=NULL;
417   PetscInt       i,numVertices,numEdges,*vtype;
418   PetscInt       *edgelist;
419   Junction       junctions=NULL;
420   Pipe           pipes=NULL;
421   PetscBool      washdist=PETSC_TRUE;
422 
423   PetscFunctionBegin;
424   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
425 
426   ierr = PetscCalloc1(1,&wash);CHKERRQ(ierr);
427   wash->comm = comm;
428   *wash_ptr  = wash;
429   wash->Q0   = 0.477432; /* RESERVOIR */
430   wash->H0   = 150.0;
431   wash->HL   = 143.488;  /* VALVE */
432   wash->QL   = 0.0;
433   wash->nnodes_loc = 0;
434 
435   numVertices = 0;
436   numEdges    = 0;
437   edgelist    = NULL;
438 
439   /* proc[0] creates a sequential wash and edgelist */
440   ierr = PetscPrintf(PETSC_COMM_WORLD,"Setup pipesCase %D\n",pipesCase);CHKERRQ(ierr);
441 
442   /* Set global number of pipes, edges, and junctions */
443   /*-------------------------------------------------*/
444   switch (pipesCase) {
445   case 0:
446     /* pipeCase 0: */
447     /* =================================================
448     (RESERVOIR) v0 --E0--> v1--E1--> v2 --E2-->v3 (VALVE)
449     ====================================================  */
450     npipes = 3;
451     ierr = PetscOptionsGetInt(NULL,NULL, "-npipes", &npipes, NULL);CHKERRQ(ierr);
452     wash->nedge   = npipes;
453     wash->nvertex = npipes + 1;
454 
455     /* Set local edges and vertices -- proc[0] sets entire network, then distributes */
456     numVertices = 0;
457     numEdges    = 0;
458     edgelist    = NULL;
459     if (!rank) {
460       numVertices = wash->nvertex;
461       numEdges    = wash->nedge;
462 
463       ierr = PetscCalloc1(2*numEdges,&edgelist);CHKERRQ(ierr);
464       for (i=0; i<numEdges; i++) {
465         edgelist[2*i] = i; edgelist[2*i+1] = i+1;
466       }
467 
468       /* Add network components */
469       /*------------------------*/
470       ierr = PetscCalloc2(numVertices,&junctions,numEdges,&pipes);CHKERRQ(ierr);
471 
472       /* vertex */
473       for (i=0; i<numVertices; i++) {
474         junctions[i].id = i;
475         junctions[i].type = JUNCTION;
476       }
477 
478       junctions[0].type             = RESERVOIR;
479       junctions[numVertices-1].type = VALVE;
480     }
481     break;
482   case 1:
483     /* pipeCase 1: */
484     /* ==========================
485                 v2 (VALVE)
486                 ^
487                 |
488                E2
489                 |
490     v0 --E0--> v3--E1--> v1
491   (RESERVOIR)            (RESERVOIR)
492     =============================  */
493     npipes = 3;
494     wash->nedge   = npipes;
495     wash->nvertex = npipes + 1;
496 
497     /* Set local edges and vertices -- proc[0] sets entire network, then distributes */
498     if (!rank) {
499       numVertices = wash->nvertex;
500       numEdges    = wash->nedge;
501 
502       ierr = PetscCalloc1(2*numEdges,&edgelist);CHKERRQ(ierr);
503       edgelist[0] = 0; edgelist[1] = 3;  /* edge[0] */
504       edgelist[2] = 3; edgelist[3] = 1;  /* edge[1] */
505       edgelist[4] = 3; edgelist[5] = 2;  /* edge[2] */
506 
507       /* Add network components */
508       /*------------------------*/
509       ierr = PetscCalloc2(numVertices,&junctions,numEdges,&pipes);CHKERRQ(ierr);
510       /* vertex */
511       for (i=0; i<numVertices; i++) {
512         junctions[i].id   = i;
513         junctions[i].type = JUNCTION;
514       }
515 
516       junctions[0].type = RESERVOIR;
517       junctions[1].type = VALVE;
518       junctions[2].type = VALVE;
519     }
520     break;
521   case 2:
522     /* pipeCase 2: */
523     /* ==========================
524     (RESERVOIR)  v2--> E2
525                        |
526             v0 --E0--> v3--E1--> v1
527     (RESERVOIR)               (VALVE)
528     =============================  */
529 
530     /* Set application parameters -- to be used in function evalutions */
531     npipes = 3;
532     wash->nedge   = npipes;
533     wash->nvertex = npipes + 1;
534 
535     /* Set local edges and vertices -- proc[0] sets entire network, then distributes */
536     if (!rank) {
537       numVertices = wash->nvertex;
538       numEdges    = wash->nedge;
539 
540       ierr = PetscCalloc1(2*numEdges,&edgelist);CHKERRQ(ierr);
541       edgelist[0] = 0; edgelist[1] = 3;  /* edge[0] */
542       edgelist[2] = 3; edgelist[3] = 1;  /* edge[1] */
543       edgelist[4] = 2; edgelist[5] = 3;  /* edge[2] */
544 
545       /* Add network components */
546       /*------------------------*/
547       ierr = PetscCalloc2(numVertices,&junctions,numEdges,&pipes);CHKERRQ(ierr);
548       /* vertex */
549       for (i=0; i<numVertices; i++) {
550         junctions[i].id = i;
551         junctions[i].type = JUNCTION;
552       }
553 
554       junctions[0].type = RESERVOIR;
555       junctions[1].type = VALVE;
556       junctions[2].type = RESERVOIR;
557     }
558     break;
559   default:
560     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"not done yet");
561   }
562 
563   /* set edge global id */
564   for (i=0; i<numEdges; i++) pipes[i].id = i;
565 
566   if (!rank) { /* set vtype for proc[0] */
567     PetscInt v;
568     ierr = PetscMalloc1(2*numEdges,&vtype);CHKERRQ(ierr);
569     for (i=0; i<2*numEdges; i++) {
570       v        = edgelist[i];
571       vtype[i] = junctions[v].type;
572     }
573     wash->vtype = vtype;
574   }
575 
576   *wash_ptr      = wash;
577   wash->nedge    = numEdges;
578   wash->nvertex  = numVertices;
579   wash->edgelist = edgelist;
580   wash->junction = junctions;
581   wash->pipe     = pipes;
582 
583   /* Distribute edgelist to other processors */
584   ierr = PetscOptionsGetBool(NULL,NULL,"-wash_distribute",&washdist,NULL);CHKERRQ(ierr);
585   if (washdist) {
586     /*
587      ierr = PetscPrintf(PETSC_COMM_WORLD," Distribute sequential wash ...\n");CHKERRQ(ierr);
588      */
589     ierr = WashNetworkDistribute(comm,wash);CHKERRQ(ierr);
590   }
591   PetscFunctionReturn(0);
592 }
593 
594 /* ------------------------------------------------------- */
main(int argc,char ** argv)595 int main(int argc,char ** argv)
596 {
597   PetscErrorCode    ierr;
598   Wash              wash;
599   Junction          junctions,junction;
600   Pipe              pipe,pipes;
601   PetscInt          KeyPipe,KeyJunction;
602   PetscInt          *edgelist = NULL,*edgelists[1],*vtype = NULL;
603   PetscInt          i,e,v,eStart,eEnd,vStart,vEnd,key;
604   PetscInt          vkey,type;
605   const PetscInt    *cone;
606   DM                networkdm;
607   PetscMPIInt       size,rank;
608   PetscReal         ftime;
609   Vec               X;
610   TS                ts;
611   PetscInt          steps=1;
612   TSConvergedReason reason;
613   PetscBool         viewpipes,monipipes=PETSC_FALSE,userJac=PETSC_TRUE,viewdm=PETSC_FALSE,viewX=PETSC_FALSE;
614   PetscBool         test=PETSC_FALSE;
615   PetscInt          pipesCase=0;
616   DMNetworkMonitor  monitor;
617   MPI_Comm          comm;
618 
619   PetscInt          nedges,nvertices; /* local num of edges and vertices */
620   PetscInt          nnodes = 6,nv,ne;
621   const PetscInt    *vtx,*edge;
622 
623   ierr = PetscInitialize(&argc,&argv,"pOption",help);if (ierr) return ierr;
624 
625   /* Read runtime options */
626   ierr = PetscOptionsGetInt(NULL,NULL, "-case", &pipesCase, NULL);CHKERRQ(ierr);
627   ierr = PetscOptionsGetBool(NULL,NULL,"-user_Jac",&userJac,NULL);CHKERRQ(ierr);
628   ierr = PetscOptionsGetBool(NULL,NULL,"-pipe_monitor",&monipipes,NULL);CHKERRQ(ierr);
629   ierr = PetscOptionsGetBool(NULL,NULL,"-viewdm",&viewdm,NULL);CHKERRQ(ierr);
630   ierr = PetscOptionsGetBool(NULL,NULL,"-viewX",&viewX,NULL);CHKERRQ(ierr);
631   ierr = PetscOptionsGetInt(NULL,NULL, "-npipenodes", &nnodes, NULL);CHKERRQ(ierr);
632   ierr = PetscOptionsGetBool(NULL,NULL,"-test",&test,NULL);CHKERRQ(ierr);
633 
634   /* Create networkdm */
635   /*------------------*/
636   ierr = DMNetworkCreate(PETSC_COMM_WORLD,&networkdm);CHKERRQ(ierr);
637   ierr = PetscObjectGetComm((PetscObject)networkdm,&comm);CHKERRQ(ierr);
638   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
639   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
640 
641   if (size == 1 && monipipes) {
642     ierr = DMNetworkMonitorCreate(networkdm,&monitor);CHKERRQ(ierr);
643   }
644 
645   /* Register the components in the network */
646   ierr = DMNetworkRegisterComponent(networkdm,"junctionstruct",sizeof(struct _p_Junction),&KeyJunction);CHKERRQ(ierr);
647   ierr = DMNetworkRegisterComponent(networkdm,"pipestruct",sizeof(struct _p_Pipe),&KeyPipe);CHKERRQ(ierr);
648 
649   /* Create a distributed wash network (user-specific) */
650   ierr = WashNetworkCreate(comm,pipesCase,&wash);CHKERRQ(ierr);
651   nedges      = wash->nedge;
652   nvertices   = wash->nvertex; /* local num of vertices, excluding ghosts */
653   edgelist    = wash->edgelist;
654   vtype       = wash->vtype;
655   junctions   = wash->junction;
656   pipes       = wash->pipe;
657 
658   /* Set up the network layout */
659   ierr = DMNetworkSetSizes(networkdm,1,&nvertices,&nedges,0,NULL);CHKERRQ(ierr);
660 
661   /* Add local edge connectivity */
662   edgelists[0] = edgelist;
663   ierr = DMNetworkSetEdgeList(networkdm,edgelists,NULL);CHKERRQ(ierr);
664   ierr = DMNetworkLayoutSetUp(networkdm);CHKERRQ(ierr);
665 
666   ierr = DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);CHKERRQ(ierr);
667   ierr = DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);CHKERRQ(ierr);
668   /* ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] eStart/End: %d - %d; vStart/End: %d - %d\n",rank,eStart,eEnd,vStart,vEnd);CHKERRQ(ierr); */
669 
670   /* Test DMNetworkGetSubnetworkInfo() */
671   if (test) {
672     ierr = DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edge);CHKERRQ(ierr);
673     if (ne != eEnd - eStart || nv != vEnd - vStart) SETERRQ2(PetscObjectComm((PetscObject)networkdm),PETSC_ERR_ARG_WRONG,"ne %D or nv %D is incorrect",ne,nv);
674   }
675 
676   if (rank) { /* junctions[] and pipes[] for proc[0] are allocated in WashNetworkCreate() */
677     /* vEnd - vStart = nvertices + num of ghost vertices! */
678     ierr = PetscCalloc2(vEnd - vStart,&junctions,nedges,&pipes);CHKERRQ(ierr);
679   }
680 
681   /* Add Pipe component to all local edges */
682   for (e = eStart; e < eEnd; e++) {
683     if (test) {
684       if (e != edge[e]) SETERRQ2(PetscObjectComm((PetscObject)networkdm),PETSC_ERR_ARG_WRONG,"e %D != edge %D from DMNetworkGetSubnetworkInfo()",e,edge[e]);
685     }
686 
687     pipes[e-eStart].nnodes = nnodes;
688     ierr = DMNetworkAddComponent(networkdm,e,KeyPipe,&pipes[e-eStart]);CHKERRQ(ierr);
689 
690     /* Add number of variables to each edge */
691     ierr = DMNetworkAddNumVariables(networkdm,e,2*pipes[e-eStart].nnodes);CHKERRQ(ierr);
692 
693     if (size == 1 && monipipes) { /* Add monitor -- show Q_{pipes[e-eStart].id}? */
694       pipes[e-eStart].length = 600.0;
695       ierr = DMNetworkMonitorAdd(monitor, "Pipe Q", e, pipes[e-eStart].nnodes, 0, 2, 0.0,pipes[e-eStart].length, -0.8, 0.8, PETSC_TRUE);CHKERRQ(ierr);
696       ierr = DMNetworkMonitorAdd(monitor, "Pipe H", e, pipes[e-eStart].nnodes, 1, 2, 0.0,pipes[e-eStart].length, -400.0, 800.0, PETSC_TRUE);CHKERRQ(ierr);
697     }
698   }
699 
700   /* Add Junction component to all local vertices, including ghost vertices! */
701   for (v = vStart; v < vEnd; v++) {
702     if (test) {
703       if (v != vtx[v-vStart]) SETERRQ2(PetscObjectComm((PetscObject)networkdm),PETSC_ERR_ARG_WRONG,"v %D != vtx %D from DMNetworkGetSubnetworkInfo()",v,vtx[v-vStart]);
704     }
705 
706     ierr = DMNetworkAddComponent(networkdm,v,KeyJunction,&junctions[v-vStart]);CHKERRQ(ierr);
707 
708     /* Add number of variables to vertex */
709     ierr = DMNetworkAddNumVariables(networkdm,v,2);CHKERRQ(ierr);
710   }
711 
712   if (size > 1) {  /* must be called before DMSetUp()???. Other partitioners do not work yet??? -- cause crash in proc[0]! */
713     DM               plexdm;
714     PetscPartitioner part;
715     ierr = DMNetworkGetPlex(networkdm,&plexdm);CHKERRQ(ierr);
716     ierr = DMPlexGetPartitioner(plexdm, &part);CHKERRQ(ierr);
717     ierr = PetscPartitionerSetType(part,PETSCPARTITIONERSIMPLE);CHKERRQ(ierr);
718     ierr = PetscOptionsSetValue(NULL,"-dm_plex_csr_via_mat","true");CHKERRQ(ierr); /* for parmetis */
719   }
720 
721   /* Set up DM for use */
722   ierr = DMSetUp(networkdm);CHKERRQ(ierr);
723   if (viewdm) {
724     ierr = PetscPrintf(PETSC_COMM_WORLD,"\nOriginal networkdm, DMView:\n");CHKERRQ(ierr);
725     ierr = DMView(networkdm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
726   }
727 
728   /* Set user physical parameters to the components */
729   for (e = eStart; e < eEnd; e++) {
730     ierr = DMNetworkGetConnectedVertices(networkdm,e,&cone);CHKERRQ(ierr);
731     /* vfrom */
732     ierr = DMNetworkGetComponent(networkdm,cone[0],0,&vkey,(void**)&junction);CHKERRQ(ierr);
733     junction->type = (VertexType)vtype[2*e];
734 
735     /* vto */
736     ierr = DMNetworkGetComponent(networkdm,cone[1],0,&vkey,(void**)&junction);CHKERRQ(ierr);
737     junction->type = (VertexType)vtype[2*e+1];
738   }
739 
740   ierr = WashNetworkCleanUp(wash);CHKERRQ(ierr);
741 
742   /* Network partitioning and distribution of data */
743   ierr = DMNetworkDistribute(&networkdm,0);CHKERRQ(ierr);
744   if (viewdm) {
745     PetscPrintf(PETSC_COMM_WORLD,"\nAfter DMNetworkDistribute, DMView:\n");CHKERRQ(ierr);
746     ierr = DMView(networkdm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
747   }
748 
749   /* create vectors */
750   ierr = DMCreateGlobalVector(networkdm,&X);CHKERRQ(ierr);
751   ierr = DMCreateLocalVector(networkdm,&wash->localX);CHKERRQ(ierr);
752   ierr = DMCreateLocalVector(networkdm,&wash->localXdot);CHKERRQ(ierr);
753 
754   /* PipeSetUp -- each process only sets its own pipes */
755   /*---------------------------------------------------*/
756   ierr = DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);CHKERRQ(ierr);
757 
758   userJac = PETSC_TRUE;
759   ierr = DMNetworkHasJacobian(networkdm,userJac,userJac);CHKERRQ(ierr);
760   ierr = DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);CHKERRQ(ierr);
761   for (e=eStart; e<eEnd; e++) { /* each edge has only one component, pipe */
762     ierr = DMNetworkGetComponent(networkdm,e,0,&type,(void**)&pipe);CHKERRQ(ierr);
763 
764     wash->nnodes_loc += pipe->nnodes; /* local total num of nodes, will be used by PipesView() */
765     ierr = PipeSetParameters(pipe,
766                              600.0,          /* length */
767                              0.5,            /* diameter */
768                              1200.0,         /* a */
769                              0.018);CHKERRQ(ierr);    /* friction */
770     ierr = PipeSetUp(pipe);CHKERRQ(ierr);
771 
772     if (userJac) {
773       /* Create Jacobian matrix structures for a Pipe */
774       Mat            *J;
775       ierr = PipeCreateJacobian(pipe,NULL,&J);CHKERRQ(ierr);
776       ierr = DMNetworkEdgeSetMatrix(networkdm,e,J);CHKERRQ(ierr);
777     }
778   }
779 
780   if (userJac) {
781     ierr = DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);CHKERRQ(ierr);
782     for (v=vStart; v<vEnd; v++) {
783       Mat            *J;
784       ierr = JunctionCreateJacobian(networkdm,v,NULL,&J);CHKERRQ(ierr);
785       ierr = DMNetworkVertexSetMatrix(networkdm,v,J);CHKERRQ(ierr);
786 
787       ierr = DMNetworkGetComponent(networkdm,v,0,&vkey,(void**)&junction);CHKERRQ(ierr);
788       junction->jacobian = J;
789     }
790   }
791 
792   /* Test DMNetworkGetSubnetworkInfo() */
793   if (test) {
794     ierr = DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);CHKERRQ(ierr);
795     ierr = DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edge);CHKERRQ(ierr);
796     if (ne != eEnd - eStart || nv != vEnd - vStart) SETERRQ2(PetscObjectComm((PetscObject)networkdm),PETSC_ERR_ARG_WRONG,"ne %D or nv %D is incorrect",ne,nv);
797 
798     for (e = eStart; e < eEnd; e++) {
799       if (e != edge[e]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"e %D != edge %D from DMNetworkGetSubnetworkInfo()",e,edge[e]);
800     }
801     for (v = vStart; v < vEnd; v++) {
802       if (v != vtx[v-vStart]) SETERRQ2(PetscObjectComm((PetscObject)networkdm),PETSC_ERR_ARG_WRONG,"v %D != vtx %D from DMNetworkGetSubnetworkInfo()",v,vtx[v-vStart]);
803     }
804   }
805 
806   /* Setup solver                                           */
807   /*--------------------------------------------------------*/
808   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
809 
810   ierr = TSSetDM(ts,(DM)networkdm);CHKERRQ(ierr);
811   ierr = TSSetIFunction(ts,NULL,WASHIFunction,wash);CHKERRQ(ierr);
812 
813   ierr = TSSetMaxSteps(ts,steps);CHKERRQ(ierr);
814   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
815   ierr = TSSetTimeStep(ts,0.1);CHKERRQ(ierr);
816   ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr);
817   if (size == 1 && monipipes) {
818     ierr = TSMonitorSet(ts, TSDMNetworkMonitor, monitor, NULL);CHKERRQ(ierr);
819   }
820   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
821 
822   ierr = WASHSetInitialSolution(networkdm,X,wash);CHKERRQ(ierr);
823 
824   ierr = TSSolve(ts,X);CHKERRQ(ierr);
825 
826   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
827   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
828   ierr = TSGetConvergedReason(ts,&reason);CHKERRQ(ierr);
829   ierr = PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps);CHKERRQ(ierr);
830   if (viewX) {
831     ierr = VecView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
832   }
833 
834   viewpipes = PETSC_FALSE;
835   ierr = PetscOptionsGetBool(NULL,NULL, "-Jac_view", &viewpipes,NULL);CHKERRQ(ierr);
836   if (viewpipes) {
837     SNES snes;
838     Mat  Jac;
839     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
840     ierr = SNESGetJacobian(snes,&Jac,NULL,NULL,NULL);CHKERRQ(ierr);
841     ierr = MatView(Jac,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);
842   }
843 
844   /* View solution q and h */
845   /* --------------------- */
846   viewpipes = PETSC_FALSE;
847   ierr = PetscOptionsGetBool(NULL,NULL, "-pipe_view", &viewpipes,NULL);CHKERRQ(ierr);
848   if (viewpipes) {
849     ierr = PipesView(X,networkdm,wash);CHKERRQ(ierr);
850   }
851 
852   /* Free spaces */
853   /* ----------- */
854   ierr = TSDestroy(&ts);CHKERRQ(ierr);
855   ierr = VecDestroy(&X);CHKERRQ(ierr);
856   ierr = VecDestroy(&wash->localX);CHKERRQ(ierr);
857   ierr = VecDestroy(&wash->localXdot);CHKERRQ(ierr);
858 
859   /* Destroy objects from each pipe that are created in PipeSetUp() */
860   ierr = DMNetworkGetEdgeRange(networkdm,&eStart, &eEnd);CHKERRQ(ierr);
861   for (i = eStart; i < eEnd; i++) {
862     ierr = DMNetworkGetComponent(networkdm,i,0,&key,(void**)&pipe);CHKERRQ(ierr);
863     ierr = PipeDestroy(&pipe);CHKERRQ(ierr);
864   }
865   if (userJac) {
866     for (v=vStart; v<vEnd; v++) {
867       ierr = DMNetworkGetComponent(networkdm,v,0,&vkey,(void**)&junction);CHKERRQ(ierr);
868       ierr = JunctionDestroyJacobian(networkdm,v,junction);CHKERRQ(ierr);
869     }
870   }
871 
872   if (size == 1 && monipipes) {
873     ierr = DMNetworkMonitorDestroy(&monitor);CHKERRQ(ierr);
874   }
875   ierr = DMDestroy(&networkdm);CHKERRQ(ierr);
876   ierr = PetscFree(wash);CHKERRQ(ierr);
877 
878   if (rank) {
879     ierr = PetscFree2(junctions,pipes);CHKERRQ(ierr);
880   }
881   ierr = PetscFinalize();
882   return ierr;
883 }
884 
885 /*TEST
886 
887    build:
888      depends: pipeInterface.c pipeImpls.c
889 
890    test:
891       args: -ts_monitor -case 1 -ts_max_steps 1 -options_left no -viewX -test
892       localrunfiles: pOption
893       output_file: output/pipes1_1.out
894 
895    test:
896       suffix: 2
897       nsize: 2
898       requires: mumps
899       args: -ts_monitor -case 1 -ts_max_steps 1 -petscpartitioner_type simple -options_left no -viewX -test
900       localrunfiles: pOption
901       output_file: output/pipes1_2.out
902 
903    test:
904       suffix: 3
905       nsize: 2
906       requires: mumps
907       args: -ts_monitor -case 0 -ts_max_steps 1 -petscpartitioner_type simple -options_left no -viewX -test
908       localrunfiles: pOption
909       output_file: output/pipes1_3.out
910 
911    test:
912       suffix: 4
913       args: -ts_monitor -case 2 -ts_max_steps 1 -options_left no -viewX -test
914       localrunfiles: pOption
915       output_file: output/pipes1_4.out
916 
917    test:
918       suffix: 5
919       nsize: 3
920       requires: mumps
921       args: -ts_monitor -case 2 -ts_max_steps 10 -petscpartitioner_type simple -options_left no -viewX -test
922       localrunfiles: pOption
923       output_file: output/pipes1_5.out
924 
925    test:
926       suffix: 6
927       nsize: 2
928       requires: mumps
929       args: -ts_monitor -case 1 -ts_max_steps 1 -petscpartitioner_type simple -options_left no -wash_distribute 0 -viewX -test
930       localrunfiles: pOption
931       output_file: output/pipes1_6.out
932 
933    test:
934       suffix: 7
935       nsize: 2
936       requires: mumps
937       args: -ts_monitor -case 2 -ts_max_steps 1 -petscpartitioner_type simple -options_left no -wash_distribute 0 -viewX -test
938       localrunfiles: pOption
939       output_file: output/pipes1_7.out
940 
941    test:
942       suffix: 8
943       nsize: 2
944       requires: mumps parmetis
945       args: -ts_monitor -case 2 -ts_max_steps 1 -petscpartitioner_type parmetis -options_left no -wash_distribute 1
946       localrunfiles: pOption
947       output_file: output/pipes1_8.out
948 
949 TEST*/
950