142dc13f1SHong Zhang static char help[] = "This example demonstrates DMNetwork. It is used for testing parallel generation of dmnetwork, then redistribute. \n\\n"; 242dc13f1SHong Zhang /* 342dc13f1SHong Zhang Example: mpiexec -n <np> ./pipes -ts_max_steps 10 442dc13f1SHong Zhang */ 542dc13f1SHong Zhang 642dc13f1SHong Zhang #include "wash.h" 742dc13f1SHong Zhang 842dc13f1SHong Zhang /* 942dc13f1SHong Zhang WashNetworkDistribute - proc[0] distributes sequential wash object 1042dc13f1SHong Zhang Input Parameters: 1142dc13f1SHong Zhang . comm - MPI communicator 1242dc13f1SHong Zhang . wash - wash context with all network data in proc[0] 1342dc13f1SHong Zhang 1442dc13f1SHong Zhang Output Parameter: 1542dc13f1SHong Zhang . wash - wash context with nedge, nvertex and edgelist distributed 1642dc13f1SHong Zhang 1742dc13f1SHong Zhang Note: The routine is used for testing parallel generation of dmnetwork, then redistribute. 1842dc13f1SHong Zhang */ 1942dc13f1SHong Zhang PetscErrorCode WashNetworkDistribute(MPI_Comm comm,Wash wash) 2042dc13f1SHong Zhang { 2142dc13f1SHong Zhang PetscErrorCode ierr; 2242dc13f1SHong Zhang PetscMPIInt rank,size,tag=0; 2342dc13f1SHong Zhang PetscInt i,e,v,numEdges,numVertices,nedges,*eowners=NULL,estart,eend,*vtype=NULL,nvertices; 2442dc13f1SHong Zhang PetscInt *edgelist = wash->edgelist,*nvtx=NULL,*vtxDone=NULL; 2542dc13f1SHong Zhang 2642dc13f1SHong Zhang PetscFunctionBegin; 2742dc13f1SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 2842dc13f1SHong Zhang if (size == 1) PetscFunctionReturn(0); 2942dc13f1SHong Zhang 3042dc13f1SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 3142dc13f1SHong Zhang numEdges = wash->nedge; 3242dc13f1SHong Zhang numVertices = wash->nvertex; 3342dc13f1SHong Zhang 3442dc13f1SHong Zhang /* (1) all processes get global and local number of edges */ 3542dc13f1SHong Zhang ierr = MPI_Bcast(&numEdges,1,MPIU_INT,0,comm);CHKERRMPI(ierr); 3642dc13f1SHong Zhang nedges = numEdges/size; /* local nedges */ 37dd400576SPatrick Sanan if (rank == 0) { 3842dc13f1SHong Zhang nedges += numEdges - size*(numEdges/size); 3942dc13f1SHong Zhang } 4042dc13f1SHong Zhang wash->Nedge = numEdges; 4142dc13f1SHong Zhang wash->nedge = nedges; 4242dc13f1SHong Zhang /* ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] nedges %d, numEdges %d\n",rank,nedges,numEdges);CHKERRQ(ierr); */ 4342dc13f1SHong Zhang 4442dc13f1SHong Zhang ierr = PetscCalloc3(size+1,&eowners,size,&nvtx,numVertices,&vtxDone);CHKERRQ(ierr); 4542dc13f1SHong Zhang ierr = MPI_Allgather(&nedges,1,MPIU_INT,eowners+1,1,MPIU_INT,PETSC_COMM_WORLD);CHKERRMPI(ierr); 4642dc13f1SHong Zhang eowners[0] = 0; 4742dc13f1SHong Zhang for (i=2; i<=size; i++) { 4842dc13f1SHong Zhang eowners[i] += eowners[i-1]; 4942dc13f1SHong Zhang } 5042dc13f1SHong Zhang 5142dc13f1SHong Zhang estart = eowners[rank]; 5242dc13f1SHong Zhang eend = eowners[rank+1]; 5342dc13f1SHong Zhang /* ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] own lists row %d - %d\n",rank,estart,eend);CHKERRQ(ierr); */ 5442dc13f1SHong Zhang 5542dc13f1SHong Zhang /* (2) distribute row block edgelist to all processors */ 56dd400576SPatrick Sanan if (rank == 0) { 5742dc13f1SHong Zhang vtype = wash->vtype; 5842dc13f1SHong Zhang for (i=1; i<size; i++) { 5942dc13f1SHong Zhang /* proc[0] sends edgelist to proc[i] */ 6042dc13f1SHong Zhang ierr = MPI_Send(edgelist+2*eowners[i],2*(eowners[i+1]-eowners[i]),MPIU_INT,i,tag,comm);CHKERRMPI(ierr); 6142dc13f1SHong Zhang 6242dc13f1SHong Zhang /* proc[0] sends vtype to proc[i] */ 6342dc13f1SHong Zhang ierr = MPI_Send(vtype+2*eowners[i],2*(eowners[i+1]-eowners[i]),MPIU_INT,i,tag,comm);CHKERRMPI(ierr); 6442dc13f1SHong Zhang } 6542dc13f1SHong Zhang } else { 6642dc13f1SHong Zhang MPI_Status status; 6742dc13f1SHong Zhang ierr = PetscMalloc1(2*(eend-estart),&vtype);CHKERRQ(ierr); 6842dc13f1SHong Zhang ierr = PetscMalloc1(2*(eend-estart),&edgelist);CHKERRQ(ierr); 6942dc13f1SHong Zhang 7042dc13f1SHong Zhang ierr = MPI_Recv(edgelist,2*(eend-estart),MPIU_INT,0,tag,comm,&status);CHKERRMPI(ierr); 7142dc13f1SHong Zhang ierr = MPI_Recv(vtype,2*(eend-estart),MPIU_INT,0,tag,comm,&status);CHKERRMPI(ierr); 7242dc13f1SHong Zhang } 7342dc13f1SHong Zhang 7442dc13f1SHong Zhang wash->edgelist = edgelist; 7542dc13f1SHong Zhang 7642dc13f1SHong Zhang /* (3) all processes get global and local number of vertices, without ghost vertices */ 77dd400576SPatrick Sanan if (rank == 0) { 7842dc13f1SHong Zhang for (i=0; i<size; i++) { 7942dc13f1SHong Zhang for (e=eowners[i]; e<eowners[i+1]; e++) { 8042dc13f1SHong Zhang v = edgelist[2*e]; 8142dc13f1SHong Zhang if (!vtxDone[v]) { 8242dc13f1SHong Zhang nvtx[i]++; vtxDone[v] = 1; 8342dc13f1SHong Zhang } 8442dc13f1SHong Zhang v = edgelist[2*e+1]; 8542dc13f1SHong Zhang if (!vtxDone[v]) { 8642dc13f1SHong Zhang nvtx[i]++; vtxDone[v] = 1; 8742dc13f1SHong Zhang } 8842dc13f1SHong Zhang } 8942dc13f1SHong Zhang } 9042dc13f1SHong Zhang } 9142dc13f1SHong Zhang ierr = MPI_Bcast(&numVertices,1,MPIU_INT,0,PETSC_COMM_WORLD);CHKERRMPI(ierr); 9242dc13f1SHong Zhang ierr = MPI_Scatter(nvtx,1,MPIU_INT,&nvertices,1,MPIU_INT,0,PETSC_COMM_WORLD);CHKERRMPI(ierr); 9342dc13f1SHong Zhang ierr = PetscFree3(eowners,nvtx,vtxDone);CHKERRQ(ierr); 9442dc13f1SHong Zhang 9542dc13f1SHong Zhang wash->Nvertex = numVertices; 9642dc13f1SHong Zhang wash->nvertex = nvertices; 9742dc13f1SHong Zhang wash->vtype = vtype; 9842dc13f1SHong Zhang PetscFunctionReturn(0); 9942dc13f1SHong Zhang } 10042dc13f1SHong Zhang 10142dc13f1SHong Zhang PetscErrorCode WASHIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void* ctx) 10242dc13f1SHong Zhang { 10342dc13f1SHong Zhang PetscErrorCode ierr; 10442dc13f1SHong Zhang Wash wash=(Wash)ctx; 10542dc13f1SHong Zhang DM networkdm; 10642dc13f1SHong Zhang Vec localX,localXdot,localF, localXold; 10742dc13f1SHong Zhang const PetscInt *cone; 10842dc13f1SHong Zhang PetscInt vfrom,vto,offsetfrom,offsetto,varoffset; 10942dc13f1SHong Zhang PetscInt v,vStart,vEnd,e,eStart,eEnd; 11042dc13f1SHong Zhang PetscInt nend,type; 11142dc13f1SHong Zhang PetscBool ghost; 11242dc13f1SHong Zhang PetscScalar *farr,*juncf, *pipef; 11342dc13f1SHong Zhang PetscReal dt; 11442dc13f1SHong Zhang Pipe pipe; 11542dc13f1SHong Zhang PipeField *pipex,*pipexdot,*juncx; 11642dc13f1SHong Zhang Junction junction; 11742dc13f1SHong Zhang DMDALocalInfo info; 11842dc13f1SHong Zhang const PetscScalar *xarr,*xdotarr, *xoldarr; 11942dc13f1SHong Zhang 12042dc13f1SHong Zhang PetscFunctionBegin; 12142dc13f1SHong Zhang localX = wash->localX; 12242dc13f1SHong Zhang localXdot = wash->localXdot; 12342dc13f1SHong Zhang 12442dc13f1SHong Zhang ierr = TSGetSolution(ts,&localXold);CHKERRQ(ierr); 12542dc13f1SHong Zhang ierr = TSGetDM(ts,&networkdm);CHKERRQ(ierr); 12642dc13f1SHong Zhang ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 12742dc13f1SHong Zhang ierr = DMGetLocalVector(networkdm,&localF);CHKERRQ(ierr); 12842dc13f1SHong Zhang 12942dc13f1SHong Zhang /* Set F and localF as zero */ 13042dc13f1SHong Zhang ierr = VecSet(F,0.0);CHKERRQ(ierr); 13142dc13f1SHong Zhang ierr = VecSet(localF,0.0);CHKERRQ(ierr); 13242dc13f1SHong Zhang 13342dc13f1SHong Zhang /* Update ghost values of locaX and locaXdot */ 13442dc13f1SHong Zhang ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 13542dc13f1SHong Zhang ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 13642dc13f1SHong Zhang 13742dc13f1SHong Zhang ierr = DMGlobalToLocalBegin(networkdm,Xdot,INSERT_VALUES,localXdot);CHKERRQ(ierr); 13842dc13f1SHong Zhang ierr = DMGlobalToLocalEnd(networkdm,Xdot,INSERT_VALUES,localXdot);CHKERRQ(ierr); 13942dc13f1SHong Zhang 14042dc13f1SHong Zhang ierr = VecGetArrayRead(localX,&xarr);CHKERRQ(ierr); 14142dc13f1SHong Zhang ierr = VecGetArrayRead(localXdot,&xdotarr);CHKERRQ(ierr); 14242dc13f1SHong Zhang ierr = VecGetArrayRead(localXold,&xoldarr);CHKERRQ(ierr); 14342dc13f1SHong Zhang ierr = VecGetArray(localF,&farr);CHKERRQ(ierr); 14442dc13f1SHong Zhang 14542dc13f1SHong Zhang /* junction->type == JUNCTION: 14642dc13f1SHong Zhang juncf[0] = -qJ + sum(qin); juncf[1] = qJ - sum(qout) 14742dc13f1SHong Zhang junction->type == RESERVOIR (upper stream): 14842dc13f1SHong Zhang juncf[0] = -hJ + H0; juncf[1] = qJ - sum(qout) 14942dc13f1SHong Zhang junction->type == VALVE (down stream): 15042dc13f1SHong Zhang juncf[0] = -qJ + sum(qin); juncf[1] = qJ 15142dc13f1SHong Zhang */ 15242dc13f1SHong Zhang /* Vertex/junction initialization */ 15342dc13f1SHong Zhang ierr = DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);CHKERRQ(ierr); 15442dc13f1SHong Zhang for (v=vStart; v<vEnd; v++) { 15542dc13f1SHong Zhang ierr = DMNetworkIsGhostVertex(networkdm,v,&ghost);CHKERRQ(ierr); 15642dc13f1SHong Zhang if (ghost) continue; 15742dc13f1SHong Zhang 15842dc13f1SHong Zhang ierr = DMNetworkGetComponent(networkdm,v,0,&type,(void**)&junction,NULL);CHKERRQ(ierr); 15942dc13f1SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,v,ALL_COMPONENTS,&varoffset);CHKERRQ(ierr); 16042dc13f1SHong Zhang juncx = (PipeField*)(xarr + varoffset); 16142dc13f1SHong Zhang juncf = (PetscScalar*)(farr + varoffset); 16242dc13f1SHong Zhang 16342dc13f1SHong Zhang juncf[0] = -juncx[0].q; 16442dc13f1SHong Zhang juncf[1] = juncx[0].q; 16542dc13f1SHong Zhang 16642dc13f1SHong Zhang if (junction->type == RESERVOIR) { /* upstream reservoir */ 16742dc13f1SHong Zhang juncf[0] = juncx[0].h - wash->H0; 16842dc13f1SHong Zhang } 16942dc13f1SHong Zhang } 17042dc13f1SHong Zhang 17142dc13f1SHong Zhang /* Edge/pipe */ 17242dc13f1SHong Zhang ierr = DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);CHKERRQ(ierr); 17342dc13f1SHong Zhang for (e=eStart; e<eEnd; e++) { 17442dc13f1SHong Zhang ierr = DMNetworkGetComponent(networkdm,e,0,&type,(void**)&pipe,NULL);CHKERRQ(ierr); 17542dc13f1SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,e,ALL_COMPONENTS,&varoffset);CHKERRQ(ierr); 17642dc13f1SHong Zhang pipex = (PipeField*)(xarr + varoffset); 17742dc13f1SHong Zhang pipexdot = (PipeField*)(xdotarr + varoffset); 17842dc13f1SHong Zhang pipef = (PetscScalar*)(farr + varoffset); 17942dc13f1SHong Zhang 18042dc13f1SHong Zhang /* Get some data into the pipe structure: note, some of these operations 18142dc13f1SHong Zhang * might be redundant. Will it consume too much time? */ 18242dc13f1SHong Zhang pipe->dt = dt; 18342dc13f1SHong Zhang pipe->xold = (PipeField*)(xoldarr + varoffset); 18442dc13f1SHong Zhang 18542dc13f1SHong Zhang /* Evaluate F over this edge/pipe: pipef[1], ...,pipef[2*nend] */ 18642dc13f1SHong Zhang ierr = DMDAGetLocalInfo(pipe->da,&info);CHKERRQ(ierr); 18742dc13f1SHong Zhang ierr = PipeIFunctionLocal_Lax(&info,t,pipex,pipexdot,pipef,pipe);CHKERRQ(ierr); 18842dc13f1SHong Zhang 18942dc13f1SHong Zhang /* Get boundary values from connected vertices */ 19042dc13f1SHong Zhang ierr = DMNetworkGetConnectedVertices(networkdm,e,&cone);CHKERRQ(ierr); 19142dc13f1SHong Zhang vfrom = cone[0]; /* local ordering */ 19242dc13f1SHong Zhang vto = cone[1]; 19342dc13f1SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vfrom,ALL_COMPONENTS,&offsetfrom);CHKERRQ(ierr); 19442dc13f1SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vto,ALL_COMPONENTS,&offsetto);CHKERRQ(ierr); 19542dc13f1SHong Zhang 19642dc13f1SHong Zhang /* Evaluate upstream boundary */ 19742dc13f1SHong Zhang ierr = DMNetworkGetComponent(networkdm,vfrom,0,&type,(void**)&junction,NULL);CHKERRQ(ierr); 198*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(junction->type != JUNCTION && junction->type != RESERVOIR,PETSC_COMM_SELF,PETSC_ERR_SUP,"junction type is not supported"); 19942dc13f1SHong Zhang juncx = (PipeField*)(xarr + offsetfrom); 20042dc13f1SHong Zhang juncf = (PetscScalar*)(farr + offsetfrom); 20142dc13f1SHong Zhang 20242dc13f1SHong Zhang pipef[0] = pipex[0].h - juncx[0].h; 20342dc13f1SHong Zhang juncf[1] -= pipex[0].q; 20442dc13f1SHong Zhang 20542dc13f1SHong Zhang /* Evaluate downstream boundary */ 20642dc13f1SHong Zhang ierr = DMNetworkGetComponent(networkdm,vto,0,&type,(void**)&junction,NULL);CHKERRQ(ierr); 207*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(junction->type != JUNCTION && junction->type != VALVE,PETSC_COMM_SELF,PETSC_ERR_SUP,"junction type is not supported"); 20842dc13f1SHong Zhang juncx = (PipeField*)(xarr + offsetto); 20942dc13f1SHong Zhang juncf = (PetscScalar*)(farr + offsetto); 21042dc13f1SHong Zhang nend = pipe->nnodes - 1; 21142dc13f1SHong Zhang 21242dc13f1SHong Zhang pipef[2*nend + 1] = pipex[nend].h - juncx[0].h; 21342dc13f1SHong Zhang juncf[0] += pipex[nend].q; 21442dc13f1SHong Zhang } 21542dc13f1SHong Zhang 21642dc13f1SHong Zhang ierr = VecRestoreArrayRead(localX,&xarr);CHKERRQ(ierr); 21742dc13f1SHong Zhang ierr = VecRestoreArrayRead(localXdot,&xdotarr);CHKERRQ(ierr); 21842dc13f1SHong Zhang ierr = VecRestoreArray(localF,&farr);CHKERRQ(ierr); 21942dc13f1SHong Zhang 22042dc13f1SHong Zhang ierr = DMLocalToGlobalBegin(networkdm,localF,ADD_VALUES,F);CHKERRQ(ierr); 22142dc13f1SHong Zhang ierr = DMLocalToGlobalEnd(networkdm,localF,ADD_VALUES,F);CHKERRQ(ierr); 22242dc13f1SHong Zhang ierr = DMRestoreLocalVector(networkdm,&localF);CHKERRQ(ierr); 22342dc13f1SHong Zhang /* 22442dc13f1SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD("F:\n");CHKERRQ(ierr); 22542dc13f1SHong Zhang ierr = VecView(F,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 22642dc13f1SHong Zhang */ 22742dc13f1SHong Zhang PetscFunctionReturn(0); 22842dc13f1SHong Zhang } 22942dc13f1SHong Zhang 23042dc13f1SHong Zhang PetscErrorCode WASHSetInitialSolution(DM networkdm,Vec X,Wash wash) 23142dc13f1SHong Zhang { 23242dc13f1SHong Zhang PetscErrorCode ierr; 23342dc13f1SHong Zhang PetscInt k,nx,vkey,vfrom,vto,offsetfrom,offsetto; 23442dc13f1SHong Zhang PetscInt type,varoffset; 23542dc13f1SHong Zhang PetscInt e,eStart,eEnd; 23642dc13f1SHong Zhang Vec localX; 23742dc13f1SHong Zhang PetscScalar *xarr; 23842dc13f1SHong Zhang Pipe pipe; 23942dc13f1SHong Zhang Junction junction; 24042dc13f1SHong Zhang const PetscInt *cone; 24142dc13f1SHong Zhang const PetscScalar *xarray; 24242dc13f1SHong Zhang 24342dc13f1SHong Zhang PetscFunctionBegin; 24442dc13f1SHong Zhang ierr = VecSet(X,0.0);CHKERRQ(ierr); 24542dc13f1SHong Zhang ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr); 24642dc13f1SHong Zhang ierr = VecGetArray(localX,&xarr);CHKERRQ(ierr); 24742dc13f1SHong Zhang 24842dc13f1SHong Zhang /* Edge */ 24942dc13f1SHong Zhang ierr = DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);CHKERRQ(ierr); 25042dc13f1SHong Zhang for (e=eStart; e<eEnd; e++) { 25142dc13f1SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,e,ALL_COMPONENTS,&varoffset);CHKERRQ(ierr); 25242dc13f1SHong Zhang ierr = DMNetworkGetComponent(networkdm,e,0,&type,(void**)&pipe,NULL);CHKERRQ(ierr); 25342dc13f1SHong Zhang 25442dc13f1SHong Zhang /* set initial values for this pipe */ 25542dc13f1SHong Zhang ierr = PipeComputeSteadyState(pipe,wash->Q0,wash->H0);CHKERRQ(ierr); 25642dc13f1SHong Zhang ierr = VecGetSize(pipe->x,&nx);CHKERRQ(ierr); 25742dc13f1SHong Zhang 25842dc13f1SHong Zhang ierr = VecGetArrayRead(pipe->x,&xarray);CHKERRQ(ierr); 25942dc13f1SHong Zhang /* copy pipe->x to xarray */ 26042dc13f1SHong Zhang for (k=0; k<nx; k++) { 26142dc13f1SHong Zhang (xarr+varoffset)[k] = xarray[k]; 26242dc13f1SHong Zhang } 26342dc13f1SHong Zhang 26442dc13f1SHong Zhang /* set boundary values into vfrom and vto */ 26542dc13f1SHong Zhang ierr = DMNetworkGetConnectedVertices(networkdm,e,&cone);CHKERRQ(ierr); 26642dc13f1SHong Zhang vfrom = cone[0]; /* local ordering */ 26742dc13f1SHong Zhang vto = cone[1]; 26842dc13f1SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vfrom,ALL_COMPONENTS,&offsetfrom);CHKERRQ(ierr); 26942dc13f1SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vto,ALL_COMPONENTS,&offsetto);CHKERRQ(ierr); 27042dc13f1SHong Zhang 27142dc13f1SHong Zhang /* if vform is a head vertex: */ 27242dc13f1SHong Zhang ierr = DMNetworkGetComponent(networkdm,vfrom,0,&vkey,(void**)&junction,NULL);CHKERRQ(ierr); 27342dc13f1SHong Zhang if (junction->type == RESERVOIR) { 27442dc13f1SHong Zhang (xarr+offsetfrom)[1] = wash->H0; /* 1st H */ 27542dc13f1SHong Zhang } 27642dc13f1SHong Zhang 27742dc13f1SHong Zhang /* if vto is an end vertex: */ 27842dc13f1SHong Zhang ierr = DMNetworkGetComponent(networkdm,vto,0,&vkey,(void**)&junction,NULL);CHKERRQ(ierr); 27942dc13f1SHong Zhang if (junction->type == VALVE) { 28042dc13f1SHong Zhang (xarr+offsetto)[0] = wash->QL; /* last Q */ 28142dc13f1SHong Zhang } 28242dc13f1SHong Zhang ierr = VecRestoreArrayRead(pipe->x,&xarray);CHKERRQ(ierr); 28342dc13f1SHong Zhang } 28442dc13f1SHong Zhang 28542dc13f1SHong Zhang ierr = VecRestoreArray(localX,&xarr);CHKERRQ(ierr); 28642dc13f1SHong Zhang ierr = DMLocalToGlobalBegin(networkdm,localX,ADD_VALUES,X);CHKERRQ(ierr); 28742dc13f1SHong Zhang ierr = DMLocalToGlobalEnd(networkdm,localX,ADD_VALUES,X);CHKERRQ(ierr); 28842dc13f1SHong Zhang ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr); 28942dc13f1SHong Zhang 29042dc13f1SHong Zhang #if 0 29142dc13f1SHong Zhang PetscInt N; 29242dc13f1SHong Zhang ierr = VecGetSize(X,&N);CHKERRQ(ierr); 29342dc13f1SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"initial solution %d:\n",N);CHKERRQ(ierr); 29442dc13f1SHong Zhang ierr = VecView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 29542dc13f1SHong Zhang #endif 29642dc13f1SHong Zhang PetscFunctionReturn(0); 29742dc13f1SHong Zhang } 29842dc13f1SHong Zhang 29942dc13f1SHong Zhang PetscErrorCode TSDMNetworkMonitor(TS ts, PetscInt step, PetscReal t, Vec x, void *context) 30042dc13f1SHong Zhang { 30142dc13f1SHong Zhang PetscErrorCode ierr; 30242dc13f1SHong Zhang DMNetworkMonitor monitor; 30342dc13f1SHong Zhang 30442dc13f1SHong Zhang PetscFunctionBegin; 30542dc13f1SHong Zhang monitor = (DMNetworkMonitor)context; 30642dc13f1SHong Zhang ierr = DMNetworkMonitorView(monitor,x);CHKERRQ(ierr); 30742dc13f1SHong Zhang PetscFunctionReturn(0); 30842dc13f1SHong Zhang } 30942dc13f1SHong Zhang 31042dc13f1SHong Zhang PetscErrorCode PipesView(DM networkdm,PetscInt KeyPipe,Vec X) 31142dc13f1SHong Zhang { 31242dc13f1SHong Zhang PetscErrorCode ierr; 31342dc13f1SHong Zhang PetscInt i,numkeys=1,*blocksize,*numselectedvariable,**selectedvariables,n; 31442dc13f1SHong Zhang IS isfrom_q,isfrom_h,isfrom; 31542dc13f1SHong Zhang Vec Xto; 31642dc13f1SHong Zhang VecScatter ctx; 31742dc13f1SHong Zhang MPI_Comm comm; 31842dc13f1SHong Zhang 31942dc13f1SHong Zhang PetscFunctionBegin; 32042dc13f1SHong Zhang ierr = PetscObjectGetComm((PetscObject)networkdm,&comm);CHKERRQ(ierr); 32142dc13f1SHong Zhang 32242dc13f1SHong Zhang /* 1. Create isfrom_q for q-variable of pipes */ 32342dc13f1SHong Zhang ierr = PetscMalloc3(numkeys,&blocksize,numkeys,&numselectedvariable,numkeys,&selectedvariables);CHKERRQ(ierr); 32442dc13f1SHong Zhang for (i=0; i<numkeys; i++) { 32542dc13f1SHong Zhang blocksize[i] = 2; 32642dc13f1SHong Zhang numselectedvariable[i] = 1; 32742dc13f1SHong Zhang ierr = PetscMalloc1(numselectedvariable[i],&selectedvariables[i]);CHKERRQ(ierr); 32842dc13f1SHong Zhang selectedvariables[i][0] = 0; /* q-variable */ 32942dc13f1SHong Zhang } 33042dc13f1SHong Zhang ierr = DMNetworkCreateIS(networkdm,numkeys,&KeyPipe,blocksize,numselectedvariable,selectedvariables,&isfrom_q);CHKERRQ(ierr); 33142dc13f1SHong Zhang 33242dc13f1SHong Zhang /* 2. Create Xto and isto */ 33342dc13f1SHong Zhang ierr = ISGetLocalSize(isfrom_q, &n);CHKERRQ(ierr); 33442dc13f1SHong Zhang ierr = VecCreate(comm,&Xto);CHKERRQ(ierr); 33542dc13f1SHong Zhang ierr = VecSetSizes(Xto,n,PETSC_DECIDE);CHKERRQ(ierr); 33642dc13f1SHong Zhang ierr = VecSetFromOptions(Xto);CHKERRQ(ierr); 33742dc13f1SHong Zhang ierr = VecSet(Xto,0.0);CHKERRQ(ierr); 33842dc13f1SHong Zhang 33942dc13f1SHong Zhang /* 3. Create scatter */ 34042dc13f1SHong Zhang ierr = VecScatterCreate(X,isfrom_q,Xto,NULL,&ctx);CHKERRQ(ierr); 34142dc13f1SHong Zhang 34242dc13f1SHong Zhang /* 4. Scatter to Xq */ 34342dc13f1SHong Zhang ierr = VecScatterBegin(ctx,X,Xto,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 34442dc13f1SHong Zhang ierr = VecScatterEnd(ctx,X,Xto,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 34542dc13f1SHong Zhang ierr = VecScatterDestroy(&ctx);CHKERRQ(ierr); 34642dc13f1SHong Zhang ierr = ISDestroy(&isfrom_q);CHKERRQ(ierr); 34742dc13f1SHong Zhang 34842dc13f1SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"Xq:\n");CHKERRQ(ierr); 34942dc13f1SHong Zhang ierr = VecView(Xto,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 35042dc13f1SHong Zhang 35142dc13f1SHong Zhang /* 5. Create isfrom_h for h-variable of pipes; Create scatter; Scatter to Xh */ 35242dc13f1SHong Zhang for (i=0; i<numkeys; i++) { 35342dc13f1SHong Zhang selectedvariables[i][0] = 1; /* h-variable */ 35442dc13f1SHong Zhang } 35542dc13f1SHong Zhang ierr = DMNetworkCreateIS(networkdm,numkeys,&KeyPipe,blocksize,numselectedvariable,selectedvariables,&isfrom_h);CHKERRQ(ierr); 35642dc13f1SHong Zhang 35742dc13f1SHong Zhang ierr = VecScatterCreate(X,isfrom_h,Xto,NULL,&ctx);CHKERRQ(ierr); 35842dc13f1SHong Zhang ierr = VecScatterBegin(ctx,X,Xto,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 35942dc13f1SHong Zhang ierr = VecScatterEnd(ctx,X,Xto,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 36042dc13f1SHong Zhang ierr = VecScatterDestroy(&ctx);CHKERRQ(ierr); 36142dc13f1SHong Zhang ierr = ISDestroy(&isfrom_h);CHKERRQ(ierr); 36242dc13f1SHong Zhang 36342dc13f1SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"Xh:\n");CHKERRQ(ierr); 36442dc13f1SHong Zhang ierr = VecView(Xto,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 36542dc13f1SHong Zhang ierr = VecDestroy(&Xto);CHKERRQ(ierr); 36642dc13f1SHong Zhang 36742dc13f1SHong Zhang /* 6. Create isfrom for all pipe variables; Create scatter; Scatter to Xpipes */ 36842dc13f1SHong Zhang for (i=0; i<numkeys; i++) { 36942dc13f1SHong Zhang blocksize[i] = -1; /* select all the variables of the i-th component */ 37042dc13f1SHong Zhang } 37142dc13f1SHong Zhang ierr = DMNetworkCreateIS(networkdm,numkeys,&KeyPipe,blocksize,NULL,NULL,&isfrom);CHKERRQ(ierr); 37242dc13f1SHong Zhang ierr = ISDestroy(&isfrom);CHKERRQ(ierr); 37342dc13f1SHong Zhang ierr = DMNetworkCreateIS(networkdm,numkeys,&KeyPipe,NULL,NULL,NULL,&isfrom);CHKERRQ(ierr); 37442dc13f1SHong Zhang 37542dc13f1SHong Zhang ierr = ISGetLocalSize(isfrom, &n);CHKERRQ(ierr); 37642dc13f1SHong Zhang ierr = VecCreate(comm,&Xto);CHKERRQ(ierr); 37742dc13f1SHong Zhang ierr = VecSetSizes(Xto,n,PETSC_DECIDE);CHKERRQ(ierr); 37842dc13f1SHong Zhang ierr = VecSetFromOptions(Xto);CHKERRQ(ierr); 37942dc13f1SHong Zhang ierr = VecSet(Xto,0.0);CHKERRQ(ierr); 38042dc13f1SHong Zhang 38142dc13f1SHong Zhang ierr = VecScatterCreate(X,isfrom,Xto,NULL,&ctx);CHKERRQ(ierr); 38242dc13f1SHong Zhang ierr = VecScatterBegin(ctx,X,Xto,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 38342dc13f1SHong Zhang ierr = VecScatterEnd(ctx,X,Xto,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 38442dc13f1SHong Zhang ierr = VecScatterDestroy(&ctx);CHKERRQ(ierr); 38542dc13f1SHong Zhang ierr = ISDestroy(&isfrom);CHKERRQ(ierr); 38642dc13f1SHong Zhang 38742dc13f1SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"Xpipes:\n");CHKERRQ(ierr); 38842dc13f1SHong Zhang ierr = VecView(Xto,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 38942dc13f1SHong Zhang 39042dc13f1SHong Zhang /* 7. Free spaces */ 39142dc13f1SHong Zhang for (i=0; i<numkeys; i++) { 39242dc13f1SHong Zhang ierr = PetscFree(selectedvariables[i]);CHKERRQ(ierr); 39342dc13f1SHong Zhang } 39442dc13f1SHong Zhang ierr = PetscFree3(blocksize,numselectedvariable,selectedvariables);CHKERRQ(ierr); 39542dc13f1SHong Zhang ierr = VecDestroy(&Xto);CHKERRQ(ierr); 39642dc13f1SHong Zhang PetscFunctionReturn(0); 39742dc13f1SHong Zhang } 39842dc13f1SHong Zhang 39942dc13f1SHong Zhang PetscErrorCode ISJunctionsView(DM networkdm,PetscInt KeyJunc) 40042dc13f1SHong Zhang { 40142dc13f1SHong Zhang PetscErrorCode ierr; 40242dc13f1SHong Zhang PetscInt numkeys=1; 40342dc13f1SHong Zhang IS isfrom; 40442dc13f1SHong Zhang MPI_Comm comm; 40542dc13f1SHong Zhang PetscMPIInt rank; 40642dc13f1SHong Zhang 40742dc13f1SHong Zhang PetscFunctionBegin; 40842dc13f1SHong Zhang ierr = PetscObjectGetComm((PetscObject)networkdm,&comm);CHKERRQ(ierr); 40942dc13f1SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 41042dc13f1SHong Zhang 41142dc13f1SHong Zhang /* Create a global isfrom for all junction variables */ 41242dc13f1SHong Zhang ierr = DMNetworkCreateIS(networkdm,numkeys,&KeyJunc,NULL,NULL,NULL,&isfrom);CHKERRQ(ierr); 41342dc13f1SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"ISJunctions:\n");CHKERRQ(ierr); 41442dc13f1SHong Zhang ierr = ISView(isfrom,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 41542dc13f1SHong Zhang ierr = ISDestroy(&isfrom);CHKERRQ(ierr); 41642dc13f1SHong Zhang 41742dc13f1SHong Zhang /* Create a local isfrom for all junction variables */ 41842dc13f1SHong Zhang ierr = DMNetworkCreateLocalIS(networkdm,numkeys,&KeyJunc,NULL,NULL,NULL,&isfrom);CHKERRQ(ierr); 41942dc13f1SHong Zhang if (!rank) { 42042dc13f1SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] ISLocalJunctions:\n",rank);CHKERRQ(ierr); 42142dc13f1SHong Zhang ierr = ISView(isfrom,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 42242dc13f1SHong Zhang } 42342dc13f1SHong Zhang ierr = ISDestroy(&isfrom);CHKERRQ(ierr); 42442dc13f1SHong Zhang PetscFunctionReturn(0); 42542dc13f1SHong Zhang } 42642dc13f1SHong Zhang 42742dc13f1SHong Zhang PetscErrorCode WashNetworkCleanUp(Wash wash) 42842dc13f1SHong Zhang { 42942dc13f1SHong Zhang PetscErrorCode ierr; 43042dc13f1SHong Zhang PetscMPIInt rank; 43142dc13f1SHong Zhang 43242dc13f1SHong Zhang PetscFunctionBegin; 43342dc13f1SHong Zhang ierr = MPI_Comm_rank(wash->comm,&rank);CHKERRMPI(ierr); 43442dc13f1SHong Zhang ierr = PetscFree(wash->edgelist);CHKERRQ(ierr); 43542dc13f1SHong Zhang ierr = PetscFree(wash->vtype);CHKERRQ(ierr); 436dd400576SPatrick Sanan if (rank == 0) { 43742dc13f1SHong Zhang ierr = PetscFree2(wash->junction,wash->pipe);CHKERRQ(ierr); 43842dc13f1SHong Zhang } 43942dc13f1SHong Zhang PetscFunctionReturn(0); 44042dc13f1SHong Zhang } 44142dc13f1SHong Zhang 44242dc13f1SHong Zhang PetscErrorCode WashNetworkCreate(MPI_Comm comm,PetscInt pipesCase,Wash *wash_ptr) 44342dc13f1SHong Zhang { 44442dc13f1SHong Zhang PetscErrorCode ierr; 44542dc13f1SHong Zhang PetscInt npipes; 44642dc13f1SHong Zhang PetscMPIInt rank; 44742dc13f1SHong Zhang Wash wash=NULL; 44842dc13f1SHong Zhang PetscInt i,numVertices,numEdges,*vtype; 44942dc13f1SHong Zhang PetscInt *edgelist; 45042dc13f1SHong Zhang Junction junctions=NULL; 45142dc13f1SHong Zhang Pipe pipes=NULL; 45242dc13f1SHong Zhang PetscBool washdist=PETSC_TRUE; 45342dc13f1SHong Zhang 45442dc13f1SHong Zhang PetscFunctionBegin; 45542dc13f1SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 45642dc13f1SHong Zhang 45742dc13f1SHong Zhang ierr = PetscCalloc1(1,&wash);CHKERRQ(ierr); 45842dc13f1SHong Zhang wash->comm = comm; 45942dc13f1SHong Zhang *wash_ptr = wash; 46042dc13f1SHong Zhang wash->Q0 = 0.477432; /* RESERVOIR */ 46142dc13f1SHong Zhang wash->H0 = 150.0; 46242dc13f1SHong Zhang wash->HL = 143.488; /* VALVE */ 46342dc13f1SHong Zhang wash->QL = 0.0; 46442dc13f1SHong Zhang wash->nnodes_loc = 0; 46542dc13f1SHong Zhang 46642dc13f1SHong Zhang numVertices = 0; 46742dc13f1SHong Zhang numEdges = 0; 46842dc13f1SHong Zhang edgelist = NULL; 46942dc13f1SHong Zhang 47042dc13f1SHong Zhang /* proc[0] creates a sequential wash and edgelist */ 47142dc13f1SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"Setup pipesCase %D\n",pipesCase);CHKERRQ(ierr); 47242dc13f1SHong Zhang 47342dc13f1SHong Zhang /* Set global number of pipes, edges, and junctions */ 47442dc13f1SHong Zhang /*-------------------------------------------------*/ 47542dc13f1SHong Zhang switch (pipesCase) { 47642dc13f1SHong Zhang case 0: 47742dc13f1SHong Zhang /* pipeCase 0: */ 47842dc13f1SHong Zhang /* ================================================= 47942dc13f1SHong Zhang (RESERVOIR) v0 --E0--> v1--E1--> v2 --E2-->v3 (VALVE) 48042dc13f1SHong Zhang ==================================================== */ 48142dc13f1SHong Zhang npipes = 3; 48242dc13f1SHong Zhang ierr = PetscOptionsGetInt(NULL,NULL, "-npipes", &npipes, NULL);CHKERRQ(ierr); 48342dc13f1SHong Zhang wash->nedge = npipes; 48442dc13f1SHong Zhang wash->nvertex = npipes + 1; 48542dc13f1SHong Zhang 48642dc13f1SHong Zhang /* Set local edges and vertices -- proc[0] sets entire network, then distributes */ 48742dc13f1SHong Zhang numVertices = 0; 48842dc13f1SHong Zhang numEdges = 0; 48942dc13f1SHong Zhang edgelist = NULL; 490dd400576SPatrick Sanan if (rank == 0) { 49142dc13f1SHong Zhang numVertices = wash->nvertex; 49242dc13f1SHong Zhang numEdges = wash->nedge; 49342dc13f1SHong Zhang 49442dc13f1SHong Zhang ierr = PetscCalloc1(2*numEdges,&edgelist);CHKERRQ(ierr); 49542dc13f1SHong Zhang for (i=0; i<numEdges; i++) { 49642dc13f1SHong Zhang edgelist[2*i] = i; edgelist[2*i+1] = i+1; 49742dc13f1SHong Zhang } 49842dc13f1SHong Zhang 49942dc13f1SHong Zhang /* Add network components */ 50042dc13f1SHong Zhang /*------------------------*/ 50142dc13f1SHong Zhang ierr = PetscCalloc2(numVertices,&junctions,numEdges,&pipes);CHKERRQ(ierr); 50242dc13f1SHong Zhang 50342dc13f1SHong Zhang /* vertex */ 50442dc13f1SHong Zhang for (i=0; i<numVertices; i++) { 50542dc13f1SHong Zhang junctions[i].id = i; 50642dc13f1SHong Zhang junctions[i].type = JUNCTION; 50742dc13f1SHong Zhang } 50842dc13f1SHong Zhang 50942dc13f1SHong Zhang junctions[0].type = RESERVOIR; 51042dc13f1SHong Zhang junctions[numVertices-1].type = VALVE; 51142dc13f1SHong Zhang } 51242dc13f1SHong Zhang break; 51342dc13f1SHong Zhang case 1: 51442dc13f1SHong Zhang /* pipeCase 1: */ 51542dc13f1SHong Zhang /* ========================== 51642dc13f1SHong Zhang v2 (VALVE) 51742dc13f1SHong Zhang ^ 51842dc13f1SHong Zhang | 51942dc13f1SHong Zhang E2 52042dc13f1SHong Zhang | 52142dc13f1SHong Zhang v0 --E0--> v3--E1--> v1 52242dc13f1SHong Zhang (RESERVOIR) (RESERVOIR) 52342dc13f1SHong Zhang ============================= */ 52442dc13f1SHong Zhang npipes = 3; 52542dc13f1SHong Zhang wash->nedge = npipes; 52642dc13f1SHong Zhang wash->nvertex = npipes + 1; 52742dc13f1SHong Zhang 52842dc13f1SHong Zhang /* Set local edges and vertices -- proc[0] sets entire network, then distributes */ 529dd400576SPatrick Sanan if (rank == 0) { 53042dc13f1SHong Zhang numVertices = wash->nvertex; 53142dc13f1SHong Zhang numEdges = wash->nedge; 53242dc13f1SHong Zhang 53342dc13f1SHong Zhang ierr = PetscCalloc1(2*numEdges,&edgelist);CHKERRQ(ierr); 53442dc13f1SHong Zhang edgelist[0] = 0; edgelist[1] = 3; /* edge[0] */ 53542dc13f1SHong Zhang edgelist[2] = 3; edgelist[3] = 1; /* edge[1] */ 53642dc13f1SHong Zhang edgelist[4] = 3; edgelist[5] = 2; /* edge[2] */ 53742dc13f1SHong Zhang 53842dc13f1SHong Zhang /* Add network components */ 53942dc13f1SHong Zhang /*------------------------*/ 54042dc13f1SHong Zhang ierr = PetscCalloc2(numVertices,&junctions,numEdges,&pipes);CHKERRQ(ierr); 54142dc13f1SHong Zhang /* vertex */ 54242dc13f1SHong Zhang for (i=0; i<numVertices; i++) { 54342dc13f1SHong Zhang junctions[i].id = i; 54442dc13f1SHong Zhang junctions[i].type = JUNCTION; 54542dc13f1SHong Zhang } 54642dc13f1SHong Zhang 54742dc13f1SHong Zhang junctions[0].type = RESERVOIR; 54842dc13f1SHong Zhang junctions[1].type = VALVE; 54942dc13f1SHong Zhang junctions[2].type = VALVE; 55042dc13f1SHong Zhang } 55142dc13f1SHong Zhang break; 55242dc13f1SHong Zhang case 2: 55342dc13f1SHong Zhang /* pipeCase 2: */ 55442dc13f1SHong Zhang /* ========================== 55542dc13f1SHong Zhang (RESERVOIR) v2--> E2 55642dc13f1SHong Zhang | 55742dc13f1SHong Zhang v0 --E0--> v3--E1--> v1 55842dc13f1SHong Zhang (RESERVOIR) (VALVE) 55942dc13f1SHong Zhang ============================= */ 56042dc13f1SHong Zhang 56142dc13f1SHong Zhang /* Set application parameters -- to be used in function evalutions */ 56242dc13f1SHong Zhang npipes = 3; 56342dc13f1SHong Zhang wash->nedge = npipes; 56442dc13f1SHong Zhang wash->nvertex = npipes + 1; 56542dc13f1SHong Zhang 56642dc13f1SHong Zhang /* Set local edges and vertices -- proc[0] sets entire network, then distributes */ 567dd400576SPatrick Sanan if (rank == 0) { 56842dc13f1SHong Zhang numVertices = wash->nvertex; 56942dc13f1SHong Zhang numEdges = wash->nedge; 57042dc13f1SHong Zhang 57142dc13f1SHong Zhang ierr = PetscCalloc1(2*numEdges,&edgelist);CHKERRQ(ierr); 57242dc13f1SHong Zhang edgelist[0] = 0; edgelist[1] = 3; /* edge[0] */ 57342dc13f1SHong Zhang edgelist[2] = 3; edgelist[3] = 1; /* edge[1] */ 57442dc13f1SHong Zhang edgelist[4] = 2; edgelist[5] = 3; /* edge[2] */ 57542dc13f1SHong Zhang 57642dc13f1SHong Zhang /* Add network components */ 57742dc13f1SHong Zhang /*------------------------*/ 57842dc13f1SHong Zhang ierr = PetscCalloc2(numVertices,&junctions,numEdges,&pipes);CHKERRQ(ierr); 57942dc13f1SHong Zhang /* vertex */ 58042dc13f1SHong Zhang for (i=0; i<numVertices; i++) { 58142dc13f1SHong Zhang junctions[i].id = i; 58242dc13f1SHong Zhang junctions[i].type = JUNCTION; 58342dc13f1SHong Zhang } 58442dc13f1SHong Zhang 58542dc13f1SHong Zhang junctions[0].type = RESERVOIR; 58642dc13f1SHong Zhang junctions[1].type = VALVE; 58742dc13f1SHong Zhang junctions[2].type = RESERVOIR; 58842dc13f1SHong Zhang } 58942dc13f1SHong Zhang break; 59042dc13f1SHong Zhang default: 59142dc13f1SHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"not done yet"); 59242dc13f1SHong Zhang } 59342dc13f1SHong Zhang 59442dc13f1SHong Zhang /* set edge global id */ 59542dc13f1SHong Zhang for (i=0; i<numEdges; i++) pipes[i].id = i; 59642dc13f1SHong Zhang 597dd400576SPatrick Sanan if (rank == 0) { /* set vtype for proc[0] */ 59842dc13f1SHong Zhang PetscInt v; 59942dc13f1SHong Zhang ierr = PetscMalloc1(2*numEdges,&vtype);CHKERRQ(ierr); 60042dc13f1SHong Zhang for (i=0; i<2*numEdges; i++) { 60142dc13f1SHong Zhang v = edgelist[i]; 60242dc13f1SHong Zhang vtype[i] = junctions[v].type; 60342dc13f1SHong Zhang } 60442dc13f1SHong Zhang wash->vtype = vtype; 60542dc13f1SHong Zhang } 60642dc13f1SHong Zhang 60742dc13f1SHong Zhang *wash_ptr = wash; 60842dc13f1SHong Zhang wash->nedge = numEdges; 60942dc13f1SHong Zhang wash->nvertex = numVertices; 61042dc13f1SHong Zhang wash->edgelist = edgelist; 61142dc13f1SHong Zhang wash->junction = junctions; 61242dc13f1SHong Zhang wash->pipe = pipes; 61342dc13f1SHong Zhang 61442dc13f1SHong Zhang /* Distribute edgelist to other processors */ 61542dc13f1SHong Zhang ierr = PetscOptionsGetBool(NULL,NULL,"-wash_distribute",&washdist,NULL);CHKERRQ(ierr); 61642dc13f1SHong Zhang if (washdist) { 61742dc13f1SHong Zhang /* 61842dc13f1SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD," Distribute sequential wash ...\n");CHKERRQ(ierr); 61942dc13f1SHong Zhang */ 62042dc13f1SHong Zhang ierr = WashNetworkDistribute(comm,wash);CHKERRQ(ierr); 62142dc13f1SHong Zhang } 62242dc13f1SHong Zhang PetscFunctionReturn(0); 62342dc13f1SHong Zhang } 62442dc13f1SHong Zhang 62542dc13f1SHong Zhang /* ------------------------------------------------------- */ 62642dc13f1SHong Zhang int main(int argc,char ** argv) 62742dc13f1SHong Zhang { 62842dc13f1SHong Zhang PetscErrorCode ierr; 62942dc13f1SHong Zhang Wash wash; 63042dc13f1SHong Zhang Junction junctions,junction; 63142dc13f1SHong Zhang Pipe pipe,pipes; 632f11a936eSBarry Smith PetscInt KeyPipe,KeyJunction,*edgelist = NULL,*vtype = NULL; 633f11a936eSBarry Smith PetscInt i,e,v,eStart,eEnd,vStart,vEnd,key,vkey,type; 634f11a936eSBarry Smith PetscInt steps=1,nedges,nnodes=6; 63542dc13f1SHong Zhang const PetscInt *cone; 63642dc13f1SHong Zhang DM networkdm; 63742dc13f1SHong Zhang PetscMPIInt size,rank; 63842dc13f1SHong Zhang PetscReal ftime; 63942dc13f1SHong Zhang Vec X; 64042dc13f1SHong Zhang TS ts; 64142dc13f1SHong Zhang TSConvergedReason reason; 64242dc13f1SHong Zhang PetscBool viewpipes,viewjuncs,monipipes=PETSC_FALSE,userJac=PETSC_TRUE,viewdm=PETSC_FALSE,viewX=PETSC_FALSE; 64342dc13f1SHong Zhang PetscInt pipesCase=0; 64442dc13f1SHong Zhang DMNetworkMonitor monitor; 64542dc13f1SHong Zhang MPI_Comm comm; 64642dc13f1SHong Zhang 64742dc13f1SHong Zhang ierr = PetscInitialize(&argc,&argv,"pOption",help);if (ierr) return ierr; 64842dc13f1SHong Zhang 64942dc13f1SHong Zhang /* Read runtime options */ 65042dc13f1SHong Zhang ierr = PetscOptionsGetInt(NULL,NULL, "-case", &pipesCase, NULL);CHKERRQ(ierr); 65142dc13f1SHong Zhang ierr = PetscOptionsGetBool(NULL,NULL,"-user_Jac",&userJac,NULL);CHKERRQ(ierr); 65242dc13f1SHong Zhang ierr = PetscOptionsGetBool(NULL,NULL,"-pipe_monitor",&monipipes,NULL);CHKERRQ(ierr); 65342dc13f1SHong Zhang ierr = PetscOptionsGetBool(NULL,NULL,"-viewdm",&viewdm,NULL);CHKERRQ(ierr); 65442dc13f1SHong Zhang ierr = PetscOptionsGetBool(NULL,NULL,"-viewX",&viewX,NULL);CHKERRQ(ierr); 65542dc13f1SHong Zhang ierr = PetscOptionsGetInt(NULL,NULL, "-npipenodes", &nnodes, NULL);CHKERRQ(ierr); 65642dc13f1SHong Zhang 65742dc13f1SHong Zhang /* Create networkdm */ 65842dc13f1SHong Zhang /*------------------*/ 65942dc13f1SHong Zhang ierr = DMNetworkCreate(PETSC_COMM_WORLD,&networkdm);CHKERRQ(ierr); 66042dc13f1SHong Zhang ierr = PetscObjectGetComm((PetscObject)networkdm,&comm);CHKERRQ(ierr); 66142dc13f1SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 66242dc13f1SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 66342dc13f1SHong Zhang 66442dc13f1SHong Zhang if (size == 1 && monipipes) { 66542dc13f1SHong Zhang ierr = DMNetworkMonitorCreate(networkdm,&monitor);CHKERRQ(ierr); 66642dc13f1SHong Zhang } 66742dc13f1SHong Zhang 66842dc13f1SHong Zhang /* Register the components in the network */ 66942dc13f1SHong Zhang ierr = DMNetworkRegisterComponent(networkdm,"junctionstruct",sizeof(struct _p_Junction),&KeyJunction);CHKERRQ(ierr); 67042dc13f1SHong Zhang ierr = DMNetworkRegisterComponent(networkdm,"pipestruct",sizeof(struct _p_Pipe),&KeyPipe);CHKERRQ(ierr); 67142dc13f1SHong Zhang 67242dc13f1SHong Zhang /* Create a distributed wash network (user-specific) */ 67342dc13f1SHong Zhang ierr = WashNetworkCreate(comm,pipesCase,&wash);CHKERRQ(ierr); 67442dc13f1SHong Zhang nedges = wash->nedge; 67542dc13f1SHong Zhang edgelist = wash->edgelist; 67642dc13f1SHong Zhang vtype = wash->vtype; 67742dc13f1SHong Zhang junctions = wash->junction; 67842dc13f1SHong Zhang pipes = wash->pipe; 67942dc13f1SHong Zhang 68042dc13f1SHong Zhang /* Set up the network layout */ 68142dc13f1SHong Zhang ierr = DMNetworkSetNumSubNetworks(networkdm,PETSC_DECIDE,1);CHKERRQ(ierr); 682f11a936eSBarry Smith ierr = DMNetworkAddSubnetwork(networkdm,NULL,nedges,edgelist,NULL);CHKERRQ(ierr); 68342dc13f1SHong Zhang 68442dc13f1SHong Zhang ierr = DMNetworkLayoutSetUp(networkdm);CHKERRQ(ierr); 68542dc13f1SHong Zhang 68642dc13f1SHong Zhang ierr = DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);CHKERRQ(ierr); 68742dc13f1SHong Zhang ierr = DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);CHKERRQ(ierr); 68842dc13f1SHong Zhang /* ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] eStart/End: %d - %d; vStart/End: %d - %d\n",rank,eStart,eEnd,vStart,vEnd);CHKERRQ(ierr); */ 68942dc13f1SHong Zhang 69042dc13f1SHong Zhang if (rank) { /* junctions[] and pipes[] for proc[0] are allocated in WashNetworkCreate() */ 69142dc13f1SHong Zhang /* vEnd - vStart = nvertices + number of ghost vertices! */ 69242dc13f1SHong Zhang ierr = PetscCalloc2(vEnd - vStart,&junctions,nedges,&pipes);CHKERRQ(ierr); 69342dc13f1SHong Zhang } 69442dc13f1SHong Zhang 69542dc13f1SHong Zhang /* Add Pipe component and number of variables to all local edges */ 69642dc13f1SHong Zhang for (e = eStart; e < eEnd; e++) { 69742dc13f1SHong Zhang pipes[e-eStart].nnodes = nnodes; 69842dc13f1SHong Zhang ierr = DMNetworkAddComponent(networkdm,e,KeyPipe,&pipes[e-eStart],2*pipes[e-eStart].nnodes);CHKERRQ(ierr); 69942dc13f1SHong Zhang 70042dc13f1SHong Zhang if (size == 1 && monipipes) { /* Add monitor -- show Q_{pipes[e-eStart].id}? */ 70142dc13f1SHong Zhang pipes[e-eStart].length = 600.0; 70242dc13f1SHong Zhang 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); 70342dc13f1SHong Zhang 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); 70442dc13f1SHong Zhang } 70542dc13f1SHong Zhang } 70642dc13f1SHong Zhang 707eac198afSGetnet /* Add Junction component and number of variables to all local vertices */ 70842dc13f1SHong Zhang for (v = vStart; v < vEnd; v++) { 70942dc13f1SHong Zhang ierr = DMNetworkAddComponent(networkdm,v,KeyJunction,&junctions[v-vStart],2);CHKERRQ(ierr); 71042dc13f1SHong Zhang } 71142dc13f1SHong Zhang 71242dc13f1SHong Zhang if (size > 1) { /* must be called before DMSetUp()???. Other partitioners do not work yet??? -- cause crash in proc[0]! */ 71342dc13f1SHong Zhang DM plexdm; 71442dc13f1SHong Zhang PetscPartitioner part; 71542dc13f1SHong Zhang ierr = DMNetworkGetPlex(networkdm,&plexdm);CHKERRQ(ierr); 71642dc13f1SHong Zhang ierr = DMPlexGetPartitioner(plexdm, &part);CHKERRQ(ierr); 71742dc13f1SHong Zhang ierr = PetscPartitionerSetType(part,PETSCPARTITIONERSIMPLE);CHKERRQ(ierr); 7185a107427SMatthew G. Knepley ierr = PetscOptionsSetValue(NULL,"-dm_plex_csr_alg","mat");CHKERRQ(ierr); /* for parmetis */ 71942dc13f1SHong Zhang } 72042dc13f1SHong Zhang 72142dc13f1SHong Zhang /* Set up DM for use */ 72242dc13f1SHong Zhang ierr = DMSetUp(networkdm);CHKERRQ(ierr); 72342dc13f1SHong Zhang if (viewdm) { 72442dc13f1SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\nOriginal networkdm, DMView:\n");CHKERRQ(ierr); 72542dc13f1SHong Zhang ierr = DMView(networkdm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 72642dc13f1SHong Zhang } 72742dc13f1SHong Zhang 72842dc13f1SHong Zhang /* Set user physical parameters to the components */ 72942dc13f1SHong Zhang for (e = eStart; e < eEnd; e++) { 73042dc13f1SHong Zhang ierr = DMNetworkGetConnectedVertices(networkdm,e,&cone);CHKERRQ(ierr); 73142dc13f1SHong Zhang /* vfrom */ 73242dc13f1SHong Zhang ierr = DMNetworkGetComponent(networkdm,cone[0],0,&vkey,(void**)&junction,NULL);CHKERRQ(ierr); 73342dc13f1SHong Zhang junction->type = (VertexType)vtype[2*e]; 73442dc13f1SHong Zhang 73542dc13f1SHong Zhang /* vto */ 73642dc13f1SHong Zhang ierr = DMNetworkGetComponent(networkdm,cone[1],0,&vkey,(void**)&junction,NULL);CHKERRQ(ierr); 73742dc13f1SHong Zhang junction->type = (VertexType)vtype[2*e+1]; 73842dc13f1SHong Zhang } 73942dc13f1SHong Zhang 74042dc13f1SHong Zhang ierr = WashNetworkCleanUp(wash);CHKERRQ(ierr); 74142dc13f1SHong Zhang 74242dc13f1SHong Zhang /* Network partitioning and distribution of data */ 74342dc13f1SHong Zhang ierr = DMNetworkDistribute(&networkdm,0);CHKERRQ(ierr); 74442dc13f1SHong Zhang if (viewdm) { 74542dc13f1SHong Zhang PetscPrintf(PETSC_COMM_WORLD,"\nAfter DMNetworkDistribute, DMView:\n");CHKERRQ(ierr); 74642dc13f1SHong Zhang ierr = DMView(networkdm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 74742dc13f1SHong Zhang } 74842dc13f1SHong Zhang 74942dc13f1SHong Zhang /* create vectors */ 75042dc13f1SHong Zhang ierr = DMCreateGlobalVector(networkdm,&X);CHKERRQ(ierr); 75142dc13f1SHong Zhang ierr = DMCreateLocalVector(networkdm,&wash->localX);CHKERRQ(ierr); 75242dc13f1SHong Zhang ierr = DMCreateLocalVector(networkdm,&wash->localXdot);CHKERRQ(ierr); 75342dc13f1SHong Zhang 75442dc13f1SHong Zhang /* PipeSetUp -- each process only sets its own pipes */ 75542dc13f1SHong Zhang /*---------------------------------------------------*/ 75642dc13f1SHong Zhang ierr = DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);CHKERRQ(ierr); 75742dc13f1SHong Zhang 75842dc13f1SHong Zhang userJac = PETSC_TRUE; 75942dc13f1SHong Zhang ierr = DMNetworkHasJacobian(networkdm,userJac,userJac);CHKERRQ(ierr); 76042dc13f1SHong Zhang ierr = DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);CHKERRQ(ierr); 76142dc13f1SHong Zhang for (e=eStart; e<eEnd; e++) { /* each edge has only one component, pipe */ 76242dc13f1SHong Zhang ierr = DMNetworkGetComponent(networkdm,e,0,&type,(void**)&pipe,NULL);CHKERRQ(ierr); 76342dc13f1SHong Zhang 76442dc13f1SHong Zhang wash->nnodes_loc += pipe->nnodes; /* local total number of nodes, will be used by PipesView() */ 76542dc13f1SHong Zhang ierr = PipeSetParameters(pipe, 76642dc13f1SHong Zhang 600.0, /* length */ 76742dc13f1SHong Zhang 0.5, /* diameter */ 76842dc13f1SHong Zhang 1200.0, /* a */ 76942dc13f1SHong Zhang 0.018);CHKERRQ(ierr); /* friction */ 77042dc13f1SHong Zhang ierr = PipeSetUp(pipe);CHKERRQ(ierr); 77142dc13f1SHong Zhang 77242dc13f1SHong Zhang if (userJac) { 77342dc13f1SHong Zhang /* Create Jacobian matrix structures for a Pipe */ 77442dc13f1SHong Zhang Mat *J; 77542dc13f1SHong Zhang ierr = PipeCreateJacobian(pipe,NULL,&J);CHKERRQ(ierr); 77642dc13f1SHong Zhang ierr = DMNetworkEdgeSetMatrix(networkdm,e,J);CHKERRQ(ierr); 77742dc13f1SHong Zhang } 77842dc13f1SHong Zhang } 77942dc13f1SHong Zhang 78042dc13f1SHong Zhang if (userJac) { 78142dc13f1SHong Zhang ierr = DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);CHKERRQ(ierr); 78242dc13f1SHong Zhang for (v=vStart; v<vEnd; v++) { 78342dc13f1SHong Zhang Mat *J; 78442dc13f1SHong Zhang ierr = JunctionCreateJacobian(networkdm,v,NULL,&J);CHKERRQ(ierr); 78542dc13f1SHong Zhang ierr = DMNetworkVertexSetMatrix(networkdm,v,J);CHKERRQ(ierr); 78642dc13f1SHong Zhang 78742dc13f1SHong Zhang ierr = DMNetworkGetComponent(networkdm,v,0,&vkey,(void**)&junction,NULL);CHKERRQ(ierr); 78842dc13f1SHong Zhang junction->jacobian = J; 78942dc13f1SHong Zhang } 79042dc13f1SHong Zhang } 79142dc13f1SHong Zhang 79242dc13f1SHong Zhang /* Setup solver */ 79342dc13f1SHong Zhang /*--------------------------------------------------------*/ 79442dc13f1SHong Zhang ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 79542dc13f1SHong Zhang 79642dc13f1SHong Zhang ierr = TSSetDM(ts,(DM)networkdm);CHKERRQ(ierr); 79742dc13f1SHong Zhang ierr = TSSetIFunction(ts,NULL,WASHIFunction,wash);CHKERRQ(ierr); 79842dc13f1SHong Zhang 79942dc13f1SHong Zhang ierr = TSSetMaxSteps(ts,steps);CHKERRQ(ierr); 80042dc13f1SHong Zhang ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 80142dc13f1SHong Zhang ierr = TSSetTimeStep(ts,0.1);CHKERRQ(ierr); 80242dc13f1SHong Zhang ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr); 80342dc13f1SHong Zhang if (size == 1 && monipipes) { 80442dc13f1SHong Zhang ierr = TSMonitorSet(ts, TSDMNetworkMonitor, monitor, NULL);CHKERRQ(ierr); 80542dc13f1SHong Zhang } 80642dc13f1SHong Zhang ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 80742dc13f1SHong Zhang 80842dc13f1SHong Zhang ierr = WASHSetInitialSolution(networkdm,X,wash);CHKERRQ(ierr); 80942dc13f1SHong Zhang 81042dc13f1SHong Zhang ierr = TSSolve(ts,X);CHKERRQ(ierr); 81142dc13f1SHong Zhang 81242dc13f1SHong Zhang ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 81342dc13f1SHong Zhang ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 81442dc13f1SHong Zhang ierr = TSGetConvergedReason(ts,&reason);CHKERRQ(ierr); 81542dc13f1SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps);CHKERRQ(ierr); 81642dc13f1SHong Zhang if (viewX) { 81742dc13f1SHong Zhang ierr = VecView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 81842dc13f1SHong Zhang } 81942dc13f1SHong Zhang 82042dc13f1SHong Zhang viewpipes = PETSC_FALSE; 82142dc13f1SHong Zhang ierr = PetscOptionsGetBool(NULL,NULL, "-Jac_view", &viewpipes,NULL);CHKERRQ(ierr); 82242dc13f1SHong Zhang if (viewpipes) { 82342dc13f1SHong Zhang SNES snes; 82442dc13f1SHong Zhang Mat Jac; 82542dc13f1SHong Zhang ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 82642dc13f1SHong Zhang ierr = SNESGetJacobian(snes,&Jac,NULL,NULL,NULL);CHKERRQ(ierr); 82742dc13f1SHong Zhang ierr = MatView(Jac,PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr); 82842dc13f1SHong Zhang } 82942dc13f1SHong Zhang 83042dc13f1SHong Zhang /* View solutions */ 83142dc13f1SHong Zhang /* -------------- */ 83242dc13f1SHong Zhang viewpipes = PETSC_FALSE; 83342dc13f1SHong Zhang ierr = PetscOptionsGetBool(NULL,NULL, "-pipe_view", &viewpipes,NULL);CHKERRQ(ierr); 83442dc13f1SHong Zhang if (viewpipes) { 83542dc13f1SHong Zhang ierr = PipesView(networkdm,KeyPipe,X);CHKERRQ(ierr); 83642dc13f1SHong Zhang } 83742dc13f1SHong Zhang 83842dc13f1SHong Zhang /* Test IS */ 83942dc13f1SHong Zhang viewjuncs = PETSC_FALSE; 84042dc13f1SHong Zhang ierr = PetscOptionsGetBool(NULL,NULL, "-isjunc_view", &viewjuncs,NULL);CHKERRQ(ierr); 84142dc13f1SHong Zhang if (viewjuncs) { 84242dc13f1SHong Zhang ierr = ISJunctionsView(networkdm,KeyJunction);CHKERRQ(ierr); 84342dc13f1SHong Zhang } 84442dc13f1SHong Zhang 84542dc13f1SHong Zhang /* Free spaces */ 84642dc13f1SHong Zhang /* ----------- */ 84742dc13f1SHong Zhang ierr = TSDestroy(&ts);CHKERRQ(ierr); 84842dc13f1SHong Zhang ierr = VecDestroy(&X);CHKERRQ(ierr); 84942dc13f1SHong Zhang ierr = VecDestroy(&wash->localX);CHKERRQ(ierr); 85042dc13f1SHong Zhang ierr = VecDestroy(&wash->localXdot);CHKERRQ(ierr); 85142dc13f1SHong Zhang 85242dc13f1SHong Zhang /* Destroy objects from each pipe that are created in PipeSetUp() */ 85342dc13f1SHong Zhang ierr = DMNetworkGetEdgeRange(networkdm,&eStart, &eEnd);CHKERRQ(ierr); 85442dc13f1SHong Zhang for (i = eStart; i < eEnd; i++) { 85542dc13f1SHong Zhang ierr = DMNetworkGetComponent(networkdm,i,0,&key,(void**)&pipe,NULL);CHKERRQ(ierr); 85642dc13f1SHong Zhang ierr = PipeDestroy(&pipe);CHKERRQ(ierr); 85742dc13f1SHong Zhang } 85842dc13f1SHong Zhang if (userJac) { 85942dc13f1SHong Zhang for (v=vStart; v<vEnd; v++) { 86042dc13f1SHong Zhang ierr = DMNetworkGetComponent(networkdm,v,0,&vkey,(void**)&junction,NULL);CHKERRQ(ierr); 86142dc13f1SHong Zhang ierr = JunctionDestroyJacobian(networkdm,v,junction);CHKERRQ(ierr); 86242dc13f1SHong Zhang } 86342dc13f1SHong Zhang } 86442dc13f1SHong Zhang 86542dc13f1SHong Zhang if (size == 1 && monipipes) { 86642dc13f1SHong Zhang ierr = DMNetworkMonitorDestroy(&monitor);CHKERRQ(ierr); 86742dc13f1SHong Zhang } 86842dc13f1SHong Zhang ierr = DMDestroy(&networkdm);CHKERRQ(ierr); 86942dc13f1SHong Zhang ierr = PetscFree(wash);CHKERRQ(ierr); 87042dc13f1SHong Zhang 87142dc13f1SHong Zhang if (rank) { 87242dc13f1SHong Zhang ierr = PetscFree2(junctions,pipes);CHKERRQ(ierr); 87342dc13f1SHong Zhang } 87442dc13f1SHong Zhang ierr = PetscFinalize(); 87542dc13f1SHong Zhang return ierr; 87642dc13f1SHong Zhang } 87742dc13f1SHong Zhang 87842dc13f1SHong Zhang /*TEST 87942dc13f1SHong Zhang 88042dc13f1SHong Zhang build: 88142dc13f1SHong Zhang depends: pipeInterface.c pipeImpls.c 88242dc13f1SHong Zhang requires: mumps 88342dc13f1SHong Zhang 88442dc13f1SHong Zhang test: 88542dc13f1SHong Zhang args: -ts_monitor -case 1 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -options_left no -viewX 88642dc13f1SHong Zhang localrunfiles: pOption 88742dc13f1SHong Zhang output_file: output/pipes_1.out 88842dc13f1SHong Zhang 88942dc13f1SHong Zhang test: 89042dc13f1SHong Zhang suffix: 2 89142dc13f1SHong Zhang nsize: 2 89242dc13f1SHong Zhang args: -ts_monitor -case 1 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -petscpartitioner_type simple -options_left no -viewX 89342dc13f1SHong Zhang localrunfiles: pOption 89442dc13f1SHong Zhang output_file: output/pipes_2.out 89542dc13f1SHong Zhang 89642dc13f1SHong Zhang test: 89742dc13f1SHong Zhang suffix: 3 89842dc13f1SHong Zhang nsize: 2 89942dc13f1SHong Zhang args: -ts_monitor -case 0 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -petscpartitioner_type simple -options_left no -viewX 90042dc13f1SHong Zhang localrunfiles: pOption 90142dc13f1SHong Zhang output_file: output/pipes_3.out 90242dc13f1SHong Zhang 90342dc13f1SHong Zhang test: 90442dc13f1SHong Zhang suffix: 4 90542dc13f1SHong Zhang args: -ts_monitor -case 2 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -options_left no -viewX 90642dc13f1SHong Zhang localrunfiles: pOption 90742dc13f1SHong Zhang output_file: output/pipes_4.out 90842dc13f1SHong Zhang 90942dc13f1SHong Zhang test: 91042dc13f1SHong Zhang suffix: 5 91142dc13f1SHong Zhang nsize: 3 91242dc13f1SHong Zhang args: -ts_monitor -case 2 -ts_max_steps 10 -pc_factor_mat_solver_type mumps -petscpartitioner_type simple -options_left no -viewX 91342dc13f1SHong Zhang localrunfiles: pOption 91442dc13f1SHong Zhang output_file: output/pipes_5.out 91542dc13f1SHong Zhang 91642dc13f1SHong Zhang test: 91742dc13f1SHong Zhang suffix: 6 91842dc13f1SHong Zhang nsize: 2 91942dc13f1SHong Zhang args: -ts_monitor -case 1 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -petscpartitioner_type simple -options_left no -wash_distribute 0 -viewX 92042dc13f1SHong Zhang localrunfiles: pOption 92142dc13f1SHong Zhang output_file: output/pipes_6.out 92242dc13f1SHong Zhang 92342dc13f1SHong Zhang test: 92442dc13f1SHong Zhang suffix: 7 92542dc13f1SHong Zhang nsize: 2 92642dc13f1SHong Zhang args: -ts_monitor -case 2 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -petscpartitioner_type simple -options_left no -wash_distribute 0 -viewX 92742dc13f1SHong Zhang localrunfiles: pOption 92842dc13f1SHong Zhang output_file: output/pipes_7.out 92942dc13f1SHong Zhang 93042dc13f1SHong Zhang test: 93142dc13f1SHong Zhang suffix: 8 93242dc13f1SHong Zhang nsize: 2 93342dc13f1SHong Zhang requires: parmetis 93442dc13f1SHong Zhang args: -ts_monitor -case 2 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -petscpartitioner_type parmetis -options_left no -wash_distribute 1 93542dc13f1SHong Zhang localrunfiles: pOption 93642dc13f1SHong Zhang output_file: output/pipes_8.out 93742dc13f1SHong Zhang 93842dc13f1SHong Zhang test: 93942dc13f1SHong Zhang suffix: 9 94042dc13f1SHong Zhang nsize: 2 94142dc13f1SHong Zhang args: -case 0 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -petscpartitioner_type simple -options_left no -wash_distribute 0 -pipe_view 94242dc13f1SHong Zhang localrunfiles: pOption 94342dc13f1SHong Zhang output_file: output/pipes_9.out 94442dc13f1SHong Zhang 94542dc13f1SHong Zhang test: 94642dc13f1SHong Zhang suffix: 10 94742dc13f1SHong Zhang nsize: 2 94842dc13f1SHong Zhang args: -case 0 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -petscpartitioner_type simple -options_left no -wash_distribute 0 -isjunc_view 94942dc13f1SHong Zhang localrunfiles: pOption 95042dc13f1SHong Zhang output_file: output/pipes_10.out 95142dc13f1SHong Zhang 95242dc13f1SHong Zhang TEST*/ 953