xref: /petsc/src/dm/impls/network/networkcreate.c (revision a81b7fe49b8f97d21ac9ed59c77b9cfc2b7a3c03)
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 
27*a81b7fe4SHong Zhang static PetscErrorCode VecArrayPrint_private(PetscViewer viewer,PetscInt n,const PetscScalar *xv)
28*a81b7fe4SHong Zhang {
29*a81b7fe4SHong Zhang   PetscErrorCode ierr;
30*a81b7fe4SHong Zhang   PetscInt       i;
31*a81b7fe4SHong Zhang 
32*a81b7fe4SHong Zhang   PetscFunctionBegin;
33*a81b7fe4SHong Zhang   for (i=0; i<n; i++) {
34*a81b7fe4SHong Zhang #if defined(PETSC_USE_COMPLEX)
35*a81b7fe4SHong Zhang     if (PetscImaginaryPart(xv[i]) > 0.0) {
36*a81b7fe4SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"    %g + %g i\n",(double)PetscRealPart(xv[i]),(double)PetscImaginaryPart(xv[i]));CHKERRQ(ierr);
37*a81b7fe4SHong Zhang     } else if (PetscImaginaryPart(xv[i]) < 0.0) {
38*a81b7fe4SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"    %g - %g i\n",(double)PetscRealPart(xv[i]),-(double)PetscImaginaryPart(xv[i]));CHKERRQ(ierr);
39*a81b7fe4SHong Zhang     } else {
40*a81b7fe4SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"    %g\n",(double)PetscRealPart(xv[i]));CHKERRQ(ierr);
41*a81b7fe4SHong Zhang     }
42*a81b7fe4SHong Zhang #else
43*a81b7fe4SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"    %g\n",(double)xv[i]);CHKERRQ(ierr);
44*a81b7fe4SHong Zhang #endif
45*a81b7fe4SHong Zhang   }
46*a81b7fe4SHong Zhang   PetscFunctionReturn(0);
47*a81b7fe4SHong Zhang }
48*a81b7fe4SHong Zhang 
494dc485aaSHong Zhang static PetscErrorCode VecView_Network_Seq(DM networkdm,Vec X,PetscViewer viewer)
504dc485aaSHong Zhang {
514dc485aaSHong Zhang   PetscErrorCode    ierr;
52*a81b7fe4SHong Zhang   PetscInt          e,v,Start,End,offset,nvar;
534dc485aaSHong Zhang   const PetscScalar *xv;
544dc485aaSHong Zhang 
554dc485aaSHong Zhang   PetscFunctionBegin;
564dc485aaSHong Zhang   ierr = VecGetArrayRead(X,&xv);CHKERRQ(ierr);
574dc485aaSHong Zhang 
584dc485aaSHong Zhang   /* iterate over edges */
594dc485aaSHong Zhang   ierr = DMNetworkGetEdgeRange(networkdm,&Start,&End);CHKERRQ(ierr);
604dc485aaSHong Zhang   for (e=Start; e<End; e++) {
614dc485aaSHong Zhang     ierr = DMNetworkGetVariableOffset(networkdm,e,&offset);CHKERRQ(ierr);
624dc485aaSHong Zhang     ierr = DMNetworkGetNumVariables(networkdm,e,&nvar);CHKERRQ(ierr);
634dc485aaSHong Zhang     if (!nvar) continue;
64*a81b7fe4SHong Zhang 
654dc485aaSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Edge %D:\n",e-Start);CHKERRQ(ierr);
66*a81b7fe4SHong Zhang     ierr = VecArrayPrint_private(viewer,nvar,xv+offset);CHKERRQ(ierr);
674dc485aaSHong Zhang   }
684dc485aaSHong Zhang 
694dc485aaSHong Zhang   /* iterate over vertices */
704dc485aaSHong Zhang   ierr = DMNetworkGetVertexRange(networkdm,&Start,&End);CHKERRQ(ierr);
714dc485aaSHong Zhang   for (v=Start; v<End; v++) {
724dc485aaSHong Zhang     ierr = DMNetworkGetVariableOffset(networkdm,v,&offset);CHKERRQ(ierr);
734dc485aaSHong Zhang     ierr = DMNetworkGetNumVariables(networkdm,v,&nvar);CHKERRQ(ierr);
744dc485aaSHong Zhang     if (!nvar) continue;
75*a81b7fe4SHong Zhang 
764dc485aaSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Vertex %D:\n",v-Start);CHKERRQ(ierr);
77*a81b7fe4SHong Zhang     ierr = VecArrayPrint_private(viewer,nvar,xv+offset);CHKERRQ(ierr);
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 
844062a5e5SHong Zhang static PetscErrorCode VecView_Network_MPI(DM networkdm,Vec X,PetscViewer viewer)
854062a5e5SHong Zhang {
864062a5e5SHong Zhang   PetscErrorCode    ierr;
87e8c2e809SHong Zhang   PetscInt          i,e,v,eStart,eEnd,vStart,vEnd,offset,nvar,len_loc,len,k;
884062a5e5SHong Zhang   const PetscScalar *xv;
894062a5e5SHong Zhang   MPI_Comm          comm;
904062a5e5SHong Zhang   PetscMPIInt       size,rank,tag = ((PetscObject)viewer)->tag;
914062a5e5SHong Zhang   Vec               localX;
924062a5e5SHong Zhang   PetscBool         ghostvtex;
934062a5e5SHong Zhang   PetscScalar       *values;
94*a81b7fe4SHong Zhang   PetscInt          j,n,ne,nv;
954062a5e5SHong Zhang   MPI_Status        status;
964062a5e5SHong Zhang 
974062a5e5SHong Zhang   PetscFunctionBegin;
984062a5e5SHong Zhang   ierr = PetscObjectGetComm((PetscObject)networkdm,&comm);CHKERRQ(ierr);
994062a5e5SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1004062a5e5SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1014062a5e5SHong Zhang 
1024062a5e5SHong Zhang   ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr);
1034062a5e5SHong Zhang   ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
1044062a5e5SHong Zhang   ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
1054062a5e5SHong Zhang   ierr = VecGetArrayRead(localX,&xv);CHKERRQ(ierr);
1064062a5e5SHong Zhang 
1074062a5e5SHong Zhang   ierr = VecGetLocalSize(localX,&len_loc);CHKERRQ(ierr);
1084062a5e5SHong Zhang 
109e8c2e809SHong Zhang   ierr = DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);CHKERRQ(ierr);
1104062a5e5SHong Zhang   ierr = DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);CHKERRQ(ierr);
111e8c2e809SHong Zhang   len_loc += 2*(1 + eEnd-eStart + vEnd-vStart);
1124062a5e5SHong Zhang 
1134062a5e5SHong Zhang   /* values = [nedges, offset, nvar, xe; nvertices, offsets, nvars, xv] */
1144062a5e5SHong Zhang   ierr = MPI_Allreduce(&len_loc,&len,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr);
1154062a5e5SHong Zhang   ierr = PetscCalloc1(len,&values);CHKERRQ(ierr);
1164062a5e5SHong Zhang 
1174062a5e5SHong Zhang   if (!rank) {
1184062a5e5SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Process [%d]\n",rank);CHKERRQ(ierr);
1194062a5e5SHong Zhang   }
1204062a5e5SHong Zhang 
1214062a5e5SHong Zhang   /* iterate over edges */
1224062a5e5SHong Zhang   k = 2;
123e8c2e809SHong Zhang   for (e=eStart; e<eEnd; e++) {
1244062a5e5SHong Zhang     ierr = DMNetworkGetVariableOffset(networkdm,e,&offset);CHKERRQ(ierr);
1254062a5e5SHong Zhang     ierr = DMNetworkGetNumVariables(networkdm,e,&nvar);CHKERRQ(ierr);
1264062a5e5SHong Zhang     if (!nvar) continue;
1274062a5e5SHong Zhang 
1284062a5e5SHong Zhang     values[0]  += 1; /* nedges */
1294062a5e5SHong Zhang     values[k++] = offset;
1304062a5e5SHong Zhang     values[k++] = nvar;
1314062a5e5SHong Zhang 
1324062a5e5SHong Zhang     if (!rank) { /* print its own entries */
133e8c2e809SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  Edge %D:\n",e-eStart);CHKERRQ(ierr);
134*a81b7fe4SHong Zhang       ierr = VecArrayPrint_private(viewer,nvar,xv+offset);CHKERRQ(ierr);
135*a81b7fe4SHong Zhang     } else {
1364062a5e5SHong Zhang       for (i=offset; i< offset+nvar; i++) {
1374062a5e5SHong Zhang         values[k++] = xv[i];
1384062a5e5SHong Zhang       }
1394062a5e5SHong Zhang     }
1404062a5e5SHong Zhang   }
1414062a5e5SHong Zhang 
1424062a5e5SHong Zhang   /* iterate over vertices */
1434062a5e5SHong Zhang   for (v=vStart; v<vEnd; v++) {
1444062a5e5SHong Zhang     ierr = DMNetworkIsGhostVertex(networkdm,v,&ghostvtex);CHKERRQ(ierr);
1454062a5e5SHong Zhang     if (ghostvtex) continue;
1464062a5e5SHong Zhang     ierr = DMNetworkGetVariableOffset(networkdm,v,&offset);CHKERRQ(ierr);
1474062a5e5SHong Zhang     ierr = DMNetworkGetNumVariables(networkdm,v,&nvar);CHKERRQ(ierr);
1484062a5e5SHong Zhang     if (!nvar) continue;
1494062a5e5SHong Zhang 
1504062a5e5SHong Zhang     values[1]  += 1; /* nvertices */
1514062a5e5SHong Zhang     values[k++] = offset;
1524062a5e5SHong Zhang     values[k++] = nvar;
1534062a5e5SHong Zhang 
1544062a5e5SHong Zhang     if (!rank) {
155e8c2e809SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  Vertex %D:\n",v-vStart);CHKERRQ(ierr);
1564062a5e5SHong Zhang     }
1574062a5e5SHong Zhang 
1584062a5e5SHong Zhang     if (!rank) {
159*a81b7fe4SHong Zhang       ierr = VecArrayPrint_private(viewer,nvar,xv+offset);CHKERRQ(ierr);
1604062a5e5SHong Zhang     } else {
161*a81b7fe4SHong Zhang       for (i=offset; i< offset+nvar; i++) {
1624062a5e5SHong Zhang         values[k++] = xv[i];
1634062a5e5SHong Zhang       }
1644062a5e5SHong Zhang     }
1654062a5e5SHong Zhang   }
1664062a5e5SHong Zhang 
1674062a5e5SHong Zhang   if (!rank) {
1684062a5e5SHong Zhang     /* proc[0] receives and prints messages */
1694062a5e5SHong Zhang     for (j=1; j<size; j++) {
1704062a5e5SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Process [%d]\n",j);CHKERRQ(ierr);
1714062a5e5SHong Zhang 
1724062a5e5SHong Zhang       ierr = MPI_Recv(values,(PetscMPIInt)len,MPIU_SCALAR,j,tag,comm,&status);CHKERRQ(ierr);
1734062a5e5SHong Zhang       ierr = MPI_Get_count(&status,MPIU_SCALAR,&n);CHKERRQ(ierr);
1744062a5e5SHong Zhang 
1754062a5e5SHong Zhang       ne = (PetscInt)values[0];
1764062a5e5SHong Zhang       nv = (PetscInt)values[1];
1774062a5e5SHong Zhang 
1784062a5e5SHong Zhang       /* print received edges */
1794062a5e5SHong Zhang       k  = 2;
1804062a5e5SHong Zhang       for (i=0; i<ne; i++) {
181*a81b7fe4SHong Zhang         offset = (PetscInt)values[k++];
182*a81b7fe4SHong Zhang         nvar   = (PetscInt)values[k++];
1834062a5e5SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  Edge %D:\n",i);CHKERRQ(ierr);
184*a81b7fe4SHong Zhang         ierr = VecArrayPrint_private(viewer,nvar,values+k);CHKERRQ(ierr);
185*a81b7fe4SHong Zhang         k   += nvar;
1864062a5e5SHong Zhang       }
1874062a5e5SHong Zhang 
1884062a5e5SHong Zhang       /* print received vertices */
1894062a5e5SHong Zhang       for (i=0; i<nv; i++) {
1904062a5e5SHong Zhang         offset = (PetscInt)values[k++];
1914062a5e5SHong Zhang         nvar   = (PetscInt)values[k++];
1924062a5e5SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  Vertex %D:\n",i);CHKERRQ(ierr);
193*a81b7fe4SHong Zhang         ierr = VecArrayPrint_private(viewer,nvar,values+k);CHKERRQ(ierr);
194*a81b7fe4SHong Zhang         k   += nvar;
1954062a5e5SHong Zhang       }
1964062a5e5SHong Zhang     }
1974062a5e5SHong Zhang   } else {
1984062a5e5SHong Zhang     /* sends values to proc[0] */
1994062a5e5SHong Zhang     ierr = MPI_Send((void*)values,k,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr);
2004062a5e5SHong Zhang   }
2014062a5e5SHong Zhang 
2024062a5e5SHong Zhang   ierr = PetscFree(values);CHKERRQ(ierr);
2034062a5e5SHong Zhang   ierr = VecRestoreArrayRead(localX,&xv);CHKERRQ(ierr);
2044062a5e5SHong Zhang   ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr);
2054062a5e5SHong Zhang   PetscFunctionReturn(0);
2064062a5e5SHong Zhang }
2074062a5e5SHong Zhang 
2084dc485aaSHong Zhang PetscErrorCode VecView_Network(Vec v,PetscViewer viewer)
2094dc485aaSHong Zhang {
2104dc485aaSHong Zhang   DM             dm;
2114dc485aaSHong Zhang   PetscErrorCode ierr;
2124dc485aaSHong Zhang   PetscBool      isseq;
2134dc485aaSHong Zhang 
2144dc485aaSHong Zhang   PetscFunctionBegin;
2154dc485aaSHong Zhang   ierr = VecGetDM(v,&dm);CHKERRQ(ierr);
2164dc485aaSHong Zhang   if (!dm) SETERRQ(PetscObjectComm((PetscObject)v),PETSC_ERR_ARG_WRONG,"Vector not generated from a DM");
2174dc485aaSHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)v,VECSEQ,&isseq);CHKERRQ(ierr);
2184dc485aaSHong Zhang   if (isseq) {
2194dc485aaSHong Zhang     ierr = VecView_Network_Seq(dm,v,viewer);CHKERRQ(ierr);
2204dc485aaSHong Zhang   } else {
2214062a5e5SHong Zhang     ierr = VecView_Network_MPI(dm,v,viewer);CHKERRQ(ierr);
2224dc485aaSHong Zhang   }
2234dc485aaSHong Zhang   PetscFunctionReturn(0);
2244dc485aaSHong Zhang }
2255f2c45f1SShri Abhyankar 
2265f2c45f1SShri Abhyankar static PetscErrorCode DMCreateGlobalVector_Network(DM dm,Vec *vec)
2275f2c45f1SShri Abhyankar {
2285f2c45f1SShri Abhyankar   PetscErrorCode ierr;
2295f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
2305f2c45f1SShri Abhyankar 
2315f2c45f1SShri Abhyankar   PetscFunctionBegin;
2325f2c45f1SShri Abhyankar   ierr = DMCreateGlobalVector(network->plex,vec);CHKERRQ(ierr);
2334dc485aaSHong Zhang   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Network);CHKERRQ(ierr);
2345f2c45f1SShri Abhyankar   ierr = VecSetDM(*vec,dm);CHKERRQ(ierr);
2355f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2365f2c45f1SShri Abhyankar }
2375f2c45f1SShri Abhyankar 
2385f2c45f1SShri Abhyankar static PetscErrorCode DMCreateLocalVector_Network(DM dm,Vec *vec)
2395f2c45f1SShri Abhyankar {
2405f2c45f1SShri Abhyankar   PetscErrorCode ierr;
2415f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
2425f2c45f1SShri Abhyankar 
2435f2c45f1SShri Abhyankar   PetscFunctionBegin;
2445f2c45f1SShri Abhyankar   ierr = DMCreateLocalVector(network->plex,vec);CHKERRQ(ierr);
2455f2c45f1SShri Abhyankar   ierr = VecSetDM(*vec,dm);CHKERRQ(ierr);
2465f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2475f2c45f1SShri Abhyankar }
2485f2c45f1SShri Abhyankar 
2495f2c45f1SShri Abhyankar PetscErrorCode DMInitialize_Network(DM dm)
2505f2c45f1SShri Abhyankar {
2515f2c45f1SShri Abhyankar 
2525f2c45f1SShri Abhyankar   PetscFunctionBegin;
2535f2c45f1SShri Abhyankar 
2545f2c45f1SShri Abhyankar   dm->ops->view                            = NULL;
2555f2c45f1SShri Abhyankar   dm->ops->setfromoptions                  = DMSetFromOptions_Network;
2568415c774SShri Abhyankar   dm->ops->clone                           = DMClone_Network;
2575f2c45f1SShri Abhyankar   dm->ops->setup                           = DMSetUp_Network;
2585f2c45f1SShri Abhyankar   dm->ops->createglobalvector              = DMCreateGlobalVector_Network;
2595f2c45f1SShri Abhyankar   dm->ops->createlocalvector               = DMCreateLocalVector_Network;
2605f2c45f1SShri Abhyankar   dm->ops->getlocaltoglobalmapping         = NULL;
2615f2c45f1SShri Abhyankar   dm->ops->createfieldis                   = NULL;
2625f2c45f1SShri Abhyankar   dm->ops->createcoordinatedm              = NULL;
2635f2c45f1SShri Abhyankar   dm->ops->getcoloring                     = 0;
2645f2c45f1SShri Abhyankar   dm->ops->creatematrix                    = DMCreateMatrix_Network;
2655f2c45f1SShri Abhyankar   dm->ops->createinterpolation             = 0;
2665f2c45f1SShri Abhyankar   dm->ops->getaggregates                   = 0;
2675f2c45f1SShri Abhyankar   dm->ops->getinjection                    = 0;
2685f2c45f1SShri Abhyankar   dm->ops->refine                          = 0;
2695f2c45f1SShri Abhyankar   dm->ops->coarsen                         = 0;
2705f2c45f1SShri Abhyankar   dm->ops->refinehierarchy                 = 0;
2715f2c45f1SShri Abhyankar   dm->ops->coarsenhierarchy                = 0;
2725f2c45f1SShri Abhyankar   dm->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Network;
2735f2c45f1SShri Abhyankar   dm->ops->globaltolocalend                = DMGlobalToLocalEnd_Network;
2745f2c45f1SShri Abhyankar   dm->ops->localtoglobalbegin              = DMLocalToGlobalBegin_Network;
2755f2c45f1SShri Abhyankar   dm->ops->localtoglobalend                = DMLocalToGlobalEnd_Network;
2765f2c45f1SShri Abhyankar   dm->ops->destroy                         = DMDestroy_Network;
2775f2c45f1SShri Abhyankar   dm->ops->createsubdm                     = NULL;
2785f2c45f1SShri Abhyankar   dm->ops->locatepoints                    = NULL;
2795f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2805f2c45f1SShri Abhyankar }
2815f2c45f1SShri Abhyankar 
2828415c774SShri Abhyankar PetscErrorCode DMClone_Network(DM dm, DM *newdm)
2838415c774SShri Abhyankar {
2848415c774SShri Abhyankar   DM_Network     *network = (DM_Network *) dm->data;
2858415c774SShri Abhyankar   PetscErrorCode ierr;
2868415c774SShri Abhyankar 
2878415c774SShri Abhyankar   PetscFunctionBegin;
2888415c774SShri Abhyankar   network->refct++;
2898415c774SShri Abhyankar   (*newdm)->data = network;
2908415c774SShri Abhyankar   ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMNETWORK);CHKERRQ(ierr);
2918415c774SShri Abhyankar   ierr = DMInitialize_Network(*newdm);CHKERRQ(ierr);
2928415c774SShri Abhyankar   PetscFunctionReturn(0);
2938415c774SShri Abhyankar }
2948415c774SShri Abhyankar 
2955f2c45f1SShri Abhyankar /*MC
2965f2c45f1SShri Abhyankar   DMNETWORK = "network" - A DM object that encapsulates an unstructured network. The implementation is based on the DM object
2975f2c45f1SShri Abhyankar                           DMPlex that manages unstructured grids. Distributed networks use a non-overlapping partitioning of
2985f2c45f1SShri Abhyankar                           the edges. In the local representation, Vecs contain all unknowns in the interior and shared boundary.
2995f2c45f1SShri Abhyankar                           This is specified by a PetscSection object. Ownership in the global representation is determined by
3005f2c45f1SShri Abhyankar                           ownership of the underlying DMPlex points. This is specified by another PetscSection object.
3015f2c45f1SShri Abhyankar 
3025f2c45f1SShri Abhyankar   Level: intermediate
3035f2c45f1SShri Abhyankar 
3045f2c45f1SShri Abhyankar .seealso: DMType, DMNetworkCreate(), DMCreate(), DMSetType()
3055f2c45f1SShri Abhyankar M*/
3065f2c45f1SShri Abhyankar 
3075f2c45f1SShri Abhyankar PETSC_EXTERN PetscErrorCode DMCreate_Network(DM dm)
3085f2c45f1SShri Abhyankar {
3095f2c45f1SShri Abhyankar   DM_Network     *network;
3105f2c45f1SShri Abhyankar   PetscErrorCode ierr;
3115f2c45f1SShri Abhyankar 
3125f2c45f1SShri Abhyankar   PetscFunctionBegin;
3135f2c45f1SShri Abhyankar   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3145f2c45f1SShri Abhyankar   ierr     = PetscNewLog(dm,&network);CHKERRQ(ierr);
3155f2c45f1SShri Abhyankar   dm->data = network;
3165f2c45f1SShri Abhyankar 
3175f2c45f1SShri Abhyankar   network->refct          = 1;
3185f2c45f1SShri Abhyankar   network->NNodes         = -1;
3195f2c45f1SShri Abhyankar   network->NEdges         = -1;
3205f2c45f1SShri Abhyankar   network->nNodes         = -1;
3215f2c45f1SShri Abhyankar   network->nEdges         = -1;
3225f2c45f1SShri Abhyankar 
32313c2a604SAdrian Maldonado 
3245f2c45f1SShri Abhyankar   ierr = DMInitialize_Network(dm);CHKERRQ(ierr);
3255f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3265f2c45f1SShri Abhyankar }
3275f2c45f1SShri Abhyankar 
3285f2c45f1SShri Abhyankar /*@
3295f2c45f1SShri Abhyankar   DMNetworkCreate - Creates a DMNetwork object, which encapsulates an unstructured network.
3305f2c45f1SShri Abhyankar 
3315f2c45f1SShri Abhyankar   Collective on MPI_Comm
3325f2c45f1SShri Abhyankar 
3335f2c45f1SShri Abhyankar   Input Parameter:
3345f2c45f1SShri Abhyankar . comm - The communicator for the DMNetwork object
3355f2c45f1SShri Abhyankar 
3365f2c45f1SShri Abhyankar   Output Parameter:
3375f2c45f1SShri Abhyankar . network  - The DMNetwork object
3385f2c45f1SShri Abhyankar 
3395f2c45f1SShri Abhyankar   Level: beginner
3405f2c45f1SShri Abhyankar 
3415f2c45f1SShri Abhyankar .keywords: DMNetwork, create
3425f2c45f1SShri Abhyankar @*/
3435f2c45f1SShri Abhyankar PetscErrorCode DMNetworkCreate(MPI_Comm comm, DM *network)
3445f2c45f1SShri Abhyankar {
3455f2c45f1SShri Abhyankar   PetscErrorCode ierr;
3465f2c45f1SShri Abhyankar 
3475f2c45f1SShri Abhyankar   PetscFunctionBegin;
3485f2c45f1SShri Abhyankar   PetscValidPointer(network,2);
3495f2c45f1SShri Abhyankar   ierr = DMCreate(comm, network);CHKERRQ(ierr);
3505f2c45f1SShri Abhyankar   ierr = DMSetType(*network, DMNETWORK);CHKERRQ(ierr);
3515f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3525f2c45f1SShri Abhyankar }
353