xref: /petsc/src/dm/impls/network/networkcreate.c (revision 9f196a0264fbaf0568fead3a30c861c7ae4cf663)
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 
5ce78bad3SBarry Smith static PetscErrorCode DMSetFromOptions_Network(DM dm, PetscOptionItems PetscOptionsObject)
6d71ae5a4SJacob Faibussowitsch {
75f2c45f1SShri Abhyankar   PetscFunctionBegin;
8d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "DMNetwork Options");
9d0609cedSBarry Smith   PetscOptionsHeadEnd();
103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115f2c45f1SShri Abhyankar }
125f2c45f1SShri Abhyankar 
135f2c45f1SShri Abhyankar /* External function declarations here */
145f2c45f1SShri Abhyankar extern PetscErrorCode DMCreateMatrix_Network(DM, Mat *);
155f2c45f1SShri Abhyankar extern PetscErrorCode DMDestroy_Network(DM);
165f2c45f1SShri Abhyankar extern PetscErrorCode DMView_Network(DM, PetscViewer);
175f2c45f1SShri Abhyankar extern PetscErrorCode DMGlobalToLocalBegin_Network(DM, Vec, InsertMode, Vec);
185f2c45f1SShri Abhyankar extern PetscErrorCode DMGlobalToLocalEnd_Network(DM, Vec, InsertMode, Vec);
195f2c45f1SShri Abhyankar extern PetscErrorCode DMLocalToGlobalBegin_Network(DM, Vec, InsertMode, Vec);
205f2c45f1SShri Abhyankar extern PetscErrorCode DMLocalToGlobalEnd_Network(DM, Vec, InsertMode, Vec);
215f2c45f1SShri Abhyankar extern PetscErrorCode DMSetUp_Network(DM);
228415c774SShri Abhyankar extern PetscErrorCode DMClone_Network(DM, DM *);
235f2c45f1SShri Abhyankar 
24d71ae5a4SJacob Faibussowitsch static PetscErrorCode VecArrayPrint_private(PetscViewer viewer, PetscInt n, const PetscScalar *xv)
25d71ae5a4SJacob Faibussowitsch {
26a81b7fe4SHong Zhang   PetscInt i;
27a81b7fe4SHong Zhang 
28a81b7fe4SHong Zhang   PetscFunctionBegin;
29a81b7fe4SHong Zhang   for (i = 0; i < n; i++) {
30a81b7fe4SHong Zhang #if defined(PETSC_USE_COMPLEX)
31a81b7fe4SHong Zhang     if (PetscImaginaryPart(xv[i]) > 0.0) {
329566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    %g + %g i\n", (double)PetscRealPart(xv[i]), (double)PetscImaginaryPart(xv[i])));
33a81b7fe4SHong Zhang     } else if (PetscImaginaryPart(xv[i]) < 0.0) {
349566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    %g - %g i\n", (double)PetscRealPart(xv[i]), -(double)PetscImaginaryPart(xv[i])));
351baa6e33SBarry Smith     } else PetscCall(PetscViewerASCIIPrintf(viewer, "    %g\n", (double)PetscRealPart(xv[i])));
36a81b7fe4SHong Zhang #else
379566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    %g\n", (double)xv[i]));
38a81b7fe4SHong Zhang #endif
39a81b7fe4SHong Zhang   }
403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41a81b7fe4SHong Zhang }
42a81b7fe4SHong Zhang 
43d71ae5a4SJacob Faibussowitsch static PetscErrorCode VecView_Network_Seq(DM networkdm, Vec X, PetscViewer viewer)
44d71ae5a4SJacob Faibussowitsch {
457b6afd5bSHong Zhang   PetscInt           e, v, Start, End, offset, nvar, id;
464dc485aaSHong Zhang   const PetscScalar *xv;
474dc485aaSHong Zhang 
484dc485aaSHong Zhang   PetscFunctionBegin;
499566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &xv));
504dc485aaSHong Zhang 
514dc485aaSHong Zhang   /* iterate over edges */
529566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetEdgeRange(networkdm, &Start, &End));
534dc485aaSHong Zhang   for (e = Start; e < End; e++) {
549566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(networkdm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
554dc485aaSHong Zhang     if (!nvar) continue;
56a81b7fe4SHong Zhang 
579566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(networkdm, e, ALL_COMPONENTS, &offset));
589566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalEdgeIndex(networkdm, e, &id));
597b6afd5bSHong Zhang 
609566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Edge %" PetscInt_FMT ":\n", id));
619566063dSJacob Faibussowitsch     PetscCall(VecArrayPrint_private(viewer, nvar, xv + offset));
624dc485aaSHong Zhang   }
634dc485aaSHong Zhang 
644dc485aaSHong Zhang   /* iterate over vertices */
659566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetVertexRange(networkdm, &Start, &End));
664dc485aaSHong Zhang   for (v = Start; v < End; v++) {
679566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(networkdm, v, ALL_COMPONENTS, NULL, NULL, &nvar));
684dc485aaSHong Zhang     if (!nvar) continue;
69a81b7fe4SHong Zhang 
709566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, ALL_COMPONENTS, &offset));
719566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalVertexIndex(networkdm, v, &id));
7262aa33baSHong Zhang 
739566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Vertex %" PetscInt_FMT ":\n", id));
749566063dSJacob Faibussowitsch     PetscCall(VecArrayPrint_private(viewer, nvar, xv + offset));
754dc485aaSHong Zhang   }
769566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &xv));
783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
794dc485aaSHong Zhang }
804dc485aaSHong Zhang 
81d71ae5a4SJacob Faibussowitsch static PetscErrorCode VecView_Network_MPI(DM networkdm, Vec X, PetscViewer viewer)
82d71ae5a4SJacob Faibussowitsch {
83835f2295SStefano Zampini   PetscInt           i, e, v, eStart, eEnd, vStart, vEnd, offset, nvar, len_loc, k;
844062a5e5SHong Zhang   const PetscScalar *xv;
854062a5e5SHong Zhang   MPI_Comm           comm;
86835f2295SStefano Zampini   PetscMPIInt        size, rank, tag = ((PetscObject)viewer)->tag, len;
874062a5e5SHong Zhang   Vec                localX;
884062a5e5SHong Zhang   PetscBool          ghostvtex;
894062a5e5SHong Zhang   PetscScalar       *values;
906497c311SBarry Smith   PetscInt           ne, nv, id;
914062a5e5SHong Zhang   MPI_Status         status;
924062a5e5SHong Zhang 
934062a5e5SHong Zhang   PetscFunctionBegin;
949566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)networkdm, &comm));
959566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
969566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
974062a5e5SHong Zhang 
989566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(networkdm, &localX));
999566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX));
1009566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX));
1019566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(localX, &xv));
1024062a5e5SHong Zhang 
1039566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(localX, &len_loc));
1044062a5e5SHong Zhang 
1059566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetEdgeRange(networkdm, &eStart, &eEnd));
1069566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetVertexRange(networkdm, &vStart, &vEnd));
107e8c2e809SHong Zhang   len_loc += 2 * (1 + eEnd - eStart + vEnd - vStart);
1084062a5e5SHong Zhang 
109e85e6aecSHong Zhang   /* values = [nedges, nvertices; id, nvar, xedge; ...; id, nvars, xvertex;...], to be sent to proc[0] */
110835f2295SStefano Zampini   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &len_loc, 1, MPIU_INT, MPI_MAX, comm));
111835f2295SStefano Zampini   PetscCall(PetscMPIIntCast(len_loc, &len));
1129566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(len, &values));
1134062a5e5SHong Zhang 
11448a46eb9SPierre Jolivet   if (rank == 0) PetscCall(PetscViewerASCIIPrintf(viewer, "Process [%d]\n", rank));
1154062a5e5SHong Zhang 
1164062a5e5SHong Zhang   /* iterate over edges */
1174062a5e5SHong Zhang   k = 2;
118e8c2e809SHong Zhang   for (e = eStart; e < eEnd; e++) {
1199566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(networkdm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
1204062a5e5SHong Zhang     if (!nvar) continue;
1214062a5e5SHong Zhang 
1229566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(networkdm, e, ALL_COMPONENTS, &offset));
1239566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalEdgeIndex(networkdm, e, &id));
12462aa33baSHong Zhang 
125dd400576SPatrick Sanan     if (rank == 0) { /* print its own entries */
1269566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Edge %" PetscInt_FMT ":\n", id));
1279566063dSJacob Faibussowitsch       PetscCall(VecArrayPrint_private(viewer, nvar, xv + offset));
128a81b7fe4SHong Zhang     } else {
1297b6afd5bSHong Zhang       values[0] += 1; /* number of edges */
1307b6afd5bSHong Zhang       values[k++] = id;
1317b6afd5bSHong Zhang       values[k++] = nvar;
1327b6afd5bSHong Zhang       for (i = offset; i < offset + nvar; i++) values[k++] = xv[i];
1334062a5e5SHong Zhang     }
1344062a5e5SHong Zhang   }
1354062a5e5SHong Zhang 
1364062a5e5SHong Zhang   /* iterate over vertices */
1374062a5e5SHong Zhang   for (v = vStart; v < vEnd; v++) {
1389566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(networkdm, v, &ghostvtex));
1394062a5e5SHong Zhang     if (ghostvtex) continue;
1409566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetComponent(networkdm, v, ALL_COMPONENTS, NULL, NULL, &nvar));
1414062a5e5SHong Zhang     if (!nvar) continue;
1424062a5e5SHong Zhang 
1439566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, ALL_COMPONENTS, &offset));
1449566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalVertexIndex(networkdm, v, &id));
1454062a5e5SHong Zhang 
146dd400576SPatrick Sanan     if (rank == 0) {
1479566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Vertex %" PetscInt_FMT ":\n", id));
1489566063dSJacob Faibussowitsch       PetscCall(VecArrayPrint_private(viewer, nvar, xv + offset));
1494062a5e5SHong Zhang     } else {
1507b6afd5bSHong Zhang       values[1] += 1; /* number of vertices */
1517b6afd5bSHong Zhang       values[k++] = id;
1527b6afd5bSHong Zhang       values[k++] = nvar;
1537b6afd5bSHong Zhang       for (i = offset; i < offset + nvar; i++) values[k++] = xv[i];
1544062a5e5SHong Zhang     }
1554062a5e5SHong Zhang   }
1564062a5e5SHong Zhang 
157dd400576SPatrick Sanan   if (rank == 0) {
1584062a5e5SHong Zhang     /* proc[0] receives and prints messages */
1596497c311SBarry Smith     for (PetscMPIInt j = 1; j < size; j++) {
1606497c311SBarry Smith       PetscCall(PetscViewerASCIIPrintf(viewer, "Process [%d]\n", j));
1614062a5e5SHong Zhang 
162835f2295SStefano Zampini       PetscCallMPI(MPI_Recv(values, len, MPIU_SCALAR, j, tag, comm, &status));
1634062a5e5SHong Zhang 
164dde233f4SSatish Balay       ne = (PetscInt)PetscAbsScalar(values[0]);
165dde233f4SSatish Balay       nv = (PetscInt)PetscAbsScalar(values[1]);
1664062a5e5SHong Zhang 
1674062a5e5SHong Zhang       /* print received edges */
1684062a5e5SHong Zhang       k = 2;
1694062a5e5SHong Zhang       for (i = 0; i < ne; i++) {
170dde233f4SSatish Balay         id   = (PetscInt)PetscAbsScalar(values[k++]);
171dde233f4SSatish Balay         nvar = (PetscInt)PetscAbsScalar(values[k++]);
1729566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Edge %" PetscInt_FMT ":\n", id));
1739566063dSJacob Faibussowitsch         PetscCall(VecArrayPrint_private(viewer, nvar, values + k));
174a81b7fe4SHong Zhang         k += nvar;
1754062a5e5SHong Zhang       }
1764062a5e5SHong Zhang 
1774062a5e5SHong Zhang       /* print received vertices */
1784062a5e5SHong Zhang       for (i = 0; i < nv; i++) {
179dde233f4SSatish Balay         id   = (PetscInt)PetscAbsScalar(values[k++]);
180dde233f4SSatish Balay         nvar = (PetscInt)PetscAbsScalar(values[k++]);
1819566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Vertex %" PetscInt_FMT ":\n", id));
1829566063dSJacob Faibussowitsch         PetscCall(VecArrayPrint_private(viewer, nvar, values + k));
183a81b7fe4SHong Zhang         k += nvar;
1844062a5e5SHong Zhang       }
1854062a5e5SHong Zhang     }
1864062a5e5SHong Zhang   } else {
1874062a5e5SHong Zhang     /* sends values to proc[0] */
1886497c311SBarry Smith     PetscCallMPI(MPIU_Send((void *)values, k, MPIU_SCALAR, 0, tag, comm));
1894062a5e5SHong Zhang   }
1904062a5e5SHong Zhang 
1919566063dSJacob Faibussowitsch   PetscCall(PetscFree(values));
1929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(localX, &xv));
1939566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(networkdm, &localX));
1943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1954062a5e5SHong Zhang }
1964062a5e5SHong Zhang 
197d6acfc2dSPierre Jolivet PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecView_MPI(Vec, PetscViewer);
198a6c4b7b7SHong Zhang 
199a4e35b19SJacob Faibussowitsch static PetscErrorCode VecView_Network(Vec v, PetscViewer viewer)
200d71ae5a4SJacob Faibussowitsch {
2014dc485aaSHong Zhang   DM        dm;
2024dc485aaSHong Zhang   PetscBool isseq;
203*9f196a02SMartin Diehl   PetscBool isascii;
2044dc485aaSHong Zhang 
2054dc485aaSHong Zhang   PetscFunctionBegin;
2069566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
2075c6496baSHong Zhang   PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
208*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
2099566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)v, VECSEQ, &isseq));
210a6c4b7b7SHong Zhang 
211a6c4b7b7SHong Zhang   /* Use VecView_Network if the viewer is ASCII; use VecView_Seq/MPI for other viewer formats */
212*9f196a02SMartin Diehl   if (isascii) {
2131baa6e33SBarry Smith     if (isseq) PetscCall(VecView_Network_Seq(dm, v, viewer));
2141baa6e33SBarry Smith     else PetscCall(VecView_Network_MPI(dm, v, viewer));
2154dc485aaSHong Zhang   } else {
2161baa6e33SBarry Smith     if (isseq) PetscCall(VecView_Seq(v, viewer));
2171baa6e33SBarry Smith     else PetscCall(VecView_MPI(v, viewer));
218a6c4b7b7SHong Zhang   }
2193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2204dc485aaSHong Zhang }
2215f2c45f1SShri Abhyankar 
222d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateGlobalVector_Network(DM dm, Vec *vec)
223d71ae5a4SJacob Faibussowitsch {
2245f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network *)dm->data;
2255f2c45f1SShri Abhyankar 
2265f2c45f1SShri Abhyankar   PetscFunctionBegin;
2279566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(network->plex, vec));
2289566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Network));
2299566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*vec, dm));
2303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2315f2c45f1SShri Abhyankar }
2325f2c45f1SShri Abhyankar 
233d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateLocalVector_Network(DM dm, Vec *vec)
234d71ae5a4SJacob Faibussowitsch {
2355f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network *)dm->data;
2365f2c45f1SShri Abhyankar 
2375f2c45f1SShri Abhyankar   PetscFunctionBegin;
2389566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(network->plex, vec));
2399566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*vec, dm));
2403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2415f2c45f1SShri Abhyankar }
2425f2c45f1SShri Abhyankar 
243d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkInitializeToDefault_NonShared(DM dm)
244d71ae5a4SJacob Faibussowitsch {
245daad07d3SAidan Hamilton   DM_Network *network = (DM_Network *)dm->data;
246daad07d3SAidan Hamilton 
247daad07d3SAidan Hamilton   PetscFunctionBegin;
248daad07d3SAidan Hamilton   network->Je                 = NULL;
249daad07d3SAidan Hamilton   network->Jv                 = NULL;
250daad07d3SAidan Hamilton   network->Jvptr              = NULL;
251daad07d3SAidan Hamilton   network->userEdgeJacobian   = PETSC_FALSE;
252daad07d3SAidan Hamilton   network->userVertexJacobian = PETSC_FALSE;
253daad07d3SAidan Hamilton 
254daad07d3SAidan Hamilton   network->vertex.DofSection       = NULL;
255daad07d3SAidan Hamilton   network->vertex.GlobalDofSection = NULL;
256daad07d3SAidan Hamilton   network->vertex.mapping          = NULL;
257daad07d3SAidan Hamilton   network->vertex.sf               = NULL;
258daad07d3SAidan Hamilton 
259daad07d3SAidan Hamilton   network->edge.DofSection       = NULL;
260daad07d3SAidan Hamilton   network->edge.GlobalDofSection = NULL;
261daad07d3SAidan Hamilton   network->edge.mapping          = NULL;
262daad07d3SAidan Hamilton   network->edge.sf               = NULL;
263daad07d3SAidan Hamilton 
264daad07d3SAidan Hamilton   network->DataSection      = NULL;
265daad07d3SAidan Hamilton   network->DofSection       = NULL;
266daad07d3SAidan Hamilton   network->GlobalDofSection = NULL;
267daad07d3SAidan Hamilton   network->componentsetup   = PETSC_FALSE;
268daad07d3SAidan Hamilton 
269daad07d3SAidan Hamilton   network->plex = NULL;
270daad07d3SAidan Hamilton 
271daad07d3SAidan Hamilton   network->component  = NULL;
272daad07d3SAidan Hamilton   network->ncomponent = 0;
273daad07d3SAidan Hamilton 
274daad07d3SAidan Hamilton   network->header             = NULL;
275daad07d3SAidan Hamilton   network->cvalue             = NULL;
276daad07d3SAidan Hamilton   network->componentdataarray = NULL;
277daad07d3SAidan Hamilton 
278daad07d3SAidan Hamilton   network->max_comps_registered = DMNETWORK_MAX_COMP_REGISTERED_DEFAULT; /* return to default */
2793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
280daad07d3SAidan Hamilton }
281daad07d3SAidan Hamilton /* Default values for the parameters in DMNetwork */
282d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkInitializeToDefault(DM dm)
283d71ae5a4SJacob Faibussowitsch {
284daad07d3SAidan Hamilton   DM_Network          *network     = (DM_Network *)dm->data;
285daad07d3SAidan Hamilton   DMNetworkCloneShared cloneshared = network->cloneshared;
286daad07d3SAidan Hamilton 
287daad07d3SAidan Hamilton   PetscFunctionBegin;
288daad07d3SAidan Hamilton   PetscCall(DMNetworkInitializeToDefault_NonShared(dm));
289daad07d3SAidan Hamilton   /* Default values for shared data */
290daad07d3SAidan Hamilton   cloneshared->refct            = 1;
291daad07d3SAidan Hamilton   cloneshared->NVertices        = 0;
292daad07d3SAidan Hamilton   cloneshared->NEdges           = 0;
293daad07d3SAidan Hamilton   cloneshared->nVertices        = 0;
294daad07d3SAidan Hamilton   cloneshared->nEdges           = 0;
295daad07d3SAidan Hamilton   cloneshared->nsubnet          = 0;
296daad07d3SAidan Hamilton   cloneshared->pStart           = -1;
297daad07d3SAidan Hamilton   cloneshared->pEnd             = -1;
298daad07d3SAidan Hamilton   cloneshared->vStart           = -1;
299daad07d3SAidan Hamilton   cloneshared->vEnd             = -1;
300daad07d3SAidan Hamilton   cloneshared->eStart           = -1;
301daad07d3SAidan Hamilton   cloneshared->eEnd             = -1;
302daad07d3SAidan Hamilton   cloneshared->vltog            = NULL;
303daad07d3SAidan Hamilton   cloneshared->distributecalled = PETSC_FALSE;
304daad07d3SAidan Hamilton 
305daad07d3SAidan Hamilton   cloneshared->subnet     = NULL;
306daad07d3SAidan Hamilton   cloneshared->subnetvtx  = NULL;
307daad07d3SAidan Hamilton   cloneshared->subnetedge = NULL;
308daad07d3SAidan Hamilton   cloneshared->svtx       = NULL;
309daad07d3SAidan Hamilton   cloneshared->nsvtx      = 0;
310daad07d3SAidan Hamilton   cloneshared->Nsvtx      = 0;
311daad07d3SAidan Hamilton   cloneshared->svertices  = NULL;
312daad07d3SAidan Hamilton   cloneshared->sedgelist  = NULL;
313daad07d3SAidan Hamilton   cloneshared->svtable    = NULL;
3143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
315daad07d3SAidan Hamilton }
316daad07d3SAidan Hamilton 
317a4e35b19SJacob Faibussowitsch static PetscErrorCode DMInitialize_Network(DM dm)
318d71ae5a4SJacob Faibussowitsch {
3195f2c45f1SShri Abhyankar   PetscFunctionBegin;
3209566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(dm, 1));
321caf410d2SHong Zhang   dm->ops->view                    = DMView_Network;
3225f2c45f1SShri Abhyankar   dm->ops->setfromoptions          = DMSetFromOptions_Network;
3238415c774SShri Abhyankar   dm->ops->clone                   = DMClone_Network;
3245f2c45f1SShri Abhyankar   dm->ops->setup                   = DMSetUp_Network;
3255f2c45f1SShri Abhyankar   dm->ops->createglobalvector      = DMCreateGlobalVector_Network;
3265f2c45f1SShri Abhyankar   dm->ops->createlocalvector       = DMCreateLocalVector_Network;
3275f2c45f1SShri Abhyankar   dm->ops->getlocaltoglobalmapping = NULL;
3285f2c45f1SShri Abhyankar   dm->ops->createfieldis           = NULL;
3298afb7921SAidan Hamilton   dm->ops->createcoordinatedm      = DMCreateCoordinateDM_Network;
33099acd26cSksagiyam   dm->ops->createcellcoordinatedm  = NULL;
331ea78f98cSLisandro Dalcin   dm->ops->getcoloring             = NULL;
3325f2c45f1SShri Abhyankar   dm->ops->creatematrix            = DMCreateMatrix_Network;
333ea78f98cSLisandro Dalcin   dm->ops->createinterpolation     = NULL;
334ea78f98cSLisandro Dalcin   dm->ops->createinjection         = NULL;
335ea78f98cSLisandro Dalcin   dm->ops->refine                  = NULL;
336ea78f98cSLisandro Dalcin   dm->ops->coarsen                 = NULL;
337ea78f98cSLisandro Dalcin   dm->ops->refinehierarchy         = NULL;
338ea78f98cSLisandro Dalcin   dm->ops->coarsenhierarchy        = NULL;
3395f2c45f1SShri Abhyankar   dm->ops->globaltolocalbegin      = DMGlobalToLocalBegin_Network;
3405f2c45f1SShri Abhyankar   dm->ops->globaltolocalend        = DMGlobalToLocalEnd_Network;
3415f2c45f1SShri Abhyankar   dm->ops->localtoglobalbegin      = DMLocalToGlobalBegin_Network;
3425f2c45f1SShri Abhyankar   dm->ops->localtoglobalend        = DMLocalToGlobalEnd_Network;
3435f2c45f1SShri Abhyankar   dm->ops->destroy                 = DMDestroy_Network;
3445f2c45f1SShri Abhyankar   dm->ops->createsubdm             = NULL;
3455f2c45f1SShri Abhyankar   dm->ops->locatepoints            = NULL;
3463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3475f2c45f1SShri Abhyankar }
348daad07d3SAidan Hamilton /*
349da81f932SPierre Jolivet   copies over the subnetid and index portions of the DMNetworkComponentHeader from original dm to the newdm
350daad07d3SAidan Hamilton */
351d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMNetworkCopyHeaderTopological(DM dm, DM newdm)
352d71ae5a4SJacob Faibussowitsch {
353daad07d3SAidan Hamilton   DM_Network *network = (DM_Network *)dm->data, *newnetwork = (DM_Network *)newdm->data;
354daad07d3SAidan Hamilton   PetscInt    p, i, np, index, subnetid;
355daad07d3SAidan Hamilton 
356daad07d3SAidan Hamilton   PetscFunctionBegin;
357daad07d3SAidan Hamilton   np = network->cloneshared->pEnd - network->cloneshared->pStart;
358daad07d3SAidan Hamilton   PetscCall(PetscCalloc2(np, &newnetwork->header, np, &newnetwork->cvalue));
359daad07d3SAidan Hamilton   for (i = 0; i < np; i++) {
360daad07d3SAidan Hamilton     p = i + network->cloneshared->pStart;
361daad07d3SAidan Hamilton     PetscCall(DMNetworkGetSubnetID(dm, p, &subnetid));
362daad07d3SAidan Hamilton     PetscCall(DMNetworkGetIndex(dm, p, &index));
363daad07d3SAidan Hamilton     newnetwork->header[i].index        = index;
364daad07d3SAidan Hamilton     newnetwork->header[i].subnetid     = subnetid;
365daad07d3SAidan Hamilton     newnetwork->header[i].size         = NULL;
366daad07d3SAidan Hamilton     newnetwork->header[i].key          = NULL;
367daad07d3SAidan Hamilton     newnetwork->header[i].offset       = NULL;
368daad07d3SAidan Hamilton     newnetwork->header[i].nvar         = NULL;
369daad07d3SAidan Hamilton     newnetwork->header[i].offsetvarrel = NULL;
370daad07d3SAidan Hamilton     newnetwork->header[i].ndata        = 0;
371daad07d3SAidan Hamilton     newnetwork->header[i].maxcomps     = DMNETWORK_MAX_COMP_AT_POINT_DEFAULT;
372daad07d3SAidan Hamilton     newnetwork->header[i].hsize        = sizeof(struct _p_DMNetworkComponentHeader) / sizeof(sizeof(DMNetworkComponentGenericDataType));
373daad07d3SAidan Hamilton   }
3743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
375daad07d3SAidan Hamilton }
3765f2c45f1SShri Abhyankar 
377d71ae5a4SJacob Faibussowitsch PetscErrorCode DMClone_Network(DM dm, DM *newdm)
378d71ae5a4SJacob Faibussowitsch {
379daad07d3SAidan Hamilton   DM_Network *network = (DM_Network *)dm->data, *newnetwork = NULL;
3808415c774SShri Abhyankar 
3818415c774SShri Abhyankar   PetscFunctionBegin;
382daad07d3SAidan Hamilton   network->cloneshared->refct++;
3834dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&newnetwork));
384daad07d3SAidan Hamilton   (*newdm)->data = newnetwork;
385daad07d3SAidan Hamilton   PetscCall(DMNetworkInitializeToDefault_NonShared(*newdm));
386daad07d3SAidan Hamilton   newnetwork->cloneshared = network->cloneshared; /* Share all data that can be cloneshared */
3879b626cebSHong Zhang 
3889b626cebSHong Zhang   PetscCheck(network->plex, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_NULL, "Must call DMNetworkLayoutSetUp() first");
389daad07d3SAidan Hamilton   PetscCall(DMClone(network->plex, &newnetwork->plex));
390daad07d3SAidan Hamilton   PetscCall(DMNetworkCopyHeaderTopological(dm, *newdm));
391daad07d3SAidan Hamilton   PetscCall(DMNetworkInitializeNonTopological(*newdm)); /* initialize all non-topological data to the state after DMNetworkLayoutSetUp as been called */
3929566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, DMNETWORK));
3939566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Network(*newdm));
3943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3958415c774SShri Abhyankar }
3969b626cebSHong Zhang 
3978afb7921SAidan Hamilton /* Developer Note: Be aware that the plex inside of the network does not have a coordinate plex.
3988afb7921SAidan Hamilton */
3998afb7921SAidan Hamilton PetscErrorCode DMCreateCoordinateDM_Network(DM dm, DM *cdm)
4008afb7921SAidan Hamilton {
4018afb7921SAidan Hamilton   DM_Network *newnetwork = NULL;
4028afb7921SAidan Hamilton   PetscInt    Nf;
403dd4c3f67SMatthew G. Knepley   const char *prefix;
4048afb7921SAidan Hamilton 
4058afb7921SAidan Hamilton   PetscFunctionBegin;
4068afb7921SAidan Hamilton   PetscCall(DMClone(dm, cdm));
4078afb7921SAidan Hamilton   newnetwork = (DM_Network *)(*cdm)->data;
4088afb7921SAidan Hamilton   PetscCall(DMGetNumFields(newnetwork->plex, &Nf));
4098afb7921SAidan Hamilton   PetscCall(DMSetNumFields(*cdm, Nf)); /* consistency with the coordinate plex */
410dd4c3f67SMatthew G. Knepley   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix));
411dd4c3f67SMatthew G. Knepley   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*cdm, prefix));
412dd4c3f67SMatthew G. Knepley   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)*cdm, "cdm_"));
4133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4148afb7921SAidan Hamilton }
4158415c774SShri Abhyankar 
4165f2c45f1SShri Abhyankar /*MC
4175f2c45f1SShri Abhyankar   DMNETWORK = "network" - A DM object that encapsulates an unstructured network. The implementation is based on the DM object
4185f2c45f1SShri Abhyankar                           DMPlex that manages unstructured grids. Distributed networks use a non-overlapping partitioning of
4195f2c45f1SShri Abhyankar                           the edges. In the local representation, Vecs contain all unknowns in the interior and shared boundary.
4205f2c45f1SShri Abhyankar                           This is specified by a PetscSection object. Ownership in the global representation is determined by
4215f2c45f1SShri Abhyankar                           ownership of the underlying DMPlex points. This is specified by another PetscSection object.
4225f2c45f1SShri Abhyankar 
4235f2c45f1SShri Abhyankar   Level: intermediate
4245f2c45f1SShri Abhyankar 
425db781477SPatrick Sanan .seealso: `DMType`, `DMNetworkCreate()`, `DMCreate()`, `DMSetType()`
4265f2c45f1SShri Abhyankar M*/
4275f2c45f1SShri Abhyankar 
428d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Network(DM dm)
429d71ae5a4SJacob Faibussowitsch {
4305f2c45f1SShri Abhyankar   DM_Network *network;
4315f2c45f1SShri Abhyankar 
4325f2c45f1SShri Abhyankar   PetscFunctionBegin;
4335f2c45f1SShri Abhyankar   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4344dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&network));
4354dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&network->cloneshared));
4365f2c45f1SShri Abhyankar   dm->data = network;
4375f2c45f1SShri Abhyankar 
438daad07d3SAidan Hamilton   PetscCall(DMNetworkInitializeToDefault(dm));
4399566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Network(dm));
4403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4415f2c45f1SShri Abhyankar }
4425f2c45f1SShri Abhyankar 
4435f2c45f1SShri Abhyankar /*@
4445f2c45f1SShri Abhyankar   DMNetworkCreate - Creates a DMNetwork object, which encapsulates an unstructured network.
4455f2c45f1SShri Abhyankar 
446d083f849SBarry Smith   Collective
4475f2c45f1SShri Abhyankar 
4485f2c45f1SShri Abhyankar   Input Parameter:
4495f2c45f1SShri Abhyankar . comm - The communicator for the DMNetwork object
4505f2c45f1SShri Abhyankar 
4515f2c45f1SShri Abhyankar   Output Parameter:
4525f2c45f1SShri Abhyankar . network - The DMNetwork object
4535f2c45f1SShri Abhyankar 
4545f2c45f1SShri Abhyankar   Level: beginner
4555f2c45f1SShri Abhyankar 
456a4e35b19SJacob Faibussowitsch .seealso: `DMCreate()`
4575f2c45f1SShri Abhyankar @*/
458d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkCreate(MPI_Comm comm, DM *network)
459d71ae5a4SJacob Faibussowitsch {
4605f2c45f1SShri Abhyankar   PetscFunctionBegin;
4614f572ea9SToby Isaac   PetscAssertPointer(network, 2);
4629566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, network));
4639566063dSJacob Faibussowitsch   PetscCall(DMSetType(*network, DMNETWORK));
4643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4655f2c45f1SShri Abhyankar }
466