Actual source code: pipes1.c

petsc-3.7.2 2016-06-05
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

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

 32:   PetscObjectGetComm((PetscObject)ts,&comm);
 33:   MPI_Comm_rank(comm,&rank);
 34:   MPI_Comm_size(comm,&size);
 35: 
 36:   VecSet(F,0.0);
 37: 
 38:   localX    = wash->localX;
 39:   localXdot = wash->localXdot;
 40: 
 41:   TSGetDM(ts,&networkdm);
 42:   DMGetLocalVector(networkdm,&localF);
 43:   VecSet(localF,0.0);

 45:   /* update ghost values of locaX and locaXdot */
 46:   DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);
 47:   DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);
 48: 
 49:   DMGlobalToLocalBegin(networkdm,Xdot,INSERT_VALUES,localXdot);
 50:   DMGlobalToLocalEnd(networkdm,Xdot,INSERT_VALUES,localXdot);

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

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

 72:   /* Edge */
 73:   DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);
 74:   for (e=eStart; e<eEnd; e++) {
 75:     DMNetworkGetComponentTypeOffset(networkdm,e,0,&type,&pipeoffset);
 76:     DMNetworkGetVariableOffset(networkdm,e,&varoffset);
 77:     pipe     = (Pipe)(nwarr + pipeoffset);
 78:     pipex    = (PipeField*)(xarr + varoffset);
 79:     pipexdot = (PipeField*)(xdotarr + varoffset);
 80:     pipef    = (PipeField*)(farr + varoffset);
 81: 
 82:     /* Get boundary values H0 and QL from connected vertices */
 83:     DMNetworkGetConnectedNodes(networkdm,e,&cone);
 84:     vfrom = cone[0]; /* local ordering */
 85:     vto   = cone[1];
 86:     DMNetworkGetVariableOffset(networkdm,vfrom,&offsetfrom);
 87:     DMNetworkGetVariableOffset(networkdm,vto,&offsetto);
 88:     if (pipe->boundary.Q0 == PIPE_CHARACTERISTIC) {
 89:       pipe->boundary.H0 = (xarr+offsetfrom)[1]; /* h_from */
 90:     } else {
 91:       pipe->boundary.Q0 = (xarr+offsetfrom)[0]; /* q_from */
 92:     }
 93:     if (pipe->boundary.HL == PIPE_CHARACTERISTIC) {
 94:       pipe->boundary.QL = (xarr+offsetto)[0];   /* q_to */
 95:     } else {
 96:       pipe->boundary.HL = (xarr+offsetto)[1];   /* h_to */
 97:     }

 99:     /* Evaluate PipeIFunctionLocal() */
100:     DMDAGetLocalInfo(pipe->da,&info);
101:     PipeIFunctionLocal(&info, t, pipex, pipexdot, pipef, pipe);
102: 
103:     /* Set F at vfrom */
104:     vf = (PetscScalar*)(farr+offsetfrom);
105:     if (pipe->boundary.Q0 == PIPE_CHARACTERISTIC) {
106:       vf[0] -= pipex[0].q; /* q_vfrom - q[0] */
107:     } else {
108:       vf[1] -= pipex[0].h; /* h_vfrom - h[0] */
109:     }

111:     /* Set F at vto */
112:     vf = (PetscScalar*)(farr+offsetto);
113:     if (pipe->boundary.HL == PIPE_CHARACTERISTIC) {
114:       vf[1] -= pipex[pipe->nnodes-1].h; /* h_vto - h[last] */
115:     } else {
116:       vf[0] -= pipex[pipe->nnodes-1].q; /* q_vto - q[last] */
117:     }
118:   }
119: 
120:   /* Set F at boundary vertices */
121:   for (v=vStart; v<vEnd; v++) {
122:     DMNetworkGetComponentTypeOffset(networkdm,v,0,&type,&voffset);
123:     DMNetworkGetVariableOffset(networkdm,v,&varoffset);
124:     junction = (Junction)(nwarr + voffset);
125:     juncf = (PetscScalar *)(farr + varoffset);
126:     if (junction->isEnd == -1) {
127:       juncf[1] -= wash->H0;
128:       } else if (junction->isEnd == 1) {
129:       juncf[0] -= wash->QL;
130:     }
131:   }

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

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

146: PetscErrorCode WASHSetInitialSolution(DM networkdm,Vec X,Wash wash)
147: {
149:   PetscInt       k,nx,vkey,vOffset,vfrom,vto,offsetfrom,offsetto;
150:   PetscInt       type,varoffset;
151:   PetscInt       e,eStart,eEnd,pipeoffset;
152:   Vec            localX;
153:   PetscScalar    *xarr;
154:   Pipe           pipe;
155:   Junction       junction;
156:   const PetscInt *cone;
157:   const PetscScalar *xarray;
158:   DMNetworkComponentGenericDataType *nwarr;
159: 
161:   VecSet(X,0.0);
162:   DMGetLocalVector(networkdm,&localX);
163:   VecGetArray(localX,&xarr);

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

168:   /* Edge */
169:   DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);
170:   for (e=eStart; e<eEnd; e++) {
171:     DMNetworkGetVariableOffset(networkdm,e,&varoffset);
172:     DMNetworkGetComponentTypeOffset(networkdm,e,0,&type,&pipeoffset);
173: 
174:     /* get from and to vertices */
175:     DMNetworkGetConnectedNodes(networkdm,e,&cone);
176:     vfrom = cone[0]; /* local ordering */
177:     vto   = cone[1];
178: 
179:     DMNetworkGetVariableOffset(networkdm,vfrom,&offsetfrom);
180:     DMNetworkGetVariableOffset(networkdm,vto,&offsetto);

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

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

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

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

207:     /* if vform is a head vertex: */
208:     DMNetworkGetComponentTypeOffset(networkdm,vfrom,0,&vkey,&vOffset);
209:     junction = (Junction)(nwarr+vOffset);
210:     if (junction->isEnd == -1) { /* head junction */
211:       (xarr+offsetfrom)[0] = 0.0;      /* 1st Q -- not used */
212:       (xarr+offsetfrom)[1] = wash->H0; /* 1st H */
213:     }
214: 
215:     /* if vto is an end vertex: */
216:     DMNetworkGetComponentTypeOffset(networkdm,vto,0,&vkey,&vOffset);
217:     junction = (Junction)(nwarr+vOffset);
218:     if (junction->isEnd == 1) { /* end junction */
219:       (xarr+offsetto)[0] = wash->QL; /* last Q */
220:       (xarr+offsetto)[1] = 0.0;      /* last H -- not used */
221:     }
222:     VecRestoreArrayRead(pipe->x,&xarray);
223:   }

225:   VecRestoreArray(localX,&xarr);
226:   DMLocalToGlobalBegin(networkdm,localX,ADD_VALUES,X);
227:   DMLocalToGlobalEnd(networkdm,localX,ADD_VALUES,X);
228:   DMRestoreLocalVector(networkdm,&localX);

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

241: PetscErrorCode TSDMNetworkMonitor(TS ts, PetscInt step, PetscReal t, Vec x, void *context)
242: {
243:   PetscErrorCode     ierr;
244:   DMNetworkMonitor   monitor;

247:   monitor = (DMNetworkMonitor)context;
248:   DMNetworkMonitorView(monitor,x);
249:   return(0);
250: }

254: PetscErrorCode PipesView(Vec X,DM networkdm,Wash wash)
255: {
256:   PetscErrorCode       ierr;
257:   Pipe                 pipe;
258:   DMNetworkComponentGenericDataType *nwarr;
259:   PetscInt             pipeOffset,key,Start,End;
260:   PetscMPIInt          rank;
261:   PetscInt             nx,nnodes,nidx,*idx1,*idx2,*idx1_h,*idx2_h,idx_start,i,k,k1,xstart,j1;
262:   Vec                  Xq,Xh,localX;
263:   IS                   is1_q,is2_q,is1_h,is2_h;
264:   VecScatter           ctx_q,ctx_h;

267:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

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

273:   VecCreate(PETSC_COMM_WORLD,&Xq);
274:   if (rank == 0) { /* all entries of Xq are in proc[0] */
275:     VecSetSizes(Xq,nx,PETSC_DECIDE);
276:   } else {
277:     VecSetSizes(Xq,0,PETSC_DECIDE);
278:   }
279:   VecSetFromOptions(Xq);
280:   VecSet(Xq,0.0);
281:   VecDuplicate(Xq,&Xh);

283:   DMGetLocalVector(networkdm,&localX);

285:   /* set idx1 and idx2 */
286:   PetscCalloc4(nidx,&idx1,nidx,&idx2,nidx,&idx1_h,nidx,&idx2_h);
287: 
288:   DMNetworkGetComponentDataArray(networkdm,&nwarr);
289:   DMNetworkGetEdgeRange(networkdm,&Start, &End);

291:   VecGetOwnershipRange(X,&xstart,NULL);
292:   k1 = 0;
293:   j1 = 0;
294:   for (i = Start; i < End; i++) {
295:     DMNetworkGetComponentTypeOffset(networkdm,i,0,&key,&pipeOffset);
296:     pipe = (Pipe)(nwarr+pipeOffset);
297:     nnodes = pipe->nnodes;
298:     idx_start = pipe->id*nnodes;
299:     for (k=0; k<nnodes; k++) {
300:       idx1[k1] = xstart + j1*2*nnodes + 2*k;
301:       idx2[k1] = idx_start + k;

303:       idx1_h[k1] = xstart + j1*2*nnodes + 2*k + 1;
304:       idx2_h[k1] = idx_start + k;
305:       k1++;
306:     }
307:     j1++;
308:   }

310:   ISCreateGeneral(PETSC_COMM_SELF,nidx,idx1,PETSC_COPY_VALUES,&is1_q);
311:   ISCreateGeneral(PETSC_COMM_SELF,nidx,idx2,PETSC_COPY_VALUES,&is2_q);
312:   VecScatterCreate(X,is1_q,Xq,is2_q,&ctx_q);
313:   VecScatterBegin(ctx_q,X,Xq,INSERT_VALUES,SCATTER_FORWARD);
314:   VecScatterEnd(ctx_q,X,Xq,INSERT_VALUES,SCATTER_FORWARD);

316:   ISCreateGeneral(PETSC_COMM_SELF,nidx,idx1_h,PETSC_COPY_VALUES,&is1_h);
317:   ISCreateGeneral(PETSC_COMM_SELF,nidx,idx2_h,PETSC_COPY_VALUES,&is2_h);
318:   VecScatterCreate(X,is1_h,Xh,is2_h,&ctx_h);
319:   VecScatterBegin(ctx_h,X,Xh,INSERT_VALUES,SCATTER_FORWARD);
320:   VecScatterEnd(ctx_h,X,Xh,INSERT_VALUES,SCATTER_FORWARD);
321: 
322:   if (!rank) printf("Xq: \n");
323:   VecView(Xq,PETSC_VIEWER_STDOUT_WORLD);
324:   if (!rank) printf("Xh: \n");
325:   VecView(Xh,PETSC_VIEWER_STDOUT_WORLD);

327: 
328:   VecScatterDestroy(&ctx_q);
329:   PetscFree4(idx1,idx2,idx1_h,idx2_h);
330:   ISDestroy(&is1_q);
331:   ISDestroy(&is2_q);

333:   VecScatterDestroy(&ctx_h);
334:   ISDestroy(&is1_h);
335:   ISDestroy(&is2_h);
336: 
337:   VecDestroy(&Xq);
338:   VecDestroy(&Xh);
339:   DMRestoreLocalVector(networkdm,&localX);
340:   return(0);
341: }

345: PetscErrorCode WashNetworkCleanUp(Wash wash,int *edgelist)
346: {
348:   PetscMPIInt    rank;

351:   MPI_Comm_rank(wash->comm,&rank);
352:   if (!rank) {
353:     PetscFree(edgelist);
354:     PetscFree2(wash->junction,wash->pipe);
355:   }
356:   return(0);
357: }

361: PetscErrorCode WashNetworkCreate(MPI_Comm comm,PetscInt pipesCase,Wash *wash_ptr,int **elist)
362: {
364:   PetscInt       nnodes,npipes;
365:   PetscMPIInt    rank;
366:   Wash           wash;
367:   PetscInt       i,numVertices,numEdges;
368:   int            *edgelist;
369:   Junction       junctions=NULL;
370:   Pipe           pipes=NULL;

373:   MPI_Comm_rank(comm,&rank);

375:   PetscCalloc1(1,&wash);
376:   wash->comm = comm;
377:   *wash_ptr  = wash;
378:   wash->Q0   = 0.477432; /* copied from initial soluiton */
379:   wash->H0   = 150.0;
380:   wash->HL   = 143.488; /* copied from initial soluiton */
381:   wash->nnodes_loc = 0;

383:   numVertices = 0;
384:   numEdges    = 0;
385:   edgelist    = NULL;

387:   if (!rank) {
388:     PetscPrintf(PETSC_COMM_SELF,"Setup pipesCase %D\n",pipesCase);
389:   }
390:   nnodes = 6;
391:   PetscOptionsGetInt(NULL,PETSC_NULL, "-npipenodes", &nnodes, PETSC_NULL);

393:   /* Set global number of pipes, edges, and junctions */
394:   /*-------------------------------------------------*/
395:   switch (pipesCase) {
396:   case 0:
397:     /* pipeCase 0: */
398:     /* =============================
399:     v0 --E0--> v1--E1--> v2 --E2-->v3
400:     ================================  */
401:     npipes = 3;
402:     PetscOptionsGetInt(NULL,PETSC_NULL, "-npipes", &npipes, PETSC_NULL);
403:     wash->nedge   = npipes;
404:     wash->nvertex = npipes + 1;
405: 
406:     /* Set local edges and vertices -- proc[0] sets entire network, then distributes */
407:     numVertices = 0;
408:     numEdges    = 0;
409:     edgelist    = NULL;
410:     if (!rank) {
411:       numVertices = wash->nvertex;
412:       numEdges    = wash->nedge;

414:       PetscCalloc1(2*numEdges,&edgelist);
415:       for (i=0; i<numEdges; i++) {
416:         edgelist[2*i] = i; edgelist[2*i+1] = i+1;
417:       }

419:       /* Add network components */
420:       /*------------------------*/
421:       PetscCalloc2(numVertices,&junctions,numEdges,&pipes);
422:       /* vertex */
423:       for (i=0; i<numVertices; i++) {
424:         junctions[i].id = i;
425:         junctions[i].isEnd = 0;
426:         junctions[i].nedges_in = 1; junctions[i].nedges_out = 1;
427: 
428:         /* Set GPS data */
429:         junctions[i].latitude  = 0.0;
430:         junctions[i].longitude = 0.0;
431:       }
432:       junctions[0].isEnd                  = -1;
433:       junctions[0].nedges_in              =  0;
434:       junctions[numVertices-1].isEnd      =  1;
435:       junctions[numVertices-1].nedges_out =  0;

437:       /* edge and pipe */
438:       for (i=0; i<numEdges; i++) {
439:         pipes[i].id   = i;
440:         pipes[i].nnodes = nnodes;
441:       }
442:     }
443:     break;
444:   case 1:
445:     /* pipeCase 1: */
446:     /* ==========================
447:                 v2
448:                 ^
449:                 |
450:                E2
451:                 |
452:     v0 --E0--> v3--E1--> v1
453:     =============================  */
454:     npipes = 3;
455:     wash->nedge   = npipes;
456:     wash->nvertex = npipes + 1;

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

463:       PetscCalloc1(2*numEdges,&edgelist);
464:       edgelist[0] = 0; edgelist[1] = 3;  /* edge[0] */
465:       edgelist[2] = 3; edgelist[3] = 1;  /* edge[1] */
466:       edgelist[4] = 3; edgelist[5] = 2;  /* edge[2] */

468:       /* Add network components */
469:       /*------------------------*/
470:       PetscCalloc2(numVertices,&junctions,numEdges,&pipes);
471:       /* vertex */
472:       for (i=0; i<numVertices; i++) {
473:         junctions[i].id = i;
474: 
475:         /* Set GPS data */
476:         junctions[i].latitude  = 0.0;
477:         junctions[i].longitude = 0.0;
478:       }
479:       junctions[0].isEnd = -1; junctions[0].nedges_in = 0; junctions[0].nedges_out = 1;
480:       junctions[1].isEnd =  1; junctions[1].nedges_in = 1; junctions[1].nedges_out = 0;
481:       junctions[2].isEnd =  1; junctions[2].nedges_in = 1; junctions[2].nedges_out = 0;
482:       junctions[3].isEnd =  0; junctions[3].nedges_in = 1; junctions[3].nedges_out = 2;

484:       /* edge and pipe */
485:       for (i=0; i<numEdges; i++) {
486:         pipes[i].id     = i;
487:         pipes[i].nnodes = nnodes;
488:       }
489:     }
490:     break;
491:   case 2:
492:     /* pipeCase 2: */
493:     /* ==========================
494:          v2--> E2
495:                 |
496:     v0 --E0--> v3--E1--> v1
497:     =============================  */

499:     /* Set application parameters -- to be used in function evalutions */
500:     npipes = 3;
501:     wash->nedge   = npipes;
502:     wash->nvertex = npipes + 1;

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

509:       PetscCalloc1(2*numEdges,&edgelist);
510:       edgelist[0] = 0; edgelist[1] = 3;  /* edge[0] */
511:       edgelist[2] = 3; edgelist[3] = 1;  /* edge[1] */
512:       edgelist[4] = 2; edgelist[5] = 3;  /* edge[2] */

514:       /* Add network components */
515:       /*------------------------*/
516:       PetscCalloc2(numVertices,&junctions,numEdges,&pipes);
517:       /* vertex */
518:       for (i=0; i<numVertices; i++) {
519:         junctions[i].id = i;
520: 
521:         /* Set GPS data */
522:         junctions[i].latitude  = 0.0;
523:         junctions[i].longitude = 0.0;
524:       }
525:       junctions[0].isEnd = -1; junctions[0].nedges_in = 0; junctions[0].nedges_out = 1;
526:       junctions[1].isEnd =  1; junctions[1].nedges_in = 1; junctions[1].nedges_out = 0;
527:       junctions[2].isEnd = -1; junctions[2].nedges_in = 0; junctions[2].nedges_out = 1;
528:       junctions[3].isEnd =  0; junctions[3].nedges_in = 2; junctions[3].nedges_out = 1;

530:       /* edge and pipe */
531:       for (i=0; i<numEdges; i++) {
532:         pipes[i].id     = i;
533:         pipes[i].nnodes = nnodes;
534:       }
535:     }
536:     break;
537:   default:
538:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"not done yet");
539:   }
540: 
541:   *wash_ptr      = wash;
542:   wash->nedge    = numEdges;
543:   wash->nvertex  = numVertices;
544:   *elist         = edgelist;
545:   wash->junction = junctions;
546:   wash->pipe     = pipes;
547:   return(0);
548: }
549: /* ------------------------------------------------------- */
552: int main(int argc,char ** argv)
553: {
554:   PetscErrorCode    ierr;
555:   Wash              wash;
556:   Junction          junctions,junction;
557:   Pipe              pipe,pipes;
558:   PetscInt          numEdges,numVertices,KeyPipe,KeyJunction;
559:   int               *edgelist = NULL;
560:   PetscInt          i,e,v,eStart,eEnd,vStart,vEnd,pipeOffset,key,frombType,tobType;
561:   PetscInt          vfrom,vto,vkey,fromOffset,toOffset,type,varoffset,pipeoffset;
562:   PetscInt          from_nedge_in,from_nedge_out,to_nedge_in;
563:   const PetscInt    *cone;
564:   DM                networkdm;
565:   PetscMPIInt       size,rank;
566:   PetscReal         ftime = 2500.0;
567:   Vec               X;
568:   TS                ts;
569:   PetscInt          maxsteps=-1,steps;
570:   TSConvergedReason reason;
571:   PetscBool         viewpipes;
572:   PetscInt          pipesCase;
573:   DMNetworkMonitor  monitor;
574:   DMNetworkComponentGenericDataType *nwarr;

576:   PetscInitialize(&argc,&argv,(char*)0,help);
577:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
578:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

580:   /* Create and setup network */
581:   /*--------------------------*/
582:   DMNetworkCreate(PETSC_COMM_WORLD,&networkdm);
583:   if (size == 1) {
584:     DMNetworkMonitorCreate(networkdm,&monitor);
585:   }
586:   /* Register the components in the network */
587:   DMNetworkRegisterComponent(networkdm,"junctionstruct",sizeof(struct _p_Junction),&KeyJunction);
588:   DMNetworkRegisterComponent(networkdm,"pipestruct",sizeof(struct _p_Pipe),&KeyPipe);

590:   /* Set global number of pipes, edges, and vertices */
591:   pipesCase = 2;
592:   PetscOptionsGetInt(NULL,PETSC_NULL, "-case", &pipesCase, PETSC_NULL);

594:   WashNetworkCreate(PETSC_COMM_WORLD,pipesCase,&wash,&edgelist);
595:   numEdges    = wash->nedge;
596:   numVertices = wash->nvertex;
597:   junctions    = wash->junction;
598:   pipes       = wash->pipe;

600:   /* Set number of vertices and edges */
601:   DMNetworkSetSizes(networkdm,numVertices,numEdges,PETSC_DETERMINE,PETSC_DETERMINE);
602:   /* Add edge connectivity */
603:   DMNetworkSetEdgeList(networkdm,edgelist);
604:   /* Set up the network layout */
605:   DMNetworkLayoutSetUp(networkdm);

607:   /* Add EDGEDATA component to all edges -- currently networkdm is a sequential network */
608:   DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);
609:   DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);
610: 
611:   for (e = eStart; e < eEnd; e++) {
612:     /* Add Pipe component to all edges -- create pipe here */
613:     DMNetworkAddComponent(networkdm,e,KeyPipe,&pipes[e-eStart]);
614: 
615:     /* Add number of variables to each edge */
616:     DMNetworkAddNumVariables(networkdm,e,2*pipes[e-eStart].nnodes);

618:     if (size == 1) { /* Add monitor -- show Q_{pipes[e-eStart].id}? */
619:       DMNetworkMonitorAdd(monitor, "Pipe Q", e, pipes[e-eStart].nnodes, 0, 2, -0.8, 0.8, PETSC_TRUE);
620:       DMNetworkMonitorAdd(monitor, "Pipe H", e, pipes[e-eStart].nnodes, 1, 2, -400.0, 800.0, PETSC_TRUE);
621:     }
622:   }

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

628:     /* Add number of variables to vertex */
629:     DMNetworkAddNumVariables(networkdm,v,2);
630:   }

632:   /* Set up DM for use */
633:   DMSetUp(networkdm);
634:   WashNetworkCleanUp(wash,edgelist);

636:   /* Network partitioning and distribution of data */
637:   if (size > 1) {
638:     DM distnetworkdm;
639:     DMNetworkDistribute(networkdm,0,&distnetworkdm);
640:     DMDestroy(&networkdm);
641:     networkdm = distnetworkdm;
642:   }
643: 
644:   /* PipeSetUp -- each process only sets its own pipes */
645:   /*---------------------------------------------------*/
646:   DMNetworkGetComponentDataArray(networkdm,&nwarr);

648:   DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);
649:   for (e=eStart; e<eEnd; e++) { /* each edge has only one component, pipe */
650:     DMNetworkGetComponentTypeOffset(networkdm,e,0,&type,&pipeoffset);
651:     DMNetworkGetVariableOffset(networkdm,e,&varoffset);
652:     pipe = (Pipe)(nwarr + pipeoffset);

654:     /* Setup conntected vertices */
655:     DMNetworkGetConnectedNodes(networkdm,e,&cone);
656:     vfrom = cone[0]; /* local ordering */
657:     vto   = cone[1];

659:     /* vfrom */
660:     DMNetworkGetComponentTypeOffset(networkdm,vfrom,0,&vkey,&fromOffset);
661:     junction = (Junction)(nwarr+fromOffset);
662:     from_nedge_in  = junction->nedges_in;
663:     from_nedge_out = junction->nedges_out;
664: 
665:     /* vto */
666:     DMNetworkGetComponentTypeOffset(networkdm,vto,0,&vkey,&toOffset);
667:     junction    = (Junction)(nwarr+toOffset);
668:     to_nedge_in = junction->nedges_in;
669: 
670:     pipe->comm = PETSC_COMM_SELF; /* must be set here, otherwise crashes in my mac??? */
671:     wash->nnodes_loc += pipe->nnodes; /* local total num of nodes, will be used by PipesView() */
672:     PipeSetParameters(pipe,
673:                              600.0,          /* length */
674:                              pipe->nnodes,   /* nnodes -- rm from PipeSetParameters */
675:                              0.5,            /* diameter */
676:                              1200.0,         /* a */
677:                              0.018);    /* friction */

679:     /* set boundary conditions for this pipe */
680:     if (from_nedge_in <= 1 && from_nedge_out > 0) {
681:       frombType = 0;
682:     } else {
683:       frombType = 1;
684:     }

686:     if (to_nedge_in == 1) {
687:       tobType = 0;
688:     } else {
689:       tobType = 1;
690:     }

692:     if (frombType == 0) {
693:       pipe->boundary.Q0 = PIPE_CHARACTERISTIC; /* will be obtained from characteristic */
694:       pipe->boundary.H0 = wash->H0;
695:     } else {
696:       pipe->boundary.Q0 = wash->Q0;
697:       pipe->boundary.H0 = PIPE_CHARACTERISTIC; /* will be obtained from characteristic */
698:     }
699:     if (tobType == 0) {
700:       pipe->boundary.QL = wash->QL;
701:       pipe->boundary.HL = PIPE_CHARACTERISTIC; /* will be obtained from characteristic */
702:     } else {
703:       pipe->boundary.QL = PIPE_CHARACTERISTIC; /* will be obtained from characteristic */
704:       pipe->boundary.HL = wash->HL;
705:     }
706: 
707:     PipeSetUp(pipe);
708:   }
709: 
710:   /* create vectors */
711:   DMCreateGlobalVector(networkdm,&X);
712:   DMCreateLocalVector(networkdm,&wash->localX);
713:   DMCreateLocalVector(networkdm,&wash->localXdot);

715:   /* Setup solver                                           */
716:   /*--------------------------------------------------------*/
717:   TSCreate(PETSC_COMM_WORLD,&ts);

719:   TSSetDM(ts,(DM)networkdm);
720:   TSSetIFunction(ts,NULL,WASHIFunction,wash);
721: 
722:   TSSetDuration(ts,maxsteps,ftime);
723:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
724:   TSSetInitialTimeStep(ts,0.0,0.1);
725:   TSSetType(ts,TSBEULER);
726:   if (size == 1) {
727:     TSMonitorSet(ts, TSDMNetworkMonitor, monitor, PETSC_NULL);
728:   }
729:   TSSetFromOptions(ts);

731:   WASHSetInitialSolution(networkdm,X,wash);

733:   TSSolve(ts,X);

735:   TSGetSolveTime(ts,&ftime);
736:   TSGetTimeStepNumber(ts,&steps);
737:   TSGetConvergedReason(ts,&reason);
738:   PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps);
739:   /* VecView(X,PETSC_VIEWER_STDOUT_WORLD); */
740: 
741:   /* View solution q and h */
742:   /* --------------------- */
743:   viewpipes = PETSC_FALSE;
744:   PetscOptionsGetBool(NULL,NULL, "-pipe_view", &viewpipes,NULL);
745:   if (viewpipes) {
746:     PipesView(X,networkdm,wash);
747:   }

749:   /* Free spaces */
750:   /* ----------- */
751:   TSDestroy(&ts);
752:   VecDestroy(&X);
753:   VecDestroy(&wash->localX);
754:   VecDestroy(&wash->localXdot);
755: 
756:   /* Destroy objects from each pipe that are created in PipeSetUp() */
757:   DMNetworkGetEdgeRange(networkdm,&eStart, &eEnd);
758:   for (i = eStart; i < eEnd; i++) {
759:     DMNetworkGetComponentTypeOffset(networkdm,i,0,&key,&pipeOffset);
760:     pipe = (Pipe)(nwarr+pipeOffset);
761:     DMDestroy(&(pipe->da));
762:     VecDestroy(&pipe->x);
763:   }
764:   if (size == 1) {
765:     DMNetworkMonitorDestroy(&monitor);
766:   }
767:   DMDestroy(&networkdm);
768:   PetscFree(wash);
769:   PetscFinalize();
770:   return 0;
771: }