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 */ 19d71ae5a4SJacob Faibussowitsch PetscErrorCode WashNetworkDistribute(MPI_Comm comm, Wash wash) 20d71ae5a4SJacob Faibussowitsch { 2142dc13f1SHong Zhang PetscMPIInt rank, size, tag = 0; 2242dc13f1SHong Zhang PetscInt i, e, v, numEdges, numVertices, nedges, *eowners = NULL, estart, eend, *vtype = NULL, nvertices; 2342dc13f1SHong Zhang PetscInt *edgelist = wash->edgelist, *nvtx = NULL, *vtxDone = NULL; 2442dc13f1SHong Zhang 2542dc13f1SHong Zhang PetscFunctionBegin; 269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 273ba16761SJacob Faibussowitsch if (size == 1) PetscFunctionReturn(PETSC_SUCCESS); 2842dc13f1SHong Zhang 299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 3042dc13f1SHong Zhang numEdges = wash->nedge; 3142dc13f1SHong Zhang numVertices = wash->nvertex; 3242dc13f1SHong Zhang 3342dc13f1SHong Zhang /* (1) all processes get global and local number of edges */ 349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&numEdges, 1, MPIU_INT, 0, comm)); 3542dc13f1SHong Zhang nedges = numEdges / size; /* local nedges */ 36ad540459SPierre Jolivet if (rank == 0) nedges += numEdges - size * (numEdges / size); 3742dc13f1SHong Zhang wash->Nedge = numEdges; 3842dc13f1SHong Zhang wash->nedge = nedges; 399566063dSJacob Faibussowitsch /* PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d] nedges %d, numEdges %d\n",rank,nedges,numEdges)); */ 4042dc13f1SHong Zhang 419566063dSJacob Faibussowitsch PetscCall(PetscCalloc3(size + 1, &eowners, size, &nvtx, numVertices, &vtxDone)); 429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nedges, 1, MPIU_INT, eowners + 1, 1, MPIU_INT, PETSC_COMM_WORLD)); 4342dc13f1SHong Zhang eowners[0] = 0; 44ad540459SPierre Jolivet for (i = 2; i <= size; i++) eowners[i] += eowners[i - 1]; 4542dc13f1SHong Zhang 4642dc13f1SHong Zhang estart = eowners[rank]; 4742dc13f1SHong Zhang eend = eowners[rank + 1]; 489566063dSJacob Faibussowitsch /* PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d] own lists row %d - %d\n",rank,estart,eend)); */ 4942dc13f1SHong Zhang 5042dc13f1SHong Zhang /* (2) distribute row block edgelist to all processors */ 51dd400576SPatrick Sanan if (rank == 0) { 5242dc13f1SHong Zhang vtype = wash->vtype; 5342dc13f1SHong Zhang for (i = 1; i < size; i++) { 5442dc13f1SHong Zhang /* proc[0] sends edgelist to proc[i] */ 559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(edgelist + 2 * eowners[i], 2 * (eowners[i + 1] - eowners[i]), MPIU_INT, i, tag, comm)); 5642dc13f1SHong Zhang 5742dc13f1SHong Zhang /* proc[0] sends vtype to proc[i] */ 589566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send(vtype + 2 * eowners[i], 2 * (eowners[i + 1] - eowners[i]), MPIU_INT, i, tag, comm)); 5942dc13f1SHong Zhang } 6042dc13f1SHong Zhang } else { 6142dc13f1SHong Zhang MPI_Status status; 629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * (eend - estart), &vtype)); 639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * (eend - estart), &edgelist)); 6442dc13f1SHong Zhang 659566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(edgelist, 2 * (eend - estart), MPIU_INT, 0, tag, comm, &status)); 669566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(vtype, 2 * (eend - estart), MPIU_INT, 0, tag, comm, &status)); 6742dc13f1SHong Zhang } 6842dc13f1SHong Zhang 6942dc13f1SHong Zhang wash->edgelist = edgelist; 7042dc13f1SHong Zhang 7142dc13f1SHong Zhang /* (3) all processes get global and local number of vertices, without ghost vertices */ 72dd400576SPatrick Sanan if (rank == 0) { 7342dc13f1SHong Zhang for (i = 0; i < size; i++) { 7442dc13f1SHong Zhang for (e = eowners[i]; e < eowners[i + 1]; e++) { 7542dc13f1SHong Zhang v = edgelist[2 * e]; 7642dc13f1SHong Zhang if (!vtxDone[v]) { 779371c9d4SSatish Balay nvtx[i]++; 789371c9d4SSatish Balay vtxDone[v] = 1; 7942dc13f1SHong Zhang } 8042dc13f1SHong Zhang v = edgelist[2 * e + 1]; 8142dc13f1SHong Zhang if (!vtxDone[v]) { 829371c9d4SSatish Balay nvtx[i]++; 839371c9d4SSatish Balay vtxDone[v] = 1; 8442dc13f1SHong Zhang } 8542dc13f1SHong Zhang } 8642dc13f1SHong Zhang } 8742dc13f1SHong Zhang } 889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&numVertices, 1, MPIU_INT, 0, PETSC_COMM_WORLD)); 899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scatter(nvtx, 1, MPIU_INT, &nvertices, 1, MPIU_INT, 0, PETSC_COMM_WORLD)); 909566063dSJacob Faibussowitsch PetscCall(PetscFree3(eowners, nvtx, vtxDone)); 9142dc13f1SHong Zhang 9242dc13f1SHong Zhang wash->Nvertex = numVertices; 9342dc13f1SHong Zhang wash->nvertex = nvertices; 9442dc13f1SHong Zhang wash->vtype = vtype; 953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9642dc13f1SHong Zhang } 9742dc13f1SHong Zhang 98d71ae5a4SJacob Faibussowitsch PetscErrorCode WASHIFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) 99d71ae5a4SJacob Faibussowitsch { 10042dc13f1SHong Zhang Wash wash = (Wash)ctx; 10142dc13f1SHong Zhang DM networkdm; 10242dc13f1SHong Zhang Vec localX, localXdot, localF, localXold; 10342dc13f1SHong Zhang const PetscInt *cone; 10442dc13f1SHong Zhang PetscInt vfrom, vto, offsetfrom, offsetto, varoffset; 10542dc13f1SHong Zhang PetscInt v, vStart, vEnd, e, eStart, eEnd; 10642dc13f1SHong Zhang PetscInt nend, type; 10742dc13f1SHong Zhang PetscBool ghost; 10842dc13f1SHong Zhang PetscScalar *farr, *juncf, *pipef; 10942dc13f1SHong Zhang PetscReal dt; 11042dc13f1SHong Zhang Pipe pipe; 11142dc13f1SHong Zhang PipeField *pipex, *pipexdot, *juncx; 11242dc13f1SHong Zhang Junction junction; 11342dc13f1SHong Zhang DMDALocalInfo info; 11442dc13f1SHong Zhang const PetscScalar *xarr, *xdotarr, *xoldarr; 11542dc13f1SHong Zhang 11642dc13f1SHong Zhang PetscFunctionBegin; 11742dc13f1SHong Zhang localX = wash->localX; 11842dc13f1SHong Zhang localXdot = wash->localXdot; 11942dc13f1SHong Zhang 1209566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &localXold)); 1219566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &networkdm)); 1229566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 1239566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(networkdm, &localF)); 12442dc13f1SHong Zhang 12542dc13f1SHong Zhang /* Set F and localF as zero */ 1269566063dSJacob Faibussowitsch PetscCall(VecSet(F, 0.0)); 1279566063dSJacob Faibussowitsch PetscCall(VecSet(localF, 0.0)); 12842dc13f1SHong Zhang 12942dc13f1SHong Zhang /* Update ghost values of locaX and locaXdot */ 1309566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX)); 1319566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX)); 13242dc13f1SHong Zhang 1339566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(networkdm, Xdot, INSERT_VALUES, localXdot)); 1349566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(networkdm, Xdot, INSERT_VALUES, localXdot)); 13542dc13f1SHong Zhang 1369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(localX, &xarr)); 1379566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(localXdot, &xdotarr)); 1389566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(localXold, &xoldarr)); 1399566063dSJacob Faibussowitsch PetscCall(VecGetArray(localF, &farr)); 14042dc13f1SHong Zhang 14142dc13f1SHong Zhang /* junction->type == JUNCTION: 14242dc13f1SHong Zhang juncf[0] = -qJ + sum(qin); juncf[1] = qJ - sum(qout) 14342dc13f1SHong Zhang junction->type == RESERVOIR (upper stream): 14442dc13f1SHong Zhang juncf[0] = -hJ + H0; juncf[1] = qJ - sum(qout) 14542dc13f1SHong Zhang junction->type == VALVE (down stream): 14642dc13f1SHong Zhang juncf[0] = -qJ + sum(qin); juncf[1] = qJ 14742dc13f1SHong Zhang */ 14842dc13f1SHong Zhang /* Vertex/junction initialization */ 1499566063dSJacob Faibussowitsch PetscCall(DMNetworkGetVertexRange(networkdm, &vStart, &vEnd)); 15042dc13f1SHong Zhang for (v = vStart; v < vEnd; v++) { 1519566063dSJacob Faibussowitsch PetscCall(DMNetworkIsGhostVertex(networkdm, v, &ghost)); 15242dc13f1SHong Zhang if (ghost) continue; 15342dc13f1SHong Zhang 1549566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, v, 0, &type, (void **)&junction, NULL)); 1559566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, ALL_COMPONENTS, &varoffset)); 15642dc13f1SHong Zhang juncx = (PipeField *)(xarr + varoffset); 15742dc13f1SHong Zhang juncf = (PetscScalar *)(farr + varoffset); 15842dc13f1SHong Zhang 15942dc13f1SHong Zhang juncf[0] = -juncx[0].q; 16042dc13f1SHong Zhang juncf[1] = juncx[0].q; 16142dc13f1SHong Zhang 16242dc13f1SHong Zhang if (junction->type == RESERVOIR) { /* upstream reservoir */ 16342dc13f1SHong Zhang juncf[0] = juncx[0].h - wash->H0; 16442dc13f1SHong Zhang } 16542dc13f1SHong Zhang } 16642dc13f1SHong Zhang 16742dc13f1SHong Zhang /* Edge/pipe */ 1689566063dSJacob Faibussowitsch PetscCall(DMNetworkGetEdgeRange(networkdm, &eStart, &eEnd)); 16942dc13f1SHong Zhang for (e = eStart; e < eEnd; e++) { 1709566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, e, 0, &type, (void **)&pipe, NULL)); 1719566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, e, ALL_COMPONENTS, &varoffset)); 17242dc13f1SHong Zhang pipex = (PipeField *)(xarr + varoffset); 17342dc13f1SHong Zhang pipexdot = (PipeField *)(xdotarr + varoffset); 17442dc13f1SHong Zhang pipef = (PetscScalar *)(farr + varoffset); 17542dc13f1SHong Zhang 17642dc13f1SHong Zhang /* Get some data into the pipe structure: note, some of these operations 17742dc13f1SHong Zhang * might be redundant. Will it consume too much time? */ 17842dc13f1SHong Zhang pipe->dt = dt; 17942dc13f1SHong Zhang pipe->xold = (PipeField *)(xoldarr + varoffset); 18042dc13f1SHong Zhang 18142dc13f1SHong Zhang /* Evaluate F over this edge/pipe: pipef[1], ...,pipef[2*nend] */ 1829566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(pipe->da, &info)); 1839566063dSJacob Faibussowitsch PetscCall(PipeIFunctionLocal_Lax(&info, t, pipex, pipexdot, pipef, pipe)); 18442dc13f1SHong Zhang 18542dc13f1SHong Zhang /* Get boundary values from connected vertices */ 1869566063dSJacob Faibussowitsch PetscCall(DMNetworkGetConnectedVertices(networkdm, e, &cone)); 18742dc13f1SHong Zhang vfrom = cone[0]; /* local ordering */ 18842dc13f1SHong Zhang vto = cone[1]; 1899566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vfrom, ALL_COMPONENTS, &offsetfrom)); 1909566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vto, ALL_COMPONENTS, &offsetto)); 19142dc13f1SHong Zhang 19242dc13f1SHong Zhang /* Evaluate upstream boundary */ 1939566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vfrom, 0, &type, (void **)&junction, NULL)); 194cad9d221SBarry Smith PetscCheck(junction->type == JUNCTION || junction->type == RESERVOIR, PETSC_COMM_SELF, PETSC_ERR_SUP, "junction type is not supported"); 19542dc13f1SHong Zhang juncx = (PipeField *)(xarr + offsetfrom); 19642dc13f1SHong Zhang juncf = (PetscScalar *)(farr + offsetfrom); 19742dc13f1SHong Zhang 19842dc13f1SHong Zhang pipef[0] = pipex[0].h - juncx[0].h; 19942dc13f1SHong Zhang juncf[1] -= pipex[0].q; 20042dc13f1SHong Zhang 20142dc13f1SHong Zhang /* Evaluate downstream boundary */ 2029566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vto, 0, &type, (void **)&junction, NULL)); 203cad9d221SBarry Smith PetscCheck(junction->type == JUNCTION || junction->type == VALVE, PETSC_COMM_SELF, PETSC_ERR_SUP, "junction type is not supported"); 20442dc13f1SHong Zhang juncx = (PipeField *)(xarr + offsetto); 20542dc13f1SHong Zhang juncf = (PetscScalar *)(farr + offsetto); 20642dc13f1SHong Zhang nend = pipe->nnodes - 1; 20742dc13f1SHong Zhang 20842dc13f1SHong Zhang pipef[2 * nend + 1] = pipex[nend].h - juncx[0].h; 20942dc13f1SHong Zhang juncf[0] += pipex[nend].q; 21042dc13f1SHong Zhang } 21142dc13f1SHong Zhang 2129566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(localX, &xarr)); 2139566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(localXdot, &xdotarr)); 2149566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localF, &farr)); 21542dc13f1SHong Zhang 2169566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(networkdm, localF, ADD_VALUES, F)); 2179566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(networkdm, localF, ADD_VALUES, F)); 2189566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(networkdm, &localF)); 21942dc13f1SHong Zhang /* 2209566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD("F:\n")); 2219566063dSJacob Faibussowitsch PetscCall(VecView(F,PETSC_VIEWER_STDOUT_WORLD)); 22242dc13f1SHong Zhang */ 2233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22442dc13f1SHong Zhang } 22542dc13f1SHong Zhang 226d71ae5a4SJacob Faibussowitsch PetscErrorCode WASHSetInitialSolution(DM networkdm, Vec X, Wash wash) 227d71ae5a4SJacob Faibussowitsch { 22842dc13f1SHong Zhang PetscInt k, nx, vkey, vfrom, vto, offsetfrom, offsetto; 22942dc13f1SHong Zhang PetscInt type, varoffset; 23042dc13f1SHong Zhang PetscInt e, eStart, eEnd; 23142dc13f1SHong Zhang Vec localX; 23242dc13f1SHong Zhang PetscScalar *xarr; 23342dc13f1SHong Zhang Pipe pipe; 23442dc13f1SHong Zhang Junction junction; 23542dc13f1SHong Zhang const PetscInt *cone; 23642dc13f1SHong Zhang const PetscScalar *xarray; 23742dc13f1SHong Zhang 23842dc13f1SHong Zhang PetscFunctionBegin; 2399566063dSJacob Faibussowitsch PetscCall(VecSet(X, 0.0)); 2409566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(networkdm, &localX)); 2419566063dSJacob Faibussowitsch PetscCall(VecGetArray(localX, &xarr)); 24242dc13f1SHong Zhang 24342dc13f1SHong Zhang /* Edge */ 2449566063dSJacob Faibussowitsch PetscCall(DMNetworkGetEdgeRange(networkdm, &eStart, &eEnd)); 24542dc13f1SHong Zhang for (e = eStart; e < eEnd; e++) { 2469566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, e, ALL_COMPONENTS, &varoffset)); 2479566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, e, 0, &type, (void **)&pipe, NULL)); 24842dc13f1SHong Zhang 24942dc13f1SHong Zhang /* set initial values for this pipe */ 2509566063dSJacob Faibussowitsch PetscCall(PipeComputeSteadyState(pipe, wash->Q0, wash->H0)); 2519566063dSJacob Faibussowitsch PetscCall(VecGetSize(pipe->x, &nx)); 25242dc13f1SHong Zhang 2539566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pipe->x, &xarray)); 25442dc13f1SHong Zhang /* copy pipe->x to xarray */ 255ad540459SPierre Jolivet for (k = 0; k < nx; k++) (xarr + varoffset)[k] = xarray[k]; 25642dc13f1SHong Zhang 25742dc13f1SHong Zhang /* set boundary values into vfrom and vto */ 2589566063dSJacob Faibussowitsch PetscCall(DMNetworkGetConnectedVertices(networkdm, e, &cone)); 25942dc13f1SHong Zhang vfrom = cone[0]; /* local ordering */ 26042dc13f1SHong Zhang vto = cone[1]; 2619566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vfrom, ALL_COMPONENTS, &offsetfrom)); 2629566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vto, ALL_COMPONENTS, &offsetto)); 26342dc13f1SHong Zhang 26442dc13f1SHong Zhang /* if vform is a head vertex: */ 2659566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vfrom, 0, &vkey, (void **)&junction, NULL)); 2669371c9d4SSatish Balay if (junction->type == RESERVOIR) { (xarr + offsetfrom)[1] = wash->H0; /* 1st H */ } 26742dc13f1SHong Zhang 26842dc13f1SHong Zhang /* if vto is an end vertex: */ 2699566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vto, 0, &vkey, (void **)&junction, NULL)); 2709371c9d4SSatish Balay if (junction->type == VALVE) { (xarr + offsetto)[0] = wash->QL; /* last Q */ } 2719566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pipe->x, &xarray)); 27242dc13f1SHong Zhang } 27342dc13f1SHong Zhang 2749566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localX, &xarr)); 2759566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(networkdm, localX, ADD_VALUES, X)); 2769566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(networkdm, localX, ADD_VALUES, X)); 2779566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(networkdm, &localX)); 27842dc13f1SHong Zhang 27942dc13f1SHong Zhang #if 0 28042dc13f1SHong Zhang PetscInt N; 2819566063dSJacob Faibussowitsch PetscCall(VecGetSize(X,&N)); 2829566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"initial solution %d:\n",N)); 2839566063dSJacob Faibussowitsch PetscCall(VecView(X,PETSC_VIEWER_STDOUT_WORLD)); 28442dc13f1SHong Zhang #endif 2853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 28642dc13f1SHong Zhang } 28742dc13f1SHong Zhang 288d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDMNetworkMonitor(TS ts, PetscInt step, PetscReal t, Vec x, void *context) 289d71ae5a4SJacob Faibussowitsch { 29042dc13f1SHong Zhang DMNetworkMonitor monitor; 29142dc13f1SHong Zhang 29242dc13f1SHong Zhang PetscFunctionBegin; 29342dc13f1SHong Zhang monitor = (DMNetworkMonitor)context; 2949566063dSJacob Faibussowitsch PetscCall(DMNetworkMonitorView(monitor, x)); 2953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 29642dc13f1SHong Zhang } 29742dc13f1SHong Zhang 298d71ae5a4SJacob Faibussowitsch PetscErrorCode PipesView(DM networkdm, PetscInt KeyPipe, Vec X) 299d71ae5a4SJacob Faibussowitsch { 30042dc13f1SHong Zhang PetscInt i, numkeys = 1, *blocksize, *numselectedvariable, **selectedvariables, n; 30142dc13f1SHong Zhang IS isfrom_q, isfrom_h, isfrom; 30242dc13f1SHong Zhang Vec Xto; 30342dc13f1SHong Zhang VecScatter ctx; 30442dc13f1SHong Zhang MPI_Comm comm; 30542dc13f1SHong Zhang 30642dc13f1SHong Zhang PetscFunctionBegin; 3079566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)networkdm, &comm)); 30842dc13f1SHong Zhang 30942dc13f1SHong Zhang /* 1. Create isfrom_q for q-variable of pipes */ 3109566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(numkeys, &blocksize, numkeys, &numselectedvariable, numkeys, &selectedvariables)); 31142dc13f1SHong Zhang for (i = 0; i < numkeys; i++) { 31242dc13f1SHong Zhang blocksize[i] = 2; 31342dc13f1SHong Zhang numselectedvariable[i] = 1; 3149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numselectedvariable[i], &selectedvariables[i])); 31542dc13f1SHong Zhang selectedvariables[i][0] = 0; /* q-variable */ 31642dc13f1SHong Zhang } 3179566063dSJacob Faibussowitsch PetscCall(DMNetworkCreateIS(networkdm, numkeys, &KeyPipe, blocksize, numselectedvariable, selectedvariables, &isfrom_q)); 31842dc13f1SHong Zhang 31942dc13f1SHong Zhang /* 2. Create Xto and isto */ 3209566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isfrom_q, &n)); 3219566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &Xto)); 3229566063dSJacob Faibussowitsch PetscCall(VecSetSizes(Xto, n, PETSC_DECIDE)); 3239566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(Xto)); 3249566063dSJacob Faibussowitsch PetscCall(VecSet(Xto, 0.0)); 32542dc13f1SHong Zhang 32642dc13f1SHong Zhang /* 3. Create scatter */ 3279566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(X, isfrom_q, Xto, NULL, &ctx)); 32842dc13f1SHong Zhang 32942dc13f1SHong Zhang /* 4. Scatter to Xq */ 3309566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx, X, Xto, INSERT_VALUES, SCATTER_FORWARD)); 3319566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx, X, Xto, INSERT_VALUES, SCATTER_FORWARD)); 3329566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ctx)); 3339566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isfrom_q)); 33442dc13f1SHong Zhang 3359566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Xq:\n")); 3369566063dSJacob Faibussowitsch PetscCall(VecView(Xto, PETSC_VIEWER_STDOUT_WORLD)); 33742dc13f1SHong Zhang 33842dc13f1SHong Zhang /* 5. Create isfrom_h for h-variable of pipes; Create scatter; Scatter to Xh */ 3399371c9d4SSatish Balay for (i = 0; i < numkeys; i++) { selectedvariables[i][0] = 1; /* h-variable */ } 3409566063dSJacob Faibussowitsch PetscCall(DMNetworkCreateIS(networkdm, numkeys, &KeyPipe, blocksize, numselectedvariable, selectedvariables, &isfrom_h)); 34142dc13f1SHong Zhang 3429566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(X, isfrom_h, Xto, NULL, &ctx)); 3439566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx, X, Xto, INSERT_VALUES, SCATTER_FORWARD)); 3449566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx, X, Xto, INSERT_VALUES, SCATTER_FORWARD)); 3459566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ctx)); 3469566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isfrom_h)); 34742dc13f1SHong Zhang 3489566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Xh:\n")); 3499566063dSJacob Faibussowitsch PetscCall(VecView(Xto, PETSC_VIEWER_STDOUT_WORLD)); 3509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Xto)); 35142dc13f1SHong Zhang 35242dc13f1SHong Zhang /* 6. Create isfrom for all pipe variables; Create scatter; Scatter to Xpipes */ 3539371c9d4SSatish Balay for (i = 0; i < numkeys; i++) { blocksize[i] = -1; /* select all the variables of the i-th component */ } 3549566063dSJacob Faibussowitsch PetscCall(DMNetworkCreateIS(networkdm, numkeys, &KeyPipe, blocksize, NULL, NULL, &isfrom)); 3559566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isfrom)); 3569566063dSJacob Faibussowitsch PetscCall(DMNetworkCreateIS(networkdm, numkeys, &KeyPipe, NULL, NULL, NULL, &isfrom)); 35742dc13f1SHong Zhang 3589566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isfrom, &n)); 3599566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &Xto)); 3609566063dSJacob Faibussowitsch PetscCall(VecSetSizes(Xto, n, PETSC_DECIDE)); 3619566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(Xto)); 3629566063dSJacob Faibussowitsch PetscCall(VecSet(Xto, 0.0)); 36342dc13f1SHong Zhang 3649566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(X, isfrom, Xto, NULL, &ctx)); 3659566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(ctx, X, Xto, INSERT_VALUES, SCATTER_FORWARD)); 3669566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(ctx, X, Xto, INSERT_VALUES, SCATTER_FORWARD)); 3679566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&ctx)); 3689566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isfrom)); 36942dc13f1SHong Zhang 3709566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Xpipes:\n")); 3719566063dSJacob Faibussowitsch PetscCall(VecView(Xto, PETSC_VIEWER_STDOUT_WORLD)); 37242dc13f1SHong Zhang 37342dc13f1SHong Zhang /* 7. Free spaces */ 37448a46eb9SPierre Jolivet for (i = 0; i < numkeys; i++) PetscCall(PetscFree(selectedvariables[i])); 3759566063dSJacob Faibussowitsch PetscCall(PetscFree3(blocksize, numselectedvariable, selectedvariables)); 3769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Xto)); 3773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 37842dc13f1SHong Zhang } 37942dc13f1SHong Zhang 380d71ae5a4SJacob Faibussowitsch PetscErrorCode ISJunctionsView(DM networkdm, PetscInt KeyJunc) 381d71ae5a4SJacob Faibussowitsch { 38242dc13f1SHong Zhang PetscInt numkeys = 1; 38342dc13f1SHong Zhang IS isfrom; 38442dc13f1SHong Zhang MPI_Comm comm; 38542dc13f1SHong Zhang PetscMPIInt rank; 38642dc13f1SHong Zhang 38742dc13f1SHong Zhang PetscFunctionBegin; 3889566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)networkdm, &comm)); 3899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 39042dc13f1SHong Zhang 39142dc13f1SHong Zhang /* Create a global isfrom for all junction variables */ 3929566063dSJacob Faibussowitsch PetscCall(DMNetworkCreateIS(networkdm, numkeys, &KeyJunc, NULL, NULL, NULL, &isfrom)); 3939566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "ISJunctions:\n")); 3949566063dSJacob Faibussowitsch PetscCall(ISView(isfrom, PETSC_VIEWER_STDOUT_WORLD)); 3959566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isfrom)); 39642dc13f1SHong Zhang 39742dc13f1SHong Zhang /* Create a local isfrom for all junction variables */ 3989566063dSJacob Faibussowitsch PetscCall(DMNetworkCreateLocalIS(networkdm, numkeys, &KeyJunc, NULL, NULL, NULL, &isfrom)); 399c5853193SPierre Jolivet if (rank == 0) { 4009566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] ISLocalJunctions:\n", rank)); 4019566063dSJacob Faibussowitsch PetscCall(ISView(isfrom, PETSC_VIEWER_STDOUT_SELF)); 40242dc13f1SHong Zhang } 4039566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isfrom)); 4043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 40542dc13f1SHong Zhang } 40642dc13f1SHong Zhang 407d71ae5a4SJacob Faibussowitsch PetscErrorCode WashNetworkCleanUp(Wash wash) 408d71ae5a4SJacob Faibussowitsch { 40942dc13f1SHong Zhang PetscMPIInt rank; 41042dc13f1SHong Zhang 41142dc13f1SHong Zhang PetscFunctionBegin; 4129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(wash->comm, &rank)); 4139566063dSJacob Faibussowitsch PetscCall(PetscFree(wash->edgelist)); 4149566063dSJacob Faibussowitsch PetscCall(PetscFree(wash->vtype)); 41548a46eb9SPierre Jolivet if (rank == 0) PetscCall(PetscFree2(wash->junction, wash->pipe)); 4163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 41742dc13f1SHong Zhang } 41842dc13f1SHong Zhang 419d71ae5a4SJacob Faibussowitsch PetscErrorCode WashNetworkCreate(MPI_Comm comm, PetscInt pipesCase, Wash *wash_ptr) 420d71ae5a4SJacob Faibussowitsch { 42142dc13f1SHong Zhang PetscInt npipes; 42242dc13f1SHong Zhang PetscMPIInt rank; 42342dc13f1SHong Zhang Wash wash = NULL; 42442dc13f1SHong Zhang PetscInt i, numVertices, numEdges, *vtype; 42542dc13f1SHong Zhang PetscInt *edgelist; 42642dc13f1SHong Zhang Junction junctions = NULL; 42742dc13f1SHong Zhang Pipe pipes = NULL; 42842dc13f1SHong Zhang PetscBool washdist = PETSC_TRUE; 42942dc13f1SHong Zhang 43042dc13f1SHong Zhang PetscFunctionBegin; 4319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 43242dc13f1SHong Zhang 4339566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(1, &wash)); 43442dc13f1SHong Zhang wash->comm = comm; 43542dc13f1SHong Zhang *wash_ptr = wash; 43642dc13f1SHong Zhang wash->Q0 = 0.477432; /* RESERVOIR */ 43742dc13f1SHong Zhang wash->H0 = 150.0; 43842dc13f1SHong Zhang wash->HL = 143.488; /* VALVE */ 43942dc13f1SHong Zhang wash->QL = 0.0; 44042dc13f1SHong Zhang wash->nnodes_loc = 0; 44142dc13f1SHong Zhang 44242dc13f1SHong Zhang numVertices = 0; 44342dc13f1SHong Zhang numEdges = 0; 44442dc13f1SHong Zhang edgelist = NULL; 44542dc13f1SHong Zhang 44642dc13f1SHong Zhang /* proc[0] creates a sequential wash and edgelist */ 44763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Setup pipesCase %" PetscInt_FMT "\n", pipesCase)); 44842dc13f1SHong Zhang 44942dc13f1SHong Zhang /* Set global number of pipes, edges, and junctions */ 45042dc13f1SHong Zhang /*-------------------------------------------------*/ 45142dc13f1SHong Zhang switch (pipesCase) { 45242dc13f1SHong Zhang case 0: 45342dc13f1SHong Zhang /* pipeCase 0: */ 45442dc13f1SHong Zhang /* ================================================= 45542dc13f1SHong Zhang (RESERVOIR) v0 --E0--> v1--E1--> v2 --E2-->v3 (VALVE) 45642dc13f1SHong Zhang ==================================================== */ 45742dc13f1SHong Zhang npipes = 3; 4589566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-npipes", &npipes, NULL)); 45942dc13f1SHong Zhang wash->nedge = npipes; 46042dc13f1SHong Zhang wash->nvertex = npipes + 1; 46142dc13f1SHong Zhang 46242dc13f1SHong Zhang /* Set local edges and vertices -- proc[0] sets entire network, then distributes */ 46342dc13f1SHong Zhang numVertices = 0; 46442dc13f1SHong Zhang numEdges = 0; 46542dc13f1SHong Zhang edgelist = NULL; 466dd400576SPatrick Sanan if (rank == 0) { 46742dc13f1SHong Zhang numVertices = wash->nvertex; 46842dc13f1SHong Zhang numEdges = wash->nedge; 46942dc13f1SHong Zhang 4709566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(2 * numEdges, &edgelist)); 47142dc13f1SHong Zhang for (i = 0; i < numEdges; i++) { 4729371c9d4SSatish Balay edgelist[2 * i] = i; 4739371c9d4SSatish Balay edgelist[2 * i + 1] = i + 1; 47442dc13f1SHong Zhang } 47542dc13f1SHong Zhang 47642dc13f1SHong Zhang /* Add network components */ 47742dc13f1SHong Zhang /*------------------------*/ 4789566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(numVertices, &junctions, numEdges, &pipes)); 47942dc13f1SHong Zhang 48042dc13f1SHong Zhang /* vertex */ 48142dc13f1SHong Zhang for (i = 0; i < numVertices; i++) { 48242dc13f1SHong Zhang junctions[i].id = i; 48342dc13f1SHong Zhang junctions[i].type = JUNCTION; 48442dc13f1SHong Zhang } 48542dc13f1SHong Zhang 48642dc13f1SHong Zhang junctions[0].type = RESERVOIR; 48742dc13f1SHong Zhang junctions[numVertices - 1].type = VALVE; 48842dc13f1SHong Zhang } 48942dc13f1SHong Zhang break; 49042dc13f1SHong Zhang case 1: 49142dc13f1SHong Zhang /* pipeCase 1: */ 49242dc13f1SHong Zhang /* ========================== 49342dc13f1SHong Zhang v2 (VALVE) 49442dc13f1SHong Zhang ^ 49542dc13f1SHong Zhang | 49642dc13f1SHong Zhang E2 49742dc13f1SHong Zhang | 49842dc13f1SHong Zhang v0 --E0--> v3--E1--> v1 49942dc13f1SHong Zhang (RESERVOIR) (RESERVOIR) 50042dc13f1SHong Zhang ============================= */ 50142dc13f1SHong Zhang npipes = 3; 50242dc13f1SHong Zhang wash->nedge = npipes; 50342dc13f1SHong Zhang wash->nvertex = npipes + 1; 50442dc13f1SHong Zhang 50542dc13f1SHong Zhang /* Set local edges and vertices -- proc[0] sets entire network, then distributes */ 506dd400576SPatrick Sanan if (rank == 0) { 50742dc13f1SHong Zhang numVertices = wash->nvertex; 50842dc13f1SHong Zhang numEdges = wash->nedge; 50942dc13f1SHong Zhang 5109566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(2 * numEdges, &edgelist)); 5119371c9d4SSatish Balay edgelist[0] = 0; 5129371c9d4SSatish Balay edgelist[1] = 3; /* edge[0] */ 5139371c9d4SSatish Balay edgelist[2] = 3; 5149371c9d4SSatish Balay edgelist[3] = 1; /* edge[1] */ 5159371c9d4SSatish Balay edgelist[4] = 3; 5169371c9d4SSatish Balay edgelist[5] = 2; /* edge[2] */ 51742dc13f1SHong Zhang 51842dc13f1SHong Zhang /* Add network components */ 51942dc13f1SHong Zhang /*------------------------*/ 5209566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(numVertices, &junctions, numEdges, &pipes)); 52142dc13f1SHong Zhang /* vertex */ 52242dc13f1SHong Zhang for (i = 0; i < numVertices; i++) { 52342dc13f1SHong Zhang junctions[i].id = i; 52442dc13f1SHong Zhang junctions[i].type = JUNCTION; 52542dc13f1SHong Zhang } 52642dc13f1SHong Zhang 52742dc13f1SHong Zhang junctions[0].type = RESERVOIR; 52842dc13f1SHong Zhang junctions[1].type = VALVE; 52942dc13f1SHong Zhang junctions[2].type = VALVE; 53042dc13f1SHong Zhang } 53142dc13f1SHong Zhang break; 53242dc13f1SHong Zhang case 2: 53342dc13f1SHong Zhang /* pipeCase 2: */ 53442dc13f1SHong Zhang /* ========================== 53542dc13f1SHong Zhang (RESERVOIR) v2--> E2 53642dc13f1SHong Zhang | 53742dc13f1SHong Zhang v0 --E0--> v3--E1--> v1 53842dc13f1SHong Zhang (RESERVOIR) (VALVE) 53942dc13f1SHong Zhang ============================= */ 54042dc13f1SHong Zhang 541d5b43468SJose E. Roman /* Set application parameters -- to be used in function evaluations */ 54242dc13f1SHong Zhang npipes = 3; 54342dc13f1SHong Zhang wash->nedge = npipes; 54442dc13f1SHong Zhang wash->nvertex = npipes + 1; 54542dc13f1SHong Zhang 54642dc13f1SHong Zhang /* Set local edges and vertices -- proc[0] sets entire network, then distributes */ 547dd400576SPatrick Sanan if (rank == 0) { 54842dc13f1SHong Zhang numVertices = wash->nvertex; 54942dc13f1SHong Zhang numEdges = wash->nedge; 55042dc13f1SHong Zhang 5519566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(2 * numEdges, &edgelist)); 5529371c9d4SSatish Balay edgelist[0] = 0; 5539371c9d4SSatish Balay edgelist[1] = 3; /* edge[0] */ 5549371c9d4SSatish Balay edgelist[2] = 3; 5559371c9d4SSatish Balay edgelist[3] = 1; /* edge[1] */ 5569371c9d4SSatish Balay edgelist[4] = 2; 5579371c9d4SSatish Balay edgelist[5] = 3; /* edge[2] */ 55842dc13f1SHong Zhang 55942dc13f1SHong Zhang /* Add network components */ 56042dc13f1SHong Zhang /*------------------------*/ 5619566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(numVertices, &junctions, numEdges, &pipes)); 56242dc13f1SHong Zhang /* vertex */ 56342dc13f1SHong Zhang for (i = 0; i < numVertices; i++) { 56442dc13f1SHong Zhang junctions[i].id = i; 56542dc13f1SHong Zhang junctions[i].type = JUNCTION; 56642dc13f1SHong Zhang } 56742dc13f1SHong Zhang 56842dc13f1SHong Zhang junctions[0].type = RESERVOIR; 56942dc13f1SHong Zhang junctions[1].type = VALVE; 57042dc13f1SHong Zhang junctions[2].type = RESERVOIR; 57142dc13f1SHong Zhang } 57242dc13f1SHong Zhang break; 573d71ae5a4SJacob Faibussowitsch default: 574d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "not done yet"); 57542dc13f1SHong Zhang } 57642dc13f1SHong Zhang 57742dc13f1SHong Zhang /* set edge global id */ 57842dc13f1SHong Zhang for (i = 0; i < numEdges; i++) pipes[i].id = i; 57942dc13f1SHong Zhang 580dd400576SPatrick Sanan if (rank == 0) { /* set vtype for proc[0] */ 58142dc13f1SHong Zhang PetscInt v; 5829566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * numEdges, &vtype)); 58342dc13f1SHong Zhang for (i = 0; i < 2 * numEdges; i++) { 58442dc13f1SHong Zhang v = edgelist[i]; 58542dc13f1SHong Zhang vtype[i] = junctions[v].type; 58642dc13f1SHong Zhang } 58742dc13f1SHong Zhang wash->vtype = vtype; 58842dc13f1SHong Zhang } 58942dc13f1SHong Zhang 59042dc13f1SHong Zhang *wash_ptr = wash; 59142dc13f1SHong Zhang wash->nedge = numEdges; 59242dc13f1SHong Zhang wash->nvertex = numVertices; 59342dc13f1SHong Zhang wash->edgelist = edgelist; 59442dc13f1SHong Zhang wash->junction = junctions; 59542dc13f1SHong Zhang wash->pipe = pipes; 59642dc13f1SHong Zhang 59742dc13f1SHong Zhang /* Distribute edgelist to other processors */ 5989566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-wash_distribute", &washdist, NULL)); 59942dc13f1SHong Zhang if (washdist) { 60042dc13f1SHong Zhang /* 6019566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Distribute sequential wash ...\n")); 60242dc13f1SHong Zhang */ 6039566063dSJacob Faibussowitsch PetscCall(WashNetworkDistribute(comm, wash)); 60442dc13f1SHong Zhang } 6053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 60642dc13f1SHong Zhang } 60742dc13f1SHong Zhang 60842dc13f1SHong Zhang /* ------------------------------------------------------- */ 609d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 610d71ae5a4SJacob Faibussowitsch { 61142dc13f1SHong Zhang Wash wash; 61242dc13f1SHong Zhang Junction junctions, junction; 61342dc13f1SHong Zhang Pipe pipe, pipes; 614f11a936eSBarry Smith PetscInt KeyPipe, KeyJunction, *edgelist = NULL, *vtype = NULL; 615f11a936eSBarry Smith PetscInt i, e, v, eStart, eEnd, vStart, vEnd, key, vkey, type; 616f11a936eSBarry Smith PetscInt steps = 1, nedges, nnodes = 6; 61742dc13f1SHong Zhang const PetscInt *cone; 61842dc13f1SHong Zhang DM networkdm; 61942dc13f1SHong Zhang PetscMPIInt size, rank; 62042dc13f1SHong Zhang PetscReal ftime; 62142dc13f1SHong Zhang Vec X; 62242dc13f1SHong Zhang TS ts; 62342dc13f1SHong Zhang TSConvergedReason reason; 62442dc13f1SHong Zhang PetscBool viewpipes, viewjuncs, monipipes = PETSC_FALSE, userJac = PETSC_TRUE, viewdm = PETSC_FALSE, viewX = PETSC_FALSE; 62542dc13f1SHong Zhang PetscInt pipesCase = 0; 62642dc13f1SHong Zhang DMNetworkMonitor monitor; 62742dc13f1SHong Zhang MPI_Comm comm; 62842dc13f1SHong Zhang 629327415f7SBarry Smith PetscFunctionBeginUser; 6309566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, "pOption", help)); 63142dc13f1SHong Zhang 63242dc13f1SHong Zhang /* Read runtime options */ 6339566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-case", &pipesCase, NULL)); 6349566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-user_Jac", &userJac, NULL)); 6359566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-pipe_monitor", &monipipes, NULL)); 6369566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewdm", &viewdm, NULL)); 6379566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewX", &viewX, NULL)); 6389566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-npipenodes", &nnodes, NULL)); 63942dc13f1SHong Zhang 64042dc13f1SHong Zhang /* Create networkdm */ 64142dc13f1SHong Zhang /*------------------*/ 6429566063dSJacob Faibussowitsch PetscCall(DMNetworkCreate(PETSC_COMM_WORLD, &networkdm)); 6439566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)networkdm, &comm)); 6449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 6459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 64642dc13f1SHong Zhang 64748a46eb9SPierre Jolivet if (size == 1 && monipipes) PetscCall(DMNetworkMonitorCreate(networkdm, &monitor)); 64842dc13f1SHong Zhang 64942dc13f1SHong Zhang /* Register the components in the network */ 6509566063dSJacob Faibussowitsch PetscCall(DMNetworkRegisterComponent(networkdm, "junctionstruct", sizeof(struct _p_Junction), &KeyJunction)); 6519566063dSJacob Faibussowitsch PetscCall(DMNetworkRegisterComponent(networkdm, "pipestruct", sizeof(struct _p_Pipe), &KeyPipe)); 65242dc13f1SHong Zhang 65342dc13f1SHong Zhang /* Create a distributed wash network (user-specific) */ 6549566063dSJacob Faibussowitsch PetscCall(WashNetworkCreate(comm, pipesCase, &wash)); 65542dc13f1SHong Zhang nedges = wash->nedge; 65642dc13f1SHong Zhang edgelist = wash->edgelist; 65742dc13f1SHong Zhang vtype = wash->vtype; 65842dc13f1SHong Zhang junctions = wash->junction; 65942dc13f1SHong Zhang pipes = wash->pipe; 66042dc13f1SHong Zhang 66142dc13f1SHong Zhang /* Set up the network layout */ 6629566063dSJacob Faibussowitsch PetscCall(DMNetworkSetNumSubNetworks(networkdm, PETSC_DECIDE, 1)); 6639566063dSJacob Faibussowitsch PetscCall(DMNetworkAddSubnetwork(networkdm, NULL, nedges, edgelist, NULL)); 66442dc13f1SHong Zhang 6659566063dSJacob Faibussowitsch PetscCall(DMNetworkLayoutSetUp(networkdm)); 66642dc13f1SHong Zhang 6679566063dSJacob Faibussowitsch PetscCall(DMNetworkGetEdgeRange(networkdm, &eStart, &eEnd)); 6689566063dSJacob Faibussowitsch PetscCall(DMNetworkGetVertexRange(networkdm, &vStart, &vEnd)); 6699566063dSJacob Faibussowitsch /* PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d] eStart/End: %d - %d; vStart/End: %d - %d\n",rank,eStart,eEnd,vStart,vEnd)); */ 67042dc13f1SHong Zhang 67142dc13f1SHong Zhang if (rank) { /* junctions[] and pipes[] for proc[0] are allocated in WashNetworkCreate() */ 67242dc13f1SHong Zhang /* vEnd - vStart = nvertices + number of ghost vertices! */ 6739566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(vEnd - vStart, &junctions, nedges, &pipes)); 67442dc13f1SHong Zhang } 67542dc13f1SHong Zhang 67642dc13f1SHong Zhang /* Add Pipe component and number of variables to all local edges */ 67742dc13f1SHong Zhang for (e = eStart; e < eEnd; e++) { 67842dc13f1SHong Zhang pipes[e - eStart].nnodes = nnodes; 6799566063dSJacob Faibussowitsch PetscCall(DMNetworkAddComponent(networkdm, e, KeyPipe, &pipes[e - eStart], 2 * pipes[e - eStart].nnodes)); 68042dc13f1SHong Zhang 68142dc13f1SHong Zhang if (size == 1 && monipipes) { /* Add monitor -- show Q_{pipes[e-eStart].id}? */ 68242dc13f1SHong Zhang pipes[e - eStart].length = 600.0; 6839566063dSJacob Faibussowitsch PetscCall(DMNetworkMonitorAdd(monitor, "Pipe Q", e, pipes[e - eStart].nnodes, 0, 2, 0.0, pipes[e - eStart].length, -0.8, 0.8, PETSC_TRUE)); 6849566063dSJacob Faibussowitsch PetscCall(DMNetworkMonitorAdd(monitor, "Pipe H", e, pipes[e - eStart].nnodes, 1, 2, 0.0, pipes[e - eStart].length, -400.0, 800.0, PETSC_TRUE)); 68542dc13f1SHong Zhang } 68642dc13f1SHong Zhang } 68742dc13f1SHong Zhang 688eac198afSGetnet /* Add Junction component and number of variables to all local vertices */ 68948a46eb9SPierre Jolivet for (v = vStart; v < vEnd; v++) PetscCall(DMNetworkAddComponent(networkdm, v, KeyJunction, &junctions[v - vStart], 2)); 69042dc13f1SHong Zhang 69142dc13f1SHong Zhang if (size > 1) { /* must be called before DMSetUp()???. Other partitioners do not work yet??? -- cause crash in proc[0]! */ 69242dc13f1SHong Zhang DM plexdm; 69342dc13f1SHong Zhang PetscPartitioner part; 6949566063dSJacob Faibussowitsch PetscCall(DMNetworkGetPlex(networkdm, &plexdm)); 6959566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(plexdm, &part)); 6969566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSIMPLE)); 6979566063dSJacob Faibussowitsch PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_csr_alg", "mat")); /* for parmetis */ 69842dc13f1SHong Zhang } 69942dc13f1SHong Zhang 70042dc13f1SHong Zhang /* Set up DM for use */ 7019566063dSJacob Faibussowitsch PetscCall(DMSetUp(networkdm)); 70242dc13f1SHong Zhang if (viewdm) { 7039566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nOriginal networkdm, DMView:\n")); 7049566063dSJacob Faibussowitsch PetscCall(DMView(networkdm, PETSC_VIEWER_STDOUT_WORLD)); 70542dc13f1SHong Zhang } 70642dc13f1SHong Zhang 70742dc13f1SHong Zhang /* Set user physical parameters to the components */ 70842dc13f1SHong Zhang for (e = eStart; e < eEnd; e++) { 7099566063dSJacob Faibussowitsch PetscCall(DMNetworkGetConnectedVertices(networkdm, e, &cone)); 71042dc13f1SHong Zhang /* vfrom */ 7119566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, cone[0], 0, &vkey, (void **)&junction, NULL)); 71242dc13f1SHong Zhang junction->type = (VertexType)vtype[2 * e]; 71342dc13f1SHong Zhang 71442dc13f1SHong Zhang /* vto */ 7159566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, cone[1], 0, &vkey, (void **)&junction, NULL)); 71642dc13f1SHong Zhang junction->type = (VertexType)vtype[2 * e + 1]; 71742dc13f1SHong Zhang } 71842dc13f1SHong Zhang 7199566063dSJacob Faibussowitsch PetscCall(WashNetworkCleanUp(wash)); 72042dc13f1SHong Zhang 72142dc13f1SHong Zhang /* Network partitioning and distribution of data */ 7229566063dSJacob Faibussowitsch PetscCall(DMNetworkDistribute(&networkdm, 0)); 72342dc13f1SHong Zhang if (viewdm) { 7249566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nAfter DMNetworkDistribute, DMView:\n")); 7259566063dSJacob Faibussowitsch PetscCall(DMView(networkdm, PETSC_VIEWER_STDOUT_WORLD)); 72642dc13f1SHong Zhang } 72742dc13f1SHong Zhang 72842dc13f1SHong Zhang /* create vectors */ 7299566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(networkdm, &X)); 7309566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(networkdm, &wash->localX)); 7319566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(networkdm, &wash->localXdot)); 73242dc13f1SHong Zhang 73342dc13f1SHong Zhang /* PipeSetUp -- each process only sets its own pipes */ 73442dc13f1SHong Zhang /*---------------------------------------------------*/ 7359566063dSJacob Faibussowitsch PetscCall(DMNetworkGetVertexRange(networkdm, &vStart, &vEnd)); 73642dc13f1SHong Zhang 73742dc13f1SHong Zhang userJac = PETSC_TRUE; 7389566063dSJacob Faibussowitsch PetscCall(DMNetworkHasJacobian(networkdm, userJac, userJac)); 7399566063dSJacob Faibussowitsch PetscCall(DMNetworkGetEdgeRange(networkdm, &eStart, &eEnd)); 74042dc13f1SHong Zhang for (e = eStart; e < eEnd; e++) { /* each edge has only one component, pipe */ 7419566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, e, 0, &type, (void **)&pipe, NULL)); 74242dc13f1SHong Zhang 74342dc13f1SHong Zhang wash->nnodes_loc += pipe->nnodes; /* local total number of nodes, will be used by PipesView() */ 7449371c9d4SSatish Balay PetscCall(PipeSetParameters(pipe, 600.0, /* length */ 74542dc13f1SHong Zhang 0.5, /* diameter */ 74642dc13f1SHong Zhang 1200.0, /* a */ 747b122ec5aSJacob Faibussowitsch 0.018)); /* friction */ 7489566063dSJacob Faibussowitsch PetscCall(PipeSetUp(pipe)); 74942dc13f1SHong Zhang 75042dc13f1SHong Zhang if (userJac) { 751*0b4b7b1cSBarry Smith /* Create Jacobian matrix nonzero structures for a Pipe */ 75242dc13f1SHong Zhang Mat *J; 7539566063dSJacob Faibussowitsch PetscCall(PipeCreateJacobian(pipe, NULL, &J)); 7549566063dSJacob Faibussowitsch PetscCall(DMNetworkEdgeSetMatrix(networkdm, e, J)); 75542dc13f1SHong Zhang } 75642dc13f1SHong Zhang } 75742dc13f1SHong Zhang 75842dc13f1SHong Zhang if (userJac) { 7599566063dSJacob Faibussowitsch PetscCall(DMNetworkGetVertexRange(networkdm, &vStart, &vEnd)); 76042dc13f1SHong Zhang for (v = vStart; v < vEnd; v++) { 76142dc13f1SHong Zhang Mat *J; 7629566063dSJacob Faibussowitsch PetscCall(JunctionCreateJacobian(networkdm, v, NULL, &J)); 7639566063dSJacob Faibussowitsch PetscCall(DMNetworkVertexSetMatrix(networkdm, v, J)); 76442dc13f1SHong Zhang 7659566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, v, 0, &vkey, (void **)&junction, NULL)); 76642dc13f1SHong Zhang junction->jacobian = J; 76742dc13f1SHong Zhang } 76842dc13f1SHong Zhang } 76942dc13f1SHong Zhang 77042dc13f1SHong Zhang /* Setup solver */ 77142dc13f1SHong Zhang /*--------------------------------------------------------*/ 7729566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 77342dc13f1SHong Zhang 7749566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, (DM)networkdm)); 7759566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, WASHIFunction, wash)); 77642dc13f1SHong Zhang 7779566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, steps)); 7789566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 7799566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.1)); 7809566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBEULER)); 7819566063dSJacob Faibussowitsch if (size == 1 && monipipes) PetscCall(TSMonitorSet(ts, TSDMNetworkMonitor, monitor, NULL)); 7829566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 78342dc13f1SHong Zhang 7849566063dSJacob Faibussowitsch PetscCall(WASHSetInitialSolution(networkdm, X, wash)); 78542dc13f1SHong Zhang 7869566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X)); 78742dc13f1SHong Zhang 7889566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 7899566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 7909566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts, &reason)); 79163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double)ftime, steps)); 7921baa6e33SBarry Smith if (viewX) PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD)); 79342dc13f1SHong Zhang 79442dc13f1SHong Zhang viewpipes = PETSC_FALSE; 7959566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-Jac_view", &viewpipes, NULL)); 79642dc13f1SHong Zhang if (viewpipes) { 79742dc13f1SHong Zhang SNES snes; 79842dc13f1SHong Zhang Mat Jac; 7999566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 8009566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, &Jac, NULL, NULL, NULL)); 8019566063dSJacob Faibussowitsch PetscCall(MatView(Jac, PETSC_VIEWER_DRAW_WORLD)); 80242dc13f1SHong Zhang } 80342dc13f1SHong Zhang 80442dc13f1SHong Zhang /* View solutions */ 80542dc13f1SHong Zhang /* -------------- */ 80642dc13f1SHong Zhang viewpipes = PETSC_FALSE; 8079566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-pipe_view", &viewpipes, NULL)); 8081baa6e33SBarry Smith if (viewpipes) PetscCall(PipesView(networkdm, KeyPipe, X)); 80942dc13f1SHong Zhang 81042dc13f1SHong Zhang /* Test IS */ 81142dc13f1SHong Zhang viewjuncs = PETSC_FALSE; 8129566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-isjunc_view", &viewjuncs, NULL)); 8131baa6e33SBarry Smith if (viewjuncs) PetscCall(ISJunctionsView(networkdm, KeyJunction)); 81442dc13f1SHong Zhang 81542dc13f1SHong Zhang /* Free spaces */ 81642dc13f1SHong Zhang /* ----------- */ 8179566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 8189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 8199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&wash->localX)); 8209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&wash->localXdot)); 82142dc13f1SHong Zhang 82242dc13f1SHong Zhang /* Destroy objects from each pipe that are created in PipeSetUp() */ 8239566063dSJacob Faibussowitsch PetscCall(DMNetworkGetEdgeRange(networkdm, &eStart, &eEnd)); 82442dc13f1SHong Zhang for (i = eStart; i < eEnd; i++) { 8259566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, i, 0, &key, (void **)&pipe, NULL)); 8269566063dSJacob Faibussowitsch PetscCall(PipeDestroy(&pipe)); 82742dc13f1SHong Zhang } 82842dc13f1SHong Zhang if (userJac) { 82942dc13f1SHong Zhang for (v = vStart; v < vEnd; v++) { 8309566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, v, 0, &vkey, (void **)&junction, NULL)); 8319566063dSJacob Faibussowitsch PetscCall(JunctionDestroyJacobian(networkdm, v, junction)); 83242dc13f1SHong Zhang } 83342dc13f1SHong Zhang } 83442dc13f1SHong Zhang 83548a46eb9SPierre Jolivet if (size == 1 && monipipes) PetscCall(DMNetworkMonitorDestroy(&monitor)); 8369566063dSJacob Faibussowitsch PetscCall(DMDestroy(&networkdm)); 8379566063dSJacob Faibussowitsch PetscCall(PetscFree(wash)); 83842dc13f1SHong Zhang 8391baa6e33SBarry Smith if (rank) PetscCall(PetscFree2(junctions, pipes)); 8409566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 841b122ec5aSJacob Faibussowitsch return 0; 84242dc13f1SHong Zhang } 84342dc13f1SHong Zhang 84442dc13f1SHong Zhang /*TEST 84542dc13f1SHong Zhang 84642dc13f1SHong Zhang build: 84742dc13f1SHong Zhang depends: pipeInterface.c pipeImpls.c 84842dc13f1SHong Zhang requires: mumps 84942dc13f1SHong Zhang 85042dc13f1SHong Zhang test: 85142dc13f1SHong Zhang args: -ts_monitor -case 1 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -options_left no -viewX 85242dc13f1SHong Zhang localrunfiles: pOption 85342dc13f1SHong Zhang output_file: output/pipes_1.out 85442dc13f1SHong Zhang 85542dc13f1SHong Zhang test: 85642dc13f1SHong Zhang suffix: 2 85742dc13f1SHong Zhang nsize: 2 85842dc13f1SHong Zhang args: -ts_monitor -case 1 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -petscpartitioner_type simple -options_left no -viewX 85942dc13f1SHong Zhang localrunfiles: pOption 86042dc13f1SHong Zhang output_file: output/pipes_2.out 86142dc13f1SHong Zhang 86242dc13f1SHong Zhang test: 86342dc13f1SHong Zhang suffix: 3 86442dc13f1SHong Zhang nsize: 2 86542dc13f1SHong Zhang args: -ts_monitor -case 0 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -petscpartitioner_type simple -options_left no -viewX 86642dc13f1SHong Zhang localrunfiles: pOption 86742dc13f1SHong Zhang output_file: output/pipes_3.out 86842dc13f1SHong Zhang 86942dc13f1SHong Zhang test: 87042dc13f1SHong Zhang suffix: 4 87142dc13f1SHong Zhang args: -ts_monitor -case 2 -ts_max_steps 1 -pc_factor_mat_solver_type mumps -options_left no -viewX 87242dc13f1SHong Zhang localrunfiles: pOption 87342dc13f1SHong Zhang output_file: output/pipes_4.out 87442dc13f1SHong Zhang 87542dc13f1SHong Zhang test: 87642dc13f1SHong Zhang suffix: 5 87742dc13f1SHong Zhang nsize: 3 87842dc13f1SHong Zhang args: -ts_monitor -case 2 -ts_max_steps 10 -pc_factor_mat_solver_type mumps -petscpartitioner_type simple -options_left no -viewX 87942dc13f1SHong Zhang localrunfiles: pOption 88042dc13f1SHong Zhang output_file: output/pipes_5.out 88142dc13f1SHong Zhang 88242dc13f1SHong Zhang test: 88342dc13f1SHong Zhang suffix: 6 88442dc13f1SHong Zhang nsize: 2 88542dc13f1SHong 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 88642dc13f1SHong Zhang localrunfiles: pOption 88742dc13f1SHong Zhang output_file: output/pipes_6.out 88842dc13f1SHong Zhang 88942dc13f1SHong Zhang test: 89042dc13f1SHong Zhang suffix: 7 89142dc13f1SHong Zhang nsize: 2 89242dc13f1SHong 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 89342dc13f1SHong Zhang localrunfiles: pOption 89442dc13f1SHong Zhang output_file: output/pipes_7.out 89542dc13f1SHong Zhang 89642dc13f1SHong Zhang test: 89742dc13f1SHong Zhang suffix: 8 89842dc13f1SHong Zhang nsize: 2 89942dc13f1SHong Zhang requires: parmetis 90042dc13f1SHong 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 90142dc13f1SHong Zhang localrunfiles: pOption 90242dc13f1SHong Zhang output_file: output/pipes_8.out 90342dc13f1SHong Zhang 90442dc13f1SHong Zhang test: 90542dc13f1SHong Zhang suffix: 9 90642dc13f1SHong Zhang nsize: 2 90742dc13f1SHong 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 90842dc13f1SHong Zhang localrunfiles: pOption 90942dc13f1SHong Zhang output_file: output/pipes_9.out 91042dc13f1SHong Zhang 91142dc13f1SHong Zhang test: 91242dc13f1SHong Zhang suffix: 10 91342dc13f1SHong Zhang nsize: 2 91442dc13f1SHong 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 91542dc13f1SHong Zhang localrunfiles: pOption 91642dc13f1SHong Zhang output_file: output/pipes_10.out 91742dc13f1SHong Zhang 91842dc13f1SHong Zhang TEST*/ 919