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