xref: /petsc/src/dm/impls/network/networkcreate.c (revision d0609ced746bc51b019815ca91d747429db24893)
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);
9*d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"DMNetwork Options");
10*d0609cedSBarry 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])));
36a81b7fe4SHong Zhang     } else {
379566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"    %g\n",(double)PetscRealPart(xv[i])));
38a81b7fe4SHong Zhang     }
39a81b7fe4SHong Zhang #else
409566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"    %g\n",(double)xv[i]));
41a81b7fe4SHong Zhang #endif
42a81b7fe4SHong Zhang   }
43a81b7fe4SHong Zhang   PetscFunctionReturn(0);
44a81b7fe4SHong Zhang }
45a81b7fe4SHong Zhang 
464dc485aaSHong Zhang static PetscErrorCode VecView_Network_Seq(DM networkdm,Vec X,PetscViewer viewer)
474dc485aaSHong Zhang {
487b6afd5bSHong Zhang   PetscInt          e,v,Start,End,offset,nvar,id;
494dc485aaSHong Zhang   const PetscScalar *xv;
504dc485aaSHong Zhang 
514dc485aaSHong Zhang   PetscFunctionBegin;
529566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&xv));
534dc485aaSHong Zhang 
544dc485aaSHong Zhang   /* iterate over edges */
559566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetEdgeRange(networkdm,&Start,&End));
564dc485aaSHong Zhang   for (e=Start; e<End; e++) {
579566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(networkdm,e,ALL_COMPONENTS,NULL,NULL,&nvar));
584dc485aaSHong Zhang     if (!nvar) continue;
59a81b7fe4SHong Zhang 
609566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(networkdm,e,ALL_COMPONENTS,&offset));
619566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalEdgeIndex(networkdm,e,&id));
627b6afd5bSHong Zhang 
639566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  Edge %" PetscInt_FMT ":\n",id));
649566063dSJacob Faibussowitsch     PetscCall(VecArrayPrint_private(viewer,nvar,xv+offset));
654dc485aaSHong Zhang   }
664dc485aaSHong Zhang 
674dc485aaSHong Zhang   /* iterate over vertices */
689566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetVertexRange(networkdm,&Start,&End));
694dc485aaSHong Zhang   for (v=Start; v<End; v++) {
709566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(networkdm,v,ALL_COMPONENTS,NULL,NULL,&nvar));
714dc485aaSHong Zhang     if (!nvar) continue;
72a81b7fe4SHong Zhang 
739566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(networkdm,v,ALL_COMPONENTS,&offset));
749566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalVertexIndex(networkdm,v,&id));
7562aa33baSHong Zhang 
769566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  Vertex %" PetscInt_FMT ":\n",id));
779566063dSJacob Faibussowitsch     PetscCall(VecArrayPrint_private(viewer,nvar,xv+offset));
784dc485aaSHong Zhang   }
799566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&xv));
814dc485aaSHong Zhang   PetscFunctionReturn(0);
824dc485aaSHong Zhang }
834dc485aaSHong Zhang 
844062a5e5SHong Zhang static PetscErrorCode VecView_Network_MPI(DM networkdm,Vec X,PetscViewer viewer)
854062a5e5SHong Zhang {
86e8c2e809SHong Zhang   PetscInt          i,e,v,eStart,eEnd,vStart,vEnd,offset,nvar,len_loc,len,k;
874062a5e5SHong Zhang   const PetscScalar *xv;
884062a5e5SHong Zhang   MPI_Comm          comm;
894062a5e5SHong Zhang   PetscMPIInt       size,rank,tag = ((PetscObject)viewer)->tag;
904062a5e5SHong Zhang   Vec               localX;
914062a5e5SHong Zhang   PetscBool         ghostvtex;
924062a5e5SHong Zhang   PetscScalar       *values;
937b6afd5bSHong Zhang   PetscInt          j,ne,nv,id;
944062a5e5SHong Zhang   MPI_Status        status;
954062a5e5SHong Zhang 
964062a5e5SHong Zhang   PetscFunctionBegin;
979566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)networkdm,&comm));
989566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
1004062a5e5SHong Zhang 
1019566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(networkdm,&localX));
1029566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX));
1039566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX));
1049566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(localX,&xv));
1054062a5e5SHong Zhang 
1069566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(localX,&len_loc));
1074062a5e5SHong Zhang 
1089566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd));
1099566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetVertexRange(networkdm,&vStart,&vEnd));
110e8c2e809SHong Zhang   len_loc += 2*(1 + eEnd-eStart + vEnd-vStart);
1114062a5e5SHong Zhang 
112e85e6aecSHong Zhang   /* values = [nedges, nvertices; id, nvar, xedge; ...; id, nvars, xvertex;...], to be sent to proc[0] */
1139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&len_loc,&len,1,MPIU_INT,MPI_MAX,comm));
1149566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(len,&values));
1154062a5e5SHong Zhang 
116dd400576SPatrick Sanan   if (rank == 0) {
1179566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"Process [%d]\n",rank));
1184062a5e5SHong Zhang   }
1194062a5e5SHong Zhang 
1204062a5e5SHong Zhang   /* iterate over edges */
1214062a5e5SHong Zhang   k = 2;
122e8c2e809SHong Zhang   for (e=eStart; e<eEnd; e++) {
1239566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(networkdm,e,ALL_COMPONENTS,NULL,NULL,&nvar));
1244062a5e5SHong Zhang     if (!nvar) continue;
1254062a5e5SHong Zhang 
1269566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(networkdm,e,ALL_COMPONENTS,&offset));
1279566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalEdgeIndex(networkdm,e,&id));
12862aa33baSHong Zhang 
129dd400576SPatrick Sanan     if (rank == 0) { /* print its own entries */
1309566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  Edge %" PetscInt_FMT ":\n",id));
1319566063dSJacob Faibussowitsch       PetscCall(VecArrayPrint_private(viewer,nvar,xv+offset));
132a81b7fe4SHong Zhang     } else {
1337b6afd5bSHong Zhang       values[0]  += 1; /* number of edges */
1347b6afd5bSHong Zhang       values[k++] = id;
1357b6afd5bSHong Zhang       values[k++] = nvar;
1367b6afd5bSHong Zhang       for (i=offset; i< offset+nvar; i++) values[k++] = xv[i];
1374062a5e5SHong Zhang     }
1384062a5e5SHong Zhang   }
1394062a5e5SHong Zhang 
1404062a5e5SHong Zhang   /* iterate over vertices */
1414062a5e5SHong Zhang   for (v=vStart; v<vEnd; v++) {
1429566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(networkdm,v,&ghostvtex));
1434062a5e5SHong Zhang     if (ghostvtex) continue;
1449566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(networkdm,v,ALL_COMPONENTS,NULL,NULL,&nvar));
1454062a5e5SHong Zhang     if (!nvar) continue;
1464062a5e5SHong Zhang 
1479566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(networkdm,v,ALL_COMPONENTS,&offset));
1489566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalVertexIndex(networkdm,v,&id));
1494062a5e5SHong Zhang 
150dd400576SPatrick Sanan     if (rank == 0) {
1519566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  Vertex %" PetscInt_FMT ":\n",id));
1529566063dSJacob Faibussowitsch       PetscCall(VecArrayPrint_private(viewer,nvar,xv+offset));
1534062a5e5SHong Zhang     } else {
1547b6afd5bSHong Zhang       values[1]  += 1; /* number of vertices */
1557b6afd5bSHong Zhang       values[k++] = id;
1567b6afd5bSHong Zhang       values[k++] = nvar;
1577b6afd5bSHong Zhang       for (i=offset; i< offset+nvar; i++) values[k++] = xv[i];
1584062a5e5SHong Zhang     }
1594062a5e5SHong Zhang   }
1604062a5e5SHong Zhang 
161dd400576SPatrick Sanan   if (rank == 0) {
1624062a5e5SHong Zhang     /* proc[0] receives and prints messages */
1634062a5e5SHong Zhang     for (j=1; j<size; j++) {
1649566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"Process [%d]\n",j));
1654062a5e5SHong Zhang 
1669566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(values,(PetscMPIInt)len,MPIU_SCALAR,j,tag,comm,&status));
1674062a5e5SHong Zhang 
168dde233f4SSatish Balay       ne = (PetscInt)PetscAbsScalar(values[0]);
169dde233f4SSatish Balay       nv = (PetscInt)PetscAbsScalar(values[1]);
1704062a5e5SHong Zhang 
1714062a5e5SHong Zhang       /* print received edges */
1724062a5e5SHong Zhang       k = 2;
1734062a5e5SHong Zhang       for (i=0; i<ne; i++) {
174dde233f4SSatish Balay         id   = (PetscInt)PetscAbsScalar(values[k++]);
175dde233f4SSatish Balay         nvar = (PetscInt)PetscAbsScalar(values[k++]);
1769566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  Edge %" PetscInt_FMT ":\n",id));
1779566063dSJacob Faibussowitsch         PetscCall(VecArrayPrint_private(viewer,nvar,values+k));
178a81b7fe4SHong Zhang         k   += nvar;
1794062a5e5SHong Zhang       }
1804062a5e5SHong Zhang 
1814062a5e5SHong Zhang       /* print received vertices */
1824062a5e5SHong Zhang       for (i=0; i<nv; i++) {
183dde233f4SSatish Balay         id   = (PetscInt)PetscAbsScalar(values[k++]);
184dde233f4SSatish Balay         nvar = (PetscInt)PetscAbsScalar(values[k++]);
1859566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  Vertex %" PetscInt_FMT ":\n",id));
1869566063dSJacob Faibussowitsch         PetscCall(VecArrayPrint_private(viewer,nvar,values+k));
187a81b7fe4SHong Zhang         k   += nvar;
1884062a5e5SHong Zhang       }
1894062a5e5SHong Zhang     }
1904062a5e5SHong Zhang   } else {
1914062a5e5SHong Zhang     /* sends values to proc[0] */
1929566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send((void*)values,k,MPIU_SCALAR,0,tag,comm));
1934062a5e5SHong Zhang   }
1944062a5e5SHong Zhang 
1959566063dSJacob Faibussowitsch   PetscCall(PetscFree(values));
1969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(localX,&xv));
1979566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(networkdm,&localX));
1984062a5e5SHong Zhang   PetscFunctionReturn(0);
1994062a5e5SHong Zhang }
2004062a5e5SHong Zhang 
20195af8c53SHong Zhang PETSC_EXTERN PetscErrorCode VecView_MPI(Vec,PetscViewer);
202a6c4b7b7SHong Zhang 
2034dc485aaSHong Zhang PetscErrorCode VecView_Network(Vec v,PetscViewer viewer)
2044dc485aaSHong Zhang {
2054dc485aaSHong Zhang   DM             dm;
2064dc485aaSHong Zhang   PetscBool      isseq;
207a6c4b7b7SHong Zhang   PetscBool      iascii;
2084dc485aaSHong Zhang 
2094dc485aaSHong Zhang   PetscFunctionBegin;
2109566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v,&dm));
2115c6496baSHong Zhang   PetscCheck(dm,PetscObjectComm((PetscObject)v),PETSC_ERR_ARG_WRONG,"Vector not generated from a DM");
2129566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
2139566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)v,VECSEQ,&isseq));
214a6c4b7b7SHong Zhang 
215a6c4b7b7SHong Zhang   /* Use VecView_Network if the viewer is ASCII; use VecView_Seq/MPI for other viewer formats */
216a6c4b7b7SHong Zhang   if (iascii) {
2174dc485aaSHong Zhang     if (isseq) {
2189566063dSJacob Faibussowitsch       PetscCall(VecView_Network_Seq(dm,v,viewer));
2194dc485aaSHong Zhang     } else {
2209566063dSJacob Faibussowitsch       PetscCall(VecView_Network_MPI(dm,v,viewer));
2214dc485aaSHong Zhang     }
222a6c4b7b7SHong Zhang   } else {
223a6c4b7b7SHong Zhang     if (isseq) {
2249566063dSJacob Faibussowitsch       PetscCall(VecView_Seq(v,viewer));
225a6c4b7b7SHong Zhang     } else {
2269566063dSJacob Faibussowitsch       PetscCall(VecView_MPI(v,viewer));
227a6c4b7b7SHong Zhang     }
228a6c4b7b7SHong Zhang   }
2294dc485aaSHong Zhang   PetscFunctionReturn(0);
2304dc485aaSHong Zhang }
2315f2c45f1SShri Abhyankar 
2325f2c45f1SShri Abhyankar static PetscErrorCode DMCreateGlobalVector_Network(DM dm,Vec *vec)
2335f2c45f1SShri Abhyankar {
2345f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
2355f2c45f1SShri Abhyankar 
2365f2c45f1SShri Abhyankar   PetscFunctionBegin;
2379566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(network->plex,vec));
2389566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Network));
2399566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*vec,dm));
2405f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2415f2c45f1SShri Abhyankar }
2425f2c45f1SShri Abhyankar 
2435f2c45f1SShri Abhyankar static PetscErrorCode DMCreateLocalVector_Network(DM dm,Vec *vec)
2445f2c45f1SShri Abhyankar {
2455f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
2465f2c45f1SShri Abhyankar 
2475f2c45f1SShri Abhyankar   PetscFunctionBegin;
2489566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(network->plex,vec));
2499566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*vec,dm));
2505f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2515f2c45f1SShri Abhyankar }
2525f2c45f1SShri Abhyankar 
2535f2c45f1SShri Abhyankar PetscErrorCode DMInitialize_Network(DM dm)
2545f2c45f1SShri Abhyankar {
2555f2c45f1SShri Abhyankar   PetscFunctionBegin;
2569566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(dm,1));
257caf410d2SHong Zhang   dm->ops->view                            = DMView_Network;
2585f2c45f1SShri Abhyankar   dm->ops->setfromoptions                  = DMSetFromOptions_Network;
2598415c774SShri Abhyankar   dm->ops->clone                           = DMClone_Network;
2605f2c45f1SShri Abhyankar   dm->ops->setup                           = DMSetUp_Network;
2615f2c45f1SShri Abhyankar   dm->ops->createglobalvector              = DMCreateGlobalVector_Network;
2625f2c45f1SShri Abhyankar   dm->ops->createlocalvector               = DMCreateLocalVector_Network;
2635f2c45f1SShri Abhyankar   dm->ops->getlocaltoglobalmapping         = NULL;
2645f2c45f1SShri Abhyankar   dm->ops->createfieldis                   = NULL;
2655f2c45f1SShri Abhyankar   dm->ops->createcoordinatedm              = NULL;
266ea78f98cSLisandro Dalcin   dm->ops->getcoloring                     = NULL;
2675f2c45f1SShri Abhyankar   dm->ops->creatematrix                    = DMCreateMatrix_Network;
268ea78f98cSLisandro Dalcin   dm->ops->createinterpolation             = NULL;
269ea78f98cSLisandro Dalcin   dm->ops->createinjection                 = NULL;
270ea78f98cSLisandro Dalcin   dm->ops->refine                          = NULL;
271ea78f98cSLisandro Dalcin   dm->ops->coarsen                         = NULL;
272ea78f98cSLisandro Dalcin   dm->ops->refinehierarchy                 = NULL;
273ea78f98cSLisandro Dalcin   dm->ops->coarsenhierarchy                = NULL;
2745f2c45f1SShri Abhyankar   dm->ops->globaltolocalbegin              = DMGlobalToLocalBegin_Network;
2755f2c45f1SShri Abhyankar   dm->ops->globaltolocalend                = DMGlobalToLocalEnd_Network;
2765f2c45f1SShri Abhyankar   dm->ops->localtoglobalbegin              = DMLocalToGlobalBegin_Network;
2775f2c45f1SShri Abhyankar   dm->ops->localtoglobalend                = DMLocalToGlobalEnd_Network;
2785f2c45f1SShri Abhyankar   dm->ops->destroy                         = DMDestroy_Network;
2795f2c45f1SShri Abhyankar   dm->ops->createsubdm                     = NULL;
2805f2c45f1SShri Abhyankar   dm->ops->locatepoints                    = NULL;
2815f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2825f2c45f1SShri Abhyankar }
2835f2c45f1SShri Abhyankar 
2848415c774SShri Abhyankar PetscErrorCode DMClone_Network(DM dm, DM *newdm)
2858415c774SShri Abhyankar {
2868415c774SShri Abhyankar   DM_Network     *network = (DM_Network *) dm->data;
2878415c774SShri Abhyankar 
2888415c774SShri Abhyankar   PetscFunctionBegin;
2898415c774SShri Abhyankar   network->refct++;
2908415c774SShri Abhyankar   (*newdm)->data = network;
2919566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject) *newdm, DMNETWORK));
2929566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Network(*newdm));
2938415c774SShri Abhyankar   PetscFunctionReturn(0);
2948415c774SShri Abhyankar }
2958415c774SShri Abhyankar 
2965f2c45f1SShri Abhyankar /*MC
2975f2c45f1SShri Abhyankar   DMNETWORK = "network" - A DM object that encapsulates an unstructured network. The implementation is based on the DM object
2985f2c45f1SShri Abhyankar                           DMPlex that manages unstructured grids. Distributed networks use a non-overlapping partitioning of
2995f2c45f1SShri Abhyankar                           the edges. In the local representation, Vecs contain all unknowns in the interior and shared boundary.
3005f2c45f1SShri Abhyankar                           This is specified by a PetscSection object. Ownership in the global representation is determined by
3015f2c45f1SShri Abhyankar                           ownership of the underlying DMPlex points. This is specified by another PetscSection object.
3025f2c45f1SShri Abhyankar 
3035f2c45f1SShri Abhyankar   Level: intermediate
3045f2c45f1SShri Abhyankar 
3055f2c45f1SShri Abhyankar .seealso: DMType, DMNetworkCreate(), DMCreate(), DMSetType()
3065f2c45f1SShri Abhyankar M*/
3075f2c45f1SShri Abhyankar 
3085f2c45f1SShri Abhyankar PETSC_EXTERN PetscErrorCode DMCreate_Network(DM dm)
3095f2c45f1SShri Abhyankar {
3105f2c45f1SShri Abhyankar   DM_Network     *network;
3115f2c45f1SShri Abhyankar 
3125f2c45f1SShri Abhyankar   PetscFunctionBegin;
3135f2c45f1SShri Abhyankar   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3149566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(dm,&network));
3155f2c45f1SShri Abhyankar   dm->data = network;
3165f2c45f1SShri Abhyankar 
3175f2c45f1SShri Abhyankar   network->refct     = 1;
318e2aaf10cSShri Abhyankar   network->NVertices = 0;
319e2aaf10cSShri Abhyankar   network->NEdges    = 0;
320e2aaf10cSShri Abhyankar   network->nVertices = 0;
321e2aaf10cSShri Abhyankar   network->nEdges    = 0;
322e2aaf10cSShri Abhyankar   network->nsubnet   = 0;
3235f2c45f1SShri Abhyankar 
32454dfd506SHong Zhang   network->max_comps_registered = 20;
32554dfd506SHong Zhang   network->component            = NULL;
32654dfd506SHong Zhang   network->header               = NULL;
32754dfd506SHong Zhang   network->cvalue               = NULL;
32813c2a604SAdrian Maldonado 
3299566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Network(dm));
3305f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3315f2c45f1SShri Abhyankar }
3325f2c45f1SShri Abhyankar 
3335f2c45f1SShri Abhyankar /*@
3345f2c45f1SShri Abhyankar   DMNetworkCreate - Creates a DMNetwork object, which encapsulates an unstructured network.
3355f2c45f1SShri Abhyankar 
336d083f849SBarry Smith   Collective
3375f2c45f1SShri Abhyankar 
3385f2c45f1SShri Abhyankar   Input Parameter:
3395f2c45f1SShri Abhyankar . comm - The communicator for the DMNetwork object
3405f2c45f1SShri Abhyankar 
3415f2c45f1SShri Abhyankar   Output Parameter:
3425f2c45f1SShri Abhyankar . network  - The DMNetwork object
3435f2c45f1SShri Abhyankar 
3445f2c45f1SShri Abhyankar   Level: beginner
3455f2c45f1SShri Abhyankar 
3465f2c45f1SShri Abhyankar @*/
3475f2c45f1SShri Abhyankar PetscErrorCode DMNetworkCreate(MPI_Comm comm, DM *network)
3485f2c45f1SShri Abhyankar {
3495f2c45f1SShri Abhyankar   PetscFunctionBegin;
3505f2c45f1SShri Abhyankar   PetscValidPointer(network,2);
3519566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, network));
3529566063dSJacob Faibussowitsch   PetscCall(DMSetType(*network, DMNETWORK));
3535f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3545f2c45f1SShri Abhyankar }
355