15f2c45f1SShri Abhyankar #define PETSCDM_DLL 2af0996ceSBarry Smith #include <petsc/private/dmnetworkimpl.h> /*I "petscdmnetwork.h" I*/ 35f2c45f1SShri Abhyankar #include <petscdmda.h> 45f2c45f1SShri Abhyankar 54416b707SBarry Smith PetscErrorCode DMSetFromOptions_Network(PetscOptionItems *PetscOptionsObject,DM dm) 65f2c45f1SShri Abhyankar { 75f2c45f1SShri Abhyankar PetscErrorCode ierr; 85f2c45f1SShri Abhyankar 95f2c45f1SShri Abhyankar PetscFunctionBegin; 105f2c45f1SShri Abhyankar PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 111a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"DMNetwork Options");CHKERRQ(ierr); 125f2c45f1SShri Abhyankar ierr = PetscOptionsTail();CHKERRQ(ierr); 135f2c45f1SShri Abhyankar PetscFunctionReturn(0); 145f2c45f1SShri Abhyankar } 155f2c45f1SShri Abhyankar 165f2c45f1SShri Abhyankar /* External function declarations here */ 175f2c45f1SShri Abhyankar extern PetscErrorCode DMCreateMatrix_Network(DM, Mat*); 185f2c45f1SShri Abhyankar extern PetscErrorCode DMDestroy_Network(DM); 195f2c45f1SShri Abhyankar extern PetscErrorCode DMView_Network(DM, PetscViewer); 205f2c45f1SShri Abhyankar extern PetscErrorCode DMGlobalToLocalBegin_Network(DM, Vec, InsertMode, Vec); 215f2c45f1SShri Abhyankar extern PetscErrorCode DMGlobalToLocalEnd_Network(DM, Vec, InsertMode, Vec); 225f2c45f1SShri Abhyankar extern PetscErrorCode DMLocalToGlobalBegin_Network(DM, Vec, InsertMode, Vec); 235f2c45f1SShri Abhyankar extern PetscErrorCode DMLocalToGlobalEnd_Network(DM, Vec, InsertMode, Vec); 245f2c45f1SShri Abhyankar extern PetscErrorCode DMSetUp_Network(DM); 258415c774SShri Abhyankar extern PetscErrorCode DMClone_Network(DM, DM*); 265f2c45f1SShri Abhyankar 274dc485aaSHong Zhang static PetscErrorCode VecView_Network_Seq(DM networkdm,Vec X,PetscViewer viewer) 284dc485aaSHong Zhang { 294dc485aaSHong Zhang PetscErrorCode ierr; 304dc485aaSHong Zhang PetscInt i,e,v,Start,End,offset,nvar; 314dc485aaSHong Zhang const PetscScalar *xv; 324dc485aaSHong Zhang 334dc485aaSHong Zhang PetscFunctionBegin; 344dc485aaSHong Zhang ierr = VecGetArrayRead(X,&xv);CHKERRQ(ierr); 354dc485aaSHong Zhang 364dc485aaSHong Zhang /* iterate over edges */ 374dc485aaSHong Zhang ierr = DMNetworkGetEdgeRange(networkdm,&Start,&End);CHKERRQ(ierr); 384dc485aaSHong Zhang for (e=Start; e<End; e++) { 394dc485aaSHong Zhang ierr = DMNetworkGetVariableOffset(networkdm,e,&offset);CHKERRQ(ierr); 404dc485aaSHong Zhang ierr = DMNetworkGetNumVariables(networkdm,e,&nvar);CHKERRQ(ierr); 414dc485aaSHong Zhang if (!nvar) continue; 424dc485aaSHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Edge %D:\n",e-Start);CHKERRQ(ierr); 434dc485aaSHong Zhang for (i=offset; i< offset+nvar; i++) { 444dc485aaSHong Zhang #if defined(PETSC_USE_COMPLEX) 454dc485aaSHong Zhang if (PetscImaginaryPart(xv[i]) > 0.0) { 464dc485aaSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g + %g i\n",(double)PetscRealPart(xv[i]),(double)PetscImaginaryPart(xv[i]));CHKERRQ(ierr); 474dc485aaSHong Zhang } else if (PetscImaginaryPart(xv[i]) < 0.0) { 484dc485aaSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g - %g i\n",(double)PetscRealPart(xv[i]),-(double)PetscImaginaryPart(xv[i]));CHKERRQ(ierr); 494dc485aaSHong Zhang } else { 504dc485aaSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g\n",(double)PetscRealPart(xv[i]));CHKERRQ(ierr); 514dc485aaSHong Zhang } 524dc485aaSHong Zhang #else 534dc485aaSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g\n",(double)xv[i]);CHKERRQ(ierr); 544dc485aaSHong Zhang #endif 554dc485aaSHong Zhang } 564dc485aaSHong Zhang } 574dc485aaSHong Zhang 584dc485aaSHong Zhang /* iterate over vertices */ 594dc485aaSHong Zhang ierr = DMNetworkGetVertexRange(networkdm,&Start,&End);CHKERRQ(ierr); 604dc485aaSHong Zhang for (v=Start; v<End; v++) { 614dc485aaSHong Zhang ierr = DMNetworkGetVariableOffset(networkdm,v,&offset);CHKERRQ(ierr); 624dc485aaSHong Zhang ierr = DMNetworkGetNumVariables(networkdm,v,&nvar);CHKERRQ(ierr); 634dc485aaSHong Zhang if (!nvar) continue; 644dc485aaSHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Vertex %D:\n",v-Start);CHKERRQ(ierr); 654dc485aaSHong Zhang for (i=offset; i< offset+nvar; i++) { 664dc485aaSHong Zhang #if defined(PETSC_USE_COMPLEX) 674dc485aaSHong Zhang if (PetscImaginaryPart(xv[i]) > 0.0) { 684dc485aaSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g + %g i\n",(double)PetscRealPart(xv[i]),(double)PetscImaginaryPart(xv[i]));CHKERRQ(ierr); 694dc485aaSHong Zhang } else if (PetscImaginaryPart(xv[i]) < 0.0) { 704dc485aaSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g - %g i\n",(double)PetscRealPart(xv[i]),-(double)PetscImaginaryPart(xv[i]));CHKERRQ(ierr); 714dc485aaSHong Zhang } else { 724dc485aaSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g\n",(double)PetscRealPart(xv[i]));CHKERRQ(ierr); 734dc485aaSHong Zhang } 744dc485aaSHong Zhang #else 754dc485aaSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g\n",(double)xv[i]);CHKERRQ(ierr); 764dc485aaSHong Zhang #endif 774dc485aaSHong Zhang } 784dc485aaSHong Zhang } 794dc485aaSHong Zhang ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 804dc485aaSHong Zhang ierr = VecRestoreArrayRead(X,&xv);CHKERRQ(ierr); 814dc485aaSHong Zhang PetscFunctionReturn(0); 824dc485aaSHong Zhang } 834dc485aaSHong Zhang 84*4062a5e5SHong Zhang static PetscErrorCode VecView_Network_MPI(DM networkdm,Vec X,PetscViewer viewer) 85*4062a5e5SHong Zhang { 86*4062a5e5SHong Zhang PetscErrorCode ierr; 87*4062a5e5SHong Zhang PetscInt i,e,v,Start,End,vStart,vEnd,offset,nvar,len_loc,len,k; 88*4062a5e5SHong Zhang const PetscScalar *xv; 89*4062a5e5SHong Zhang MPI_Comm comm; 90*4062a5e5SHong Zhang PetscMPIInt size,rank,tag = ((PetscObject)viewer)->tag; 91*4062a5e5SHong Zhang Vec localX; 92*4062a5e5SHong Zhang PetscBool ghostvtex; 93*4062a5e5SHong Zhang PetscScalar *values; 94*4062a5e5SHong Zhang PetscInt j,n,ne,nv,ii; 95*4062a5e5SHong Zhang MPI_Status status; 96*4062a5e5SHong Zhang 97*4062a5e5SHong Zhang PetscFunctionBegin; 98*4062a5e5SHong Zhang ierr = PetscObjectGetComm((PetscObject)networkdm,&comm);CHKERRQ(ierr); 99*4062a5e5SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 100*4062a5e5SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 101*4062a5e5SHong Zhang 102*4062a5e5SHong Zhang ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr); 103*4062a5e5SHong Zhang ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 104*4062a5e5SHong Zhang ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 105*4062a5e5SHong Zhang ierr = VecGetArrayRead(localX,&xv);CHKERRQ(ierr); 106*4062a5e5SHong Zhang 107*4062a5e5SHong Zhang ierr = VecGetLocalSize(localX,&len_loc);CHKERRQ(ierr); 108*4062a5e5SHong Zhang 109*4062a5e5SHong Zhang ierr = DMNetworkGetEdgeRange(networkdm,&Start,&End);CHKERRQ(ierr); 110*4062a5e5SHong Zhang ierr = DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);CHKERRQ(ierr); 111*4062a5e5SHong Zhang len_loc += 2*(1 + End-Start + vEnd-vStart); 112*4062a5e5SHong Zhang 113*4062a5e5SHong Zhang /* values = [nedges, offset, nvar, xe; nvertices, offsets, nvars, xv] */ 114*4062a5e5SHong Zhang ierr = MPI_Allreduce(&len_loc,&len,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr); 115*4062a5e5SHong Zhang ierr = PetscCalloc1(len,&values);CHKERRQ(ierr); 116*4062a5e5SHong Zhang 117*4062a5e5SHong Zhang if (!rank) { 118*4062a5e5SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Process [%d]\n",rank);CHKERRQ(ierr); 119*4062a5e5SHong Zhang } 120*4062a5e5SHong Zhang 121*4062a5e5SHong Zhang /* iterate over edges */ 122*4062a5e5SHong Zhang k = 2; 123*4062a5e5SHong Zhang for (e=Start; e<End; e++) { 124*4062a5e5SHong Zhang ierr = DMNetworkGetVariableOffset(networkdm,e,&offset);CHKERRQ(ierr); 125*4062a5e5SHong Zhang ierr = DMNetworkGetNumVariables(networkdm,e,&nvar);CHKERRQ(ierr); 126*4062a5e5SHong Zhang if (!nvar) continue; 127*4062a5e5SHong Zhang 128*4062a5e5SHong Zhang values[0] += 1; /* nedges */ 129*4062a5e5SHong Zhang values[k++] = offset; 130*4062a5e5SHong Zhang values[k++] = nvar; 131*4062a5e5SHong Zhang 132*4062a5e5SHong Zhang if (!rank) { /* print its own entries */ 133*4062a5e5SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Edge %D:\n",e-Start);CHKERRQ(ierr); 134*4062a5e5SHong Zhang } 135*4062a5e5SHong Zhang 136*4062a5e5SHong Zhang for (i=offset; i< offset+nvar; i++) { 137*4062a5e5SHong Zhang if (!rank) { 138*4062a5e5SHong Zhang #if defined(PETSC_USE_COMPLEX) 139*4062a5e5SHong Zhang if (PetscImaginaryPart(xv[i]) > 0.0) { 140*4062a5e5SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g + %g i\n",(double)PetscRealPart(xv[i]),(double)PetscImaginaryPart(xv[i]));CHKERRQ(ierr); 141*4062a5e5SHong Zhang } else if (PetscImaginaryPart(xv[i]) < 0.0) { 142*4062a5e5SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g - %g i\n",(double)PetscRealPart(xv[i]),-(double)PetscImaginaryPart(xv[i]));CHKERRQ(ierr); 143*4062a5e5SHong Zhang } else { 144*4062a5e5SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g\n",(double)PetscRealPart(xv[i]));CHKERRQ(ierr); 145*4062a5e5SHong Zhang } 146*4062a5e5SHong Zhang #else 147*4062a5e5SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g\n",(double)xv[i]);CHKERRQ(ierr); 148*4062a5e5SHong Zhang #endif 149*4062a5e5SHong Zhang } else { 150*4062a5e5SHong Zhang values[k++] = xv[i]; 151*4062a5e5SHong Zhang } 152*4062a5e5SHong Zhang } 153*4062a5e5SHong Zhang } 154*4062a5e5SHong Zhang 155*4062a5e5SHong Zhang /* iterate over vertices */ 156*4062a5e5SHong Zhang for (v=vStart; v<vEnd; v++) { 157*4062a5e5SHong Zhang ierr = DMNetworkIsGhostVertex(networkdm,v,&ghostvtex);CHKERRQ(ierr); 158*4062a5e5SHong Zhang if (ghostvtex) continue; 159*4062a5e5SHong Zhang ierr = DMNetworkGetVariableOffset(networkdm,v,&offset);CHKERRQ(ierr); 160*4062a5e5SHong Zhang ierr = DMNetworkGetNumVariables(networkdm,v,&nvar);CHKERRQ(ierr); 161*4062a5e5SHong Zhang if (!nvar) continue; 162*4062a5e5SHong Zhang 163*4062a5e5SHong Zhang values[1] += 1; /* nvertices */ 164*4062a5e5SHong Zhang values[k++] = offset; 165*4062a5e5SHong Zhang values[k++] = nvar; 166*4062a5e5SHong Zhang 167*4062a5e5SHong Zhang if (!rank) { 168*4062a5e5SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Vertex %D:\n",v-Start);CHKERRQ(ierr); 169*4062a5e5SHong Zhang } 170*4062a5e5SHong Zhang 171*4062a5e5SHong Zhang for (i=offset; i< offset+nvar; i++) { 172*4062a5e5SHong Zhang if (!rank) { 173*4062a5e5SHong Zhang #if defined(PETSC_USE_COMPLEX) 174*4062a5e5SHong Zhang if (PetscImaginaryPart(xv[i]) > 0.0) { 175*4062a5e5SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g + %g i\n",(double)PetscRealPart(xv[i]),(double)PetscImaginaryPart(xv[i]));CHKERRQ(ierr); 176*4062a5e5SHong Zhang } else if (PetscImaginaryPart(xv[i]) < 0.0) { 177*4062a5e5SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g - %g i\n",(double)PetscRealPart(xv[i]),-(double)PetscImaginaryPart(xv[i]));CHKERRQ(ierr); 178*4062a5e5SHong Zhang } else { 179*4062a5e5SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g\n",(double)PetscRealPart(xv[i]));CHKERRQ(ierr); 180*4062a5e5SHong Zhang } 181*4062a5e5SHong Zhang #else 182*4062a5e5SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g\n",(double)xv[i]);CHKERRQ(ierr); 183*4062a5e5SHong Zhang #endif 184*4062a5e5SHong Zhang } else { 185*4062a5e5SHong Zhang values[k++] = xv[i]; 186*4062a5e5SHong Zhang } 187*4062a5e5SHong Zhang } 188*4062a5e5SHong Zhang } 189*4062a5e5SHong Zhang 190*4062a5e5SHong Zhang if (!rank) { 191*4062a5e5SHong Zhang /* proc[0] receives and prints messages */ 192*4062a5e5SHong Zhang for (j=1; j<size; j++) { 193*4062a5e5SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Process [%d]\n",j);CHKERRQ(ierr); 194*4062a5e5SHong Zhang 195*4062a5e5SHong Zhang ierr = MPI_Recv(values,(PetscMPIInt)len,MPIU_SCALAR,j,tag,comm,&status);CHKERRQ(ierr); 196*4062a5e5SHong Zhang ierr = MPI_Get_count(&status,MPIU_SCALAR,&n);CHKERRQ(ierr); 197*4062a5e5SHong Zhang 198*4062a5e5SHong Zhang ne = (PetscInt)values[0]; 199*4062a5e5SHong Zhang nv = (PetscInt)values[1]; 200*4062a5e5SHong Zhang 201*4062a5e5SHong Zhang /* print received edges */ 202*4062a5e5SHong Zhang k = 2; 203*4062a5e5SHong Zhang for (i=0; i<ne; i++) { 204*4062a5e5SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Edge %D:\n",i);CHKERRQ(ierr); 205*4062a5e5SHong Zhang } 206*4062a5e5SHong Zhang 207*4062a5e5SHong Zhang /* print received vertices */ 208*4062a5e5SHong Zhang for (i=0; i<nv; i++) { 209*4062a5e5SHong Zhang offset = (PetscInt)values[k++]; 210*4062a5e5SHong Zhang nvar = (PetscInt)values[k++]; 211*4062a5e5SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Vertex %D:\n",i);CHKERRQ(ierr); 212*4062a5e5SHong Zhang for (ii=0; ii<nvar; ii++) { 213*4062a5e5SHong Zhang #if defined(PETSC_USE_COMPLEX) 214*4062a5e5SHong Zhang if (PetscImaginaryPart(xv[i]) > 0.0) { 215*4062a5e5SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g + %g i\n",(double)PetscRealPart(xv[k]),(double)PetscImaginaryPart(xv[k]));CHKERRQ(ierr); 216*4062a5e5SHong Zhang k++; 217*4062a5e5SHong Zhang } else if (PetscImaginaryPart(xv[i]) < 0.0) { 218*4062a5e5SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g - %g i\n",(double)PetscRealPart(xv[k]),-(double)PetscImaginaryPart(xv[k]));CHKERRQ(ierr); 219*4062a5e5SHong Zhang k++; 220*4062a5e5SHong Zhang } else { 221*4062a5e5SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g\n",(double)PetscRealPart(xv[k++]));CHKERRQ(ierr); 222*4062a5e5SHong Zhang } 223*4062a5e5SHong Zhang #else 224*4062a5e5SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %g\n",(double)values[k++]);CHKERRQ(ierr); 225*4062a5e5SHong Zhang #endif 226*4062a5e5SHong Zhang } 227*4062a5e5SHong Zhang } 228*4062a5e5SHong Zhang } 229*4062a5e5SHong Zhang } else { 230*4062a5e5SHong Zhang /* sends values to proc[0] */ 231*4062a5e5SHong Zhang ierr = MPI_Send((void*)values,k,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr); 232*4062a5e5SHong Zhang } 233*4062a5e5SHong Zhang 234*4062a5e5SHong Zhang ierr = PetscFree(values);CHKERRQ(ierr); 235*4062a5e5SHong Zhang ierr = VecRestoreArrayRead(localX,&xv);CHKERRQ(ierr); 236*4062a5e5SHong Zhang ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr); 237*4062a5e5SHong Zhang PetscFunctionReturn(0); 238*4062a5e5SHong Zhang } 239*4062a5e5SHong Zhang 2404dc485aaSHong Zhang PetscErrorCode VecView_Network(Vec v,PetscViewer viewer) 2414dc485aaSHong Zhang { 2424dc485aaSHong Zhang DM dm; 2434dc485aaSHong Zhang PetscErrorCode ierr; 2444dc485aaSHong Zhang PetscBool isseq; 2454dc485aaSHong Zhang 2464dc485aaSHong Zhang PetscFunctionBegin; 2474dc485aaSHong Zhang ierr = VecGetDM(v,&dm);CHKERRQ(ierr); 2484dc485aaSHong Zhang if (!dm) SETERRQ(PetscObjectComm((PetscObject)v),PETSC_ERR_ARG_WRONG,"Vector not generated from a DM"); 2494dc485aaSHong Zhang ierr = PetscObjectTypeCompare((PetscObject)v,VECSEQ,&isseq);CHKERRQ(ierr); 2504dc485aaSHong Zhang if (isseq) { 2514dc485aaSHong Zhang ierr = VecView_Network_Seq(dm,v,viewer);CHKERRQ(ierr); 2524dc485aaSHong Zhang } else { 253*4062a5e5SHong Zhang ierr = VecView_Network_MPI(dm,v,viewer);CHKERRQ(ierr); 2544dc485aaSHong Zhang } 2554dc485aaSHong Zhang PetscFunctionReturn(0); 2564dc485aaSHong Zhang } 2575f2c45f1SShri Abhyankar 2585f2c45f1SShri Abhyankar static PetscErrorCode DMCreateGlobalVector_Network(DM dm,Vec *vec) 2595f2c45f1SShri Abhyankar { 2605f2c45f1SShri Abhyankar PetscErrorCode ierr; 2615f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 2625f2c45f1SShri Abhyankar 2635f2c45f1SShri Abhyankar PetscFunctionBegin; 2645f2c45f1SShri Abhyankar ierr = DMCreateGlobalVector(network->plex,vec);CHKERRQ(ierr); 2654dc485aaSHong Zhang ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Network);CHKERRQ(ierr); 2665f2c45f1SShri Abhyankar ierr = VecSetDM(*vec,dm);CHKERRQ(ierr); 2675f2c45f1SShri Abhyankar PetscFunctionReturn(0); 2685f2c45f1SShri Abhyankar } 2695f2c45f1SShri Abhyankar 2705f2c45f1SShri Abhyankar static PetscErrorCode DMCreateLocalVector_Network(DM dm,Vec *vec) 2715f2c45f1SShri Abhyankar { 2725f2c45f1SShri Abhyankar PetscErrorCode ierr; 2735f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 2745f2c45f1SShri Abhyankar 2755f2c45f1SShri Abhyankar PetscFunctionBegin; 2765f2c45f1SShri Abhyankar ierr = DMCreateLocalVector(network->plex,vec);CHKERRQ(ierr); 2775f2c45f1SShri Abhyankar ierr = VecSetDM(*vec,dm);CHKERRQ(ierr); 2785f2c45f1SShri Abhyankar PetscFunctionReturn(0); 2795f2c45f1SShri Abhyankar } 2805f2c45f1SShri Abhyankar 2815f2c45f1SShri Abhyankar PetscErrorCode DMInitialize_Network(DM dm) 2825f2c45f1SShri Abhyankar { 2835f2c45f1SShri Abhyankar 2845f2c45f1SShri Abhyankar PetscFunctionBegin; 2855f2c45f1SShri Abhyankar 2865f2c45f1SShri Abhyankar dm->ops->view = NULL; 2875f2c45f1SShri Abhyankar dm->ops->setfromoptions = DMSetFromOptions_Network; 2888415c774SShri Abhyankar dm->ops->clone = DMClone_Network; 2895f2c45f1SShri Abhyankar dm->ops->setup = DMSetUp_Network; 2905f2c45f1SShri Abhyankar dm->ops->createglobalvector = DMCreateGlobalVector_Network; 2915f2c45f1SShri Abhyankar dm->ops->createlocalvector = DMCreateLocalVector_Network; 2925f2c45f1SShri Abhyankar dm->ops->getlocaltoglobalmapping = NULL; 2935f2c45f1SShri Abhyankar dm->ops->createfieldis = NULL; 2945f2c45f1SShri Abhyankar dm->ops->createcoordinatedm = NULL; 2955f2c45f1SShri Abhyankar dm->ops->getcoloring = 0; 2965f2c45f1SShri Abhyankar dm->ops->creatematrix = DMCreateMatrix_Network; 2975f2c45f1SShri Abhyankar dm->ops->createinterpolation = 0; 2985f2c45f1SShri Abhyankar dm->ops->getaggregates = 0; 2995f2c45f1SShri Abhyankar dm->ops->getinjection = 0; 3005f2c45f1SShri Abhyankar dm->ops->refine = 0; 3015f2c45f1SShri Abhyankar dm->ops->coarsen = 0; 3025f2c45f1SShri Abhyankar dm->ops->refinehierarchy = 0; 3035f2c45f1SShri Abhyankar dm->ops->coarsenhierarchy = 0; 3045f2c45f1SShri Abhyankar dm->ops->globaltolocalbegin = DMGlobalToLocalBegin_Network; 3055f2c45f1SShri Abhyankar dm->ops->globaltolocalend = DMGlobalToLocalEnd_Network; 3065f2c45f1SShri Abhyankar dm->ops->localtoglobalbegin = DMLocalToGlobalBegin_Network; 3075f2c45f1SShri Abhyankar dm->ops->localtoglobalend = DMLocalToGlobalEnd_Network; 3085f2c45f1SShri Abhyankar dm->ops->destroy = DMDestroy_Network; 3095f2c45f1SShri Abhyankar dm->ops->createsubdm = NULL; 3105f2c45f1SShri Abhyankar dm->ops->locatepoints = NULL; 3115f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3125f2c45f1SShri Abhyankar } 3135f2c45f1SShri Abhyankar 3148415c774SShri Abhyankar PetscErrorCode DMClone_Network(DM dm, DM *newdm) 3158415c774SShri Abhyankar { 3168415c774SShri Abhyankar DM_Network *network = (DM_Network *) dm->data; 3178415c774SShri Abhyankar PetscErrorCode ierr; 3188415c774SShri Abhyankar 3198415c774SShri Abhyankar PetscFunctionBegin; 3208415c774SShri Abhyankar network->refct++; 3218415c774SShri Abhyankar (*newdm)->data = network; 3228415c774SShri Abhyankar ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMNETWORK);CHKERRQ(ierr); 3238415c774SShri Abhyankar ierr = DMInitialize_Network(*newdm);CHKERRQ(ierr); 3248415c774SShri Abhyankar PetscFunctionReturn(0); 3258415c774SShri Abhyankar } 3268415c774SShri Abhyankar 3275f2c45f1SShri Abhyankar /*MC 3285f2c45f1SShri Abhyankar DMNETWORK = "network" - A DM object that encapsulates an unstructured network. The implementation is based on the DM object 3295f2c45f1SShri Abhyankar DMPlex that manages unstructured grids. Distributed networks use a non-overlapping partitioning of 3305f2c45f1SShri Abhyankar the edges. In the local representation, Vecs contain all unknowns in the interior and shared boundary. 3315f2c45f1SShri Abhyankar This is specified by a PetscSection object. Ownership in the global representation is determined by 3325f2c45f1SShri Abhyankar ownership of the underlying DMPlex points. This is specified by another PetscSection object. 3335f2c45f1SShri Abhyankar 3345f2c45f1SShri Abhyankar Level: intermediate 3355f2c45f1SShri Abhyankar 3365f2c45f1SShri Abhyankar .seealso: DMType, DMNetworkCreate(), DMCreate(), DMSetType() 3375f2c45f1SShri Abhyankar M*/ 3385f2c45f1SShri Abhyankar 3395f2c45f1SShri Abhyankar PETSC_EXTERN PetscErrorCode DMCreate_Network(DM dm) 3405f2c45f1SShri Abhyankar { 3415f2c45f1SShri Abhyankar DM_Network *network; 3425f2c45f1SShri Abhyankar PetscErrorCode ierr; 3435f2c45f1SShri Abhyankar 3445f2c45f1SShri Abhyankar PetscFunctionBegin; 3455f2c45f1SShri Abhyankar PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3465f2c45f1SShri Abhyankar ierr = PetscNewLog(dm,&network);CHKERRQ(ierr); 3475f2c45f1SShri Abhyankar dm->data = network; 3485f2c45f1SShri Abhyankar 3495f2c45f1SShri Abhyankar network->refct = 1; 3505f2c45f1SShri Abhyankar network->NNodes = -1; 3515f2c45f1SShri Abhyankar network->NEdges = -1; 3525f2c45f1SShri Abhyankar network->nNodes = -1; 3535f2c45f1SShri Abhyankar network->nEdges = -1; 3545f2c45f1SShri Abhyankar 35513c2a604SAdrian Maldonado 3565f2c45f1SShri Abhyankar ierr = DMInitialize_Network(dm);CHKERRQ(ierr); 3575f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3585f2c45f1SShri Abhyankar } 3595f2c45f1SShri Abhyankar 3605f2c45f1SShri Abhyankar /*@ 3615f2c45f1SShri Abhyankar DMNetworkCreate - Creates a DMNetwork object, which encapsulates an unstructured network. 3625f2c45f1SShri Abhyankar 3635f2c45f1SShri Abhyankar Collective on MPI_Comm 3645f2c45f1SShri Abhyankar 3655f2c45f1SShri Abhyankar Input Parameter: 3665f2c45f1SShri Abhyankar . comm - The communicator for the DMNetwork object 3675f2c45f1SShri Abhyankar 3685f2c45f1SShri Abhyankar Output Parameter: 3695f2c45f1SShri Abhyankar . network - The DMNetwork object 3705f2c45f1SShri Abhyankar 3715f2c45f1SShri Abhyankar Level: beginner 3725f2c45f1SShri Abhyankar 3735f2c45f1SShri Abhyankar .keywords: DMNetwork, create 3745f2c45f1SShri Abhyankar @*/ 3755f2c45f1SShri Abhyankar PetscErrorCode DMNetworkCreate(MPI_Comm comm, DM *network) 3765f2c45f1SShri Abhyankar { 3775f2c45f1SShri Abhyankar PetscErrorCode ierr; 3785f2c45f1SShri Abhyankar 3795f2c45f1SShri Abhyankar PetscFunctionBegin; 3805f2c45f1SShri Abhyankar PetscValidPointer(network,2); 3815f2c45f1SShri Abhyankar ierr = DMCreate(comm, network);CHKERRQ(ierr); 3825f2c45f1SShri Abhyankar ierr = DMSetType(*network, DMNETWORK);CHKERRQ(ierr); 3835f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3845f2c45f1SShri Abhyankar } 385