142dc13f1SHong Zhang #include "wash.h" 242dc13f1SHong Zhang 342dc13f1SHong Zhang /* Subroutines for Pipe */ 442dc13f1SHong Zhang /* -------------------------------------------------------*/ 542dc13f1SHong Zhang 642dc13f1SHong Zhang /* 742dc13f1SHong Zhang PipeCreate - Create Pipe object. 842dc13f1SHong Zhang 942dc13f1SHong Zhang Input Parameters: 1042dc13f1SHong Zhang comm - MPI communicator 1142dc13f1SHong Zhang 1242dc13f1SHong Zhang Output Parameter: 1342dc13f1SHong Zhang . pipe - location to put the PIPE context 1442dc13f1SHong Zhang */ 1542dc13f1SHong Zhang PetscErrorCode PipeCreate(MPI_Comm comm,Pipe *pipe) 1642dc13f1SHong Zhang { 1742dc13f1SHong Zhang PetscErrorCode ierr; 1842dc13f1SHong Zhang 1942dc13f1SHong Zhang PetscFunctionBegin; 2042dc13f1SHong Zhang ierr = PetscNew(pipe);CHKERRQ(ierr); 2142dc13f1SHong Zhang PetscFunctionReturn(0); 2242dc13f1SHong Zhang } 2342dc13f1SHong Zhang 2442dc13f1SHong Zhang /* 2542dc13f1SHong Zhang PipeDestroy - Destroy Pipe object. 2642dc13f1SHong Zhang 2742dc13f1SHong Zhang Input Parameters: 2842dc13f1SHong Zhang pipe - Reference to pipe intended to be destroyed. 2942dc13f1SHong Zhang */ 3042dc13f1SHong Zhang PetscErrorCode PipeDestroy(Pipe *pipe) 3142dc13f1SHong Zhang { 3242dc13f1SHong Zhang PetscErrorCode ierr; 3342dc13f1SHong Zhang 3442dc13f1SHong Zhang PetscFunctionBegin; 3542dc13f1SHong Zhang if (!*pipe) PetscFunctionReturn(0); 3642dc13f1SHong Zhang 3742dc13f1SHong Zhang ierr = PipeDestroyJacobian(*pipe);CHKERRQ(ierr); 3842dc13f1SHong Zhang ierr = VecDestroy(&(*pipe)->x);CHKERRQ(ierr); 3942dc13f1SHong Zhang ierr = DMDestroy(&(*pipe)->da);CHKERRQ(ierr); 4042dc13f1SHong Zhang PetscFunctionReturn(0); 4142dc13f1SHong Zhang } 4242dc13f1SHong Zhang 4342dc13f1SHong Zhang /* 4442dc13f1SHong Zhang PipeSetParameters - Set parameters for Pipe context 4542dc13f1SHong Zhang 4642dc13f1SHong Zhang Input Parameter: 4742dc13f1SHong Zhang + pipe - PIPE object 4842dc13f1SHong Zhang . length - 4942dc13f1SHong Zhang . nnodes - 5042dc13f1SHong Zhang . D - 5142dc13f1SHong Zhang . a - 5242dc13f1SHong Zhang - fric - 5342dc13f1SHong Zhang */ 5442dc13f1SHong Zhang PetscErrorCode PipeSetParameters(Pipe pipe,PetscReal length,PetscReal D,PetscReal a,PetscReal fric) 5542dc13f1SHong Zhang { 5642dc13f1SHong Zhang PetscFunctionBegin; 5742dc13f1SHong Zhang pipe->length = length; 5842dc13f1SHong Zhang pipe->D = D; 5942dc13f1SHong Zhang pipe->a = a; 6042dc13f1SHong Zhang pipe->fric = fric; 6142dc13f1SHong Zhang PetscFunctionReturn(0); 6242dc13f1SHong Zhang } 6342dc13f1SHong Zhang 6442dc13f1SHong Zhang /* 6542dc13f1SHong Zhang PipeSetUp - Set up pipe based on set parameters. 6642dc13f1SHong Zhang */ 6742dc13f1SHong Zhang PetscErrorCode PipeSetUp(Pipe pipe) 6842dc13f1SHong Zhang { 6942dc13f1SHong Zhang DMDALocalInfo info; 7042dc13f1SHong Zhang PetscErrorCode ierr; 7142dc13f1SHong Zhang 7242dc13f1SHong Zhang PetscFunctionBegin; 7342dc13f1SHong Zhang ierr = DMDACreate1d(PETSC_COMM_SELF, DM_BOUNDARY_GHOSTED, pipe->nnodes, 2, 1, NULL, &pipe->da);CHKERRQ(ierr); 7442dc13f1SHong Zhang ierr = DMSetFromOptions(pipe->da);CHKERRQ(ierr); 7542dc13f1SHong Zhang ierr = DMSetUp(pipe->da);CHKERRQ(ierr); 7642dc13f1SHong Zhang ierr = DMDASetFieldName(pipe->da, 0, "Q");CHKERRQ(ierr); 7742dc13f1SHong Zhang ierr = DMDASetFieldName(pipe->da, 1, "H");CHKERRQ(ierr); 7842dc13f1SHong Zhang ierr = DMDASetUniformCoordinates(pipe->da, 0, pipe->length, 0, 0, 0, 0);CHKERRQ(ierr); 7942dc13f1SHong Zhang ierr = DMCreateGlobalVector(pipe->da, &(pipe->x));CHKERRQ(ierr); 8042dc13f1SHong Zhang 8142dc13f1SHong Zhang ierr = DMDAGetLocalInfo(pipe->da, &info);CHKERRQ(ierr); 8242dc13f1SHong Zhang 8342dc13f1SHong Zhang pipe->rad = pipe->D / 2; 8442dc13f1SHong Zhang pipe->A = PETSC_PI*pipe->rad*pipe->rad; 8542dc13f1SHong Zhang pipe->R = pipe->fric / (2*pipe->D*pipe->A); 8642dc13f1SHong Zhang PetscFunctionReturn(0); 8742dc13f1SHong Zhang } 8842dc13f1SHong Zhang 8942dc13f1SHong Zhang /* 9042dc13f1SHong Zhang PipeCreateJacobian - Create Jacobian matrix structures for a Pipe. 9142dc13f1SHong Zhang 9242dc13f1SHong Zhang Collective on Pipe 9342dc13f1SHong Zhang 9442dc13f1SHong Zhang Input Parameter: 9542dc13f1SHong Zhang + pipe - the Pipe object 9642dc13f1SHong Zhang - Jin - array of three constructed Jacobian matrices to be reused. Set NULL if it is not available 9742dc13f1SHong Zhang 9842dc13f1SHong Zhang Output Parameter: 9942dc13f1SHong Zhang . J - array of three empty Jacobian matrices 10042dc13f1SHong Zhang 10142dc13f1SHong Zhang Level: beginner 10242dc13f1SHong Zhang */ 10342dc13f1SHong Zhang PetscErrorCode PipeCreateJacobian(Pipe pipe,Mat *Jin,Mat *J[]) 10442dc13f1SHong Zhang { 10542dc13f1SHong Zhang PetscErrorCode ierr; 10642dc13f1SHong Zhang Mat *Jpipe; 10742dc13f1SHong Zhang PetscInt M,rows[2],cols[2],*nz; 10842dc13f1SHong Zhang PetscScalar *aa; 10942dc13f1SHong Zhang 11042dc13f1SHong Zhang PetscFunctionBegin; 11142dc13f1SHong Zhang if (Jin) { 11242dc13f1SHong Zhang *J = Jin; 11342dc13f1SHong Zhang pipe->jacobian = Jin; 11442dc13f1SHong Zhang ierr = PetscObjectReference((PetscObject)(Jin[0]));CHKERRQ(ierr); 11542dc13f1SHong Zhang PetscFunctionReturn(0); 11642dc13f1SHong Zhang } 11742dc13f1SHong Zhang ierr = PetscMalloc1(3,&Jpipe);CHKERRQ(ierr); 11842dc13f1SHong Zhang 11942dc13f1SHong Zhang /* Jacobian for this pipe */ 12042dc13f1SHong Zhang ierr = DMSetMatrixStructureOnly(pipe->da,PETSC_TRUE);CHKERRQ(ierr); 12142dc13f1SHong Zhang ierr = DMCreateMatrix(pipe->da,&Jpipe[0]);CHKERRQ(ierr); 12242dc13f1SHong Zhang ierr = DMSetMatrixStructureOnly(pipe->da,PETSC_FALSE);CHKERRQ(ierr); 12342dc13f1SHong Zhang 12442dc13f1SHong Zhang /* Jacobian for upstream vertex */ 12542dc13f1SHong Zhang ierr = MatGetSize(Jpipe[0],&M,NULL);CHKERRQ(ierr); 12642dc13f1SHong Zhang ierr = PetscCalloc2(M,&nz,4,&aa);CHKERRQ(ierr); 12742dc13f1SHong Zhang 12842dc13f1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Jpipe[1]);CHKERRQ(ierr); 12942dc13f1SHong Zhang ierr = MatSetSizes(Jpipe[1],PETSC_DECIDE,PETSC_DECIDE,M,2);CHKERRQ(ierr); 13042dc13f1SHong Zhang ierr = MatSetFromOptions(Jpipe[1]);CHKERRQ(ierr); 13142dc13f1SHong Zhang ierr = MatSetOption(Jpipe[1],MAT_STRUCTURE_ONLY,PETSC_TRUE);CHKERRQ(ierr); 13242dc13f1SHong Zhang nz[0] = 2; nz[1] = 2; 13342dc13f1SHong Zhang rows[0] = 0; rows[1] = 1; 13442dc13f1SHong Zhang cols[0] = 0; cols[1] = 1; 13542dc13f1SHong Zhang ierr = MatSeqAIJSetPreallocation(Jpipe[1],0,nz);CHKERRQ(ierr); 13642dc13f1SHong Zhang ierr = MatSetValues(Jpipe[1],2,rows,2,cols,aa,INSERT_VALUES);CHKERRQ(ierr); 13742dc13f1SHong Zhang ierr = MatAssemblyBegin(Jpipe[1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13842dc13f1SHong Zhang ierr = MatAssemblyEnd(Jpipe[1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13942dc13f1SHong Zhang 14042dc13f1SHong Zhang /* Jacobian for downstream vertex */ 14142dc13f1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Jpipe[2]);CHKERRQ(ierr); 14242dc13f1SHong Zhang ierr = MatSetSizes(Jpipe[2],PETSC_DECIDE,PETSC_DECIDE,M,2);CHKERRQ(ierr); 14342dc13f1SHong Zhang ierr = MatSetFromOptions(Jpipe[2]);CHKERRQ(ierr); 14442dc13f1SHong Zhang ierr = MatSetOption(Jpipe[2],MAT_STRUCTURE_ONLY,PETSC_TRUE);CHKERRQ(ierr); 14542dc13f1SHong Zhang nz[0] = 0; nz[1] = 0; nz[M-2] = 2; nz[M-1] = 2; 14642dc13f1SHong Zhang rows[0] = M - 2; rows[1] = M - 1; 14742dc13f1SHong Zhang ierr = MatSeqAIJSetPreallocation(Jpipe[2],0,nz);CHKERRQ(ierr); 14842dc13f1SHong Zhang ierr = MatSetValues(Jpipe[2],2,rows,2,cols,aa,INSERT_VALUES);CHKERRQ(ierr); 14942dc13f1SHong Zhang ierr = MatAssemblyBegin(Jpipe[2],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15042dc13f1SHong Zhang ierr = MatAssemblyEnd(Jpipe[2],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15142dc13f1SHong Zhang 15242dc13f1SHong Zhang ierr = PetscFree2(nz,aa);CHKERRQ(ierr); 15342dc13f1SHong Zhang 15442dc13f1SHong Zhang *J = Jpipe; 15542dc13f1SHong Zhang pipe->jacobian = Jpipe; 15642dc13f1SHong Zhang PetscFunctionReturn(0); 15742dc13f1SHong Zhang } 15842dc13f1SHong Zhang 15942dc13f1SHong Zhang PetscErrorCode PipeDestroyJacobian(Pipe pipe) 16042dc13f1SHong Zhang { 16142dc13f1SHong Zhang PetscErrorCode ierr; 16242dc13f1SHong Zhang Mat *Jpipe = pipe->jacobian; 16342dc13f1SHong Zhang PetscInt i; 16442dc13f1SHong Zhang 16542dc13f1SHong Zhang PetscFunctionBegin; 16642dc13f1SHong Zhang if (Jpipe) { 16742dc13f1SHong Zhang for (i=0; i<3; i++) { 16842dc13f1SHong Zhang ierr = MatDestroy(&Jpipe[i]);CHKERRQ(ierr); 16942dc13f1SHong Zhang } 17042dc13f1SHong Zhang } 17142dc13f1SHong Zhang ierr = PetscFree(Jpipe);CHKERRQ(ierr); 17242dc13f1SHong Zhang PetscFunctionReturn(0); 17342dc13f1SHong Zhang } 17442dc13f1SHong Zhang 17542dc13f1SHong Zhang /* 17642dc13f1SHong Zhang JunctionCreateJacobian - Create Jacobian matrices for a vertex. 17742dc13f1SHong Zhang 17842dc13f1SHong Zhang Collective on Pipe 17942dc13f1SHong Zhang 18042dc13f1SHong Zhang Input Parameter: 18142dc13f1SHong Zhang + dm - the DMNetwork object 18242dc13f1SHong Zhang . v - vertex point 18342dc13f1SHong Zhang - Jin - Jacobian patterns created by JunctionCreateJacobianSample() for reuse 18442dc13f1SHong Zhang 18542dc13f1SHong Zhang Output Parameter: 18642dc13f1SHong Zhang . J - array of Jacobian matrices (see dmnetworkimpl.h) 18742dc13f1SHong Zhang 18842dc13f1SHong Zhang Level: beginner 18942dc13f1SHong Zhang */ 19042dc13f1SHong Zhang PetscErrorCode JunctionCreateJacobian(DM dm,PetscInt v,Mat *Jin,Mat *J[]) 19142dc13f1SHong Zhang { 19242dc13f1SHong Zhang PetscErrorCode ierr; 19342dc13f1SHong Zhang Mat *Jv; 19442dc13f1SHong Zhang PetscInt nedges,e,i,M,N,*rows,*cols; 19542dc13f1SHong Zhang PetscBool isSelf; 19642dc13f1SHong Zhang const PetscInt *edges,*cone; 19742dc13f1SHong Zhang PetscScalar *zeros; 19842dc13f1SHong Zhang 19942dc13f1SHong Zhang PetscFunctionBegin; 200a5b23f4aSJose E. Roman /* Get array size of Jv */ 20142dc13f1SHong Zhang ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 202*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(nedges <= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"%d vertex, nedges %d",v,nedges); 20342dc13f1SHong Zhang 20442dc13f1SHong Zhang /* two Jacobians for each connected edge: J(v,e) and J(v,vc); adding J(v,v), total 2*nedges+1 Jacobians */ 20542dc13f1SHong Zhang ierr = PetscCalloc1(2*nedges+1,&Jv);CHKERRQ(ierr); 20642dc13f1SHong Zhang 20742dc13f1SHong Zhang /* Create dense zero block for this vertex: J[0] = Jacobian(v,v) */ 20842dc13f1SHong Zhang ierr = DMNetworkGetComponent(dm,v,-1,NULL,NULL,&M);CHKERRQ(ierr); 209*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(M !=2,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"M != 2",M); 21042dc13f1SHong Zhang ierr = PetscMalloc3(M,&rows,M,&cols,M*M,&zeros);CHKERRQ(ierr); 21142dc13f1SHong Zhang ierr = PetscArrayzero(zeros,M*M);CHKERRQ(ierr); 21242dc13f1SHong Zhang for (i=0; i<M; i++) rows[i] = i; 21342dc13f1SHong Zhang 21442dc13f1SHong Zhang for (e=0; e<nedges; e++) { 21542dc13f1SHong Zhang /* create Jv[2*e+1] = Jacobian(v,e), e: supporting edge */ 21642dc13f1SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr); 21742dc13f1SHong Zhang isSelf = (v == cone[0]) ? PETSC_TRUE:PETSC_FALSE; 21842dc13f1SHong Zhang 21942dc13f1SHong Zhang if (Jin) { 22042dc13f1SHong Zhang if (isSelf) { 22142dc13f1SHong Zhang Jv[2*e+1] = Jin[0]; 22242dc13f1SHong Zhang } else { 22342dc13f1SHong Zhang Jv[2*e+1] = Jin[1]; 22442dc13f1SHong Zhang } 22542dc13f1SHong Zhang Jv[2*e+2] = Jin[2]; 22642dc13f1SHong Zhang ierr = PetscObjectReference((PetscObject)(Jv[2*e+1]));CHKERRQ(ierr); 22742dc13f1SHong Zhang ierr = PetscObjectReference((PetscObject)(Jv[2*e+2]));CHKERRQ(ierr); 22842dc13f1SHong Zhang } else { 22942dc13f1SHong Zhang /* create J(v,e) */ 23042dc13f1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Jv[2*e+1]);CHKERRQ(ierr); 23142dc13f1SHong Zhang ierr = DMNetworkGetComponent(dm,edges[e],-1,NULL,NULL,&N);CHKERRQ(ierr); 23242dc13f1SHong Zhang ierr = MatSetSizes(Jv[2*e+1],PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 23342dc13f1SHong Zhang ierr = MatSetFromOptions(Jv[2*e+1]);CHKERRQ(ierr); 23442dc13f1SHong Zhang ierr = MatSetOption(Jv[2*e+1],MAT_STRUCTURE_ONLY,PETSC_TRUE);CHKERRQ(ierr); 23542dc13f1SHong Zhang ierr = MatSeqAIJSetPreallocation(Jv[2*e+1],2,NULL);CHKERRQ(ierr); 23642dc13f1SHong Zhang if (N) { 23742dc13f1SHong Zhang if (isSelf) { /* coupling at upstream */ 23842dc13f1SHong Zhang for (i=0; i<2; i++) cols[i] = i; 23942dc13f1SHong Zhang } else { /* coupling at downstream */ 24042dc13f1SHong Zhang cols[0] = N-2; cols[1] = N-1; 24142dc13f1SHong Zhang } 24242dc13f1SHong Zhang ierr = MatSetValues(Jv[2*e+1],2,rows,2,cols,zeros,INSERT_VALUES);CHKERRQ(ierr); 24342dc13f1SHong Zhang } 24442dc13f1SHong Zhang ierr = MatAssemblyBegin(Jv[2*e+1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 24542dc13f1SHong Zhang ierr = MatAssemblyEnd(Jv[2*e+1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 24642dc13f1SHong Zhang 24742dc13f1SHong Zhang /* create Jv[2*e+2] = Jacobian(v,vc), vc: connected vertex. 24842dc13f1SHong Zhang In WashNetwork, v and vc are not connected, thus Jacobian(v,vc) is empty */ 24942dc13f1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Jv[2*e+2]);CHKERRQ(ierr); 25042dc13f1SHong Zhang ierr = MatSetSizes(Jv[2*e+2],PETSC_DECIDE,PETSC_DECIDE,M,M);CHKERRQ(ierr); /* empty matrix, sizes can be arbitrary */ 25142dc13f1SHong Zhang ierr = MatSetFromOptions(Jv[2*e+2]);CHKERRQ(ierr); 25242dc13f1SHong Zhang ierr = MatSetOption(Jv[2*e+2],MAT_STRUCTURE_ONLY,PETSC_TRUE);CHKERRQ(ierr); 25342dc13f1SHong Zhang ierr = MatSeqAIJSetPreallocation(Jv[2*e+2],1,NULL);CHKERRQ(ierr); 25442dc13f1SHong Zhang ierr = MatAssemblyBegin(Jv[2*e+2],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 25542dc13f1SHong Zhang ierr = MatAssemblyEnd(Jv[2*e+2],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 25642dc13f1SHong Zhang } 25742dc13f1SHong Zhang } 25842dc13f1SHong Zhang ierr = PetscFree3(rows,cols,zeros);CHKERRQ(ierr); 25942dc13f1SHong Zhang 26042dc13f1SHong Zhang *J = Jv; 26142dc13f1SHong Zhang PetscFunctionReturn(0); 26242dc13f1SHong Zhang } 26342dc13f1SHong Zhang 26442dc13f1SHong Zhang PetscErrorCode JunctionDestroyJacobian(DM dm,PetscInt v,Junction junc) 26542dc13f1SHong Zhang { 26642dc13f1SHong Zhang PetscErrorCode ierr; 26742dc13f1SHong Zhang Mat *Jv=junc->jacobian; 26842dc13f1SHong Zhang const PetscInt *edges; 26942dc13f1SHong Zhang PetscInt nedges,e; 27042dc13f1SHong Zhang 27142dc13f1SHong Zhang PetscFunctionBegin; 27242dc13f1SHong Zhang if (!Jv) PetscFunctionReturn(0); 27342dc13f1SHong Zhang 27442dc13f1SHong Zhang ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 27542dc13f1SHong Zhang for (e=0; e<nedges; e++) { 27642dc13f1SHong Zhang ierr = MatDestroy(&Jv[2*e+1]);CHKERRQ(ierr); 27742dc13f1SHong Zhang ierr = MatDestroy(&Jv[2*e+2]);CHKERRQ(ierr); 27842dc13f1SHong Zhang } 27942dc13f1SHong Zhang ierr = PetscFree(Jv);CHKERRQ(ierr); 28042dc13f1SHong Zhang PetscFunctionReturn(0); 28142dc13f1SHong Zhang } 282