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 5d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetFromOptions_Network(DM dm, PetscOptionItems *PetscOptionsObject) 6d71ae5a4SJacob Faibussowitsch { 75f2c45f1SShri Abhyankar PetscFunctionBegin; 8d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "DMNetwork Options"); 9d0609cedSBarry Smith PetscOptionsHeadEnd(); 105f2c45f1SShri Abhyankar PetscFunctionReturn(0); 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 *); 238afb7921SAidan Hamilton extern PetscErrorCode DMCreateCoordinateDM_Network(DM, DM *); 245f2c45f1SShri Abhyankar 25d71ae5a4SJacob Faibussowitsch static PetscErrorCode VecArrayPrint_private(PetscViewer viewer, PetscInt n, const PetscScalar *xv) 26d71ae5a4SJacob Faibussowitsch { 27a81b7fe4SHong Zhang PetscInt i; 28a81b7fe4SHong Zhang 29a81b7fe4SHong Zhang PetscFunctionBegin; 30a81b7fe4SHong Zhang for (i = 0; i < n; i++) { 31a81b7fe4SHong Zhang #if defined(PETSC_USE_COMPLEX) 32a81b7fe4SHong Zhang if (PetscImaginaryPart(xv[i]) > 0.0) { 339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " %g + %g i\n", (double)PetscRealPart(xv[i]), (double)PetscImaginaryPart(xv[i]))); 34a81b7fe4SHong Zhang } else if (PetscImaginaryPart(xv[i]) < 0.0) { 359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " %g - %g i\n", (double)PetscRealPart(xv[i]), -(double)PetscImaginaryPart(xv[i]))); 361baa6e33SBarry Smith } else PetscCall(PetscViewerASCIIPrintf(viewer, " %g\n", (double)PetscRealPart(xv[i]))); 37a81b7fe4SHong Zhang #else 389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " %g\n", (double)xv[i])); 39a81b7fe4SHong Zhang #endif 40a81b7fe4SHong Zhang } 41a81b7fe4SHong Zhang PetscFunctionReturn(0); 42a81b7fe4SHong Zhang } 43a81b7fe4SHong Zhang 44d71ae5a4SJacob Faibussowitsch static PetscErrorCode VecView_Network_Seq(DM networkdm, Vec X, PetscViewer viewer) 45d71ae5a4SJacob Faibussowitsch { 467b6afd5bSHong Zhang PetscInt e, v, Start, End, offset, nvar, id; 474dc485aaSHong Zhang const PetscScalar *xv; 484dc485aaSHong Zhang 494dc485aaSHong Zhang PetscFunctionBegin; 509566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &xv)); 514dc485aaSHong Zhang 524dc485aaSHong Zhang /* iterate over edges */ 539566063dSJacob Faibussowitsch PetscCall(DMNetworkGetEdgeRange(networkdm, &Start, &End)); 544dc485aaSHong Zhang for (e = Start; e < End; e++) { 559566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, e, ALL_COMPONENTS, NULL, NULL, &nvar)); 564dc485aaSHong Zhang if (!nvar) continue; 57a81b7fe4SHong Zhang 589566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, e, ALL_COMPONENTS, &offset)); 599566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalEdgeIndex(networkdm, e, &id)); 607b6afd5bSHong Zhang 619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Edge %" PetscInt_FMT ":\n", id)); 629566063dSJacob Faibussowitsch PetscCall(VecArrayPrint_private(viewer, nvar, xv + offset)); 634dc485aaSHong Zhang } 644dc485aaSHong Zhang 654dc485aaSHong Zhang /* iterate over vertices */ 669566063dSJacob Faibussowitsch PetscCall(DMNetworkGetVertexRange(networkdm, &Start, &End)); 674dc485aaSHong Zhang for (v = Start; v < End; v++) { 689566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, v, ALL_COMPONENTS, NULL, NULL, &nvar)); 694dc485aaSHong Zhang if (!nvar) continue; 70a81b7fe4SHong Zhang 719566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, ALL_COMPONENTS, &offset)); 729566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVertexIndex(networkdm, v, &id)); 7362aa33baSHong Zhang 749566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Vertex %" PetscInt_FMT ":\n", id)); 759566063dSJacob Faibussowitsch PetscCall(VecArrayPrint_private(viewer, nvar, xv + offset)); 764dc485aaSHong Zhang } 779566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 789566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &xv)); 794dc485aaSHong Zhang PetscFunctionReturn(0); 804dc485aaSHong Zhang } 814dc485aaSHong Zhang 82d71ae5a4SJacob Faibussowitsch static PetscErrorCode VecView_Network_MPI(DM networkdm, Vec X, PetscViewer viewer) 83d71ae5a4SJacob Faibussowitsch { 84e8c2e809SHong Zhang PetscInt i, e, v, eStart, eEnd, vStart, vEnd, offset, nvar, len_loc, len, k; 854062a5e5SHong Zhang const PetscScalar *xv; 864062a5e5SHong Zhang MPI_Comm comm; 874062a5e5SHong Zhang PetscMPIInt size, rank, tag = ((PetscObject)viewer)->tag; 884062a5e5SHong Zhang Vec localX; 894062a5e5SHong Zhang PetscBool ghostvtex; 904062a5e5SHong Zhang PetscScalar *values; 917b6afd5bSHong Zhang PetscInt j, ne, nv, id; 924062a5e5SHong Zhang MPI_Status status; 934062a5e5SHong Zhang 944062a5e5SHong Zhang PetscFunctionBegin; 959566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)networkdm, &comm)); 969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 984062a5e5SHong Zhang 999566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(networkdm, &localX)); 1009566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX)); 1019566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX)); 1029566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(localX, &xv)); 1034062a5e5SHong Zhang 1049566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(localX, &len_loc)); 1054062a5e5SHong Zhang 1069566063dSJacob Faibussowitsch PetscCall(DMNetworkGetEdgeRange(networkdm, &eStart, &eEnd)); 1079566063dSJacob Faibussowitsch PetscCall(DMNetworkGetVertexRange(networkdm, &vStart, &vEnd)); 108e8c2e809SHong Zhang len_loc += 2 * (1 + eEnd - eStart + vEnd - vStart); 1094062a5e5SHong Zhang 110e85e6aecSHong Zhang /* values = [nedges, nvertices; id, nvar, xedge; ...; id, nvars, xvertex;...], to be sent to proc[0] */ 1119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&len_loc, &len, 1, MPIU_INT, MPI_MAX, comm)); 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 */ 1594062a5e5SHong Zhang for (j = 1; j < size; j++) { 16063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Process [%" PetscInt_FMT "]\n", j)); 1614062a5e5SHong Zhang 1629566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(values, (PetscMPIInt)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] */ 1889566063dSJacob Faibussowitsch PetscCallMPI(MPI_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)); 1944062a5e5SHong Zhang PetscFunctionReturn(0); 1954062a5e5SHong Zhang } 1964062a5e5SHong Zhang 19795af8c53SHong Zhang PETSC_EXTERN PetscErrorCode VecView_MPI(Vec, PetscViewer); 198a6c4b7b7SHong Zhang 199d71ae5a4SJacob Faibussowitsch PetscErrorCode VecView_Network(Vec v, PetscViewer viewer) 200d71ae5a4SJacob Faibussowitsch { 2014dc485aaSHong Zhang DM dm; 2024dc485aaSHong Zhang PetscBool isseq; 203a6c4b7b7SHong Zhang PetscBool iascii; 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"); 2089566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 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 */ 212a6c4b7b7SHong Zhang if (iascii) { 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 } 2194dc485aaSHong Zhang PetscFunctionReturn(0); 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)); 2305f2c45f1SShri Abhyankar PetscFunctionReturn(0); 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)); 2405f2c45f1SShri Abhyankar PetscFunctionReturn(0); 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 */ 279daad07d3SAidan Hamilton PetscFunctionReturn(0); 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; 314daad07d3SAidan Hamilton PetscFunctionReturn(0); 315daad07d3SAidan Hamilton } 316daad07d3SAidan Hamilton 317d71ae5a4SJacob Faibussowitsch 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; 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; 3455f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3465f2c45f1SShri Abhyankar } 347daad07d3SAidan Hamilton /* 348daad07d3SAidan Hamilton copies over the subnetid and index portions of the DMNetworkComponentHeader from orignal dm to the newdm 349daad07d3SAidan Hamilton */ 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 } 373daad07d3SAidan Hamilton PetscFunctionReturn(0); 374daad07d3SAidan Hamilton } 3755f2c45f1SShri Abhyankar 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 */ 386daad07d3SAidan Hamilton PetscCall(DMClone(network->plex, &newnetwork->plex)); 387daad07d3SAidan Hamilton PetscCall(DMNetworkCopyHeaderTopological(dm, *newdm)); 388daad07d3SAidan Hamilton PetscCall(DMNetworkInitializeNonTopological(*newdm)); /* initialize all non-topological data to the state after DMNetworkLayoutSetUp as been called */ 3899566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, DMNETWORK)); 3909566063dSJacob Faibussowitsch PetscCall(DMInitialize_Network(*newdm)); 3918415c774SShri Abhyankar PetscFunctionReturn(0); 3928415c774SShri Abhyankar } 3938afb7921SAidan Hamilton /* Developer Note: Be aware that the plex inside of the network does not have a coordinate plex. 3948afb7921SAidan Hamilton */ 3958afb7921SAidan Hamilton PetscErrorCode DMCreateCoordinateDM_Network(DM dm, DM *cdm) 3968afb7921SAidan Hamilton { 3978afb7921SAidan Hamilton DM_Network *newnetwork = NULL; 3988afb7921SAidan Hamilton PetscInt Nf; 399*dd4c3f67SMatthew G. Knepley const char *prefix; 4008afb7921SAidan Hamilton 4018afb7921SAidan Hamilton PetscFunctionBegin; 4028afb7921SAidan Hamilton PetscCall(DMClone(dm, cdm)); 4038afb7921SAidan Hamilton newnetwork = (DM_Network *)(*cdm)->data; 4048afb7921SAidan Hamilton PetscCall(DMGetNumFields(newnetwork->plex, &Nf)); 4058afb7921SAidan Hamilton PetscCall(DMSetNumFields(*cdm, Nf)); /* consistency with the coordinate plex */ 406*dd4c3f67SMatthew G. Knepley PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 407*dd4c3f67SMatthew G. Knepley PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*cdm, prefix)); 408*dd4c3f67SMatthew G. Knepley PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)*cdm, "cdm_")); 4098afb7921SAidan Hamilton PetscFunctionReturn(0); 4108afb7921SAidan Hamilton } 4118415c774SShri Abhyankar 4125f2c45f1SShri Abhyankar /*MC 4135f2c45f1SShri Abhyankar DMNETWORK = "network" - A DM object that encapsulates an unstructured network. The implementation is based on the DM object 4145f2c45f1SShri Abhyankar DMPlex that manages unstructured grids. Distributed networks use a non-overlapping partitioning of 4155f2c45f1SShri Abhyankar the edges. In the local representation, Vecs contain all unknowns in the interior and shared boundary. 4165f2c45f1SShri Abhyankar This is specified by a PetscSection object. Ownership in the global representation is determined by 4175f2c45f1SShri Abhyankar ownership of the underlying DMPlex points. This is specified by another PetscSection object. 4185f2c45f1SShri Abhyankar 4195f2c45f1SShri Abhyankar Level: intermediate 4205f2c45f1SShri Abhyankar 421db781477SPatrick Sanan .seealso: `DMType`, `DMNetworkCreate()`, `DMCreate()`, `DMSetType()` 4225f2c45f1SShri Abhyankar M*/ 4235f2c45f1SShri Abhyankar 424d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Network(DM dm) 425d71ae5a4SJacob Faibussowitsch { 4265f2c45f1SShri Abhyankar DM_Network *network; 4275f2c45f1SShri Abhyankar 4285f2c45f1SShri Abhyankar PetscFunctionBegin; 4295f2c45f1SShri Abhyankar PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4304dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&network)); 4314dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&network->cloneshared)); 4325f2c45f1SShri Abhyankar dm->data = network; 4335f2c45f1SShri Abhyankar 434daad07d3SAidan Hamilton PetscCall(DMNetworkInitializeToDefault(dm)); 4359566063dSJacob Faibussowitsch PetscCall(DMInitialize_Network(dm)); 4365f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4375f2c45f1SShri Abhyankar } 4385f2c45f1SShri Abhyankar 4395f2c45f1SShri Abhyankar /*@ 4405f2c45f1SShri Abhyankar DMNetworkCreate - Creates a DMNetwork object, which encapsulates an unstructured network. 4415f2c45f1SShri Abhyankar 442d083f849SBarry Smith Collective 4435f2c45f1SShri Abhyankar 4445f2c45f1SShri Abhyankar Input Parameter: 4455f2c45f1SShri Abhyankar . comm - The communicator for the DMNetwork object 4465f2c45f1SShri Abhyankar 4475f2c45f1SShri Abhyankar Output Parameter: 4485f2c45f1SShri Abhyankar . network - The DMNetwork object 4495f2c45f1SShri Abhyankar 4505f2c45f1SShri Abhyankar Level: beginner 4515f2c45f1SShri Abhyankar 4525f2c45f1SShri Abhyankar @*/ 453d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkCreate(MPI_Comm comm, DM *network) 454d71ae5a4SJacob Faibussowitsch { 4555f2c45f1SShri Abhyankar PetscFunctionBegin; 4565f2c45f1SShri Abhyankar PetscValidPointer(network, 2); 4579566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, network)); 4589566063dSJacob Faibussowitsch PetscCall(DMSetType(*network, DMNETWORK)); 4595f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4605f2c45f1SShri Abhyankar } 461