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 5*a4e35b19SJacob Faibussowitsch 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 { 83e8c2e809SHong Zhang PetscInt i, e, v, eStart, eEnd, vStart, vEnd, offset, nvar, len_loc, len, k; 844062a5e5SHong Zhang const PetscScalar *xv; 854062a5e5SHong Zhang MPI_Comm comm; 864062a5e5SHong Zhang PetscMPIInt size, rank, tag = ((PetscObject)viewer)->tag; 874062a5e5SHong Zhang Vec localX; 884062a5e5SHong Zhang PetscBool ghostvtex; 894062a5e5SHong Zhang PetscScalar *values; 907b6afd5bSHong Zhang PetscInt j, 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] */ 110712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(&len_loc, &len, 1, MPIU_INT, MPI_MAX, comm)); 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 */ 1584062a5e5SHong Zhang for (j = 1; j < size; j++) { 15963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Process [%" PetscInt_FMT "]\n", j)); 1604062a5e5SHong Zhang 1619566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(values, (PetscMPIInt)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] */ 1879566063dSJacob Faibussowitsch PetscCallMPI(MPI_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 19695af8c53SHong Zhang PETSC_EXTERN PetscErrorCode VecView_MPI(Vec, PetscViewer); 197a6c4b7b7SHong Zhang 198*a4e35b19SJacob Faibussowitsch static PetscErrorCode VecView_Network(Vec v, PetscViewer viewer) 199d71ae5a4SJacob Faibussowitsch { 2004dc485aaSHong Zhang DM dm; 2014dc485aaSHong Zhang PetscBool isseq; 202a6c4b7b7SHong Zhang PetscBool iascii; 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"); 2079566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 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 */ 211a6c4b7b7SHong Zhang if (iascii) { 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 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)); 2279566063dSJacob Faibussowitsch PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Network)); 2289566063dSJacob Faibussowitsch PetscCall(VecSetDM(*vec, dm)); 2293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2305f2c45f1SShri Abhyankar } 2315f2c45f1SShri Abhyankar 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 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 */ 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 316*a4e35b19SJacob 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; 329ea78f98cSLisandro Dalcin dm->ops->getcoloring = NULL; 3305f2c45f1SShri Abhyankar dm->ops->creatematrix = DMCreateMatrix_Network; 331ea78f98cSLisandro Dalcin dm->ops->createinterpolation = NULL; 332ea78f98cSLisandro Dalcin dm->ops->createinjection = NULL; 333ea78f98cSLisandro Dalcin dm->ops->refine = NULL; 334ea78f98cSLisandro Dalcin dm->ops->coarsen = NULL; 335ea78f98cSLisandro Dalcin dm->ops->refinehierarchy = NULL; 336ea78f98cSLisandro Dalcin dm->ops->coarsenhierarchy = NULL; 3375f2c45f1SShri Abhyankar dm->ops->globaltolocalbegin = DMGlobalToLocalBegin_Network; 3385f2c45f1SShri Abhyankar dm->ops->globaltolocalend = DMGlobalToLocalEnd_Network; 3395f2c45f1SShri Abhyankar dm->ops->localtoglobalbegin = DMLocalToGlobalBegin_Network; 3405f2c45f1SShri Abhyankar dm->ops->localtoglobalend = DMLocalToGlobalEnd_Network; 3415f2c45f1SShri Abhyankar dm->ops->destroy = DMDestroy_Network; 3425f2c45f1SShri Abhyankar dm->ops->createsubdm = NULL; 3435f2c45f1SShri Abhyankar dm->ops->locatepoints = NULL; 3443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3455f2c45f1SShri Abhyankar } 346daad07d3SAidan Hamilton /* 347da81f932SPierre Jolivet copies over the subnetid and index portions of the DMNetworkComponentHeader from original dm to the newdm 348daad07d3SAidan Hamilton */ 349d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMNetworkCopyHeaderTopological(DM dm, DM newdm) 350d71ae5a4SJacob Faibussowitsch { 351daad07d3SAidan Hamilton DM_Network *network = (DM_Network *)dm->data, *newnetwork = (DM_Network *)newdm->data; 352daad07d3SAidan Hamilton PetscInt p, i, np, index, subnetid; 353daad07d3SAidan Hamilton 354daad07d3SAidan Hamilton PetscFunctionBegin; 355daad07d3SAidan Hamilton np = network->cloneshared->pEnd - network->cloneshared->pStart; 356daad07d3SAidan Hamilton PetscCall(PetscCalloc2(np, &newnetwork->header, np, &newnetwork->cvalue)); 357daad07d3SAidan Hamilton for (i = 0; i < np; i++) { 358daad07d3SAidan Hamilton p = i + network->cloneshared->pStart; 359daad07d3SAidan Hamilton PetscCall(DMNetworkGetSubnetID(dm, p, &subnetid)); 360daad07d3SAidan Hamilton PetscCall(DMNetworkGetIndex(dm, p, &index)); 361daad07d3SAidan Hamilton newnetwork->header[i].index = index; 362daad07d3SAidan Hamilton newnetwork->header[i].subnetid = subnetid; 363daad07d3SAidan Hamilton newnetwork->header[i].size = NULL; 364daad07d3SAidan Hamilton newnetwork->header[i].key = NULL; 365daad07d3SAidan Hamilton newnetwork->header[i].offset = NULL; 366daad07d3SAidan Hamilton newnetwork->header[i].nvar = NULL; 367daad07d3SAidan Hamilton newnetwork->header[i].offsetvarrel = NULL; 368daad07d3SAidan Hamilton newnetwork->header[i].ndata = 0; 369daad07d3SAidan Hamilton newnetwork->header[i].maxcomps = DMNETWORK_MAX_COMP_AT_POINT_DEFAULT; 370daad07d3SAidan Hamilton newnetwork->header[i].hsize = sizeof(struct _p_DMNetworkComponentHeader) / sizeof(sizeof(DMNetworkComponentGenericDataType)); 371daad07d3SAidan Hamilton } 3723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 373daad07d3SAidan Hamilton } 3745f2c45f1SShri Abhyankar 375d71ae5a4SJacob Faibussowitsch PetscErrorCode DMClone_Network(DM dm, DM *newdm) 376d71ae5a4SJacob Faibussowitsch { 377daad07d3SAidan Hamilton DM_Network *network = (DM_Network *)dm->data, *newnetwork = NULL; 3788415c774SShri Abhyankar 3798415c774SShri Abhyankar PetscFunctionBegin; 380daad07d3SAidan Hamilton network->cloneshared->refct++; 3814dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&newnetwork)); 382daad07d3SAidan Hamilton (*newdm)->data = newnetwork; 383daad07d3SAidan Hamilton PetscCall(DMNetworkInitializeToDefault_NonShared(*newdm)); 384daad07d3SAidan Hamilton newnetwork->cloneshared = network->cloneshared; /* Share all data that can be cloneshared */ 3859b626cebSHong Zhang 3869b626cebSHong Zhang PetscCheck(network->plex, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_NULL, "Must call DMNetworkLayoutSetUp() first"); 387daad07d3SAidan Hamilton PetscCall(DMClone(network->plex, &newnetwork->plex)); 388daad07d3SAidan Hamilton PetscCall(DMNetworkCopyHeaderTopological(dm, *newdm)); 389daad07d3SAidan Hamilton PetscCall(DMNetworkInitializeNonTopological(*newdm)); /* initialize all non-topological data to the state after DMNetworkLayoutSetUp as been called */ 3909566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, DMNETWORK)); 3919566063dSJacob Faibussowitsch PetscCall(DMInitialize_Network(*newdm)); 3923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3938415c774SShri Abhyankar } 3949b626cebSHong Zhang 3958afb7921SAidan Hamilton /* Developer Note: Be aware that the plex inside of the network does not have a coordinate plex. 3968afb7921SAidan Hamilton */ 3978afb7921SAidan Hamilton PetscErrorCode DMCreateCoordinateDM_Network(DM dm, DM *cdm) 3988afb7921SAidan Hamilton { 3998afb7921SAidan Hamilton DM_Network *newnetwork = NULL; 4008afb7921SAidan Hamilton PetscInt Nf; 401dd4c3f67SMatthew G. Knepley const char *prefix; 4028afb7921SAidan Hamilton 4038afb7921SAidan Hamilton PetscFunctionBegin; 4048afb7921SAidan Hamilton PetscCall(DMClone(dm, cdm)); 4058afb7921SAidan Hamilton newnetwork = (DM_Network *)(*cdm)->data; 4068afb7921SAidan Hamilton PetscCall(DMGetNumFields(newnetwork->plex, &Nf)); 4078afb7921SAidan Hamilton PetscCall(DMSetNumFields(*cdm, Nf)); /* consistency with the coordinate plex */ 408dd4c3f67SMatthew G. Knepley PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 409dd4c3f67SMatthew G. Knepley PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*cdm, prefix)); 410dd4c3f67SMatthew G. Knepley PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)*cdm, "cdm_")); 4113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4128afb7921SAidan Hamilton } 4138415c774SShri Abhyankar 4145f2c45f1SShri Abhyankar /*MC 4155f2c45f1SShri Abhyankar DMNETWORK = "network" - A DM object that encapsulates an unstructured network. The implementation is based on the DM object 4165f2c45f1SShri Abhyankar DMPlex that manages unstructured grids. Distributed networks use a non-overlapping partitioning of 4175f2c45f1SShri Abhyankar the edges. In the local representation, Vecs contain all unknowns in the interior and shared boundary. 4185f2c45f1SShri Abhyankar This is specified by a PetscSection object. Ownership in the global representation is determined by 4195f2c45f1SShri Abhyankar ownership of the underlying DMPlex points. This is specified by another PetscSection object. 4205f2c45f1SShri Abhyankar 4215f2c45f1SShri Abhyankar Level: intermediate 4225f2c45f1SShri Abhyankar 423db781477SPatrick Sanan .seealso: `DMType`, `DMNetworkCreate()`, `DMCreate()`, `DMSetType()` 4245f2c45f1SShri Abhyankar M*/ 4255f2c45f1SShri Abhyankar 426d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Network(DM dm) 427d71ae5a4SJacob Faibussowitsch { 4285f2c45f1SShri Abhyankar DM_Network *network; 4295f2c45f1SShri Abhyankar 4305f2c45f1SShri Abhyankar PetscFunctionBegin; 4315f2c45f1SShri Abhyankar PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4324dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&network)); 4334dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&network->cloneshared)); 4345f2c45f1SShri Abhyankar dm->data = network; 4355f2c45f1SShri Abhyankar 436daad07d3SAidan Hamilton PetscCall(DMNetworkInitializeToDefault(dm)); 4379566063dSJacob Faibussowitsch PetscCall(DMInitialize_Network(dm)); 4383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4395f2c45f1SShri Abhyankar } 4405f2c45f1SShri Abhyankar 4415f2c45f1SShri Abhyankar /*@ 4425f2c45f1SShri Abhyankar DMNetworkCreate - Creates a DMNetwork object, which encapsulates an unstructured network. 4435f2c45f1SShri Abhyankar 444d083f849SBarry Smith Collective 4455f2c45f1SShri Abhyankar 4465f2c45f1SShri Abhyankar Input Parameter: 4475f2c45f1SShri Abhyankar . comm - The communicator for the DMNetwork object 4485f2c45f1SShri Abhyankar 4495f2c45f1SShri Abhyankar Output Parameter: 4505f2c45f1SShri Abhyankar . network - The DMNetwork object 4515f2c45f1SShri Abhyankar 4525f2c45f1SShri Abhyankar Level: beginner 4535f2c45f1SShri Abhyankar 454*a4e35b19SJacob Faibussowitsch .seealso: `DMCreate()` 4555f2c45f1SShri Abhyankar @*/ 456d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkCreate(MPI_Comm comm, DM *network) 457d71ae5a4SJacob Faibussowitsch { 4585f2c45f1SShri Abhyankar PetscFunctionBegin; 4594f572ea9SToby Isaac PetscAssertPointer(network, 2); 4609566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, network)); 4619566063dSJacob Faibussowitsch PetscCall(DMSetType(*network, DMNETWORK)); 4623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4635f2c45f1SShri Abhyankar } 464