Actual source code: pipes1.c

petsc-3.8.3 2017-12-09
Report Typos and Errors
  1:      static char help[] = "This example demonstrates the use of DMNetwork \n\\n";

  3: /*
  4:   Example: mpiexec -n <np> ./pipes1 -ts_max_steps 10
  5: */

  7: #include "wash.h"
  8:  #include <petscdmnetwork.h>

 10: PetscErrorCode WASHIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void* ctx)
 11: {
 13:   Wash           wash=(Wash)ctx;
 14:   DM             networkdm;
 15:   Vec            localX,localXdot,localF;
 16:   const PetscInt *cone;
 17:   PetscInt       vfrom,vto,offsetfrom,offsetto,type,varoffset,voffset;
 18:   PetscInt       v,vStart,vEnd,e,eStart,eEnd,pipeoffset;
 19:   PetscBool      ghost;
 20:   PetscScalar    *farr,*vf,*juncx,*juncf;
 21:   Pipe           pipe;
 22:   PipeField      *pipex,*pipexdot,*pipef;
 23:   DMDALocalInfo  info;
 24:   Junction       junction;
 25:   MPI_Comm       comm;
 26:   PetscMPIInt    rank,size;
 27:   const PetscScalar *xarr,*xdotarr;
 28:   DMNetworkComponentGenericDataType *nwarr;

 31:   PetscObjectGetComm((PetscObject)ts,&comm);
 32:   MPI_Comm_rank(comm,&rank);
 33:   MPI_Comm_size(comm,&size);

 35:   VecSet(F,0.0);

 37:   localX    = wash->localX;
 38:   localXdot = wash->localXdot;

 40:   TSGetDM(ts,&networkdm);
 41:   DMGetLocalVector(networkdm,&localF);
 42:   VecSet(localF,0.0);

 44:   /* update ghost values of locaX and locaXdot */
 45:   DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);
 46:   DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);

 48:   DMGlobalToLocalBegin(networkdm,Xdot,INSERT_VALUES,localXdot);
 49:   DMGlobalToLocalEnd(networkdm,Xdot,INSERT_VALUES,localXdot);

 51:   VecGetArrayRead(localX,&xarr);
 52:   VecGetArrayRead(localXdot,&xdotarr);
 53:   VecGetArray(localF,&farr);

 55:   /* Initialize localF = localX at non-ghost vertices */
 56:   DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);
 57:   for (v=vStart; v<vEnd; v++) {
 58:     DMNetworkIsGhostVertex(networkdm,v,&ghost);
 59:     if (!ghost) {
 60:       DMNetworkGetVariableOffset(networkdm,v,&varoffset);
 61:       juncx  = (PetscScalar*)(xarr+varoffset);
 62:       juncf  = (PetscScalar*)(farr+varoffset);
 63:       juncf[0] = juncx[0];
 64:       juncf[1] = juncx[1];
 65:     }
 66:   }

 68:   /* Get component(application) data array */
 69:   DMNetworkGetComponentDataArray(networkdm,&nwarr);

 71:   /* Edge */
 72:   DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);
 73:   for (e=eStart; e<eEnd; e++) {
 74:     DMNetworkGetComponentKeyOffset(networkdm,e,0,&type,&pipeoffset);
 75:     DMNetworkGetVariableOffset(networkdm,e,&varoffset);
 76:     pipe     = (Pipe)(nwarr + pipeoffset);
 77:     pipex    = (PipeField*)(xarr + varoffset);
 78:     pipexdot = (PipeField*)(xdotarr + varoffset);
 79:     pipef    = (PipeField*)(farr + varoffset);

 81:     /* Get boundary values H0 and QL from connected vertices */
 82:     DMNetworkGetConnectedVertices(networkdm,e,&cone);
 83:     vfrom = cone[0]; /* local ordering */
 84:     vto   = cone[1];
 85:     DMNetworkGetVariableOffset(networkdm,vfrom,&offsetfrom);
 86:     DMNetworkGetVariableOffset(networkdm,vto,&offsetto);
 87:     if (pipe->boundary.Q0 == PIPE_CHARACTERISTIC) {
 88:       pipe->boundary.H0 = (xarr+offsetfrom)[1]; /* h_from */
 89:     } else {
 90:       pipe->boundary.Q0 = (xarr+offsetfrom)[0]; /* q_from */
 91:     }
 92:     if (pipe->boundary.HL == PIPE_CHARACTERISTIC) {
 93:       pipe->boundary.QL = (xarr+offsetto)[0];   /* q_to */
 94:     } else {
 95:       pipe->boundary.HL = (xarr+offsetto)[1];   /* h_to */
 96:     }

 98:     /* Evaluate PipeIFunctionLocal() */
 99:     DMDAGetLocalInfo(pipe->da,&info);
100:     PipeIFunctionLocal(&info, t, pipex, pipexdot, pipef, pipe);

102:     /* Set F at vfrom */
103:     vf = (PetscScalar*)(farr+offsetfrom);
104:     if (pipe->boundary.Q0 == PIPE_CHARACTERISTIC) {
105:       vf[0] -= pipex[0].q; /* q_vfrom - q[0] */
106:     } else {
107:       vf[1] -= pipex[0].h; /* h_vfrom - h[0] */
108:     }

110:     /* Set F at vto */
111:     vf = (PetscScalar*)(farr+offsetto);
112:     if (pipe->boundary.HL == PIPE_CHARACTERISTIC) {
113:       vf[1] -= pipex[pipe->nnodes-1].h; /* h_vto - h[last] */
114:     } else {
115:       vf[0] -= pipex[pipe->nnodes-1].q; /* q_vto - q[last] */
116:     }
117:   }

119:   /* Set F at boundary vertices */
120:   for (v=vStart; v<vEnd; v++) {
121:     DMNetworkGetComponentKeyOffset(networkdm,v,0,&type,&voffset);
122:     DMNetworkGetVariableOffset(networkdm,v,&varoffset);
123:     junction = (Junction)(nwarr + voffset);
124:     juncf = (PetscScalar *)(farr + varoffset);
125:     if (junction->isEnd == -1) {
126:       juncf[1] -= wash->H0;
127:       } else if (junction->isEnd == 1) {
128:       juncf[0] -= wash->QL;
129:     }
130:   }

132:   VecRestoreArrayRead(localX,&xarr);
133:   VecRestoreArrayRead(localXdot,&xdotarr);
134:   VecRestoreArray(localF,&farr);

136:   DMLocalToGlobalBegin(networkdm,localF,ADD_VALUES,F);
137:   DMLocalToGlobalEnd(networkdm,localF,ADD_VALUES,F);
138:   DMRestoreLocalVector(networkdm,&localF);
139:   /* VecView(F,PETSC_VIEWER_STDOUT_WORLD); */
140:   return(0);
141: }

143: PetscErrorCode WASHSetInitialSolution(DM networkdm,Vec X,Wash wash)
144: {
146:   PetscInt       k,nx,vkey,vOffset,vfrom,vto,offsetfrom,offsetto;
147:   PetscInt       type,varoffset;
148:   PetscInt       e,eStart,eEnd,pipeoffset;
149:   Vec            localX;
150:   PetscScalar    *xarr;
151:   Pipe           pipe;
152:   Junction       junction;
153:   const PetscInt *cone;
154:   const PetscScalar *xarray;
155:   DMNetworkComponentGenericDataType *nwarr;

158:   VecSet(X,0.0);
159:   DMGetLocalVector(networkdm,&localX);
160:   VecGetArray(localX,&xarr);

162:    /* Get component(application) data array */
163:   DMNetworkGetComponentDataArray(networkdm,&nwarr);

165:   /* Edge */
166:   DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);
167:   for (e=eStart; e<eEnd; e++) {
168:     DMNetworkGetVariableOffset(networkdm,e,&varoffset);
169:     DMNetworkGetComponentKeyOffset(networkdm,e,0,&type,&pipeoffset);

171:     /* get from and to vertices */
172:     DMNetworkGetConnectedVertices(networkdm,e,&cone);
173:     vfrom = cone[0]; /* local ordering */
174:     vto   = cone[1];

176:     DMNetworkGetVariableOffset(networkdm,vfrom,&offsetfrom);
177:     DMNetworkGetVariableOffset(networkdm,vto,&offsetto);

179:     /* set initial values for this pipe */
180:     /* Q0=0.477432; H0=150.0; needs to be updated by its succeeding pipe. Use SNESSolve()? */
181:     pipe     = (Pipe)(nwarr + pipeoffset);
182:     PipeComputeSteadyState(pipe, 0.477432, wash->H0);
183:     VecGetSize(pipe->x,&nx);

185:     VecGetArrayRead(pipe->x,&xarray);
186:     /* copy pipe->x to xarray */
187:     for (k=0; k<nx; k++) {
188:       (xarr+varoffset)[k] = xarray[k];
189:     }

191:     /* set boundary values into vfrom and vto */
192:     if (pipe->boundary.Q0 == PIPE_CHARACTERISTIC) {
193:       (xarr+offsetfrom)[0] += xarray[0];    /* Q0 -> vfrom[0] */
194:     } else {
195:       (xarr+offsetfrom)[1] += xarray[1];    /* H0 -> vfrom[1] */
196:     }

198:     if (pipe->boundary.HL == PIPE_CHARACTERISTIC) {
199:       (xarr+offsetto)[1]   += xarray[nx-1]; /* HL -> vto[1]   */
200:     } else {
201:       (xarr+offsetto)[0]   += xarray[nx-2]; /* QL -> vto[0]   */
202:     }

204:     /* if vform is a head vertex: */
205:     DMNetworkGetComponentKeyOffset(networkdm,vfrom,0,&vkey,&vOffset);
206:     junction = (Junction)(nwarr+vOffset);
207:     if (junction->isEnd == -1) { /* head junction */
208:       (xarr+offsetfrom)[0] = 0.0;      /* 1st Q -- not used */
209:       (xarr+offsetfrom)[1] = wash->H0; /* 1st H */
210:     }

212:     /* if vto is an end vertex: */
213:     DMNetworkGetComponentKeyOffset(networkdm,vto,0,&vkey,&vOffset);
214:     junction = (Junction)(nwarr+vOffset);
215:     if (junction->isEnd == 1) { /* end junction */
216:       (xarr+offsetto)[0] = wash->QL; /* last Q */
217:       (xarr+offsetto)[1] = 0.0;      /* last H -- not used */
218:     }
219:     VecRestoreArrayRead(pipe->x,&xarray);
220:   }

222:   VecRestoreArray(localX,&xarr);
223:   DMLocalToGlobalBegin(networkdm,localX,ADD_VALUES,X);
224:   DMLocalToGlobalEnd(networkdm,localX,ADD_VALUES,X);
225:   DMRestoreLocalVector(networkdm,&localX);

227: #if 0
228:   PetscInt N;
229:   VecGetSize(X,&N);
230:   printf("initial solution %d:\n",N);
231:   VecView(X,PETSC_VIEWER_STDOUT_WORLD);
232: #endif
233:   return(0);
234: }

236: PetscErrorCode TSDMNetworkMonitor(TS ts, PetscInt step, PetscReal t, Vec x, void *context)
237: {
238:   PetscErrorCode     ierr;
239:   DMNetworkMonitor   monitor;

242:   monitor = (DMNetworkMonitor)context;
243:   DMNetworkMonitorView(monitor,x);
244:   return(0);
245: }

247: PetscErrorCode PipesView(Vec X,DM networkdm,Wash wash)
248: {
249:   PetscErrorCode       ierr;
250:   Pipe                 pipe;
251:   DMNetworkComponentGenericDataType *nwarr;
252:   PetscInt             pipeOffset,key,Start,End;
253:   PetscMPIInt          rank;
254:   PetscInt             nx,nnodes,nidx,*idx1,*idx2,*idx1_h,*idx2_h,idx_start,i,k,k1,xstart,j1;
255:   Vec                  Xq,Xh,localX;
256:   IS                   is1_q,is2_q,is1_h,is2_h;
257:   VecScatter           ctx_q,ctx_h;

260:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

262:   /* get num of local and global total nnodes */
263:   nidx = wash->nnodes_loc;
264:   MPIU_Allreduce(&nidx,&nx,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD);

266:   VecCreate(PETSC_COMM_WORLD,&Xq);
267:   if (rank == 0) { /* all entries of Xq are in proc[0] */
268:     VecSetSizes(Xq,nx,PETSC_DECIDE);
269:   } else {
270:     VecSetSizes(Xq,0,PETSC_DECIDE);
271:   }
272:   VecSetFromOptions(Xq);
273:   VecSet(Xq,0.0);
274:   VecDuplicate(Xq,&Xh);

276:   DMGetLocalVector(networkdm,&localX);

278:   /* set idx1 and idx2 */
279:   PetscCalloc4(nidx,&idx1,nidx,&idx2,nidx,&idx1_h,nidx,&idx2_h);

281:   DMNetworkGetComponentDataArray(networkdm,&nwarr);
282:   DMNetworkGetEdgeRange(networkdm,&Start, &End);

284:   VecGetOwnershipRange(X,&xstart,NULL);
285:   k1 = 0;
286:   j1 = 0;
287:   for (i = Start; i < End; i++) {
288:     DMNetworkGetComponentKeyOffset(networkdm,i,0,&key,&pipeOffset);
289:     pipe = (Pipe)(nwarr+pipeOffset);
290:     nnodes = pipe->nnodes;
291:     idx_start = pipe->id*nnodes;
292:     for (k=0; k<nnodes; k++) {
293:       idx1[k1] = xstart + j1*2*nnodes + 2*k;
294:       idx2[k1] = idx_start + k;

296:       idx1_h[k1] = xstart + j1*2*nnodes + 2*k + 1;
297:       idx2_h[k1] = idx_start + k;
298:       k1++;
299:     }
300:     j1++;
301:   }

303:   ISCreateGeneral(PETSC_COMM_SELF,nidx,idx1,PETSC_COPY_VALUES,&is1_q);
304:   ISCreateGeneral(PETSC_COMM_SELF,nidx,idx2,PETSC_COPY_VALUES,&is2_q);
305:   VecScatterCreate(X,is1_q,Xq,is2_q,&ctx_q);
306:   VecScatterBegin(ctx_q,X,Xq,INSERT_VALUES,SCATTER_FORWARD);
307:   VecScatterEnd(ctx_q,X,Xq,INSERT_VALUES,SCATTER_FORWARD);

309:   ISCreateGeneral(PETSC_COMM_SELF,nidx,idx1_h,PETSC_COPY_VALUES,&is1_h);
310:   ISCreateGeneral(PETSC_COMM_SELF,nidx,idx2_h,PETSC_COPY_VALUES,&is2_h);
311:   VecScatterCreate(X,is1_h,Xh,is2_h,&ctx_h);
312:   VecScatterBegin(ctx_h,X,Xh,INSERT_VALUES,SCATTER_FORWARD);
313:   VecScatterEnd(ctx_h,X,Xh,INSERT_VALUES,SCATTER_FORWARD);

315:   if (!rank) printf("Xq: \n");
316:   VecView(Xq,PETSC_VIEWER_STDOUT_WORLD);
317:   if (!rank) printf("Xh: \n");
318:   VecView(Xh,PETSC_VIEWER_STDOUT_WORLD);

320:   VecScatterDestroy(&ctx_q);
321:   PetscFree4(idx1,idx2,idx1_h,idx2_h);
322:   ISDestroy(&is1_q);
323:   ISDestroy(&is2_q);

325:   VecScatterDestroy(&ctx_h);
326:   ISDestroy(&is1_h);
327:   ISDestroy(&is2_h);

329:   VecDestroy(&Xq);
330:   VecDestroy(&Xh);
331:   DMRestoreLocalVector(networkdm,&localX);
332:   return(0);
333: }

335: PetscErrorCode WashNetworkCleanUp(Wash wash,int *edgelist)
336: {
338:   PetscMPIInt    rank;

341:   MPI_Comm_rank(wash->comm,&rank);
342:   if (!rank) {
343:     PetscFree(edgelist);
344:     PetscFree2(wash->junction,wash->pipe);
345:   }
346:   return(0);
347: }

349: PetscErrorCode WashNetworkCreate(MPI_Comm comm,PetscInt pipesCase,Wash *wash_ptr,int **elist)
350: {
352:   PetscInt       nnodes,npipes;
353:   PetscMPIInt    rank;
354:   Wash           wash;
355:   PetscInt       i,numVertices,numEdges;
356:   int            *edgelist;
357:   Junction       junctions=NULL;
358:   Pipe           pipes=NULL;

361:   MPI_Comm_rank(comm,&rank);

363:   PetscCalloc1(1,&wash);
364:   wash->comm = comm;
365:   *wash_ptr  = wash;
366:   wash->Q0   = 0.477432; /* copied from initial soluiton */
367:   wash->H0   = 150.0;
368:   wash->HL   = 143.488; /* copied from initial soluiton */
369:   wash->nnodes_loc = 0;

371:   numVertices = 0;
372:   numEdges    = 0;
373:   edgelist    = NULL;

375:   if (!rank) {
376:     PetscPrintf(PETSC_COMM_SELF,"Setup pipesCase %D\n",pipesCase);
377:   }
378:   nnodes = 6;
379:   PetscOptionsGetInt(NULL,NULL, "-npipenodes", &nnodes, NULL);

381:   /* Set global number of pipes, edges, and junctions */
382:   /*-------------------------------------------------*/
383:   switch (pipesCase) {
384:   case 0:
385:     /* pipeCase 0: */
386:     /* =============================
387:     v0 --E0--> v1--E1--> v2 --E2-->v3
388:     ================================  */
389:     npipes = 3;
390:     PetscOptionsGetInt(NULL,NULL, "-npipes", &npipes, NULL);
391:     wash->nedge   = npipes;
392:     wash->nvertex = npipes + 1;

394:     /* Set local edges and vertices -- proc[0] sets entire network, then distributes */
395:     numVertices = 0;
396:     numEdges    = 0;
397:     edgelist    = NULL;
398:     if (!rank) {
399:       numVertices = wash->nvertex;
400:       numEdges    = wash->nedge;

402:       PetscCalloc1(2*numEdges,&edgelist);
403:       for (i=0; i<numEdges; i++) {
404:         edgelist[2*i] = i; edgelist[2*i+1] = i+1;
405:       }

407:       /* Add network components */
408:       /*------------------------*/
409:       PetscCalloc2(numVertices,&junctions,numEdges,&pipes);
410:       /* vertex */
411:       for (i=0; i<numVertices; i++) {
412:         junctions[i].id = i;
413:         junctions[i].isEnd = 0;
414:         junctions[i].nedges_in = 1; junctions[i].nedges_out = 1;

416:         /* Set GPS data */
417:         junctions[i].latitude  = 0.0;
418:         junctions[i].longitude = 0.0;
419:       }
420:       junctions[0].isEnd                  = -1;
421:       junctions[0].nedges_in              =  0;
422:       junctions[numVertices-1].isEnd      =  1;
423:       junctions[numVertices-1].nedges_out =  0;

425:       /* edge and pipe */
426:       for (i=0; i<numEdges; i++) {
427:         pipes[i].id   = i;
428:         pipes[i].nnodes = nnodes;
429:       }
430:     }
431:     break;
432:   case 1:
433:     /* pipeCase 1: */
434:     /* ==========================
435:                 v2
436:                 ^
437:                 |
438:                E2
439:                 |
440:     v0 --E0--> v3--E1--> v1
441:     =============================  */
442:     npipes = 3;
443:     wash->nedge   = npipes;
444:     wash->nvertex = npipes + 1;

446:     /* Set local edges and vertices -- proc[0] sets entire network, then distributes */
447:     if (!rank) {
448:       numVertices = wash->nvertex;
449:       numEdges    = wash->nedge;

451:       PetscCalloc1(2*numEdges,&edgelist);
452:       edgelist[0] = 0; edgelist[1] = 3;  /* edge[0] */
453:       edgelist[2] = 3; edgelist[3] = 1;  /* edge[1] */
454:       edgelist[4] = 3; edgelist[5] = 2;  /* edge[2] */

456:       /* Add network components */
457:       /*------------------------*/
458:       PetscCalloc2(numVertices,&junctions,numEdges,&pipes);
459:       /* vertex */
460:       for (i=0; i<numVertices; i++) {
461:         junctions[i].id = i;

463:         /* Set GPS data */
464:         junctions[i].latitude  = 0.0;
465:         junctions[i].longitude = 0.0;
466:       }
467:       junctions[0].isEnd = -1; junctions[0].nedges_in = 0; junctions[0].nedges_out = 1;
468:       junctions[1].isEnd =  1; junctions[1].nedges_in = 1; junctions[1].nedges_out = 0;
469:       junctions[2].isEnd =  1; junctions[2].nedges_in = 1; junctions[2].nedges_out = 0;
470:       junctions[3].isEnd =  0; junctions[3].nedges_in = 1; junctions[3].nedges_out = 2;

472:       /* edge and pipe */
473:       for (i=0; i<numEdges; i++) {
474:         pipes[i].id     = i;
475:         pipes[i].nnodes = nnodes;
476:       }
477:     }
478:     break;
479:   case 2:
480:     /* pipeCase 2: */
481:     /* ==========================
482:          v2--> E2
483:                 |
484:     v0 --E0--> v3--E1--> v1
485:     =============================  */

487:     /* Set application parameters -- to be used in function evalutions */
488:     npipes = 3;
489:     wash->nedge   = npipes;
490:     wash->nvertex = npipes + 1;

492:     /* Set local edges and vertices -- proc[0] sets entire network, then distributes */
493:     if (!rank) {
494:       numVertices = wash->nvertex;
495:       numEdges    = wash->nedge;

497:       PetscCalloc1(2*numEdges,&edgelist);
498:       edgelist[0] = 0; edgelist[1] = 3;  /* edge[0] */
499:       edgelist[2] = 3; edgelist[3] = 1;  /* edge[1] */
500:       edgelist[4] = 2; edgelist[5] = 3;  /* edge[2] */

502:       /* Add network components */
503:       /*------------------------*/
504:       PetscCalloc2(numVertices,&junctions,numEdges,&pipes);
505:       /* vertex */
506:       for (i=0; i<numVertices; i++) {
507:         junctions[i].id = i;

509:         /* Set GPS data */
510:         junctions[i].latitude  = 0.0;
511:         junctions[i].longitude = 0.0;
512:       }
513:       junctions[0].isEnd = -1; junctions[0].nedges_in = 0; junctions[0].nedges_out = 1;
514:       junctions[1].isEnd =  1; junctions[1].nedges_in = 1; junctions[1].nedges_out = 0;
515:       junctions[2].isEnd = -1; junctions[2].nedges_in = 0; junctions[2].nedges_out = 1;
516:       junctions[3].isEnd =  0; junctions[3].nedges_in = 2; junctions[3].nedges_out = 1;

518:       /* edge and pipe */
519:       for (i=0; i<numEdges; i++) {
520:         pipes[i].id     = i;
521:         pipes[i].nnodes = nnodes;
522:       }
523:     }
524:     break;
525:   default:
526:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"not done yet");
527:   }

529:   *wash_ptr      = wash;
530:   wash->nedge    = numEdges;
531:   wash->nvertex  = numVertices;
532:   *elist         = edgelist;
533:   wash->junction = junctions;
534:   wash->pipe     = pipes;
535:   return(0);
536: }
537: /* ------------------------------------------------------- */
538: int main(int argc,char ** argv)
539: {
540:   PetscErrorCode    ierr;
541:   Wash              wash;
542:   Junction          junctions,junction;
543:   Pipe              pipe,pipes;
544:   PetscInt          numEdges,numVertices,KeyPipe,KeyJunction;
545:   int               *edgelist = NULL;
546:   PetscInt          i,e,v,eStart,eEnd,vStart,vEnd,pipeOffset,key,frombType,tobType;
547:   PetscInt          vfrom,vto,vkey,fromOffset,toOffset,type,varoffset,pipeoffset;
548:   PetscInt          from_nedge_in,from_nedge_out,to_nedge_in;
549:   const PetscInt    *cone;
550:   DM                networkdm;
551:   PetscMPIInt       size,rank;
552:   PetscReal         ftime = 2500.0;
553:   Vec               X;
554:   TS                ts;
555:   PetscInt          steps;
556:   TSConvergedReason reason;
557:   PetscBool         viewpipes;
558:   PetscInt          pipesCase;
559:   DMNetworkMonitor  monitor;
560:   DMNetworkComponentGenericDataType *nwarr;

562:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
563:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
564:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

566:   /* Create and setup network */
567:   /*--------------------------*/
568:   DMNetworkCreate(PETSC_COMM_WORLD,&networkdm);
569:   if (size == 1) {
570:     DMNetworkMonitorCreate(networkdm,&monitor);
571:   }
572:   /* Register the components in the network */
573:   DMNetworkRegisterComponent(networkdm,"junctionstruct",sizeof(struct _p_Junction),&KeyJunction);
574:   DMNetworkRegisterComponent(networkdm,"pipestruct",sizeof(struct _p_Pipe),&KeyPipe);

576:   /* Set global number of pipes, edges, and vertices */
577:   pipesCase = 2;
578:   PetscOptionsGetInt(NULL,NULL, "-case", &pipesCase, NULL);

580:   WashNetworkCreate(PETSC_COMM_WORLD,pipesCase,&wash,&edgelist);
581:   numEdges    = wash->nedge;
582:   numVertices = wash->nvertex;
583:   junctions    = wash->junction;
584:   pipes       = wash->pipe;

586:   /* Set number of vertices and edges */
587:   DMNetworkSetSizes(networkdm,numVertices,numEdges,PETSC_DETERMINE,PETSC_DETERMINE);
588:   /* Add edge connectivity */
589:   DMNetworkSetEdgeList(networkdm,edgelist);
590:   /* Set up the network layout */
591:   DMNetworkLayoutSetUp(networkdm);

593:   /* Add EDGEDATA component to all edges -- currently networkdm is a sequential network */
594:   DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);
595:   DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);

597:   for (e = eStart; e < eEnd; e++) {
598:     /* Add Pipe component to all edges -- create pipe here */
599:     DMNetworkAddComponent(networkdm,e,KeyPipe,&pipes[e-eStart]);

601:     /* Add number of variables to each edge */
602:     DMNetworkAddNumVariables(networkdm,e,2*pipes[e-eStart].nnodes);

604:     if (size == 1) { /* Add monitor -- show Q_{pipes[e-eStart].id}? */
605:       DMNetworkMonitorAdd(monitor, "Pipe Q", e, pipes[e-eStart].nnodes, 0, 2, -0.8, 0.8, PETSC_TRUE);
606:       DMNetworkMonitorAdd(monitor, "Pipe H", e, pipes[e-eStart].nnodes, 1, 2, -400.0, 800.0, PETSC_TRUE);
607:     }
608:   }

610:   /* Add Junction component to all vertices */
611:   for (v = vStart; v < vEnd; v++) {
612:     DMNetworkAddComponent(networkdm,v,KeyJunction,&junctions[v-vStart]);

614:     /* Add number of variables to vertex */
615:     DMNetworkAddNumVariables(networkdm,v,2);
616:   }

618:   /* Set up DM for use */
619:   DMSetUp(networkdm);
620:   WashNetworkCleanUp(wash,edgelist);

622:   /* Network partitioning and distribution of data */
623:   DMNetworkDistribute(&networkdm,0);

625:   /* PipeSetUp -- each process only sets its own pipes */
626:   /*---------------------------------------------------*/
627:   DMNetworkGetComponentDataArray(networkdm,&nwarr);

629:   DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);
630:   for (e=eStart; e<eEnd; e++) { /* each edge has only one component, pipe */
631:     DMNetworkGetComponentKeyOffset(networkdm,e,0,&type,&pipeoffset);
632:     DMNetworkGetVariableOffset(networkdm,e,&varoffset);
633:     pipe = (Pipe)(nwarr + pipeoffset);

635:     /* Setup conntected vertices */
636:     DMNetworkGetConnectedVertices(networkdm,e,&cone);
637:     vfrom = cone[0]; /* local ordering */
638:     vto   = cone[1];

640:     /* vfrom */
641:     DMNetworkGetComponentKeyOffset(networkdm,vfrom,0,&vkey,&fromOffset);
642:     junction = (Junction)(nwarr+fromOffset);
643:     from_nedge_in  = junction->nedges_in;
644:     from_nedge_out = junction->nedges_out;

646:     /* vto */
647:     DMNetworkGetComponentKeyOffset(networkdm,vto,0,&vkey,&toOffset);
648:     junction    = (Junction)(nwarr+toOffset);
649:     to_nedge_in = junction->nedges_in;

651:     pipe->comm = PETSC_COMM_SELF; /* must be set here, otherwise crashes in my mac??? */
652:     wash->nnodes_loc += pipe->nnodes; /* local total num of nodes, will be used by PipesView() */
653:     PipeSetParameters(pipe,
654:                              600.0,          /* length */
655:                              pipe->nnodes,   /* nnodes -- rm from PipeSetParameters */
656:                              0.5,            /* diameter */
657:                              1200.0,         /* a */
658:                              0.018);    /* friction */

660:     /* set boundary conditions for this pipe */
661:     if (from_nedge_in <= 1 && from_nedge_out > 0) {
662:       frombType = 0;
663:     } else {
664:       frombType = 1;
665:     }

667:     if (to_nedge_in == 1) {
668:       tobType = 0;
669:     } else {
670:       tobType = 1;
671:     }

673:     if (frombType == 0) {
674:       pipe->boundary.Q0 = PIPE_CHARACTERISTIC; /* will be obtained from characteristic */
675:       pipe->boundary.H0 = wash->H0;
676:     } else {
677:       pipe->boundary.Q0 = wash->Q0;
678:       pipe->boundary.H0 = PIPE_CHARACTERISTIC; /* will be obtained from characteristic */
679:     }
680:     if (tobType == 0) {
681:       pipe->boundary.QL = wash->QL;
682:       pipe->boundary.HL = PIPE_CHARACTERISTIC; /* will be obtained from characteristic */
683:     } else {
684:       pipe->boundary.QL = PIPE_CHARACTERISTIC; /* will be obtained from characteristic */
685:       pipe->boundary.HL = wash->HL;
686:     }

688:     PipeSetUp(pipe);
689:   }

691:   /* create vectors */
692:   DMCreateGlobalVector(networkdm,&X);
693:   DMCreateLocalVector(networkdm,&wash->localX);
694:   DMCreateLocalVector(networkdm,&wash->localXdot);

696:   /* Setup solver                                           */
697:   /*--------------------------------------------------------*/
698:   TSCreate(PETSC_COMM_WORLD,&ts);

700:   TSSetDM(ts,(DM)networkdm);
701:   TSSetIFunction(ts,NULL,WASHIFunction,wash);

703:   TSSetMaxTime(ts,ftime);
704:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
705:   TSSetTimeStep(ts,0.1);
706:   TSSetType(ts,TSBEULER);
707:   if (size == 1) {
708:     TSMonitorSet(ts, TSDMNetworkMonitor, monitor, NULL);
709:   }
710:   TSSetFromOptions(ts);

712:   WASHSetInitialSolution(networkdm,X,wash);

714:   TSSolve(ts,X);

716:   TSGetSolveTime(ts,&ftime);
717:   TSGetStepNumber(ts,&steps);
718:   TSGetConvergedReason(ts,&reason);
719:   PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps);
720:   /* VecView(X,PETSC_VIEWER_STDOUT_WORLD); */

722:   /* View solution q and h */
723:   /* --------------------- */
724:   viewpipes = PETSC_FALSE;
725:   PetscOptionsGetBool(NULL,NULL, "-pipe_view", &viewpipes,NULL);
726:   if (viewpipes) {
727:     PipesView(X,networkdm,wash);
728:   }

730:   /* Free spaces */
731:   /* ----------- */
732:   TSDestroy(&ts);
733:   VecDestroy(&X);
734:   VecDestroy(&wash->localX);
735:   VecDestroy(&wash->localXdot);

737:   /* Destroy objects from each pipe that are created in PipeSetUp() */
738:   DMNetworkGetEdgeRange(networkdm,&eStart, &eEnd);
739:   for (i = eStart; i < eEnd; i++) {
740:     DMNetworkGetComponentKeyOffset(networkdm,i,0,&key,&pipeOffset);
741:     pipe = (Pipe)(nwarr+pipeOffset);
742:     DMDestroy(&(pipe->da));
743:     VecDestroy(&pipe->x);
744:   }
745:   if (size == 1) {
746:     DMNetworkMonitorDestroy(&monitor);
747:   }
748:   DMDestroy(&networkdm);
749:   PetscFree(wash);
750:   PetscFinalize();
751:   return ierr;
752: }