1*42dc13f1SHong Zhang #include "wash.h" 2*42dc13f1SHong Zhang 3*42dc13f1SHong Zhang /* Subroutines for Pipe */ 4*42dc13f1SHong Zhang /* -------------------------------------------------------*/ 5*42dc13f1SHong Zhang 6*42dc13f1SHong Zhang /* 7*42dc13f1SHong Zhang PipeCreate - Create Pipe object. 8*42dc13f1SHong Zhang 9*42dc13f1SHong Zhang Input Parameters: 10*42dc13f1SHong Zhang comm - MPI communicator 11*42dc13f1SHong Zhang 12*42dc13f1SHong Zhang Output Parameter: 13*42dc13f1SHong Zhang . pipe - location to put the PIPE context 14*42dc13f1SHong Zhang */ 15*42dc13f1SHong Zhang PetscErrorCode PipeCreate(MPI_Comm comm,Pipe *pipe) 16*42dc13f1SHong Zhang { 17*42dc13f1SHong Zhang PetscErrorCode ierr; 18*42dc13f1SHong Zhang 19*42dc13f1SHong Zhang PetscFunctionBegin; 20*42dc13f1SHong Zhang ierr = PetscNew(pipe);CHKERRQ(ierr); 21*42dc13f1SHong Zhang PetscFunctionReturn(0); 22*42dc13f1SHong Zhang } 23*42dc13f1SHong Zhang 24*42dc13f1SHong Zhang /* 25*42dc13f1SHong Zhang PipeDestroy - Destroy Pipe object. 26*42dc13f1SHong Zhang 27*42dc13f1SHong Zhang Input Parameters: 28*42dc13f1SHong Zhang pipe - Reference to pipe intended to be destroyed. 29*42dc13f1SHong Zhang */ 30*42dc13f1SHong Zhang PetscErrorCode PipeDestroy(Pipe *pipe) 31*42dc13f1SHong Zhang { 32*42dc13f1SHong Zhang PetscErrorCode ierr; 33*42dc13f1SHong Zhang 34*42dc13f1SHong Zhang PetscFunctionBegin; 35*42dc13f1SHong Zhang if (!*pipe) PetscFunctionReturn(0); 36*42dc13f1SHong Zhang 37*42dc13f1SHong Zhang ierr = PipeDestroyJacobian(*pipe);CHKERRQ(ierr); 38*42dc13f1SHong Zhang ierr = VecDestroy(&(*pipe)->x);CHKERRQ(ierr); 39*42dc13f1SHong Zhang ierr = DMDestroy(&(*pipe)->da);CHKERRQ(ierr); 40*42dc13f1SHong Zhang PetscFunctionReturn(0); 41*42dc13f1SHong Zhang } 42*42dc13f1SHong Zhang 43*42dc13f1SHong Zhang /* 44*42dc13f1SHong Zhang PipeSetParameters - Set parameters for Pipe context 45*42dc13f1SHong Zhang 46*42dc13f1SHong Zhang Input Parameter: 47*42dc13f1SHong Zhang + pipe - PIPE object 48*42dc13f1SHong Zhang . length - 49*42dc13f1SHong Zhang . nnodes - 50*42dc13f1SHong Zhang . D - 51*42dc13f1SHong Zhang . a - 52*42dc13f1SHong Zhang - fric - 53*42dc13f1SHong Zhang */ 54*42dc13f1SHong Zhang PetscErrorCode PipeSetParameters(Pipe pipe,PetscReal length,PetscReal D,PetscReal a,PetscReal fric) 55*42dc13f1SHong Zhang { 56*42dc13f1SHong Zhang PetscFunctionBegin; 57*42dc13f1SHong Zhang pipe->length = length; 58*42dc13f1SHong Zhang pipe->D = D; 59*42dc13f1SHong Zhang pipe->a = a; 60*42dc13f1SHong Zhang pipe->fric = fric; 61*42dc13f1SHong Zhang PetscFunctionReturn(0); 62*42dc13f1SHong Zhang } 63*42dc13f1SHong Zhang 64*42dc13f1SHong Zhang /* 65*42dc13f1SHong Zhang PipeSetUp - Set up pipe based on set parameters. 66*42dc13f1SHong Zhang */ 67*42dc13f1SHong Zhang PetscErrorCode PipeSetUp(Pipe pipe) 68*42dc13f1SHong Zhang { 69*42dc13f1SHong Zhang DMDALocalInfo info; 70*42dc13f1SHong Zhang PetscErrorCode ierr; 71*42dc13f1SHong Zhang 72*42dc13f1SHong Zhang PetscFunctionBegin; 73*42dc13f1SHong Zhang ierr = DMDACreate1d(PETSC_COMM_SELF, DM_BOUNDARY_GHOSTED, pipe->nnodes, 2, 1, NULL, &pipe->da);CHKERRQ(ierr); 74*42dc13f1SHong Zhang ierr = DMSetFromOptions(pipe->da);CHKERRQ(ierr); 75*42dc13f1SHong Zhang ierr = DMSetUp(pipe->da);CHKERRQ(ierr); 76*42dc13f1SHong Zhang ierr = DMDASetFieldName(pipe->da, 0, "Q");CHKERRQ(ierr); 77*42dc13f1SHong Zhang ierr = DMDASetFieldName(pipe->da, 1, "H");CHKERRQ(ierr); 78*42dc13f1SHong Zhang ierr = DMDASetUniformCoordinates(pipe->da, 0, pipe->length, 0, 0, 0, 0);CHKERRQ(ierr); 79*42dc13f1SHong Zhang ierr = DMCreateGlobalVector(pipe->da, &(pipe->x));CHKERRQ(ierr); 80*42dc13f1SHong Zhang 81*42dc13f1SHong Zhang ierr = DMDAGetLocalInfo(pipe->da, &info);CHKERRQ(ierr); 82*42dc13f1SHong Zhang 83*42dc13f1SHong Zhang pipe->rad = pipe->D / 2; 84*42dc13f1SHong Zhang pipe->A = PETSC_PI*pipe->rad*pipe->rad; 85*42dc13f1SHong Zhang pipe->R = pipe->fric / (2*pipe->D*pipe->A); 86*42dc13f1SHong Zhang PetscFunctionReturn(0); 87*42dc13f1SHong Zhang } 88*42dc13f1SHong Zhang 89*42dc13f1SHong Zhang /* 90*42dc13f1SHong Zhang PipeCreateJacobian - Create Jacobian matrix structures for a Pipe. 91*42dc13f1SHong Zhang 92*42dc13f1SHong Zhang Collective on Pipe 93*42dc13f1SHong Zhang 94*42dc13f1SHong Zhang Input Parameter: 95*42dc13f1SHong Zhang + pipe - the Pipe object 96*42dc13f1SHong Zhang - Jin - array of three constructed Jacobian matrices to be reused. Set NULL if it is not available 97*42dc13f1SHong Zhang 98*42dc13f1SHong Zhang Output Parameter: 99*42dc13f1SHong Zhang . J - array of three empty Jacobian matrices 100*42dc13f1SHong Zhang 101*42dc13f1SHong Zhang Level: beginner 102*42dc13f1SHong Zhang */ 103*42dc13f1SHong Zhang PetscErrorCode PipeCreateJacobian(Pipe pipe,Mat *Jin,Mat *J[]) 104*42dc13f1SHong Zhang { 105*42dc13f1SHong Zhang PetscErrorCode ierr; 106*42dc13f1SHong Zhang Mat *Jpipe; 107*42dc13f1SHong Zhang PetscInt M,rows[2],cols[2],*nz; 108*42dc13f1SHong Zhang PetscScalar *aa; 109*42dc13f1SHong Zhang 110*42dc13f1SHong Zhang PetscFunctionBegin; 111*42dc13f1SHong Zhang if (Jin) { 112*42dc13f1SHong Zhang *J = Jin; 113*42dc13f1SHong Zhang pipe->jacobian = Jin; 114*42dc13f1SHong Zhang ierr = PetscObjectReference((PetscObject)(Jin[0]));CHKERRQ(ierr); 115*42dc13f1SHong Zhang PetscFunctionReturn(0); 116*42dc13f1SHong Zhang } 117*42dc13f1SHong Zhang ierr = PetscMalloc1(3,&Jpipe);CHKERRQ(ierr); 118*42dc13f1SHong Zhang 119*42dc13f1SHong Zhang /* Jacobian for this pipe */ 120*42dc13f1SHong Zhang ierr = DMSetMatrixStructureOnly(pipe->da,PETSC_TRUE);CHKERRQ(ierr); 121*42dc13f1SHong Zhang ierr = DMCreateMatrix(pipe->da,&Jpipe[0]);CHKERRQ(ierr); 122*42dc13f1SHong Zhang ierr = DMSetMatrixStructureOnly(pipe->da,PETSC_FALSE);CHKERRQ(ierr); 123*42dc13f1SHong Zhang 124*42dc13f1SHong Zhang /* Jacobian for upstream vertex */ 125*42dc13f1SHong Zhang ierr = MatGetSize(Jpipe[0],&M,NULL);CHKERRQ(ierr); 126*42dc13f1SHong Zhang ierr = PetscCalloc2(M,&nz,4,&aa);CHKERRQ(ierr); 127*42dc13f1SHong Zhang 128*42dc13f1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Jpipe[1]);CHKERRQ(ierr); 129*42dc13f1SHong Zhang ierr = MatSetSizes(Jpipe[1],PETSC_DECIDE,PETSC_DECIDE,M,2);CHKERRQ(ierr); 130*42dc13f1SHong Zhang ierr = MatSetFromOptions(Jpipe[1]);CHKERRQ(ierr); 131*42dc13f1SHong Zhang ierr = MatSetOption(Jpipe[1],MAT_STRUCTURE_ONLY,PETSC_TRUE);CHKERRQ(ierr); 132*42dc13f1SHong Zhang nz[0] = 2; nz[1] = 2; 133*42dc13f1SHong Zhang rows[0] = 0; rows[1] = 1; 134*42dc13f1SHong Zhang cols[0] = 0; cols[1] = 1; 135*42dc13f1SHong Zhang ierr = MatSeqAIJSetPreallocation(Jpipe[1],0,nz);CHKERRQ(ierr); 136*42dc13f1SHong Zhang ierr = MatSetValues(Jpipe[1],2,rows,2,cols,aa,INSERT_VALUES);CHKERRQ(ierr); 137*42dc13f1SHong Zhang ierr = MatAssemblyBegin(Jpipe[1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 138*42dc13f1SHong Zhang ierr = MatAssemblyEnd(Jpipe[1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 139*42dc13f1SHong Zhang 140*42dc13f1SHong Zhang /* Jacobian for downstream vertex */ 141*42dc13f1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Jpipe[2]);CHKERRQ(ierr); 142*42dc13f1SHong Zhang ierr = MatSetSizes(Jpipe[2],PETSC_DECIDE,PETSC_DECIDE,M,2);CHKERRQ(ierr); 143*42dc13f1SHong Zhang ierr = MatSetFromOptions(Jpipe[2]);CHKERRQ(ierr); 144*42dc13f1SHong Zhang ierr = MatSetOption(Jpipe[2],MAT_STRUCTURE_ONLY,PETSC_TRUE);CHKERRQ(ierr); 145*42dc13f1SHong Zhang nz[0] = 0; nz[1] = 0; nz[M-2] = 2; nz[M-1] = 2; 146*42dc13f1SHong Zhang rows[0] = M - 2; rows[1] = M - 1; 147*42dc13f1SHong Zhang ierr = MatSeqAIJSetPreallocation(Jpipe[2],0,nz);CHKERRQ(ierr); 148*42dc13f1SHong Zhang ierr = MatSetValues(Jpipe[2],2,rows,2,cols,aa,INSERT_VALUES);CHKERRQ(ierr); 149*42dc13f1SHong Zhang ierr = MatAssemblyBegin(Jpipe[2],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 150*42dc13f1SHong Zhang ierr = MatAssemblyEnd(Jpipe[2],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 151*42dc13f1SHong Zhang 152*42dc13f1SHong Zhang ierr = PetscFree2(nz,aa);CHKERRQ(ierr); 153*42dc13f1SHong Zhang 154*42dc13f1SHong Zhang *J = Jpipe; 155*42dc13f1SHong Zhang pipe->jacobian = Jpipe; 156*42dc13f1SHong Zhang PetscFunctionReturn(0); 157*42dc13f1SHong Zhang } 158*42dc13f1SHong Zhang 159*42dc13f1SHong Zhang PetscErrorCode PipeDestroyJacobian(Pipe pipe) 160*42dc13f1SHong Zhang { 161*42dc13f1SHong Zhang PetscErrorCode ierr; 162*42dc13f1SHong Zhang Mat *Jpipe = pipe->jacobian; 163*42dc13f1SHong Zhang PetscInt i; 164*42dc13f1SHong Zhang 165*42dc13f1SHong Zhang PetscFunctionBegin; 166*42dc13f1SHong Zhang if (Jpipe) { 167*42dc13f1SHong Zhang for (i=0; i<3; i++) { 168*42dc13f1SHong Zhang ierr = MatDestroy(&Jpipe[i]);CHKERRQ(ierr); 169*42dc13f1SHong Zhang } 170*42dc13f1SHong Zhang } 171*42dc13f1SHong Zhang ierr = PetscFree(Jpipe);CHKERRQ(ierr); 172*42dc13f1SHong Zhang PetscFunctionReturn(0); 173*42dc13f1SHong Zhang } 174*42dc13f1SHong Zhang 175*42dc13f1SHong Zhang /* 176*42dc13f1SHong Zhang JunctionCreateJacobian - Create Jacobian matrices for a vertex. 177*42dc13f1SHong Zhang 178*42dc13f1SHong Zhang Collective on Pipe 179*42dc13f1SHong Zhang 180*42dc13f1SHong Zhang Input Parameter: 181*42dc13f1SHong Zhang + dm - the DMNetwork object 182*42dc13f1SHong Zhang . v - vertex point 183*42dc13f1SHong Zhang - Jin - Jacobian patterns created by JunctionCreateJacobianSample() for reuse 184*42dc13f1SHong Zhang 185*42dc13f1SHong Zhang Output Parameter: 186*42dc13f1SHong Zhang . J - array of Jacobian matrices (see dmnetworkimpl.h) 187*42dc13f1SHong Zhang 188*42dc13f1SHong Zhang Level: beginner 189*42dc13f1SHong Zhang */ 190*42dc13f1SHong Zhang PetscErrorCode JunctionCreateJacobian(DM dm,PetscInt v,Mat *Jin,Mat *J[]) 191*42dc13f1SHong Zhang { 192*42dc13f1SHong Zhang PetscErrorCode ierr; 193*42dc13f1SHong Zhang Mat *Jv; 194*42dc13f1SHong Zhang PetscInt nedges,e,i,M,N,*rows,*cols; 195*42dc13f1SHong Zhang PetscBool isSelf; 196*42dc13f1SHong Zhang const PetscInt *edges,*cone; 197*42dc13f1SHong Zhang PetscScalar *zeros; 198*42dc13f1SHong Zhang 199*42dc13f1SHong Zhang PetscFunctionBegin; 200*42dc13f1SHong Zhang /* Get arrary size of Jv */ 201*42dc13f1SHong Zhang ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 202*42dc13f1SHong Zhang if (nedges <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"%d vertex, nedges %d\n",v,nedges); 203*42dc13f1SHong Zhang 204*42dc13f1SHong Zhang /* two Jacobians for each connected edge: J(v,e) and J(v,vc); adding J(v,v), total 2*nedges+1 Jacobians */ 205*42dc13f1SHong Zhang ierr = PetscCalloc1(2*nedges+1,&Jv);CHKERRQ(ierr); 206*42dc13f1SHong Zhang 207*42dc13f1SHong Zhang /* Create dense zero block for this vertex: J[0] = Jacobian(v,v) */ 208*42dc13f1SHong Zhang ierr = DMNetworkGetComponent(dm,v,-1,NULL,NULL,&M);CHKERRQ(ierr); 209*42dc13f1SHong Zhang if (M !=2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"M != 2",M); 210*42dc13f1SHong Zhang ierr = PetscMalloc3(M,&rows,M,&cols,M*M,&zeros);CHKERRQ(ierr); 211*42dc13f1SHong Zhang ierr = PetscArrayzero(zeros,M*M);CHKERRQ(ierr); 212*42dc13f1SHong Zhang for (i=0; i<M; i++) rows[i] = i; 213*42dc13f1SHong Zhang 214*42dc13f1SHong Zhang for (e=0; e<nedges; e++) { 215*42dc13f1SHong Zhang /* create Jv[2*e+1] = Jacobian(v,e), e: supporting edge */ 216*42dc13f1SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr); 217*42dc13f1SHong Zhang isSelf = (v == cone[0]) ? PETSC_TRUE:PETSC_FALSE; 218*42dc13f1SHong Zhang 219*42dc13f1SHong Zhang if (Jin) { 220*42dc13f1SHong Zhang if (isSelf) { 221*42dc13f1SHong Zhang Jv[2*e+1] = Jin[0]; 222*42dc13f1SHong Zhang } else { 223*42dc13f1SHong Zhang Jv[2*e+1] = Jin[1]; 224*42dc13f1SHong Zhang } 225*42dc13f1SHong Zhang Jv[2*e+2] = Jin[2]; 226*42dc13f1SHong Zhang ierr = PetscObjectReference((PetscObject)(Jv[2*e+1]));CHKERRQ(ierr); 227*42dc13f1SHong Zhang ierr = PetscObjectReference((PetscObject)(Jv[2*e+2]));CHKERRQ(ierr); 228*42dc13f1SHong Zhang } else { 229*42dc13f1SHong Zhang /* create J(v,e) */ 230*42dc13f1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Jv[2*e+1]);CHKERRQ(ierr); 231*42dc13f1SHong Zhang ierr = DMNetworkGetComponent(dm,edges[e],-1,NULL,NULL,&N);CHKERRQ(ierr); 232*42dc13f1SHong Zhang ierr = MatSetSizes(Jv[2*e+1],PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 233*42dc13f1SHong Zhang ierr = MatSetFromOptions(Jv[2*e+1]);CHKERRQ(ierr); 234*42dc13f1SHong Zhang ierr = MatSetOption(Jv[2*e+1],MAT_STRUCTURE_ONLY,PETSC_TRUE);CHKERRQ(ierr); 235*42dc13f1SHong Zhang ierr = MatSeqAIJSetPreallocation(Jv[2*e+1],2,NULL);CHKERRQ(ierr); 236*42dc13f1SHong Zhang if (N) { 237*42dc13f1SHong Zhang if (isSelf) { /* coupling at upstream */ 238*42dc13f1SHong Zhang for (i=0; i<2; i++) cols[i] = i; 239*42dc13f1SHong Zhang } else { /* coupling at downstream */ 240*42dc13f1SHong Zhang cols[0] = N-2; cols[1] = N-1; 241*42dc13f1SHong Zhang } 242*42dc13f1SHong Zhang ierr = MatSetValues(Jv[2*e+1],2,rows,2,cols,zeros,INSERT_VALUES);CHKERRQ(ierr); 243*42dc13f1SHong Zhang } 244*42dc13f1SHong Zhang ierr = MatAssemblyBegin(Jv[2*e+1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 245*42dc13f1SHong Zhang ierr = MatAssemblyEnd(Jv[2*e+1],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 246*42dc13f1SHong Zhang 247*42dc13f1SHong Zhang /* create Jv[2*e+2] = Jacobian(v,vc), vc: connected vertex. 248*42dc13f1SHong Zhang In WashNetwork, v and vc are not connected, thus Jacobian(v,vc) is empty */ 249*42dc13f1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Jv[2*e+2]);CHKERRQ(ierr); 250*42dc13f1SHong Zhang ierr = MatSetSizes(Jv[2*e+2],PETSC_DECIDE,PETSC_DECIDE,M,M);CHKERRQ(ierr); /* empty matrix, sizes can be arbitrary */ 251*42dc13f1SHong Zhang ierr = MatSetFromOptions(Jv[2*e+2]);CHKERRQ(ierr); 252*42dc13f1SHong Zhang ierr = MatSetOption(Jv[2*e+2],MAT_STRUCTURE_ONLY,PETSC_TRUE);CHKERRQ(ierr); 253*42dc13f1SHong Zhang ierr = MatSeqAIJSetPreallocation(Jv[2*e+2],1,NULL);CHKERRQ(ierr); 254*42dc13f1SHong Zhang ierr = MatAssemblyBegin(Jv[2*e+2],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 255*42dc13f1SHong Zhang ierr = MatAssemblyEnd(Jv[2*e+2],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 256*42dc13f1SHong Zhang } 257*42dc13f1SHong Zhang } 258*42dc13f1SHong Zhang ierr = PetscFree3(rows,cols,zeros);CHKERRQ(ierr); 259*42dc13f1SHong Zhang 260*42dc13f1SHong Zhang *J = Jv; 261*42dc13f1SHong Zhang PetscFunctionReturn(0); 262*42dc13f1SHong Zhang } 263*42dc13f1SHong Zhang 264*42dc13f1SHong Zhang PetscErrorCode JunctionDestroyJacobian(DM dm,PetscInt v,Junction junc) 265*42dc13f1SHong Zhang { 266*42dc13f1SHong Zhang PetscErrorCode ierr; 267*42dc13f1SHong Zhang Mat *Jv=junc->jacobian; 268*42dc13f1SHong Zhang const PetscInt *edges; 269*42dc13f1SHong Zhang PetscInt nedges,e; 270*42dc13f1SHong Zhang 271*42dc13f1SHong Zhang PetscFunctionBegin; 272*42dc13f1SHong Zhang if (!Jv) PetscFunctionReturn(0); 273*42dc13f1SHong Zhang 274*42dc13f1SHong Zhang ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 275*42dc13f1SHong Zhang for (e=0; e<nedges; e++) { 276*42dc13f1SHong Zhang ierr = MatDestroy(&Jv[2*e+1]);CHKERRQ(ierr); 277*42dc13f1SHong Zhang ierr = MatDestroy(&Jv[2*e+2]);CHKERRQ(ierr); 278*42dc13f1SHong Zhang } 279*42dc13f1SHong Zhang ierr = PetscFree(Jv);CHKERRQ(ierr); 280*42dc13f1SHong Zhang PetscFunctionReturn(0); 281*42dc13f1SHong Zhang } 282