xref: /petsc/src/dm/impls/network/networkcreate.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558) !
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 
59371c9d4SSatish Balay PetscErrorCode DMSetFromOptions_Network(DM dm, PetscOptionItems *PetscOptionsObject) {
65f2c45f1SShri Abhyankar   PetscFunctionBegin;
7d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "DMNetwork Options");
8d0609cedSBarry Smith   PetscOptionsHeadEnd();
95f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
105f2c45f1SShri Abhyankar }
115f2c45f1SShri Abhyankar 
125f2c45f1SShri Abhyankar /* External function declarations here */
135f2c45f1SShri Abhyankar extern PetscErrorCode DMCreateMatrix_Network(DM, Mat *);
145f2c45f1SShri Abhyankar extern PetscErrorCode DMDestroy_Network(DM);
155f2c45f1SShri Abhyankar extern PetscErrorCode DMView_Network(DM, PetscViewer);
165f2c45f1SShri Abhyankar extern PetscErrorCode DMGlobalToLocalBegin_Network(DM, Vec, InsertMode, Vec);
175f2c45f1SShri Abhyankar extern PetscErrorCode DMGlobalToLocalEnd_Network(DM, Vec, InsertMode, Vec);
185f2c45f1SShri Abhyankar extern PetscErrorCode DMLocalToGlobalBegin_Network(DM, Vec, InsertMode, Vec);
195f2c45f1SShri Abhyankar extern PetscErrorCode DMLocalToGlobalEnd_Network(DM, Vec, InsertMode, Vec);
205f2c45f1SShri Abhyankar extern PetscErrorCode DMSetUp_Network(DM);
218415c774SShri Abhyankar extern PetscErrorCode DMClone_Network(DM, DM *);
225f2c45f1SShri Abhyankar 
239371c9d4SSatish Balay static PetscErrorCode VecArrayPrint_private(PetscViewer viewer, PetscInt n, const PetscScalar *xv) {
24a81b7fe4SHong Zhang   PetscInt i;
25a81b7fe4SHong Zhang 
26a81b7fe4SHong Zhang   PetscFunctionBegin;
27a81b7fe4SHong Zhang   for (i = 0; i < n; i++) {
28a81b7fe4SHong Zhang #if defined(PETSC_USE_COMPLEX)
29a81b7fe4SHong Zhang     if (PetscImaginaryPart(xv[i]) > 0.0) {
309566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    %g + %g i\n", (double)PetscRealPart(xv[i]), (double)PetscImaginaryPart(xv[i])));
31a81b7fe4SHong Zhang     } else if (PetscImaginaryPart(xv[i]) < 0.0) {
329566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    %g - %g i\n", (double)PetscRealPart(xv[i]), -(double)PetscImaginaryPart(xv[i])));
331baa6e33SBarry Smith     } else PetscCall(PetscViewerASCIIPrintf(viewer, "    %g\n", (double)PetscRealPart(xv[i])));
34a81b7fe4SHong Zhang #else
359566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    %g\n", (double)xv[i]));
36a81b7fe4SHong Zhang #endif
37a81b7fe4SHong Zhang   }
38a81b7fe4SHong Zhang   PetscFunctionReturn(0);
39a81b7fe4SHong Zhang }
40a81b7fe4SHong Zhang 
419371c9d4SSatish Balay static PetscErrorCode VecView_Network_Seq(DM networkdm, Vec X, PetscViewer viewer) {
427b6afd5bSHong Zhang   PetscInt           e, v, Start, End, offset, nvar, id;
434dc485aaSHong Zhang   const PetscScalar *xv;
444dc485aaSHong Zhang 
454dc485aaSHong Zhang   PetscFunctionBegin;
469566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &xv));
474dc485aaSHong Zhang 
484dc485aaSHong Zhang   /* iterate over edges */
499566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetEdgeRange(networkdm, &Start, &End));
504dc485aaSHong Zhang   for (e = Start; e < End; e++) {
519566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(networkdm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
524dc485aaSHong Zhang     if (!nvar) continue;
53a81b7fe4SHong Zhang 
549566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(networkdm, e, ALL_COMPONENTS, &offset));
559566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalEdgeIndex(networkdm, e, &id));
567b6afd5bSHong Zhang 
579566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Edge %" PetscInt_FMT ":\n", id));
589566063dSJacob Faibussowitsch     PetscCall(VecArrayPrint_private(viewer, nvar, xv + offset));
594dc485aaSHong Zhang   }
604dc485aaSHong Zhang 
614dc485aaSHong Zhang   /* iterate over vertices */
629566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetVertexRange(networkdm, &Start, &End));
634dc485aaSHong Zhang   for (v = Start; v < End; v++) {
649566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(networkdm, v, ALL_COMPONENTS, NULL, NULL, &nvar));
654dc485aaSHong Zhang     if (!nvar) continue;
66a81b7fe4SHong Zhang 
679566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, ALL_COMPONENTS, &offset));
689566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalVertexIndex(networkdm, v, &id));
6962aa33baSHong Zhang 
709566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Vertex %" PetscInt_FMT ":\n", id));
719566063dSJacob Faibussowitsch     PetscCall(VecArrayPrint_private(viewer, nvar, xv + offset));
724dc485aaSHong Zhang   }
739566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &xv));
754dc485aaSHong Zhang   PetscFunctionReturn(0);
764dc485aaSHong Zhang }
774dc485aaSHong Zhang 
789371c9d4SSatish Balay static PetscErrorCode VecView_Network_MPI(DM networkdm, Vec X, PetscViewer viewer) {
79e8c2e809SHong Zhang   PetscInt           i, e, v, eStart, eEnd, vStart, vEnd, offset, nvar, len_loc, len, k;
804062a5e5SHong Zhang   const PetscScalar *xv;
814062a5e5SHong Zhang   MPI_Comm           comm;
824062a5e5SHong Zhang   PetscMPIInt        size, rank, tag = ((PetscObject)viewer)->tag;
834062a5e5SHong Zhang   Vec                localX;
844062a5e5SHong Zhang   PetscBool          ghostvtex;
854062a5e5SHong Zhang   PetscScalar       *values;
867b6afd5bSHong Zhang   PetscInt           j, ne, nv, id;
874062a5e5SHong Zhang   MPI_Status         status;
884062a5e5SHong Zhang 
894062a5e5SHong Zhang   PetscFunctionBegin;
909566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)networkdm, &comm));
919566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
929566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
934062a5e5SHong Zhang 
949566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(networkdm, &localX));
959566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX));
969566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX));
979566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(localX, &xv));
984062a5e5SHong Zhang 
999566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(localX, &len_loc));
1004062a5e5SHong Zhang 
1019566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetEdgeRange(networkdm, &eStart, &eEnd));
1029566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetVertexRange(networkdm, &vStart, &vEnd));
103e8c2e809SHong Zhang   len_loc += 2 * (1 + eEnd - eStart + vEnd - vStart);
1044062a5e5SHong Zhang 
105e85e6aecSHong Zhang   /* values = [nedges, nvertices; id, nvar, xedge; ...; id, nvars, xvertex;...], to be sent to proc[0] */
1069566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&len_loc, &len, 1, MPIU_INT, MPI_MAX, comm));
1079566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(len, &values));
1084062a5e5SHong Zhang 
109*48a46eb9SPierre Jolivet   if (rank == 0) PetscCall(PetscViewerASCIIPrintf(viewer, "Process [%d]\n", rank));
1104062a5e5SHong Zhang 
1114062a5e5SHong Zhang   /* iterate over edges */
1124062a5e5SHong Zhang   k = 2;
113e8c2e809SHong Zhang   for (e = eStart; e < eEnd; e++) {
1149566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(networkdm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
1154062a5e5SHong Zhang     if (!nvar) continue;
1164062a5e5SHong Zhang 
1179566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(networkdm, e, ALL_COMPONENTS, &offset));
1189566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalEdgeIndex(networkdm, e, &id));
11962aa33baSHong Zhang 
120dd400576SPatrick Sanan     if (rank == 0) { /* print its own entries */
1219566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Edge %" PetscInt_FMT ":\n", id));
1229566063dSJacob Faibussowitsch       PetscCall(VecArrayPrint_private(viewer, nvar, xv + offset));
123a81b7fe4SHong Zhang     } else {
1247b6afd5bSHong Zhang       values[0] += 1; /* number of edges */
1257b6afd5bSHong Zhang       values[k++] = id;
1267b6afd5bSHong Zhang       values[k++] = nvar;
1277b6afd5bSHong Zhang       for (i = offset; i < offset + nvar; i++) values[k++] = xv[i];
1284062a5e5SHong Zhang     }
1294062a5e5SHong Zhang   }
1304062a5e5SHong Zhang 
1314062a5e5SHong Zhang   /* iterate over vertices */
1324062a5e5SHong Zhang   for (v = vStart; v < vEnd; v++) {
1339566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(networkdm, v, &ghostvtex));
1344062a5e5SHong Zhang     if (ghostvtex) continue;
1359566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(networkdm, v, ALL_COMPONENTS, NULL, NULL, &nvar));
1364062a5e5SHong Zhang     if (!nvar) continue;
1374062a5e5SHong Zhang 
1389566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, ALL_COMPONENTS, &offset));
1399566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalVertexIndex(networkdm, v, &id));
1404062a5e5SHong Zhang 
141dd400576SPatrick Sanan     if (rank == 0) {
1429566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Vertex %" PetscInt_FMT ":\n", id));
1439566063dSJacob Faibussowitsch       PetscCall(VecArrayPrint_private(viewer, nvar, xv + offset));
1444062a5e5SHong Zhang     } else {
1457b6afd5bSHong Zhang       values[1] += 1; /* number of vertices */
1467b6afd5bSHong Zhang       values[k++] = id;
1477b6afd5bSHong Zhang       values[k++] = nvar;
1487b6afd5bSHong Zhang       for (i = offset; i < offset + nvar; i++) values[k++] = xv[i];
1494062a5e5SHong Zhang     }
1504062a5e5SHong Zhang   }
1514062a5e5SHong Zhang 
152dd400576SPatrick Sanan   if (rank == 0) {
1534062a5e5SHong Zhang     /* proc[0] receives and prints messages */
1544062a5e5SHong Zhang     for (j = 1; j < size; j++) {
15563a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Process [%" PetscInt_FMT "]\n", j));
1564062a5e5SHong Zhang 
1579566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(values, (PetscMPIInt)len, MPIU_SCALAR, j, tag, comm, &status));
1584062a5e5SHong Zhang 
159dde233f4SSatish Balay       ne = (PetscInt)PetscAbsScalar(values[0]);
160dde233f4SSatish Balay       nv = (PetscInt)PetscAbsScalar(values[1]);
1614062a5e5SHong Zhang 
1624062a5e5SHong Zhang       /* print received edges */
1634062a5e5SHong Zhang       k = 2;
1644062a5e5SHong Zhang       for (i = 0; i < ne; i++) {
165dde233f4SSatish Balay         id   = (PetscInt)PetscAbsScalar(values[k++]);
166dde233f4SSatish Balay         nvar = (PetscInt)PetscAbsScalar(values[k++]);
1679566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Edge %" PetscInt_FMT ":\n", id));
1689566063dSJacob Faibussowitsch         PetscCall(VecArrayPrint_private(viewer, nvar, values + k));
169a81b7fe4SHong Zhang         k += nvar;
1704062a5e5SHong Zhang       }
1714062a5e5SHong Zhang 
1724062a5e5SHong Zhang       /* print received vertices */
1734062a5e5SHong Zhang       for (i = 0; i < nv; i++) {
174dde233f4SSatish Balay         id   = (PetscInt)PetscAbsScalar(values[k++]);
175dde233f4SSatish Balay         nvar = (PetscInt)PetscAbsScalar(values[k++]);
1769566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Vertex %" PetscInt_FMT ":\n", id));
1779566063dSJacob Faibussowitsch         PetscCall(VecArrayPrint_private(viewer, nvar, values + k));
178a81b7fe4SHong Zhang         k += nvar;
1794062a5e5SHong Zhang       }
1804062a5e5SHong Zhang     }
1814062a5e5SHong Zhang   } else {
1824062a5e5SHong Zhang     /* sends values to proc[0] */
1839566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Send((void *)values, k, MPIU_SCALAR, 0, tag, comm));
1844062a5e5SHong Zhang   }
1854062a5e5SHong Zhang 
1869566063dSJacob Faibussowitsch   PetscCall(PetscFree(values));
1879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(localX, &xv));
1889566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(networkdm, &localX));
1894062a5e5SHong Zhang   PetscFunctionReturn(0);
1904062a5e5SHong Zhang }
1914062a5e5SHong Zhang 
19295af8c53SHong Zhang PETSC_EXTERN PetscErrorCode VecView_MPI(Vec, PetscViewer);
193a6c4b7b7SHong Zhang 
1949371c9d4SSatish Balay PetscErrorCode VecView_Network(Vec v, PetscViewer viewer) {
1954dc485aaSHong Zhang   DM        dm;
1964dc485aaSHong Zhang   PetscBool isseq;
197a6c4b7b7SHong Zhang   PetscBool iascii;
1984dc485aaSHong Zhang 
1994dc485aaSHong Zhang   PetscFunctionBegin;
2009566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
2015c6496baSHong Zhang   PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
2029566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
2039566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)v, VECSEQ, &isseq));
204a6c4b7b7SHong Zhang 
205a6c4b7b7SHong Zhang   /* Use VecView_Network if the viewer is ASCII; use VecView_Seq/MPI for other viewer formats */
206a6c4b7b7SHong Zhang   if (iascii) {
2071baa6e33SBarry Smith     if (isseq) PetscCall(VecView_Network_Seq(dm, v, viewer));
2081baa6e33SBarry Smith     else PetscCall(VecView_Network_MPI(dm, v, viewer));
2094dc485aaSHong Zhang   } else {
2101baa6e33SBarry Smith     if (isseq) PetscCall(VecView_Seq(v, viewer));
2111baa6e33SBarry Smith     else PetscCall(VecView_MPI(v, viewer));
212a6c4b7b7SHong Zhang   }
2134dc485aaSHong Zhang   PetscFunctionReturn(0);
2144dc485aaSHong Zhang }
2155f2c45f1SShri Abhyankar 
2169371c9d4SSatish Balay static PetscErrorCode DMCreateGlobalVector_Network(DM dm, Vec *vec) {
2175f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network *)dm->data;
2185f2c45f1SShri Abhyankar 
2195f2c45f1SShri Abhyankar   PetscFunctionBegin;
2209566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(network->plex, vec));
2219566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Network));
2229566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*vec, dm));
2235f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2245f2c45f1SShri Abhyankar }
2255f2c45f1SShri Abhyankar 
2269371c9d4SSatish Balay static PetscErrorCode DMCreateLocalVector_Network(DM dm, Vec *vec) {
2275f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network *)dm->data;
2285f2c45f1SShri Abhyankar 
2295f2c45f1SShri Abhyankar   PetscFunctionBegin;
2309566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(network->plex, vec));
2319566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*vec, dm));
2325f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2335f2c45f1SShri Abhyankar }
2345f2c45f1SShri Abhyankar 
2359371c9d4SSatish Balay PetscErrorCode DMNetworkInitializeToDefault_NonShared(DM dm) {
236daad07d3SAidan Hamilton   DM_Network *network = (DM_Network *)dm->data;
237daad07d3SAidan Hamilton 
238daad07d3SAidan Hamilton   PetscFunctionBegin;
239daad07d3SAidan Hamilton   network->Je                 = NULL;
240daad07d3SAidan Hamilton   network->Jv                 = NULL;
241daad07d3SAidan Hamilton   network->Jvptr              = NULL;
242daad07d3SAidan Hamilton   network->userEdgeJacobian   = PETSC_FALSE;
243daad07d3SAidan Hamilton   network->userVertexJacobian = PETSC_FALSE;
244daad07d3SAidan Hamilton 
245daad07d3SAidan Hamilton   network->vertex.DofSection       = NULL;
246daad07d3SAidan Hamilton   network->vertex.GlobalDofSection = NULL;
247daad07d3SAidan Hamilton   network->vertex.mapping          = NULL;
248daad07d3SAidan Hamilton   network->vertex.sf               = NULL;
249daad07d3SAidan Hamilton 
250daad07d3SAidan Hamilton   network->edge.DofSection       = NULL;
251daad07d3SAidan Hamilton   network->edge.GlobalDofSection = NULL;
252daad07d3SAidan Hamilton   network->edge.mapping          = NULL;
253daad07d3SAidan Hamilton   network->edge.sf               = NULL;
254daad07d3SAidan Hamilton 
255daad07d3SAidan Hamilton   network->DataSection      = NULL;
256daad07d3SAidan Hamilton   network->DofSection       = NULL;
257daad07d3SAidan Hamilton   network->GlobalDofSection = NULL;
258daad07d3SAidan Hamilton   network->componentsetup   = PETSC_FALSE;
259daad07d3SAidan Hamilton 
260daad07d3SAidan Hamilton   network->plex = NULL;
261daad07d3SAidan Hamilton 
262daad07d3SAidan Hamilton   network->component  = NULL;
263daad07d3SAidan Hamilton   network->ncomponent = 0;
264daad07d3SAidan Hamilton 
265daad07d3SAidan Hamilton   network->header             = NULL;
266daad07d3SAidan Hamilton   network->cvalue             = NULL;
267daad07d3SAidan Hamilton   network->componentdataarray = NULL;
268daad07d3SAidan Hamilton 
269daad07d3SAidan Hamilton   network->max_comps_registered = DMNETWORK_MAX_COMP_REGISTERED_DEFAULT; /* return to default */
270daad07d3SAidan Hamilton   PetscFunctionReturn(0);
271daad07d3SAidan Hamilton }
272daad07d3SAidan Hamilton /* Default values for the parameters in DMNetwork */
2739371c9d4SSatish Balay PetscErrorCode DMNetworkInitializeToDefault(DM dm) {
274daad07d3SAidan Hamilton   DM_Network          *network     = (DM_Network *)dm->data;
275daad07d3SAidan Hamilton   DMNetworkCloneShared cloneshared = network->cloneshared;
276daad07d3SAidan Hamilton 
277daad07d3SAidan Hamilton   PetscFunctionBegin;
278daad07d3SAidan Hamilton   PetscCall(DMNetworkInitializeToDefault_NonShared(dm));
279daad07d3SAidan Hamilton   /* Default values for shared data */
280daad07d3SAidan Hamilton   cloneshared->refct            = 1;
281daad07d3SAidan Hamilton   cloneshared->NVertices        = 0;
282daad07d3SAidan Hamilton   cloneshared->NEdges           = 0;
283daad07d3SAidan Hamilton   cloneshared->nVertices        = 0;
284daad07d3SAidan Hamilton   cloneshared->nEdges           = 0;
285daad07d3SAidan Hamilton   cloneshared->nsubnet          = 0;
286daad07d3SAidan Hamilton   cloneshared->pStart           = -1;
287daad07d3SAidan Hamilton   cloneshared->pEnd             = -1;
288daad07d3SAidan Hamilton   cloneshared->vStart           = -1;
289daad07d3SAidan Hamilton   cloneshared->vEnd             = -1;
290daad07d3SAidan Hamilton   cloneshared->eStart           = -1;
291daad07d3SAidan Hamilton   cloneshared->eEnd             = -1;
292daad07d3SAidan Hamilton   cloneshared->vltog            = NULL;
293daad07d3SAidan Hamilton   cloneshared->distributecalled = PETSC_FALSE;
294daad07d3SAidan Hamilton 
295daad07d3SAidan Hamilton   cloneshared->subnet     = NULL;
296daad07d3SAidan Hamilton   cloneshared->subnetvtx  = NULL;
297daad07d3SAidan Hamilton   cloneshared->subnetedge = NULL;
298daad07d3SAidan Hamilton   cloneshared->svtx       = NULL;
299daad07d3SAidan Hamilton   cloneshared->nsvtx      = 0;
300daad07d3SAidan Hamilton   cloneshared->Nsvtx      = 0;
301daad07d3SAidan Hamilton   cloneshared->svertices  = NULL;
302daad07d3SAidan Hamilton   cloneshared->sedgelist  = NULL;
303daad07d3SAidan Hamilton   cloneshared->svtable    = NULL;
304daad07d3SAidan Hamilton   PetscFunctionReturn(0);
305daad07d3SAidan Hamilton }
306daad07d3SAidan Hamilton 
3079371c9d4SSatish Balay PetscErrorCode DMInitialize_Network(DM dm) {
3085f2c45f1SShri Abhyankar   PetscFunctionBegin;
3099566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(dm, 1));
310caf410d2SHong Zhang   dm->ops->view                    = DMView_Network;
3115f2c45f1SShri Abhyankar   dm->ops->setfromoptions          = DMSetFromOptions_Network;
3128415c774SShri Abhyankar   dm->ops->clone                   = DMClone_Network;
3135f2c45f1SShri Abhyankar   dm->ops->setup                   = DMSetUp_Network;
3145f2c45f1SShri Abhyankar   dm->ops->createglobalvector      = DMCreateGlobalVector_Network;
3155f2c45f1SShri Abhyankar   dm->ops->createlocalvector       = DMCreateLocalVector_Network;
3165f2c45f1SShri Abhyankar   dm->ops->getlocaltoglobalmapping = NULL;
3175f2c45f1SShri Abhyankar   dm->ops->createfieldis           = NULL;
3185f2c45f1SShri Abhyankar   dm->ops->createcoordinatedm      = NULL;
319ea78f98cSLisandro Dalcin   dm->ops->getcoloring             = NULL;
3205f2c45f1SShri Abhyankar   dm->ops->creatematrix            = DMCreateMatrix_Network;
321ea78f98cSLisandro Dalcin   dm->ops->createinterpolation     = NULL;
322ea78f98cSLisandro Dalcin   dm->ops->createinjection         = NULL;
323ea78f98cSLisandro Dalcin   dm->ops->refine                  = NULL;
324ea78f98cSLisandro Dalcin   dm->ops->coarsen                 = NULL;
325ea78f98cSLisandro Dalcin   dm->ops->refinehierarchy         = NULL;
326ea78f98cSLisandro Dalcin   dm->ops->coarsenhierarchy        = NULL;
3275f2c45f1SShri Abhyankar   dm->ops->globaltolocalbegin      = DMGlobalToLocalBegin_Network;
3285f2c45f1SShri Abhyankar   dm->ops->globaltolocalend        = DMGlobalToLocalEnd_Network;
3295f2c45f1SShri Abhyankar   dm->ops->localtoglobalbegin      = DMLocalToGlobalBegin_Network;
3305f2c45f1SShri Abhyankar   dm->ops->localtoglobalend        = DMLocalToGlobalEnd_Network;
3315f2c45f1SShri Abhyankar   dm->ops->destroy                 = DMDestroy_Network;
3325f2c45f1SShri Abhyankar   dm->ops->createsubdm             = NULL;
3335f2c45f1SShri Abhyankar   dm->ops->locatepoints            = NULL;
3345f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3355f2c45f1SShri Abhyankar }
336daad07d3SAidan Hamilton /*
337daad07d3SAidan Hamilton   copies over the subnetid and index portions of the DMNetworkComponentHeader from orignal dm to the newdm
338daad07d3SAidan Hamilton */
3399371c9d4SSatish Balay static PetscErrorCode DMNetworkCopyHeaderTopological(DM dm, DM newdm) {
340daad07d3SAidan Hamilton   DM_Network *network = (DM_Network *)dm->data, *newnetwork = (DM_Network *)newdm->data;
341daad07d3SAidan Hamilton   PetscInt    p, i, np, index, subnetid;
342daad07d3SAidan Hamilton 
343daad07d3SAidan Hamilton   PetscFunctionBegin;
344daad07d3SAidan Hamilton   np = network->cloneshared->pEnd - network->cloneshared->pStart;
345daad07d3SAidan Hamilton   PetscCall(PetscCalloc2(np, &newnetwork->header, np, &newnetwork->cvalue));
346daad07d3SAidan Hamilton   for (i = 0; i < np; i++) {
347daad07d3SAidan Hamilton     p = i + network->cloneshared->pStart;
348daad07d3SAidan Hamilton     PetscCall(DMNetworkGetSubnetID(dm, p, &subnetid));
349daad07d3SAidan Hamilton     PetscCall(DMNetworkGetIndex(dm, p, &index));
350daad07d3SAidan Hamilton     newnetwork->header[i].index        = index;
351daad07d3SAidan Hamilton     newnetwork->header[i].subnetid     = subnetid;
352daad07d3SAidan Hamilton     newnetwork->header[i].size         = NULL;
353daad07d3SAidan Hamilton     newnetwork->header[i].key          = NULL;
354daad07d3SAidan Hamilton     newnetwork->header[i].offset       = NULL;
355daad07d3SAidan Hamilton     newnetwork->header[i].nvar         = NULL;
356daad07d3SAidan Hamilton     newnetwork->header[i].offsetvarrel = NULL;
357daad07d3SAidan Hamilton     newnetwork->header[i].ndata        = 0;
358daad07d3SAidan Hamilton     newnetwork->header[i].maxcomps     = DMNETWORK_MAX_COMP_AT_POINT_DEFAULT;
359daad07d3SAidan Hamilton     newnetwork->header[i].hsize        = sizeof(struct _p_DMNetworkComponentHeader) / sizeof(sizeof(DMNetworkComponentGenericDataType));
360daad07d3SAidan Hamilton   }
361daad07d3SAidan Hamilton   PetscFunctionReturn(0);
362daad07d3SAidan Hamilton }
3635f2c45f1SShri Abhyankar 
3649371c9d4SSatish Balay PetscErrorCode DMClone_Network(DM dm, DM *newdm) {
365daad07d3SAidan Hamilton   DM_Network *network = (DM_Network *)dm->data, *newnetwork = NULL;
3668415c774SShri Abhyankar 
3678415c774SShri Abhyankar   PetscFunctionBegin;
368daad07d3SAidan Hamilton   network->cloneshared->refct++;
369daad07d3SAidan Hamilton   PetscCall(PetscNewLog(*newdm, &newnetwork));
370daad07d3SAidan Hamilton   (*newdm)->data = newnetwork;
371daad07d3SAidan Hamilton   PetscCall(DMNetworkInitializeToDefault_NonShared(*newdm));
372daad07d3SAidan Hamilton   newnetwork->cloneshared = network->cloneshared; /* Share all data that can be cloneshared */
373daad07d3SAidan Hamilton   PetscCall(DMClone(network->plex, &newnetwork->plex));
374daad07d3SAidan Hamilton   PetscCall(DMNetworkCopyHeaderTopological(dm, *newdm));
375daad07d3SAidan Hamilton   PetscCall(DMNetworkInitializeNonTopological(*newdm)); /* initialize all non-topological data to the state after DMNetworkLayoutSetUp as been called */
3769566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, DMNETWORK));
3779566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Network(*newdm));
3788415c774SShri Abhyankar   PetscFunctionReturn(0);
3798415c774SShri Abhyankar }
3808415c774SShri Abhyankar 
3815f2c45f1SShri Abhyankar /*MC
3825f2c45f1SShri Abhyankar   DMNETWORK = "network" - A DM object that encapsulates an unstructured network. The implementation is based on the DM object
3835f2c45f1SShri Abhyankar                           DMPlex that manages unstructured grids. Distributed networks use a non-overlapping partitioning of
3845f2c45f1SShri Abhyankar                           the edges. In the local representation, Vecs contain all unknowns in the interior and shared boundary.
3855f2c45f1SShri Abhyankar                           This is specified by a PetscSection object. Ownership in the global representation is determined by
3865f2c45f1SShri Abhyankar                           ownership of the underlying DMPlex points. This is specified by another PetscSection object.
3875f2c45f1SShri Abhyankar 
3885f2c45f1SShri Abhyankar   Level: intermediate
3895f2c45f1SShri Abhyankar 
390db781477SPatrick Sanan .seealso: `DMType`, `DMNetworkCreate()`, `DMCreate()`, `DMSetType()`
3915f2c45f1SShri Abhyankar M*/
3925f2c45f1SShri Abhyankar 
3939371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMCreate_Network(DM dm) {
3945f2c45f1SShri Abhyankar   DM_Network *network;
3955f2c45f1SShri Abhyankar 
3965f2c45f1SShri Abhyankar   PetscFunctionBegin;
3975f2c45f1SShri Abhyankar   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3989566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(dm, &network));
399daad07d3SAidan Hamilton   PetscCall(PetscNewLog(dm, &network->cloneshared));
4005f2c45f1SShri Abhyankar   dm->data = network;
4015f2c45f1SShri Abhyankar 
402daad07d3SAidan Hamilton   PetscCall(DMNetworkInitializeToDefault(dm));
4039566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Network(dm));
4045f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4055f2c45f1SShri Abhyankar }
4065f2c45f1SShri Abhyankar 
4075f2c45f1SShri Abhyankar /*@
4085f2c45f1SShri Abhyankar   DMNetworkCreate - Creates a DMNetwork object, which encapsulates an unstructured network.
4095f2c45f1SShri Abhyankar 
410d083f849SBarry Smith   Collective
4115f2c45f1SShri Abhyankar 
4125f2c45f1SShri Abhyankar   Input Parameter:
4135f2c45f1SShri Abhyankar . comm - The communicator for the DMNetwork object
4145f2c45f1SShri Abhyankar 
4155f2c45f1SShri Abhyankar   Output Parameter:
4165f2c45f1SShri Abhyankar . network  - The DMNetwork object
4175f2c45f1SShri Abhyankar 
4185f2c45f1SShri Abhyankar   Level: beginner
4195f2c45f1SShri Abhyankar 
4205f2c45f1SShri Abhyankar @*/
4219371c9d4SSatish Balay PetscErrorCode DMNetworkCreate(MPI_Comm comm, DM *network) {
4225f2c45f1SShri Abhyankar   PetscFunctionBegin;
4235f2c45f1SShri Abhyankar   PetscValidPointer(network, 2);
4249566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, network));
4259566063dSJacob Faibussowitsch   PetscCall(DMSetType(*network, DMNETWORK));
4265f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4275f2c45f1SShri Abhyankar }
428