1af0996ceSBarry Smith #include <petsc/private/dmnetworkimpl.h> /*I "petscdmnetwork.h" I*/
295af8c53SHong Zhang #include <petsc/private/vecimpl.h>
35f2c45f1SShri Abhyankar
DMSetFromOptions_Network(DM dm,PetscOptionItems PetscOptionsObject)4ce78bad3SBarry Smith static PetscErrorCode DMSetFromOptions_Network(DM dm, PetscOptionItems PetscOptionsObject)
5d71ae5a4SJacob Faibussowitsch {
65f2c45f1SShri Abhyankar PetscFunctionBegin;
7d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "DMNetwork Options");
8d0609cedSBarry Smith PetscOptionsHeadEnd();
93ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
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
VecArrayPrint_private(PetscViewer viewer,PetscInt n,const PetscScalar * xv)23d71ae5a4SJacob Faibussowitsch static PetscErrorCode VecArrayPrint_private(PetscViewer viewer, PetscInt n, const PetscScalar *xv)
24d71ae5a4SJacob Faibussowitsch {
25a81b7fe4SHong Zhang PetscInt i;
26a81b7fe4SHong Zhang
27a81b7fe4SHong Zhang PetscFunctionBegin;
28a81b7fe4SHong Zhang for (i = 0; i < n; i++) {
29a81b7fe4SHong Zhang #if defined(PETSC_USE_COMPLEX)
30a81b7fe4SHong Zhang if (PetscImaginaryPart(xv[i]) > 0.0) {
319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " %g + %g i\n", (double)PetscRealPart(xv[i]), (double)PetscImaginaryPart(xv[i])));
32a81b7fe4SHong Zhang } else if (PetscImaginaryPart(xv[i]) < 0.0) {
339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " %g - %g i\n", (double)PetscRealPart(xv[i]), -(double)PetscImaginaryPart(xv[i])));
341baa6e33SBarry Smith } else PetscCall(PetscViewerASCIIPrintf(viewer, " %g\n", (double)PetscRealPart(xv[i])));
35a81b7fe4SHong Zhang #else
369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " %g\n", (double)xv[i]));
37a81b7fe4SHong Zhang #endif
38a81b7fe4SHong Zhang }
393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
40a81b7fe4SHong Zhang }
41a81b7fe4SHong Zhang
VecView_Network_Seq(DM networkdm,Vec X,PetscViewer viewer)42d71ae5a4SJacob Faibussowitsch static PetscErrorCode VecView_Network_Seq(DM networkdm, Vec X, PetscViewer viewer)
43d71ae5a4SJacob Faibussowitsch {
447b6afd5bSHong Zhang PetscInt e, v, Start, End, offset, nvar, id;
454dc485aaSHong Zhang const PetscScalar *xv;
464dc485aaSHong Zhang
474dc485aaSHong Zhang PetscFunctionBegin;
489566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &xv));
494dc485aaSHong Zhang
504dc485aaSHong Zhang /* iterate over edges */
519566063dSJacob Faibussowitsch PetscCall(DMNetworkGetEdgeRange(networkdm, &Start, &End));
524dc485aaSHong Zhang for (e = Start; e < End; e++) {
539566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
544dc485aaSHong Zhang if (!nvar) continue;
55a81b7fe4SHong Zhang
569566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, e, ALL_COMPONENTS, &offset));
579566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalEdgeIndex(networkdm, e, &id));
587b6afd5bSHong Zhang
599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Edge %" PetscInt_FMT ":\n", id));
609566063dSJacob Faibussowitsch PetscCall(VecArrayPrint_private(viewer, nvar, xv + offset));
614dc485aaSHong Zhang }
624dc485aaSHong Zhang
634dc485aaSHong Zhang /* iterate over vertices */
649566063dSJacob Faibussowitsch PetscCall(DMNetworkGetVertexRange(networkdm, &Start, &End));
654dc485aaSHong Zhang for (v = Start; v < End; v++) {
669566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, v, ALL_COMPONENTS, NULL, NULL, &nvar));
674dc485aaSHong Zhang if (!nvar) continue;
68a81b7fe4SHong Zhang
699566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, ALL_COMPONENTS, &offset));
709566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVertexIndex(networkdm, v, &id));
7162aa33baSHong Zhang
729566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Vertex %" PetscInt_FMT ":\n", id));
739566063dSJacob Faibussowitsch PetscCall(VecArrayPrint_private(viewer, nvar, xv + offset));
744dc485aaSHong Zhang }
759566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer));
769566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &xv));
773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
784dc485aaSHong Zhang }
794dc485aaSHong Zhang
VecView_Network_MPI(DM networkdm,Vec X,PetscViewer viewer)80d71ae5a4SJacob Faibussowitsch static PetscErrorCode VecView_Network_MPI(DM networkdm, Vec X, PetscViewer viewer)
81d71ae5a4SJacob Faibussowitsch {
82835f2295SStefano Zampini PetscInt i, e, v, eStart, eEnd, vStart, vEnd, offset, nvar, len_loc, k;
834062a5e5SHong Zhang const PetscScalar *xv;
844062a5e5SHong Zhang MPI_Comm comm;
85835f2295SStefano Zampini PetscMPIInt size, rank, tag = ((PetscObject)viewer)->tag, len;
864062a5e5SHong Zhang Vec localX;
874062a5e5SHong Zhang PetscBool ghostvtex;
884062a5e5SHong Zhang PetscScalar *values;
896497c311SBarry Smith PetscInt ne, nv, id;
904062a5e5SHong Zhang MPI_Status status;
914062a5e5SHong Zhang
924062a5e5SHong Zhang PetscFunctionBegin;
939566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)networkdm, &comm));
949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size));
959566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank));
964062a5e5SHong Zhang
979566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(networkdm, &localX));
989566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX));
999566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX));
1009566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(localX, &xv));
1014062a5e5SHong Zhang
1029566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(localX, &len_loc));
1034062a5e5SHong Zhang
1049566063dSJacob Faibussowitsch PetscCall(DMNetworkGetEdgeRange(networkdm, &eStart, &eEnd));
1059566063dSJacob Faibussowitsch PetscCall(DMNetworkGetVertexRange(networkdm, &vStart, &vEnd));
106e8c2e809SHong Zhang len_loc += 2 * (1 + eEnd - eStart + vEnd - vStart);
1074062a5e5SHong Zhang
108e85e6aecSHong Zhang /* values = [nedges, nvertices; id, nvar, xedge; ...; id, nvars, xvertex;...], to be sent to proc[0] */
109835f2295SStefano Zampini PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &len_loc, 1, MPIU_INT, MPI_MAX, comm));
110835f2295SStefano Zampini PetscCall(PetscMPIIntCast(len_loc, &len));
1119566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(len, &values));
1124062a5e5SHong Zhang
11348a46eb9SPierre Jolivet if (rank == 0) PetscCall(PetscViewerASCIIPrintf(viewer, "Process [%d]\n", rank));
1144062a5e5SHong Zhang
1154062a5e5SHong Zhang /* iterate over edges */
1164062a5e5SHong Zhang k = 2;
117e8c2e809SHong Zhang for (e = eStart; e < eEnd; e++) {
1189566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, e, ALL_COMPONENTS, NULL, NULL, &nvar));
1194062a5e5SHong Zhang if (!nvar) continue;
1204062a5e5SHong Zhang
1219566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, e, ALL_COMPONENTS, &offset));
1229566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalEdgeIndex(networkdm, e, &id));
12362aa33baSHong Zhang
124dd400576SPatrick Sanan if (rank == 0) { /* print its own entries */
1259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Edge %" PetscInt_FMT ":\n", id));
1269566063dSJacob Faibussowitsch PetscCall(VecArrayPrint_private(viewer, nvar, xv + offset));
127a81b7fe4SHong Zhang } else {
1287b6afd5bSHong Zhang values[0] += 1; /* number of edges */
1297b6afd5bSHong Zhang values[k++] = id;
1307b6afd5bSHong Zhang values[k++] = nvar;
1317b6afd5bSHong Zhang for (i = offset; i < offset + nvar; i++) values[k++] = xv[i];
1324062a5e5SHong Zhang }
1334062a5e5SHong Zhang }
1344062a5e5SHong Zhang
1354062a5e5SHong Zhang /* iterate over vertices */
1364062a5e5SHong Zhang for (v = vStart; v < vEnd; v++) {
1379566063dSJacob Faibussowitsch PetscCall(DMNetworkIsGhostVertex(networkdm, v, &ghostvtex));
1384062a5e5SHong Zhang if (ghostvtex) continue;
1399566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, v, ALL_COMPONENTS, NULL, NULL, &nvar));
1404062a5e5SHong Zhang if (!nvar) continue;
1414062a5e5SHong Zhang
1429566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, ALL_COMPONENTS, &offset));
1439566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVertexIndex(networkdm, v, &id));
1444062a5e5SHong Zhang
145dd400576SPatrick Sanan if (rank == 0) {
1469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Vertex %" PetscInt_FMT ":\n", id));
1479566063dSJacob Faibussowitsch PetscCall(VecArrayPrint_private(viewer, nvar, xv + offset));
1484062a5e5SHong Zhang } else {
1497b6afd5bSHong Zhang values[1] += 1; /* number of vertices */
1507b6afd5bSHong Zhang values[k++] = id;
1517b6afd5bSHong Zhang values[k++] = nvar;
1527b6afd5bSHong Zhang for (i = offset; i < offset + nvar; i++) values[k++] = xv[i];
1534062a5e5SHong Zhang }
1544062a5e5SHong Zhang }
1554062a5e5SHong Zhang
156dd400576SPatrick Sanan if (rank == 0) {
1574062a5e5SHong Zhang /* proc[0] receives and prints messages */
1586497c311SBarry Smith for (PetscMPIInt j = 1; j < size; j++) {
1596497c311SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "Process [%d]\n", j));
1604062a5e5SHong Zhang
161835f2295SStefano Zampini PetscCallMPI(MPI_Recv(values, len, MPIU_SCALAR, j, tag, comm, &status));
1624062a5e5SHong Zhang
163dde233f4SSatish Balay ne = (PetscInt)PetscAbsScalar(values[0]);
164dde233f4SSatish Balay nv = (PetscInt)PetscAbsScalar(values[1]);
1654062a5e5SHong Zhang
1664062a5e5SHong Zhang /* print received edges */
1674062a5e5SHong Zhang k = 2;
1684062a5e5SHong Zhang for (i = 0; i < ne; i++) {
169dde233f4SSatish Balay id = (PetscInt)PetscAbsScalar(values[k++]);
170dde233f4SSatish Balay nvar = (PetscInt)PetscAbsScalar(values[k++]);
1719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Edge %" PetscInt_FMT ":\n", id));
1729566063dSJacob Faibussowitsch PetscCall(VecArrayPrint_private(viewer, nvar, values + k));
173a81b7fe4SHong Zhang k += nvar;
1744062a5e5SHong Zhang }
1754062a5e5SHong Zhang
1764062a5e5SHong Zhang /* print received vertices */
1774062a5e5SHong Zhang for (i = 0; i < nv; i++) {
178dde233f4SSatish Balay id = (PetscInt)PetscAbsScalar(values[k++]);
179dde233f4SSatish Balay nvar = (PetscInt)PetscAbsScalar(values[k++]);
1809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Vertex %" PetscInt_FMT ":\n", id));
1819566063dSJacob Faibussowitsch PetscCall(VecArrayPrint_private(viewer, nvar, values + k));
182a81b7fe4SHong Zhang k += nvar;
1834062a5e5SHong Zhang }
1844062a5e5SHong Zhang }
1854062a5e5SHong Zhang } else {
1864062a5e5SHong Zhang /* sends values to proc[0] */
1876497c311SBarry Smith PetscCallMPI(MPIU_Send((void *)values, k, MPIU_SCALAR, 0, tag, comm));
1884062a5e5SHong Zhang }
1894062a5e5SHong Zhang
1909566063dSJacob Faibussowitsch PetscCall(PetscFree(values));
1919566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(localX, &xv));
1929566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(networkdm, &localX));
1933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1944062a5e5SHong Zhang }
1954062a5e5SHong Zhang
196d6acfc2dSPierre Jolivet PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecView_MPI(Vec, PetscViewer);
197a6c4b7b7SHong Zhang
VecView_Network(Vec v,PetscViewer viewer)198a4e35b19SJacob Faibussowitsch static PetscErrorCode VecView_Network(Vec v, PetscViewer viewer)
199d71ae5a4SJacob Faibussowitsch {
2004dc485aaSHong Zhang DM dm;
2014dc485aaSHong Zhang PetscBool isseq;
2029f196a02SMartin Diehl PetscBool isascii;
2034dc485aaSHong Zhang
2044dc485aaSHong Zhang PetscFunctionBegin;
2059566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm));
2065c6496baSHong Zhang PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
2079f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
2089566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)v, VECSEQ, &isseq));
209a6c4b7b7SHong Zhang
210a6c4b7b7SHong Zhang /* Use VecView_Network if the viewer is ASCII; use VecView_Seq/MPI for other viewer formats */
2119f196a02SMartin Diehl if (isascii) {
2121baa6e33SBarry Smith if (isseq) PetscCall(VecView_Network_Seq(dm, v, viewer));
2131baa6e33SBarry Smith else PetscCall(VecView_Network_MPI(dm, v, viewer));
2144dc485aaSHong Zhang } else {
2151baa6e33SBarry Smith if (isseq) PetscCall(VecView_Seq(v, viewer));
2161baa6e33SBarry Smith else PetscCall(VecView_MPI(v, viewer));
217a6c4b7b7SHong Zhang }
2183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2194dc485aaSHong Zhang }
2205f2c45f1SShri Abhyankar
DMCreateGlobalVector_Network(DM dm,Vec * vec)221d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateGlobalVector_Network(DM dm, Vec *vec)
222d71ae5a4SJacob Faibussowitsch {
2235f2c45f1SShri Abhyankar DM_Network *network = (DM_Network *)dm->data;
2245f2c45f1SShri Abhyankar
2255f2c45f1SShri Abhyankar PetscFunctionBegin;
2269566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(network->plex, vec));
227*57d50842SBarry Smith PetscCall(VecSetOperation(*vec, VECOP_VIEW, (PetscErrorCodeFn *)VecView_Network));
2289566063dSJacob Faibussowitsch PetscCall(VecSetDM(*vec, dm));
2293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2305f2c45f1SShri Abhyankar }
2315f2c45f1SShri Abhyankar
DMCreateLocalVector_Network(DM dm,Vec * vec)232d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateLocalVector_Network(DM dm, Vec *vec)
233d71ae5a4SJacob Faibussowitsch {
2345f2c45f1SShri Abhyankar DM_Network *network = (DM_Network *)dm->data;
2355f2c45f1SShri Abhyankar
2365f2c45f1SShri Abhyankar PetscFunctionBegin;
2379566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(network->plex, vec));
2389566063dSJacob Faibussowitsch PetscCall(VecSetDM(*vec, dm));
2393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2405f2c45f1SShri Abhyankar }
2415f2c45f1SShri Abhyankar
DMNetworkInitializeToDefault_NonShared(DM dm)242d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkInitializeToDefault_NonShared(DM dm)
243d71ae5a4SJacob Faibussowitsch {
244daad07d3SAidan Hamilton DM_Network *network = (DM_Network *)dm->data;
245daad07d3SAidan Hamilton
246daad07d3SAidan Hamilton PetscFunctionBegin;
247daad07d3SAidan Hamilton network->Je = NULL;
248daad07d3SAidan Hamilton network->Jv = NULL;
249daad07d3SAidan Hamilton network->Jvptr = NULL;
250daad07d3SAidan Hamilton network->userEdgeJacobian = PETSC_FALSE;
251daad07d3SAidan Hamilton network->userVertexJacobian = PETSC_FALSE;
252daad07d3SAidan Hamilton
253daad07d3SAidan Hamilton network->vertex.DofSection = NULL;
254daad07d3SAidan Hamilton network->vertex.GlobalDofSection = NULL;
255daad07d3SAidan Hamilton network->vertex.mapping = NULL;
256daad07d3SAidan Hamilton network->vertex.sf = NULL;
257daad07d3SAidan Hamilton
258daad07d3SAidan Hamilton network->edge.DofSection = NULL;
259daad07d3SAidan Hamilton network->edge.GlobalDofSection = NULL;
260daad07d3SAidan Hamilton network->edge.mapping = NULL;
261daad07d3SAidan Hamilton network->edge.sf = NULL;
262daad07d3SAidan Hamilton
263daad07d3SAidan Hamilton network->DataSection = NULL;
264daad07d3SAidan Hamilton network->DofSection = NULL;
265daad07d3SAidan Hamilton network->GlobalDofSection = NULL;
266daad07d3SAidan Hamilton network->componentsetup = PETSC_FALSE;
267daad07d3SAidan Hamilton
268daad07d3SAidan Hamilton network->plex = NULL;
269daad07d3SAidan Hamilton
270daad07d3SAidan Hamilton network->component = NULL;
271daad07d3SAidan Hamilton network->ncomponent = 0;
272daad07d3SAidan Hamilton
273daad07d3SAidan Hamilton network->header = NULL;
274daad07d3SAidan Hamilton network->cvalue = NULL;
275daad07d3SAidan Hamilton network->componentdataarray = NULL;
276daad07d3SAidan Hamilton
277daad07d3SAidan Hamilton network->max_comps_registered = DMNETWORK_MAX_COMP_REGISTERED_DEFAULT; /* return to default */
2783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
279daad07d3SAidan Hamilton }
280daad07d3SAidan Hamilton /* Default values for the parameters in DMNetwork */
DMNetworkInitializeToDefault(DM dm)281d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkInitializeToDefault(DM dm)
282d71ae5a4SJacob Faibussowitsch {
283daad07d3SAidan Hamilton DM_Network *network = (DM_Network *)dm->data;
284daad07d3SAidan Hamilton DMNetworkCloneShared cloneshared = network->cloneshared;
285daad07d3SAidan Hamilton
286daad07d3SAidan Hamilton PetscFunctionBegin;
287daad07d3SAidan Hamilton PetscCall(DMNetworkInitializeToDefault_NonShared(dm));
288daad07d3SAidan Hamilton /* Default values for shared data */
289daad07d3SAidan Hamilton cloneshared->refct = 1;
290daad07d3SAidan Hamilton cloneshared->NVertices = 0;
291daad07d3SAidan Hamilton cloneshared->NEdges = 0;
292daad07d3SAidan Hamilton cloneshared->nVertices = 0;
293daad07d3SAidan Hamilton cloneshared->nEdges = 0;
294daad07d3SAidan Hamilton cloneshared->nsubnet = 0;
295daad07d3SAidan Hamilton cloneshared->pStart = -1;
296daad07d3SAidan Hamilton cloneshared->pEnd = -1;
297daad07d3SAidan Hamilton cloneshared->vStart = -1;
298daad07d3SAidan Hamilton cloneshared->vEnd = -1;
299daad07d3SAidan Hamilton cloneshared->eStart = -1;
300daad07d3SAidan Hamilton cloneshared->eEnd = -1;
301daad07d3SAidan Hamilton cloneshared->vltog = NULL;
302daad07d3SAidan Hamilton cloneshared->distributecalled = PETSC_FALSE;
303daad07d3SAidan Hamilton
304daad07d3SAidan Hamilton cloneshared->subnet = NULL;
305daad07d3SAidan Hamilton cloneshared->subnetvtx = NULL;
306daad07d3SAidan Hamilton cloneshared->subnetedge = NULL;
307daad07d3SAidan Hamilton cloneshared->svtx = NULL;
308daad07d3SAidan Hamilton cloneshared->nsvtx = 0;
309daad07d3SAidan Hamilton cloneshared->Nsvtx = 0;
310daad07d3SAidan Hamilton cloneshared->svertices = NULL;
311daad07d3SAidan Hamilton cloneshared->sedgelist = NULL;
312daad07d3SAidan Hamilton cloneshared->svtable = NULL;
3133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
314daad07d3SAidan Hamilton }
315daad07d3SAidan Hamilton
DMInitialize_Network(DM dm)316a4e35b19SJacob Faibussowitsch static PetscErrorCode DMInitialize_Network(DM dm)
317d71ae5a4SJacob Faibussowitsch {
3185f2c45f1SShri Abhyankar PetscFunctionBegin;
3199566063dSJacob Faibussowitsch PetscCall(DMSetDimension(dm, 1));
320caf410d2SHong Zhang dm->ops->view = DMView_Network;
3215f2c45f1SShri Abhyankar dm->ops->setfromoptions = DMSetFromOptions_Network;
3228415c774SShri Abhyankar dm->ops->clone = DMClone_Network;
3235f2c45f1SShri Abhyankar dm->ops->setup = DMSetUp_Network;
3245f2c45f1SShri Abhyankar dm->ops->createglobalvector = DMCreateGlobalVector_Network;
3255f2c45f1SShri Abhyankar dm->ops->createlocalvector = DMCreateLocalVector_Network;
3265f2c45f1SShri Abhyankar dm->ops->getlocaltoglobalmapping = NULL;
3275f2c45f1SShri Abhyankar dm->ops->createfieldis = NULL;
3288afb7921SAidan Hamilton dm->ops->createcoordinatedm = DMCreateCoordinateDM_Network;
32999acd26cSksagiyam dm->ops->createcellcoordinatedm = NULL;
330ea78f98cSLisandro Dalcin dm->ops->getcoloring = NULL;
3315f2c45f1SShri Abhyankar dm->ops->creatematrix = DMCreateMatrix_Network;
332ea78f98cSLisandro Dalcin dm->ops->createinterpolation = NULL;
333ea78f98cSLisandro Dalcin dm->ops->createinjection = NULL;
334ea78f98cSLisandro Dalcin dm->ops->refine = NULL;
335ea78f98cSLisandro Dalcin dm->ops->coarsen = NULL;
336ea78f98cSLisandro Dalcin dm->ops->refinehierarchy = NULL;
337ea78f98cSLisandro Dalcin dm->ops->coarsenhierarchy = NULL;
3385f2c45f1SShri Abhyankar dm->ops->globaltolocalbegin = DMGlobalToLocalBegin_Network;
3395f2c45f1SShri Abhyankar dm->ops->globaltolocalend = DMGlobalToLocalEnd_Network;
3405f2c45f1SShri Abhyankar dm->ops->localtoglobalbegin = DMLocalToGlobalBegin_Network;
3415f2c45f1SShri Abhyankar dm->ops->localtoglobalend = DMLocalToGlobalEnd_Network;
3425f2c45f1SShri Abhyankar dm->ops->destroy = DMDestroy_Network;
3435f2c45f1SShri Abhyankar dm->ops->createsubdm = NULL;
3445f2c45f1SShri Abhyankar dm->ops->locatepoints = NULL;
3453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3465f2c45f1SShri Abhyankar }
347daad07d3SAidan Hamilton /*
348da81f932SPierre Jolivet copies over the subnetid and index portions of the DMNetworkComponentHeader from original dm to the newdm
349daad07d3SAidan Hamilton */
DMNetworkCopyHeaderTopological(DM dm,DM newdm)350d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMNetworkCopyHeaderTopological(DM dm, DM newdm)
351d71ae5a4SJacob Faibussowitsch {
352daad07d3SAidan Hamilton DM_Network *network = (DM_Network *)dm->data, *newnetwork = (DM_Network *)newdm->data;
353daad07d3SAidan Hamilton PetscInt p, i, np, index, subnetid;
354daad07d3SAidan Hamilton
355daad07d3SAidan Hamilton PetscFunctionBegin;
356daad07d3SAidan Hamilton np = network->cloneshared->pEnd - network->cloneshared->pStart;
357daad07d3SAidan Hamilton PetscCall(PetscCalloc2(np, &newnetwork->header, np, &newnetwork->cvalue));
358daad07d3SAidan Hamilton for (i = 0; i < np; i++) {
359daad07d3SAidan Hamilton p = i + network->cloneshared->pStart;
360daad07d3SAidan Hamilton PetscCall(DMNetworkGetSubnetID(dm, p, &subnetid));
361daad07d3SAidan Hamilton PetscCall(DMNetworkGetIndex(dm, p, &index));
362daad07d3SAidan Hamilton newnetwork->header[i].index = index;
363daad07d3SAidan Hamilton newnetwork->header[i].subnetid = subnetid;
364daad07d3SAidan Hamilton newnetwork->header[i].size = NULL;
365daad07d3SAidan Hamilton newnetwork->header[i].key = NULL;
366daad07d3SAidan Hamilton newnetwork->header[i].offset = NULL;
367daad07d3SAidan Hamilton newnetwork->header[i].nvar = NULL;
368daad07d3SAidan Hamilton newnetwork->header[i].offsetvarrel = NULL;
369daad07d3SAidan Hamilton newnetwork->header[i].ndata = 0;
370daad07d3SAidan Hamilton newnetwork->header[i].maxcomps = DMNETWORK_MAX_COMP_AT_POINT_DEFAULT;
371daad07d3SAidan Hamilton newnetwork->header[i].hsize = sizeof(struct _p_DMNetworkComponentHeader) / sizeof(sizeof(DMNetworkComponentGenericDataType));
372daad07d3SAidan Hamilton }
3733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
374daad07d3SAidan Hamilton }
3755f2c45f1SShri Abhyankar
DMClone_Network(DM dm,DM * newdm)376d71ae5a4SJacob Faibussowitsch PetscErrorCode DMClone_Network(DM dm, DM *newdm)
377d71ae5a4SJacob Faibussowitsch {
378daad07d3SAidan Hamilton DM_Network *network = (DM_Network *)dm->data, *newnetwork = NULL;
3798415c774SShri Abhyankar
3808415c774SShri Abhyankar PetscFunctionBegin;
381daad07d3SAidan Hamilton network->cloneshared->refct++;
3824dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&newnetwork));
383daad07d3SAidan Hamilton (*newdm)->data = newnetwork;
384daad07d3SAidan Hamilton PetscCall(DMNetworkInitializeToDefault_NonShared(*newdm));
385daad07d3SAidan Hamilton newnetwork->cloneshared = network->cloneshared; /* Share all data that can be cloneshared */
3869b626cebSHong Zhang
3879b626cebSHong Zhang PetscCheck(network->plex, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_NULL, "Must call DMNetworkLayoutSetUp() first");
388daad07d3SAidan Hamilton PetscCall(DMClone(network->plex, &newnetwork->plex));
389daad07d3SAidan Hamilton PetscCall(DMNetworkCopyHeaderTopological(dm, *newdm));
390daad07d3SAidan Hamilton PetscCall(DMNetworkInitializeNonTopological(*newdm)); /* initialize all non-topological data to the state after DMNetworkLayoutSetUp as been called */
3919566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, DMNETWORK));
3929566063dSJacob Faibussowitsch PetscCall(DMInitialize_Network(*newdm));
3933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3948415c774SShri Abhyankar }
3959b626cebSHong Zhang
3968afb7921SAidan Hamilton /* Developer Note: Be aware that the plex inside of the network does not have a coordinate plex.
3978afb7921SAidan Hamilton */
DMCreateCoordinateDM_Network(DM dm,DM * cdm)3988afb7921SAidan Hamilton PetscErrorCode DMCreateCoordinateDM_Network(DM dm, DM *cdm)
3998afb7921SAidan Hamilton {
4008afb7921SAidan Hamilton DM_Network *newnetwork = NULL;
4018afb7921SAidan Hamilton PetscInt Nf;
402dd4c3f67SMatthew G. Knepley const char *prefix;
4038afb7921SAidan Hamilton
4048afb7921SAidan Hamilton PetscFunctionBegin;
4058afb7921SAidan Hamilton PetscCall(DMClone(dm, cdm));
4068afb7921SAidan Hamilton newnetwork = (DM_Network *)(*cdm)->data;
4078afb7921SAidan Hamilton PetscCall(DMGetNumFields(newnetwork->plex, &Nf));
4088afb7921SAidan Hamilton PetscCall(DMSetNumFields(*cdm, Nf)); /* consistency with the coordinate plex */
409dd4c3f67SMatthew G. Knepley PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix));
410dd4c3f67SMatthew G. Knepley PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*cdm, prefix));
411dd4c3f67SMatthew G. Knepley PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)*cdm, "cdm_"));
4123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4138afb7921SAidan Hamilton }
4148415c774SShri Abhyankar
4155f2c45f1SShri Abhyankar /*MC
4165f2c45f1SShri Abhyankar DMNETWORK = "network" - A DM object that encapsulates an unstructured network. The implementation is based on the DM object
4175f2c45f1SShri Abhyankar DMPlex that manages unstructured grids. Distributed networks use a non-overlapping partitioning of
4185f2c45f1SShri Abhyankar the edges. In the local representation, Vecs contain all unknowns in the interior and shared boundary.
4195f2c45f1SShri Abhyankar This is specified by a PetscSection object. Ownership in the global representation is determined by
4205f2c45f1SShri Abhyankar ownership of the underlying DMPlex points. This is specified by another PetscSection object.
4215f2c45f1SShri Abhyankar
4225f2c45f1SShri Abhyankar Level: intermediate
4235f2c45f1SShri Abhyankar
424db781477SPatrick Sanan .seealso: `DMType`, `DMNetworkCreate()`, `DMCreate()`, `DMSetType()`
4255f2c45f1SShri Abhyankar M*/
4265f2c45f1SShri Abhyankar
DMCreate_Network(DM dm)427d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Network(DM dm)
428d71ae5a4SJacob Faibussowitsch {
4295f2c45f1SShri Abhyankar DM_Network *network;
4305f2c45f1SShri Abhyankar
4315f2c45f1SShri Abhyankar PetscFunctionBegin;
4325f2c45f1SShri Abhyankar PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4334dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&network));
4344dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&network->cloneshared));
4355f2c45f1SShri Abhyankar dm->data = network;
4365f2c45f1SShri Abhyankar
437daad07d3SAidan Hamilton PetscCall(DMNetworkInitializeToDefault(dm));
4389566063dSJacob Faibussowitsch PetscCall(DMInitialize_Network(dm));
4393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4405f2c45f1SShri Abhyankar }
4415f2c45f1SShri Abhyankar
4425f2c45f1SShri Abhyankar /*@
4435f2c45f1SShri Abhyankar DMNetworkCreate - Creates a DMNetwork object, which encapsulates an unstructured network.
4445f2c45f1SShri Abhyankar
445d083f849SBarry Smith Collective
4465f2c45f1SShri Abhyankar
4475f2c45f1SShri Abhyankar Input Parameter:
4485f2c45f1SShri Abhyankar . comm - The communicator for the DMNetwork object
4495f2c45f1SShri Abhyankar
4505f2c45f1SShri Abhyankar Output Parameter:
4515f2c45f1SShri Abhyankar . network - The DMNetwork object
4525f2c45f1SShri Abhyankar
4535f2c45f1SShri Abhyankar Level: beginner
4545f2c45f1SShri Abhyankar
455a4e35b19SJacob Faibussowitsch .seealso: `DMCreate()`
4565f2c45f1SShri Abhyankar @*/
DMNetworkCreate(MPI_Comm comm,DM * network)457d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkCreate(MPI_Comm comm, DM *network)
458d71ae5a4SJacob Faibussowitsch {
4595f2c45f1SShri Abhyankar PetscFunctionBegin;
4604f572ea9SToby Isaac PetscAssertPointer(network, 2);
4619566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, network));
4629566063dSJacob Faibussowitsch PetscCall(DMSetType(*network, DMNETWORK));
4633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4645f2c45f1SShri Abhyankar }
465