15f2c45f1SShri Abhyankar #define PETSCDM_DLL 2af0996ceSBarry Smith #include <petsc/private/dmnetworkimpl.h> /*I "petscdmnetwork.h" I*/ 395af8c53SHong Zhang #include <petsc/private/vecimpl.h> 45f2c45f1SShri Abhyankar 54416b707SBarry Smith PetscErrorCode DMSetFromOptions_Network(PetscOptionItems *PetscOptionsObject,DM dm) 65f2c45f1SShri Abhyankar { 75f2c45f1SShri Abhyankar PetscFunctionBegin; 8064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 9d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"DMNetwork Options"); 10d0609cedSBarry Smith PetscOptionsHeadEnd(); 115f2c45f1SShri Abhyankar PetscFunctionReturn(0); 125f2c45f1SShri Abhyankar } 135f2c45f1SShri Abhyankar 145f2c45f1SShri Abhyankar /* External function declarations here */ 155f2c45f1SShri Abhyankar extern PetscErrorCode DMCreateMatrix_Network(DM, Mat*); 165f2c45f1SShri Abhyankar extern PetscErrorCode DMDestroy_Network(DM); 175f2c45f1SShri Abhyankar extern PetscErrorCode DMView_Network(DM, PetscViewer); 185f2c45f1SShri Abhyankar extern PetscErrorCode DMGlobalToLocalBegin_Network(DM, Vec, InsertMode, Vec); 195f2c45f1SShri Abhyankar extern PetscErrorCode DMGlobalToLocalEnd_Network(DM, Vec, InsertMode, Vec); 205f2c45f1SShri Abhyankar extern PetscErrorCode DMLocalToGlobalBegin_Network(DM, Vec, InsertMode, Vec); 215f2c45f1SShri Abhyankar extern PetscErrorCode DMLocalToGlobalEnd_Network(DM, Vec, InsertMode, Vec); 225f2c45f1SShri Abhyankar extern PetscErrorCode DMSetUp_Network(DM); 238415c774SShri Abhyankar extern PetscErrorCode DMClone_Network(DM, DM*); 245f2c45f1SShri Abhyankar 25a81b7fe4SHong Zhang static PetscErrorCode VecArrayPrint_private(PetscViewer viewer,PetscInt n,const PetscScalar *xv) 26a81b7fe4SHong Zhang { 27a81b7fe4SHong Zhang PetscInt i; 28a81b7fe4SHong Zhang 29a81b7fe4SHong Zhang PetscFunctionBegin; 30a81b7fe4SHong Zhang for (i=0; i<n; i++) { 31a81b7fe4SHong Zhang #if defined(PETSC_USE_COMPLEX) 32a81b7fe4SHong Zhang if (PetscImaginaryPart(xv[i]) > 0.0) { 339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," %g + %g i\n",(double)PetscRealPart(xv[i]),(double)PetscImaginaryPart(xv[i]))); 34a81b7fe4SHong Zhang } else if (PetscImaginaryPart(xv[i]) < 0.0) { 359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," %g - %g i\n",(double)PetscRealPart(xv[i]),-(double)PetscImaginaryPart(xv[i]))); 361baa6e33SBarry Smith } else PetscCall(PetscViewerASCIIPrintf(viewer," %g\n",(double)PetscRealPart(xv[i]))); 37a81b7fe4SHong Zhang #else 389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," %g\n",(double)xv[i])); 39a81b7fe4SHong Zhang #endif 40a81b7fe4SHong Zhang } 41a81b7fe4SHong Zhang PetscFunctionReturn(0); 42a81b7fe4SHong Zhang } 43a81b7fe4SHong Zhang 444dc485aaSHong Zhang static PetscErrorCode VecView_Network_Seq(DM networkdm,Vec X,PetscViewer viewer) 454dc485aaSHong Zhang { 467b6afd5bSHong Zhang PetscInt e,v,Start,End,offset,nvar,id; 474dc485aaSHong Zhang const PetscScalar *xv; 484dc485aaSHong Zhang 494dc485aaSHong Zhang PetscFunctionBegin; 509566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&xv)); 514dc485aaSHong Zhang 524dc485aaSHong Zhang /* iterate over edges */ 539566063dSJacob Faibussowitsch PetscCall(DMNetworkGetEdgeRange(networkdm,&Start,&End)); 544dc485aaSHong Zhang for (e=Start; e<End; e++) { 559566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm,e,ALL_COMPONENTS,NULL,NULL,&nvar)); 564dc485aaSHong Zhang if (!nvar) continue; 57a81b7fe4SHong Zhang 589566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm,e,ALL_COMPONENTS,&offset)); 599566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalEdgeIndex(networkdm,e,&id)); 607b6afd5bSHong Zhang 619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Edge %" PetscInt_FMT ":\n",id)); 629566063dSJacob Faibussowitsch PetscCall(VecArrayPrint_private(viewer,nvar,xv+offset)); 634dc485aaSHong Zhang } 644dc485aaSHong Zhang 654dc485aaSHong Zhang /* iterate over vertices */ 669566063dSJacob Faibussowitsch PetscCall(DMNetworkGetVertexRange(networkdm,&Start,&End)); 674dc485aaSHong Zhang for (v=Start; v<End; v++) { 689566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm,v,ALL_COMPONENTS,NULL,NULL,&nvar)); 694dc485aaSHong Zhang if (!nvar) continue; 70a81b7fe4SHong Zhang 719566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm,v,ALL_COMPONENTS,&offset)); 729566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVertexIndex(networkdm,v,&id)); 7362aa33baSHong Zhang 749566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Vertex %" PetscInt_FMT ":\n",id)); 759566063dSJacob Faibussowitsch PetscCall(VecArrayPrint_private(viewer,nvar,xv+offset)); 764dc485aaSHong Zhang } 779566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 789566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&xv)); 794dc485aaSHong Zhang PetscFunctionReturn(0); 804dc485aaSHong Zhang } 814dc485aaSHong Zhang 824062a5e5SHong Zhang static PetscErrorCode VecView_Network_MPI(DM networkdm,Vec X,PetscViewer viewer) 834062a5e5SHong Zhang { 84e8c2e809SHong Zhang PetscInt i,e,v,eStart,eEnd,vStart,vEnd,offset,nvar,len_loc,len,k; 854062a5e5SHong Zhang const PetscScalar *xv; 864062a5e5SHong Zhang MPI_Comm comm; 874062a5e5SHong Zhang PetscMPIInt size,rank,tag = ((PetscObject)viewer)->tag; 884062a5e5SHong Zhang Vec localX; 894062a5e5SHong Zhang PetscBool ghostvtex; 904062a5e5SHong Zhang PetscScalar *values; 917b6afd5bSHong Zhang PetscInt j,ne,nv,id; 924062a5e5SHong Zhang MPI_Status status; 934062a5e5SHong Zhang 944062a5e5SHong Zhang PetscFunctionBegin; 959566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)networkdm,&comm)); 969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 984062a5e5SHong Zhang 999566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(networkdm,&localX)); 1009566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX)); 1019566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX)); 1029566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(localX,&xv)); 1034062a5e5SHong Zhang 1049566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(localX,&len_loc)); 1054062a5e5SHong Zhang 1069566063dSJacob Faibussowitsch PetscCall(DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd)); 1079566063dSJacob Faibussowitsch PetscCall(DMNetworkGetVertexRange(networkdm,&vStart,&vEnd)); 108e8c2e809SHong Zhang len_loc += 2*(1 + eEnd-eStart + vEnd-vStart); 1094062a5e5SHong Zhang 110e85e6aecSHong Zhang /* values = [nedges, nvertices; id, nvar, xedge; ...; id, nvars, xvertex;...], to be sent to proc[0] */ 1119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&len_loc,&len,1,MPIU_INT,MPI_MAX,comm)); 1129566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(len,&values)); 1134062a5e5SHong Zhang 114dd400576SPatrick Sanan if (rank == 0) { 1159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Process [%d]\n",rank)); 1164062a5e5SHong Zhang } 1174062a5e5SHong Zhang 1184062a5e5SHong Zhang /* iterate over edges */ 1194062a5e5SHong Zhang k = 2; 120e8c2e809SHong Zhang for (e=eStart; e<eEnd; e++) { 1219566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm,e,ALL_COMPONENTS,NULL,NULL,&nvar)); 1224062a5e5SHong Zhang if (!nvar) continue; 1234062a5e5SHong Zhang 1249566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm,e,ALL_COMPONENTS,&offset)); 1259566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalEdgeIndex(networkdm,e,&id)); 12662aa33baSHong Zhang 127dd400576SPatrick Sanan if (rank == 0) { /* print its own entries */ 1289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Edge %" PetscInt_FMT ":\n",id)); 1299566063dSJacob Faibussowitsch PetscCall(VecArrayPrint_private(viewer,nvar,xv+offset)); 130a81b7fe4SHong Zhang } else { 1317b6afd5bSHong Zhang values[0] += 1; /* number of edges */ 1327b6afd5bSHong Zhang values[k++] = id; 1337b6afd5bSHong Zhang values[k++] = nvar; 1347b6afd5bSHong Zhang for (i=offset; i< offset+nvar; i++) values[k++] = xv[i]; 1354062a5e5SHong Zhang } 1364062a5e5SHong Zhang } 1374062a5e5SHong Zhang 1384062a5e5SHong Zhang /* iterate over vertices */ 1394062a5e5SHong Zhang for (v=vStart; v<vEnd; v++) { 1409566063dSJacob Faibussowitsch PetscCall(DMNetworkIsGhostVertex(networkdm,v,&ghostvtex)); 1414062a5e5SHong Zhang if (ghostvtex) continue; 1429566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm,v,ALL_COMPONENTS,NULL,NULL,&nvar)); 1434062a5e5SHong Zhang if (!nvar) continue; 1444062a5e5SHong Zhang 1459566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm,v,ALL_COMPONENTS,&offset)); 1469566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVertexIndex(networkdm,v,&id)); 1474062a5e5SHong Zhang 148dd400576SPatrick Sanan if (rank == 0) { 1499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Vertex %" PetscInt_FMT ":\n",id)); 1509566063dSJacob Faibussowitsch PetscCall(VecArrayPrint_private(viewer,nvar,xv+offset)); 1514062a5e5SHong Zhang } else { 1527b6afd5bSHong Zhang values[1] += 1; /* number of vertices */ 1537b6afd5bSHong Zhang values[k++] = id; 1547b6afd5bSHong Zhang values[k++] = nvar; 1557b6afd5bSHong Zhang for (i=offset; i< offset+nvar; i++) values[k++] = xv[i]; 1564062a5e5SHong Zhang } 1574062a5e5SHong Zhang } 1584062a5e5SHong Zhang 159dd400576SPatrick Sanan if (rank == 0) { 1604062a5e5SHong Zhang /* proc[0] receives and prints messages */ 1614062a5e5SHong Zhang for (j=1; j<size; j++) { 16263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Process [%" PetscInt_FMT "]\n",j)); 1634062a5e5SHong Zhang 1649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(values,(PetscMPIInt)len,MPIU_SCALAR,j,tag,comm,&status)); 1654062a5e5SHong Zhang 166dde233f4SSatish Balay ne = (PetscInt)PetscAbsScalar(values[0]); 167dde233f4SSatish Balay nv = (PetscInt)PetscAbsScalar(values[1]); 1684062a5e5SHong Zhang 1694062a5e5SHong Zhang /* print received edges */ 1704062a5e5SHong Zhang k = 2; 1714062a5e5SHong Zhang for (i=0; i<ne; i++) { 172dde233f4SSatish Balay id = (PetscInt)PetscAbsScalar(values[k++]); 173dde233f4SSatish Balay nvar = (PetscInt)PetscAbsScalar(values[k++]); 1749566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Edge %" PetscInt_FMT ":\n",id)); 1759566063dSJacob Faibussowitsch PetscCall(VecArrayPrint_private(viewer,nvar,values+k)); 176a81b7fe4SHong Zhang k += nvar; 1774062a5e5SHong Zhang } 1784062a5e5SHong Zhang 1794062a5e5SHong Zhang /* print received vertices */ 1804062a5e5SHong Zhang for (i=0; i<nv; i++) { 181dde233f4SSatish Balay id = (PetscInt)PetscAbsScalar(values[k++]); 182dde233f4SSatish Balay nvar = (PetscInt)PetscAbsScalar(values[k++]); 1839566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Vertex %" PetscInt_FMT ":\n",id)); 1849566063dSJacob Faibussowitsch PetscCall(VecArrayPrint_private(viewer,nvar,values+k)); 185a81b7fe4SHong Zhang k += nvar; 1864062a5e5SHong Zhang } 1874062a5e5SHong Zhang } 1884062a5e5SHong Zhang } else { 1894062a5e5SHong Zhang /* sends values to proc[0] */ 1909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void*)values,k,MPIU_SCALAR,0,tag,comm)); 1914062a5e5SHong Zhang } 1924062a5e5SHong Zhang 1939566063dSJacob Faibussowitsch PetscCall(PetscFree(values)); 1949566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(localX,&xv)); 1959566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(networkdm,&localX)); 1964062a5e5SHong Zhang PetscFunctionReturn(0); 1974062a5e5SHong Zhang } 1984062a5e5SHong Zhang 19995af8c53SHong Zhang PETSC_EXTERN PetscErrorCode VecView_MPI(Vec,PetscViewer); 200a6c4b7b7SHong Zhang 2014dc485aaSHong Zhang PetscErrorCode VecView_Network(Vec v,PetscViewer viewer) 2024dc485aaSHong Zhang { 2034dc485aaSHong Zhang DM dm; 2044dc485aaSHong Zhang PetscBool isseq; 205a6c4b7b7SHong Zhang PetscBool iascii; 2064dc485aaSHong Zhang 2074dc485aaSHong Zhang PetscFunctionBegin; 2089566063dSJacob Faibussowitsch PetscCall(VecGetDM(v,&dm)); 2095c6496baSHong Zhang PetscCheck(dm,PetscObjectComm((PetscObject)v),PETSC_ERR_ARG_WRONG,"Vector not generated from a DM"); 2109566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 2119566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)v,VECSEQ,&isseq)); 212a6c4b7b7SHong Zhang 213a6c4b7b7SHong Zhang /* Use VecView_Network if the viewer is ASCII; use VecView_Seq/MPI for other viewer formats */ 214a6c4b7b7SHong Zhang if (iascii) { 2151baa6e33SBarry Smith if (isseq) PetscCall(VecView_Network_Seq(dm,v,viewer)); 2161baa6e33SBarry Smith else PetscCall(VecView_Network_MPI(dm,v,viewer)); 2174dc485aaSHong Zhang } else { 2181baa6e33SBarry Smith if (isseq) PetscCall(VecView_Seq(v,viewer)); 2191baa6e33SBarry Smith else PetscCall(VecView_MPI(v,viewer)); 220a6c4b7b7SHong Zhang } 2214dc485aaSHong Zhang PetscFunctionReturn(0); 2224dc485aaSHong Zhang } 2235f2c45f1SShri Abhyankar 2245f2c45f1SShri Abhyankar static PetscErrorCode DMCreateGlobalVector_Network(DM dm,Vec *vec) 2255f2c45f1SShri Abhyankar { 2265f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 2275f2c45f1SShri Abhyankar 2285f2c45f1SShri Abhyankar PetscFunctionBegin; 2299566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(network->plex,vec)); 2309566063dSJacob Faibussowitsch PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Network)); 2319566063dSJacob Faibussowitsch PetscCall(VecSetDM(*vec,dm)); 2325f2c45f1SShri Abhyankar PetscFunctionReturn(0); 2335f2c45f1SShri Abhyankar } 2345f2c45f1SShri Abhyankar 2355f2c45f1SShri Abhyankar static PetscErrorCode DMCreateLocalVector_Network(DM dm,Vec *vec) 2365f2c45f1SShri Abhyankar { 2375f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 2385f2c45f1SShri Abhyankar 2395f2c45f1SShri Abhyankar PetscFunctionBegin; 2409566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(network->plex,vec)); 2419566063dSJacob Faibussowitsch PetscCall(VecSetDM(*vec,dm)); 2425f2c45f1SShri Abhyankar PetscFunctionReturn(0); 2435f2c45f1SShri Abhyankar } 2445f2c45f1SShri Abhyankar 245*daad07d3SAidan Hamilton PetscErrorCode DMNetworkInitializeToDefault_NonShared(DM dm) 246*daad07d3SAidan Hamilton { 247*daad07d3SAidan Hamilton DM_Network *network = (DM_Network*) dm->data; 248*daad07d3SAidan Hamilton 249*daad07d3SAidan Hamilton PetscFunctionBegin; 250*daad07d3SAidan Hamilton network->Je = NULL; 251*daad07d3SAidan Hamilton network->Jv = NULL; 252*daad07d3SAidan Hamilton network->Jvptr = NULL; 253*daad07d3SAidan Hamilton network->userEdgeJacobian = PETSC_FALSE; 254*daad07d3SAidan Hamilton network->userVertexJacobian = PETSC_FALSE; 255*daad07d3SAidan Hamilton 256*daad07d3SAidan Hamilton network->vertex.DofSection = NULL; 257*daad07d3SAidan Hamilton network->vertex.GlobalDofSection = NULL; 258*daad07d3SAidan Hamilton network->vertex.mapping = NULL; 259*daad07d3SAidan Hamilton network->vertex.sf = NULL; 260*daad07d3SAidan Hamilton 261*daad07d3SAidan Hamilton network->edge.DofSection = NULL; 262*daad07d3SAidan Hamilton network->edge.GlobalDofSection = NULL; 263*daad07d3SAidan Hamilton network->edge.mapping = NULL; 264*daad07d3SAidan Hamilton network->edge.sf = NULL; 265*daad07d3SAidan Hamilton 266*daad07d3SAidan Hamilton network->DataSection = NULL; 267*daad07d3SAidan Hamilton network->DofSection = NULL; 268*daad07d3SAidan Hamilton network->GlobalDofSection = NULL; 269*daad07d3SAidan Hamilton network->componentsetup = PETSC_FALSE; 270*daad07d3SAidan Hamilton 271*daad07d3SAidan Hamilton network->plex = NULL; 272*daad07d3SAidan Hamilton 273*daad07d3SAidan Hamilton network->component = NULL; 274*daad07d3SAidan Hamilton network->ncomponent = 0; 275*daad07d3SAidan Hamilton 276*daad07d3SAidan Hamilton network->header = NULL; 277*daad07d3SAidan Hamilton network->cvalue = NULL; 278*daad07d3SAidan Hamilton network->componentdataarray = NULL; 279*daad07d3SAidan Hamilton 280*daad07d3SAidan Hamilton network->max_comps_registered = DMNETWORK_MAX_COMP_REGISTERED_DEFAULT; /* return to default */ 281*daad07d3SAidan Hamilton PetscFunctionReturn(0); 282*daad07d3SAidan Hamilton } 283*daad07d3SAidan Hamilton /* Default values for the parameters in DMNetwork */ 284*daad07d3SAidan Hamilton PetscErrorCode DMNetworkInitializeToDefault(DM dm) 285*daad07d3SAidan Hamilton { 286*daad07d3SAidan Hamilton DM_Network *network = (DM_Network*) dm->data; 287*daad07d3SAidan Hamilton DMNetworkCloneShared cloneshared = network->cloneshared; 288*daad07d3SAidan Hamilton 289*daad07d3SAidan Hamilton PetscFunctionBegin; 290*daad07d3SAidan Hamilton PetscCall(DMNetworkInitializeToDefault_NonShared(dm)); 291*daad07d3SAidan Hamilton /* Default values for shared data */ 292*daad07d3SAidan Hamilton cloneshared->refct = 1; 293*daad07d3SAidan Hamilton cloneshared->NVertices = 0; 294*daad07d3SAidan Hamilton cloneshared->NEdges = 0; 295*daad07d3SAidan Hamilton cloneshared->nVertices = 0; 296*daad07d3SAidan Hamilton cloneshared->nEdges = 0; 297*daad07d3SAidan Hamilton cloneshared->nsubnet = 0; 298*daad07d3SAidan Hamilton cloneshared->pStart = -1; 299*daad07d3SAidan Hamilton cloneshared->pEnd = -1; 300*daad07d3SAidan Hamilton cloneshared->vStart = -1; 301*daad07d3SAidan Hamilton cloneshared->vEnd = -1; 302*daad07d3SAidan Hamilton cloneshared->eStart = -1; 303*daad07d3SAidan Hamilton cloneshared->eEnd = -1; 304*daad07d3SAidan Hamilton cloneshared->vltog = NULL; 305*daad07d3SAidan Hamilton cloneshared->distributecalled = PETSC_FALSE; 306*daad07d3SAidan Hamilton 307*daad07d3SAidan Hamilton cloneshared->subnet = NULL; 308*daad07d3SAidan Hamilton cloneshared->subnetvtx = NULL; 309*daad07d3SAidan Hamilton cloneshared->subnetedge = NULL; 310*daad07d3SAidan Hamilton cloneshared->svtx = NULL; 311*daad07d3SAidan Hamilton cloneshared->nsvtx = 0; 312*daad07d3SAidan Hamilton cloneshared->Nsvtx = 0; 313*daad07d3SAidan Hamilton cloneshared->svertices = NULL; 314*daad07d3SAidan Hamilton cloneshared->sedgelist = NULL; 315*daad07d3SAidan Hamilton cloneshared->svtable = NULL; 316*daad07d3SAidan Hamilton PetscFunctionReturn(0); 317*daad07d3SAidan Hamilton } 318*daad07d3SAidan Hamilton 3195f2c45f1SShri Abhyankar PetscErrorCode DMInitialize_Network(DM dm) 3205f2c45f1SShri Abhyankar { 3215f2c45f1SShri Abhyankar PetscFunctionBegin; 3229566063dSJacob Faibussowitsch PetscCall(DMSetDimension(dm,1)); 323caf410d2SHong Zhang dm->ops->view = DMView_Network; 3245f2c45f1SShri Abhyankar dm->ops->setfromoptions = DMSetFromOptions_Network; 3258415c774SShri Abhyankar dm->ops->clone = DMClone_Network; 3265f2c45f1SShri Abhyankar dm->ops->setup = DMSetUp_Network; 3275f2c45f1SShri Abhyankar dm->ops->createglobalvector = DMCreateGlobalVector_Network; 3285f2c45f1SShri Abhyankar dm->ops->createlocalvector = DMCreateLocalVector_Network; 3295f2c45f1SShri Abhyankar dm->ops->getlocaltoglobalmapping = NULL; 3305f2c45f1SShri Abhyankar dm->ops->createfieldis = NULL; 3315f2c45f1SShri Abhyankar dm->ops->createcoordinatedm = NULL; 332ea78f98cSLisandro Dalcin dm->ops->getcoloring = NULL; 3335f2c45f1SShri Abhyankar dm->ops->creatematrix = DMCreateMatrix_Network; 334ea78f98cSLisandro Dalcin dm->ops->createinterpolation = NULL; 335ea78f98cSLisandro Dalcin dm->ops->createinjection = NULL; 336ea78f98cSLisandro Dalcin dm->ops->refine = NULL; 337ea78f98cSLisandro Dalcin dm->ops->coarsen = NULL; 338ea78f98cSLisandro Dalcin dm->ops->refinehierarchy = NULL; 339ea78f98cSLisandro Dalcin dm->ops->coarsenhierarchy = NULL; 3405f2c45f1SShri Abhyankar dm->ops->globaltolocalbegin = DMGlobalToLocalBegin_Network; 3415f2c45f1SShri Abhyankar dm->ops->globaltolocalend = DMGlobalToLocalEnd_Network; 3425f2c45f1SShri Abhyankar dm->ops->localtoglobalbegin = DMLocalToGlobalBegin_Network; 3435f2c45f1SShri Abhyankar dm->ops->localtoglobalend = DMLocalToGlobalEnd_Network; 3445f2c45f1SShri Abhyankar dm->ops->destroy = DMDestroy_Network; 3455f2c45f1SShri Abhyankar dm->ops->createsubdm = NULL; 3465f2c45f1SShri Abhyankar dm->ops->locatepoints = NULL; 3475f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3485f2c45f1SShri Abhyankar } 349*daad07d3SAidan Hamilton /* 350*daad07d3SAidan Hamilton copies over the subnetid and index portions of the DMNetworkComponentHeader from orignal dm to the newdm 351*daad07d3SAidan Hamilton */ 352*daad07d3SAidan Hamilton static PetscErrorCode DMNetworkCopyHeaderTopological(DM dm, DM newdm) 353*daad07d3SAidan Hamilton { 354*daad07d3SAidan Hamilton DM_Network *network = (DM_Network *) dm->data, *newnetwork = (DM_Network *) newdm->data; 355*daad07d3SAidan Hamilton PetscInt p,i,np,index,subnetid; 356*daad07d3SAidan Hamilton 357*daad07d3SAidan Hamilton PetscFunctionBegin; 358*daad07d3SAidan Hamilton np = network->cloneshared->pEnd - network->cloneshared->pStart; 359*daad07d3SAidan Hamilton PetscCall(PetscCalloc2(np,&newnetwork->header,np,&newnetwork->cvalue)); 360*daad07d3SAidan Hamilton for (i = 0; i < np; i++ ) { 361*daad07d3SAidan Hamilton p = i + network->cloneshared->pStart; 362*daad07d3SAidan Hamilton PetscCall(DMNetworkGetSubnetID(dm,p,&subnetid)); 363*daad07d3SAidan Hamilton PetscCall(DMNetworkGetIndex(dm,p,&index)); 364*daad07d3SAidan Hamilton newnetwork->header[i].index = index; 365*daad07d3SAidan Hamilton newnetwork->header[i].subnetid = subnetid; 366*daad07d3SAidan Hamilton newnetwork->header[i].size = NULL; 367*daad07d3SAidan Hamilton newnetwork->header[i].key = NULL; 368*daad07d3SAidan Hamilton newnetwork->header[i].offset = NULL; 369*daad07d3SAidan Hamilton newnetwork->header[i].nvar = NULL; 370*daad07d3SAidan Hamilton newnetwork->header[i].offsetvarrel = NULL; 371*daad07d3SAidan Hamilton newnetwork->header[i].ndata = 0; 372*daad07d3SAidan Hamilton newnetwork->header[i].maxcomps = DMNETWORK_MAX_COMP_AT_POINT_DEFAULT; 373*daad07d3SAidan Hamilton newnetwork->header[i].hsize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(sizeof(DMNetworkComponentGenericDataType)); 374*daad07d3SAidan Hamilton } 375*daad07d3SAidan Hamilton PetscFunctionReturn(0); 376*daad07d3SAidan Hamilton } 3775f2c45f1SShri Abhyankar 3788415c774SShri Abhyankar PetscErrorCode DMClone_Network(DM dm, DM *newdm) 3798415c774SShri Abhyankar { 380*daad07d3SAidan Hamilton DM_Network *network = (DM_Network *) dm->data, *newnetwork = NULL; 3818415c774SShri Abhyankar 3828415c774SShri Abhyankar PetscFunctionBegin; 383*daad07d3SAidan Hamilton network->cloneshared->refct++; 384*daad07d3SAidan Hamilton PetscCall(PetscNewLog(*newdm,&newnetwork)); 385*daad07d3SAidan Hamilton (*newdm)->data = newnetwork; 386*daad07d3SAidan Hamilton PetscCall(DMNetworkInitializeToDefault_NonShared(*newdm)); 387*daad07d3SAidan Hamilton newnetwork->cloneshared = network->cloneshared; /* Share all data that can be cloneshared */ 388*daad07d3SAidan Hamilton PetscCall(DMClone(network->plex,&newnetwork->plex)); 389*daad07d3SAidan Hamilton PetscCall(DMNetworkCopyHeaderTopological(dm, *newdm)); 390*daad07d3SAidan Hamilton PetscCall(DMNetworkInitializeNonTopological(*newdm)); /* initialize all non-topological data to the state after DMNetworkLayoutSetUp as been called */ 3919566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject) *newdm, DMNETWORK)); 3929566063dSJacob Faibussowitsch PetscCall(DMInitialize_Network(*newdm)); 3938415c774SShri Abhyankar PetscFunctionReturn(0); 3948415c774SShri Abhyankar } 3958415c774SShri Abhyankar 3965f2c45f1SShri Abhyankar /*MC 3975f2c45f1SShri Abhyankar DMNETWORK = "network" - A DM object that encapsulates an unstructured network. The implementation is based on the DM object 3985f2c45f1SShri Abhyankar DMPlex that manages unstructured grids. Distributed networks use a non-overlapping partitioning of 3995f2c45f1SShri Abhyankar the edges. In the local representation, Vecs contain all unknowns in the interior and shared boundary. 4005f2c45f1SShri Abhyankar This is specified by a PetscSection object. Ownership in the global representation is determined by 4015f2c45f1SShri Abhyankar ownership of the underlying DMPlex points. This is specified by another PetscSection object. 4025f2c45f1SShri Abhyankar 4035f2c45f1SShri Abhyankar Level: intermediate 4045f2c45f1SShri Abhyankar 405db781477SPatrick Sanan .seealso: `DMType`, `DMNetworkCreate()`, `DMCreate()`, `DMSetType()` 4065f2c45f1SShri Abhyankar M*/ 4075f2c45f1SShri Abhyankar 4085f2c45f1SShri Abhyankar PETSC_EXTERN PetscErrorCode DMCreate_Network(DM dm) 4095f2c45f1SShri Abhyankar { 4105f2c45f1SShri Abhyankar DM_Network *network; 4115f2c45f1SShri Abhyankar 4125f2c45f1SShri Abhyankar PetscFunctionBegin; 4135f2c45f1SShri Abhyankar PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4149566063dSJacob Faibussowitsch PetscCall(PetscNewLog(dm,&network)); 415*daad07d3SAidan Hamilton PetscCall(PetscNewLog(dm,&network->cloneshared)); 4165f2c45f1SShri Abhyankar dm->data = network; 4175f2c45f1SShri Abhyankar 418*daad07d3SAidan Hamilton PetscCall(DMNetworkInitializeToDefault(dm)); 4199566063dSJacob Faibussowitsch PetscCall(DMInitialize_Network(dm)); 4205f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4215f2c45f1SShri Abhyankar } 4225f2c45f1SShri Abhyankar 4235f2c45f1SShri Abhyankar /*@ 4245f2c45f1SShri Abhyankar DMNetworkCreate - Creates a DMNetwork object, which encapsulates an unstructured network. 4255f2c45f1SShri Abhyankar 426d083f849SBarry Smith Collective 4275f2c45f1SShri Abhyankar 4285f2c45f1SShri Abhyankar Input Parameter: 4295f2c45f1SShri Abhyankar . comm - The communicator for the DMNetwork object 4305f2c45f1SShri Abhyankar 4315f2c45f1SShri Abhyankar Output Parameter: 4325f2c45f1SShri Abhyankar . network - The DMNetwork object 4335f2c45f1SShri Abhyankar 4345f2c45f1SShri Abhyankar Level: beginner 4355f2c45f1SShri Abhyankar 4365f2c45f1SShri Abhyankar @*/ 4375f2c45f1SShri Abhyankar PetscErrorCode DMNetworkCreate(MPI_Comm comm, DM *network) 4385f2c45f1SShri Abhyankar { 4395f2c45f1SShri Abhyankar PetscFunctionBegin; 4405f2c45f1SShri Abhyankar PetscValidPointer(network,2); 4419566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, network)); 4429566063dSJacob Faibussowitsch PetscCall(DMSetType(*network, DMNETWORK)); 4435f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4445f2c45f1SShri Abhyankar } 445