xref: /petsc/src/dm/impls/network/networkcreate.c (revision 7b6afd5b8febdfe9eecd254e43b0d168d52e271f)
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 
27a81b7fe4SHong Zhang static PetscErrorCode VecArrayPrint_private(PetscViewer viewer,PetscInt n,const PetscScalar *xv)
28a81b7fe4SHong Zhang {
29a81b7fe4SHong Zhang   PetscErrorCode ierr;
30a81b7fe4SHong Zhang   PetscInt       i;
31a81b7fe4SHong Zhang 
32a81b7fe4SHong Zhang   PetscFunctionBegin;
33a81b7fe4SHong Zhang   for (i=0; i<n; i++) {
34a81b7fe4SHong Zhang #if defined(PETSC_USE_COMPLEX)
35a81b7fe4SHong Zhang     if (PetscImaginaryPart(xv[i]) > 0.0) {
36a81b7fe4SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"    %g + %g i\n",(double)PetscRealPart(xv[i]),(double)PetscImaginaryPart(xv[i]));CHKERRQ(ierr);
37a81b7fe4SHong Zhang     } else if (PetscImaginaryPart(xv[i]) < 0.0) {
38a81b7fe4SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"    %g - %g i\n",(double)PetscRealPart(xv[i]),-(double)PetscImaginaryPart(xv[i]));CHKERRQ(ierr);
39a81b7fe4SHong Zhang     } else {
40a81b7fe4SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"    %g\n",(double)PetscRealPart(xv[i]));CHKERRQ(ierr);
41a81b7fe4SHong Zhang     }
42a81b7fe4SHong Zhang #else
43a81b7fe4SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"    %g\n",(double)xv[i]);CHKERRQ(ierr);
44a81b7fe4SHong Zhang #endif
45a81b7fe4SHong Zhang   }
46a81b7fe4SHong Zhang   PetscFunctionReturn(0);
47a81b7fe4SHong Zhang }
48a81b7fe4SHong Zhang 
494dc485aaSHong Zhang static PetscErrorCode VecView_Network_Seq(DM networkdm,Vec X,PetscViewer viewer)
504dc485aaSHong Zhang {
514dc485aaSHong Zhang   PetscErrorCode    ierr;
52*7b6afd5bSHong Zhang   PetscInt          e,v,Start,End,offset,nvar,id;
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 = DMNetworkGetNumVariables(networkdm,e,&nvar);CHKERRQ(ierr);
624dc485aaSHong Zhang     if (!nvar) continue;
63a81b7fe4SHong Zhang 
64*7b6afd5bSHong Zhang     ierr = DMNetworkGetVariableOffset(networkdm,e,&offset);CHKERRQ(ierr);
65*7b6afd5bSHong Zhang     ierr = DMNetworkGetEdgeVertexID(networkdm,e,&id);CHKERRQ(ierr);
66*7b6afd5bSHong Zhang 
67*7b6afd5bSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Edge %D:\n",id);CHKERRQ(ierr);
68a81b7fe4SHong Zhang     ierr = VecArrayPrint_private(viewer,nvar,xv+offset);CHKERRQ(ierr);
694dc485aaSHong Zhang   }
704dc485aaSHong Zhang 
714dc485aaSHong Zhang   /* iterate over vertices */
724dc485aaSHong Zhang   ierr = DMNetworkGetVertexRange(networkdm,&Start,&End);CHKERRQ(ierr);
734dc485aaSHong Zhang   for (v=Start; v<End; v++) {
744dc485aaSHong Zhang     ierr = DMNetworkGetNumVariables(networkdm,v,&nvar);CHKERRQ(ierr);
754dc485aaSHong Zhang     if (!nvar) continue;
76a81b7fe4SHong Zhang 
77*7b6afd5bSHong Zhang     ierr = DMNetworkGetVariableOffset(networkdm,v,&offset);CHKERRQ(ierr);
78*7b6afd5bSHong Zhang     ierr = DMNetworkGetEdgeVertexID(networkdm,v,&id);CHKERRQ(ierr);
7962aa33baSHong Zhang 
80*7b6afd5bSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Vertex %D:\n",id);CHKERRQ(ierr);
81a81b7fe4SHong Zhang     ierr = VecArrayPrint_private(viewer,nvar,xv+offset);CHKERRQ(ierr);
824dc485aaSHong Zhang   }
834dc485aaSHong Zhang   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
844dc485aaSHong Zhang   ierr = VecRestoreArrayRead(X,&xv);CHKERRQ(ierr);
854dc485aaSHong Zhang   PetscFunctionReturn(0);
864dc485aaSHong Zhang }
874dc485aaSHong Zhang 
884062a5e5SHong Zhang static PetscErrorCode VecView_Network_MPI(DM networkdm,Vec X,PetscViewer viewer)
894062a5e5SHong Zhang {
904062a5e5SHong Zhang   PetscErrorCode    ierr;
91e8c2e809SHong Zhang   PetscInt          i,e,v,eStart,eEnd,vStart,vEnd,offset,nvar,len_loc,len,k;
924062a5e5SHong Zhang   const PetscScalar *xv;
934062a5e5SHong Zhang   MPI_Comm          comm;
944062a5e5SHong Zhang   PetscMPIInt       size,rank,tag = ((PetscObject)viewer)->tag;
954062a5e5SHong Zhang   Vec               localX;
964062a5e5SHong Zhang   PetscBool         ghostvtex;
974062a5e5SHong Zhang   PetscScalar       *values;
98*7b6afd5bSHong Zhang   PetscInt          j,ne,nv,id;
994062a5e5SHong Zhang   MPI_Status        status;
1004062a5e5SHong Zhang 
1014062a5e5SHong Zhang   PetscFunctionBegin;
1024062a5e5SHong Zhang   ierr = PetscObjectGetComm((PetscObject)networkdm,&comm);CHKERRQ(ierr);
1034062a5e5SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1044062a5e5SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1054062a5e5SHong Zhang 
1064062a5e5SHong Zhang   ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr);
1074062a5e5SHong Zhang   ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
1084062a5e5SHong Zhang   ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
1094062a5e5SHong Zhang   ierr = VecGetArrayRead(localX,&xv);CHKERRQ(ierr);
1104062a5e5SHong Zhang 
1114062a5e5SHong Zhang   ierr = VecGetLocalSize(localX,&len_loc);CHKERRQ(ierr);
1124062a5e5SHong Zhang 
113e8c2e809SHong Zhang   ierr = DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);CHKERRQ(ierr);
1144062a5e5SHong Zhang   ierr = DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);CHKERRQ(ierr);
115e8c2e809SHong Zhang   len_loc += 2*(1 + eEnd-eStart + vEnd-vStart);
1164062a5e5SHong Zhang 
1174062a5e5SHong Zhang   /* values = [nedges, offset, nvar, xe; nvertices, offsets, nvars, xv] */
1184062a5e5SHong Zhang   ierr = MPI_Allreduce(&len_loc,&len,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr);
1194062a5e5SHong Zhang   ierr = PetscCalloc1(len,&values);CHKERRQ(ierr);
1204062a5e5SHong Zhang 
1214062a5e5SHong Zhang   if (!rank) {
1224062a5e5SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Process [%d]\n",rank);CHKERRQ(ierr);
1234062a5e5SHong Zhang   }
1244062a5e5SHong Zhang 
1254062a5e5SHong Zhang   /* iterate over edges */
1264062a5e5SHong Zhang   k = 2;
127e8c2e809SHong Zhang   for (e=eStart; e<eEnd; e++) {
1284062a5e5SHong Zhang     ierr = DMNetworkGetNumVariables(networkdm,e,&nvar);CHKERRQ(ierr);
1294062a5e5SHong Zhang     if (!nvar) continue;
1304062a5e5SHong Zhang 
131*7b6afd5bSHong Zhang     ierr = DMNetworkGetVariableOffset(networkdm,e,&offset);CHKERRQ(ierr);
132*7b6afd5bSHong Zhang     ierr = DMNetworkGetEdgeVertexID(networkdm,e,&id);CHKERRQ(ierr);
13362aa33baSHong Zhang 
1344062a5e5SHong Zhang     if (!rank) { /* print its own entries */
135*7b6afd5bSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  Edge %D:\n",id);CHKERRQ(ierr);
136a81b7fe4SHong Zhang       ierr = VecArrayPrint_private(viewer,nvar,xv+offset);CHKERRQ(ierr);
137a81b7fe4SHong Zhang     } else {
138*7b6afd5bSHong Zhang       values[0]  += 1; /* number of edges */
139*7b6afd5bSHong Zhang       values[k++] = id;
140*7b6afd5bSHong Zhang       values[k++] = nvar;
141*7b6afd5bSHong Zhang       for (i=offset; i< offset+nvar; i++) values[k++] = xv[i];
1424062a5e5SHong Zhang     }
1434062a5e5SHong Zhang   }
1444062a5e5SHong Zhang 
1454062a5e5SHong Zhang   /* iterate over vertices */
1464062a5e5SHong Zhang   for (v=vStart; v<vEnd; v++) {
1474062a5e5SHong Zhang     ierr = DMNetworkIsGhostVertex(networkdm,v,&ghostvtex);CHKERRQ(ierr);
1484062a5e5SHong Zhang     if (ghostvtex) continue;
1494062a5e5SHong Zhang     ierr = DMNetworkGetNumVariables(networkdm,v,&nvar);CHKERRQ(ierr);
1504062a5e5SHong Zhang     if (!nvar) continue;
1514062a5e5SHong Zhang 
152*7b6afd5bSHong Zhang     ierr = DMNetworkGetVariableOffset(networkdm,v,&offset);CHKERRQ(ierr);
153*7b6afd5bSHong Zhang     ierr = DMNetworkGetEdgeVertexID(networkdm,v,&id);CHKERRQ(ierr);
1544062a5e5SHong Zhang 
1554062a5e5SHong Zhang     if (!rank) {
156*7b6afd5bSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  Vertex %D:\n",id);CHKERRQ(ierr);
157a81b7fe4SHong Zhang       ierr = VecArrayPrint_private(viewer,nvar,xv+offset);CHKERRQ(ierr);
1584062a5e5SHong Zhang     } else {
159*7b6afd5bSHong Zhang       values[1]  += 1; /* number of vertices */
160*7b6afd5bSHong Zhang       values[k++] = id;
161*7b6afd5bSHong Zhang       values[k++] = nvar;
162*7b6afd5bSHong Zhang       for (i=offset; i< offset+nvar; i++) values[k++] = xv[i];
1634062a5e5SHong Zhang     }
1644062a5e5SHong Zhang   }
1654062a5e5SHong Zhang 
1664062a5e5SHong Zhang   if (!rank) {
1674062a5e5SHong Zhang     /* proc[0] receives and prints messages */
1684062a5e5SHong Zhang     for (j=1; j<size; j++) {
1694062a5e5SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Process [%d]\n",j);CHKERRQ(ierr);
1704062a5e5SHong Zhang 
1714062a5e5SHong Zhang       ierr = MPI_Recv(values,(PetscMPIInt)len,MPIU_SCALAR,j,tag,comm,&status);CHKERRQ(ierr);
1724062a5e5SHong Zhang 
1734062a5e5SHong Zhang       ne = (PetscInt)values[0];
1744062a5e5SHong Zhang       nv = (PetscInt)values[1];
1754062a5e5SHong Zhang 
1764062a5e5SHong Zhang       /* print received edges */
1774062a5e5SHong Zhang       k = 2;
1784062a5e5SHong Zhang       for (i=0; i<ne; i++) {
179*7b6afd5bSHong Zhang         id   = (PetscInt)values[k++];
180a81b7fe4SHong Zhang         nvar = (PetscInt)values[k++];
181*7b6afd5bSHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  Edge %D:\n",id);CHKERRQ(ierr);
182a81b7fe4SHong Zhang         ierr = VecArrayPrint_private(viewer,nvar,values+k);CHKERRQ(ierr);
183a81b7fe4SHong Zhang         k   += nvar;
1844062a5e5SHong Zhang       }
1854062a5e5SHong Zhang 
1864062a5e5SHong Zhang       /* print received vertices */
1874062a5e5SHong Zhang       for (i=0; i<nv; i++) {
188*7b6afd5bSHong Zhang         id   = (PetscInt)values[k++];
1894062a5e5SHong Zhang         nvar = (PetscInt)values[k++];
190*7b6afd5bSHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  Vertex %D:\n",id);CHKERRQ(ierr);
191a81b7fe4SHong Zhang         ierr = VecArrayPrint_private(viewer,nvar,values+k);CHKERRQ(ierr);
192a81b7fe4SHong Zhang         k   += nvar;
1934062a5e5SHong Zhang       }
1944062a5e5SHong Zhang     }
1954062a5e5SHong Zhang   } else {
1964062a5e5SHong Zhang     /* sends values to proc[0] */
1974062a5e5SHong Zhang     ierr = MPI_Send((void*)values,k,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr);
1984062a5e5SHong Zhang   }
1994062a5e5SHong Zhang 
2004062a5e5SHong Zhang   ierr = PetscFree(values);CHKERRQ(ierr);
2014062a5e5SHong Zhang   ierr = VecRestoreArrayRead(localX,&xv);CHKERRQ(ierr);
2024062a5e5SHong Zhang   ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr);
2034062a5e5SHong Zhang   PetscFunctionReturn(0);
2044062a5e5SHong Zhang }
2054062a5e5SHong Zhang 
2064dc485aaSHong Zhang PetscErrorCode VecView_Network(Vec v,PetscViewer viewer)
2074dc485aaSHong Zhang {
2084dc485aaSHong Zhang   DM             dm;
2094dc485aaSHong Zhang   PetscErrorCode ierr;
2104dc485aaSHong Zhang   PetscBool      isseq;
2114dc485aaSHong Zhang 
2124dc485aaSHong Zhang   PetscFunctionBegin;
2134dc485aaSHong Zhang   ierr = VecGetDM(v,&dm);CHKERRQ(ierr);
2144dc485aaSHong Zhang   if (!dm) SETERRQ(PetscObjectComm((PetscObject)v),PETSC_ERR_ARG_WRONG,"Vector not generated from a DM");
2154dc485aaSHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)v,VECSEQ,&isseq);CHKERRQ(ierr);
2164dc485aaSHong Zhang   if (isseq) {
2174dc485aaSHong Zhang     ierr = VecView_Network_Seq(dm,v,viewer);CHKERRQ(ierr);
2184dc485aaSHong Zhang   } else {
2194062a5e5SHong Zhang     ierr = VecView_Network_MPI(dm,v,viewer);CHKERRQ(ierr);
2204dc485aaSHong Zhang   }
2214dc485aaSHong Zhang   PetscFunctionReturn(0);
2224dc485aaSHong Zhang }
2235f2c45f1SShri Abhyankar 
2245f2c45f1SShri Abhyankar static PetscErrorCode DMCreateGlobalVector_Network(DM dm,Vec *vec)
2255f2c45f1SShri Abhyankar {
2265f2c45f1SShri Abhyankar   PetscErrorCode ierr;
2275f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
2285f2c45f1SShri Abhyankar 
2295f2c45f1SShri Abhyankar   PetscFunctionBegin;
2305f2c45f1SShri Abhyankar   ierr = DMCreateGlobalVector(network->plex,vec);CHKERRQ(ierr);
2314dc485aaSHong Zhang   ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Network);CHKERRQ(ierr);
2325f2c45f1SShri Abhyankar   ierr = VecSetDM(*vec,dm);CHKERRQ(ierr);
2335f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2345f2c45f1SShri Abhyankar }
2355f2c45f1SShri Abhyankar 
2365f2c45f1SShri Abhyankar static PetscErrorCode DMCreateLocalVector_Network(DM dm,Vec *vec)
2375f2c45f1SShri Abhyankar {
2385f2c45f1SShri Abhyankar   PetscErrorCode ierr;
2395f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
2405f2c45f1SShri Abhyankar 
2415f2c45f1SShri Abhyankar   PetscFunctionBegin;
2425f2c45f1SShri Abhyankar   ierr = DMCreateLocalVector(network->plex,vec);CHKERRQ(ierr);
2435f2c45f1SShri Abhyankar   ierr = VecSetDM(*vec,dm);CHKERRQ(ierr);
2445f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2455f2c45f1SShri Abhyankar }
2465f2c45f1SShri Abhyankar 
2475f2c45f1SShri Abhyankar PetscErrorCode DMInitialize_Network(DM dm)
2485f2c45f1SShri Abhyankar {
2495f2c45f1SShri Abhyankar 
2505f2c45f1SShri Abhyankar   PetscFunctionBegin;
2515f2c45f1SShri Abhyankar 
2525f2c45f1SShri Abhyankar   dm->ops->view                            = NULL;
2535f2c45f1SShri Abhyankar   dm->ops->setfromoptions                  = DMSetFromOptions_Network;
2548415c774SShri Abhyankar   dm->ops->clone                           = DMClone_Network;
2555f2c45f1SShri Abhyankar   dm->ops->setup                           = DMSetUp_Network;
2565f2c45f1SShri Abhyankar   dm->ops->createglobalvector              = DMCreateGlobalVector_Network;
2575f2c45f1SShri Abhyankar   dm->ops->createlocalvector               = DMCreateLocalVector_Network;
2585f2c45f1SShri Abhyankar   dm->ops->getlocaltoglobalmapping         = NULL;
2595f2c45f1SShri Abhyankar   dm->ops->createfieldis                   = NULL;
2605f2c45f1SShri Abhyankar   dm->ops->createcoordinatedm              = NULL;
2615f2c45f1SShri Abhyankar   dm->ops->getcoloring                     = 0;
2625f2c45f1SShri Abhyankar   dm->ops->creatematrix                    = DMCreateMatrix_Network;
2635f2c45f1SShri Abhyankar   dm->ops->createinterpolation             = 0;
2645f2c45f1SShri Abhyankar   dm->ops->getaggregates                   = 0;
2655f2c45f1SShri Abhyankar   dm->ops->getinjection                    = 0;
2665f2c45f1SShri Abhyankar   dm->ops->refine                          = 0;
2675f2c45f1SShri Abhyankar   dm->ops->coarsen                         = 0;
2685f2c45f1SShri Abhyankar   dm->ops->refinehierarchy                 = 0;
2695f2c45f1SShri Abhyankar   dm->ops->coarsenhierarchy                = 0;
2705f2c45f1SShri Abhyankar   dm->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Network;
2715f2c45f1SShri Abhyankar   dm->ops->globaltolocalend                = DMGlobalToLocalEnd_Network;
2725f2c45f1SShri Abhyankar   dm->ops->localtoglobalbegin              = DMLocalToGlobalBegin_Network;
2735f2c45f1SShri Abhyankar   dm->ops->localtoglobalend                = DMLocalToGlobalEnd_Network;
2745f2c45f1SShri Abhyankar   dm->ops->destroy                         = DMDestroy_Network;
2755f2c45f1SShri Abhyankar   dm->ops->createsubdm                     = NULL;
2765f2c45f1SShri Abhyankar   dm->ops->locatepoints                    = NULL;
2775f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2785f2c45f1SShri Abhyankar }
2795f2c45f1SShri Abhyankar 
2808415c774SShri Abhyankar PetscErrorCode DMClone_Network(DM dm, DM *newdm)
2818415c774SShri Abhyankar {
2828415c774SShri Abhyankar   DM_Network     *network = (DM_Network *) dm->data;
2838415c774SShri Abhyankar   PetscErrorCode ierr;
2848415c774SShri Abhyankar 
2858415c774SShri Abhyankar   PetscFunctionBegin;
2868415c774SShri Abhyankar   network->refct++;
2878415c774SShri Abhyankar   (*newdm)->data = network;
2888415c774SShri Abhyankar   ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMNETWORK);CHKERRQ(ierr);
2898415c774SShri Abhyankar   ierr = DMInitialize_Network(*newdm);CHKERRQ(ierr);
2908415c774SShri Abhyankar   PetscFunctionReturn(0);
2918415c774SShri Abhyankar }
2928415c774SShri Abhyankar 
2935f2c45f1SShri Abhyankar /*MC
2945f2c45f1SShri Abhyankar   DMNETWORK = "network" - A DM object that encapsulates an unstructured network. The implementation is based on the DM object
2955f2c45f1SShri Abhyankar                           DMPlex that manages unstructured grids. Distributed networks use a non-overlapping partitioning of
2965f2c45f1SShri Abhyankar                           the edges. In the local representation, Vecs contain all unknowns in the interior and shared boundary.
2975f2c45f1SShri Abhyankar                           This is specified by a PetscSection object. Ownership in the global representation is determined by
2985f2c45f1SShri Abhyankar                           ownership of the underlying DMPlex points. This is specified by another PetscSection object.
2995f2c45f1SShri Abhyankar 
3005f2c45f1SShri Abhyankar   Level: intermediate
3015f2c45f1SShri Abhyankar 
3025f2c45f1SShri Abhyankar .seealso: DMType, DMNetworkCreate(), DMCreate(), DMSetType()
3035f2c45f1SShri Abhyankar M*/
3045f2c45f1SShri Abhyankar 
3055f2c45f1SShri Abhyankar PETSC_EXTERN PetscErrorCode DMCreate_Network(DM dm)
3065f2c45f1SShri Abhyankar {
3075f2c45f1SShri Abhyankar   DM_Network     *network;
3085f2c45f1SShri Abhyankar   PetscErrorCode ierr;
3095f2c45f1SShri Abhyankar 
3105f2c45f1SShri Abhyankar   PetscFunctionBegin;
3115f2c45f1SShri Abhyankar   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3125f2c45f1SShri Abhyankar   ierr     = PetscNewLog(dm,&network);CHKERRQ(ierr);
3135f2c45f1SShri Abhyankar   dm->data = network;
3145f2c45f1SShri Abhyankar 
3155f2c45f1SShri Abhyankar   network->refct          = 1;
3165f2c45f1SShri Abhyankar   network->NNodes         = -1;
3175f2c45f1SShri Abhyankar   network->NEdges         = -1;
3185f2c45f1SShri Abhyankar   network->nNodes         = -1;
3195f2c45f1SShri Abhyankar   network->nEdges         = -1;
3205f2c45f1SShri Abhyankar 
32113c2a604SAdrian Maldonado 
3225f2c45f1SShri Abhyankar   ierr = DMInitialize_Network(dm);CHKERRQ(ierr);
3235f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3245f2c45f1SShri Abhyankar }
3255f2c45f1SShri Abhyankar 
3265f2c45f1SShri Abhyankar /*@
3275f2c45f1SShri Abhyankar   DMNetworkCreate - Creates a DMNetwork object, which encapsulates an unstructured network.
3285f2c45f1SShri Abhyankar 
3295f2c45f1SShri Abhyankar   Collective on MPI_Comm
3305f2c45f1SShri Abhyankar 
3315f2c45f1SShri Abhyankar   Input Parameter:
3325f2c45f1SShri Abhyankar . comm - The communicator for the DMNetwork object
3335f2c45f1SShri Abhyankar 
3345f2c45f1SShri Abhyankar   Output Parameter:
3355f2c45f1SShri Abhyankar . network  - The DMNetwork object
3365f2c45f1SShri Abhyankar 
3375f2c45f1SShri Abhyankar   Level: beginner
3385f2c45f1SShri Abhyankar 
3395f2c45f1SShri Abhyankar .keywords: DMNetwork, create
3405f2c45f1SShri Abhyankar @*/
3415f2c45f1SShri Abhyankar PetscErrorCode DMNetworkCreate(MPI_Comm comm, DM *network)
3425f2c45f1SShri Abhyankar {
3435f2c45f1SShri Abhyankar   PetscErrorCode ierr;
3445f2c45f1SShri Abhyankar 
3455f2c45f1SShri Abhyankar   PetscFunctionBegin;
3465f2c45f1SShri Abhyankar   PetscValidPointer(network,2);
3475f2c45f1SShri Abhyankar   ierr = DMCreate(comm, network);CHKERRQ(ierr);
3485f2c45f1SShri Abhyankar   ierr = DMSetType(*network, DMNETWORK);CHKERRQ(ierr);
3495f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3505f2c45f1SShri Abhyankar }
351