1af0996ceSBarry Smith #include <petsc/private/dmimpl.h> 28ac4e037SJed Brown #include <petscdmredundant.h> /*I "petscdmredundant.h" I*/ 38ac4e037SJed Brown 48ac4e037SJed Brown typedef struct { 5907376e6SBarry Smith PetscMPIInt rank; /* owner */ 68ac4e037SJed Brown PetscInt N; /* total number of dofs */ 78ac4e037SJed Brown PetscInt n; /* owned number of dofs, n=N on owner, n=0 on non-owners */ 88ac4e037SJed Brown } DM_Redundant; 98ac4e037SJed Brown 109371c9d4SSatish Balay static PetscErrorCode DMCreateMatrix_Redundant(DM dm, Mat *J) { 118ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data; 1245b6f7e9SBarry Smith ISLocalToGlobalMapping ltog; 13d2e2b695SJed Brown PetscInt i, rstart, rend, *cols; 14d2e2b695SJed Brown PetscScalar *vals; 158ac4e037SJed Brown 168ac4e037SJed Brown PetscFunctionBegin; 179566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), J)); 189566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*J, red->n, red->n, red->N, red->N)); 199566063dSJacob Faibussowitsch PetscCall(MatSetType(*J, dm->mattype)); 209566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*J, red->n, NULL)); 219566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(*J, 1, red->n, NULL)); 229566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*J, red->n, NULL, red->N - red->n, NULL)); 239566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(*J, 1, red->n, NULL, red->N - red->n, NULL)); 248ac4e037SJed Brown 259566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dm, <og)); 269566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(*J, ltog, ltog)); 279566063dSJacob Faibussowitsch PetscCall(MatSetDM(*J, dm)); 28d2e2b695SJed Brown 299566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(red->N, &cols, red->N, &vals)); 30d2e2b695SJed Brown for (i = 0; i < red->N; i++) { 31d2e2b695SJed Brown cols[i] = i; 32d2e2b695SJed Brown vals[i] = 0.0; 33d2e2b695SJed Brown } 349566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(*J, &rstart, &rend)); 3548a46eb9SPierre Jolivet for (i = rstart; i < rend; i++) PetscCall(MatSetValues(*J, 1, &i, red->N, cols, vals, INSERT_VALUES)); 369566063dSJacob Faibussowitsch PetscCall(PetscFree2(cols, vals)); 379566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY)); 389566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY)); 398ac4e037SJed Brown PetscFunctionReturn(0); 408ac4e037SJed Brown } 418ac4e037SJed Brown 429371c9d4SSatish Balay static PetscErrorCode DMDestroy_Redundant(DM dm) { 438ac4e037SJed Brown PetscFunctionBegin; 449566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMRedundantSetSize_C", NULL)); 459566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMRedundantGetSize_C", NULL)); 469566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", NULL)); 47435a35e8SMatthew G Knepley /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */ 489566063dSJacob Faibussowitsch PetscCall(PetscFree(dm->data)); 498ac4e037SJed Brown PetscFunctionReturn(0); 508ac4e037SJed Brown } 518ac4e037SJed Brown 529371c9d4SSatish Balay static PetscErrorCode DMCreateGlobalVector_Redundant(DM dm, Vec *gvec) { 538ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data; 5445b6f7e9SBarry Smith ISLocalToGlobalMapping ltog; 558ac4e037SJed Brown 568ac4e037SJed Brown PetscFunctionBegin; 578ac4e037SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 588ac4e037SJed Brown PetscValidPointer(gvec, 2); 59ea78f98cSLisandro Dalcin *gvec = NULL; 609566063dSJacob Faibussowitsch PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), gvec)); 619566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*gvec, red->n, red->N)); 629566063dSJacob Faibussowitsch PetscCall(VecSetType(*gvec, dm->vectype)); 639566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dm, <og)); 649566063dSJacob Faibussowitsch PetscCall(VecSetLocalToGlobalMapping(*gvec, ltog)); 659566063dSJacob Faibussowitsch PetscCall(VecSetDM(*gvec, dm)); 668ac4e037SJed Brown PetscFunctionReturn(0); 678ac4e037SJed Brown } 688ac4e037SJed Brown 699371c9d4SSatish Balay static PetscErrorCode DMCreateLocalVector_Redundant(DM dm, Vec *lvec) { 708ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data; 718ac4e037SJed Brown 728ac4e037SJed Brown PetscFunctionBegin; 738ac4e037SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 748ac4e037SJed Brown PetscValidPointer(lvec, 2); 75ea78f98cSLisandro Dalcin *lvec = NULL; 769566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, lvec)); 779566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*lvec, red->N, red->N)); 789566063dSJacob Faibussowitsch PetscCall(VecSetType(*lvec, dm->vectype)); 799566063dSJacob Faibussowitsch PetscCall(VecSetDM(*lvec, dm)); 808ac4e037SJed Brown PetscFunctionReturn(0); 818ac4e037SJed Brown } 828ac4e037SJed Brown 839371c9d4SSatish Balay static PetscErrorCode DMLocalToGlobalBegin_Redundant(DM dm, Vec l, InsertMode imode, Vec g) { 848ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data; 858ac4e037SJed Brown const PetscScalar *lv; 868ac4e037SJed Brown PetscScalar *gv; 878ac4e037SJed Brown PetscMPIInt rank; 888ac4e037SJed Brown 898ac4e037SJed Brown PetscFunctionBegin; 909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 919566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(l, &lv)); 929566063dSJacob Faibussowitsch PetscCall(VecGetArray(g, &gv)); 938ac4e037SJed Brown switch (imode) { 948ac4e037SJed Brown case ADD_VALUES: 959371c9d4SSatish Balay case MAX_VALUES: { 968ac4e037SJed Brown void *source; 978ac4e037SJed Brown PetscScalar *buffer; 988ac4e037SJed Brown PetscInt i; 998ac4e037SJed Brown if (rank == red->rank) { 1008ac4e037SJed Brown buffer = gv; 1018ac4e037SJed Brown source = MPI_IN_PLACE; 1029371c9d4SSatish Balay if (imode == ADD_VALUES) 1039371c9d4SSatish Balay for (i = 0; i < red->N; i++) buffer[i] = gv[i] + lv[i]; 1048ac4e037SJed Brown #if !defined(PETSC_USE_COMPLEX) 1059371c9d4SSatish Balay if (imode == MAX_VALUES) 1069371c9d4SSatish Balay for (i = 0; i < red->N; i++) buffer[i] = PetscMax(gv[i], lv[i]); 1078ac4e037SJed Brown #endif 1088865f1eaSKarl Rupp } else source = (void *)lv; 1099566063dSJacob Faibussowitsch PetscCallMPI(MPI_Reduce(source, gv, red->N, MPIU_SCALAR, (imode == ADD_VALUES) ? MPIU_SUM : MPIU_MAX, red->rank, PetscObjectComm((PetscObject)dm))); 1108ac4e037SJed Brown } break; 1119371c9d4SSatish Balay case INSERT_VALUES: PetscCall(PetscArraycpy(gv, lv, red->n)); break; 112ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "InsertMode not supported"); 1138ac4e037SJed Brown } 1149566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(l, &lv)); 1159566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(g, &gv)); 1168ac4e037SJed Brown PetscFunctionReturn(0); 1178ac4e037SJed Brown } 1188ac4e037SJed Brown 1199371c9d4SSatish Balay static PetscErrorCode DMLocalToGlobalEnd_Redundant(DM dm, Vec l, InsertMode imode, Vec g) { 1208ac4e037SJed Brown PetscFunctionBegin; 1218ac4e037SJed Brown PetscFunctionReturn(0); 1228ac4e037SJed Brown } 1238ac4e037SJed Brown 1249371c9d4SSatish Balay static PetscErrorCode DMGlobalToLocalBegin_Redundant(DM dm, Vec g, InsertMode imode, Vec l) { 1258ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data; 1268ac4e037SJed Brown const PetscScalar *gv; 1278ac4e037SJed Brown PetscScalar *lv; 1288ac4e037SJed Brown 1298ac4e037SJed Brown PetscFunctionBegin; 1309566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(g, &gv)); 1319566063dSJacob Faibussowitsch PetscCall(VecGetArray(l, &lv)); 1328ac4e037SJed Brown switch (imode) { 1338ac4e037SJed Brown case INSERT_VALUES: 1349566063dSJacob Faibussowitsch if (red->n) PetscCall(PetscArraycpy(lv, gv, red->n)); 1359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(lv, red->N, MPIU_SCALAR, red->rank, PetscObjectComm((PetscObject)dm))); 1368ac4e037SJed Brown break; 137ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "InsertMode not supported"); 1388ac4e037SJed Brown } 1399566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(g, &gv)); 1409566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(l, &lv)); 1418ac4e037SJed Brown PetscFunctionReturn(0); 1428ac4e037SJed Brown } 1438ac4e037SJed Brown 1449371c9d4SSatish Balay static PetscErrorCode DMGlobalToLocalEnd_Redundant(DM dm, Vec g, InsertMode imode, Vec l) { 1458ac4e037SJed Brown PetscFunctionBegin; 1468ac4e037SJed Brown PetscFunctionReturn(0); 1478ac4e037SJed Brown } 1488ac4e037SJed Brown 1499371c9d4SSatish Balay static PetscErrorCode DMSetUp_Redundant(DM dm) { 1508ac4e037SJed Brown PetscFunctionBegin; 1518ac4e037SJed Brown PetscFunctionReturn(0); 1528ac4e037SJed Brown } 1538ac4e037SJed Brown 1549371c9d4SSatish Balay static PetscErrorCode DMView_Redundant(DM dm, PetscViewer viewer) { 1558ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data; 1568ac4e037SJed Brown PetscBool iascii; 1578ac4e037SJed Brown 1588ac4e037SJed Brown PetscFunctionBegin; 1599566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 16048a46eb9SPierre Jolivet if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, "redundant: rank=%d N=%" PetscInt_FMT "\n", red->rank, red->N)); 1618ac4e037SJed Brown PetscFunctionReturn(0); 1628ac4e037SJed Brown } 1638ac4e037SJed Brown 1649371c9d4SSatish Balay static PetscErrorCode DMCreateColoring_Redundant(DM dm, ISColoringType ctype, ISColoring *coloring) { 16569527181SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data; 16669527181SJed Brown PetscInt i, nloc; 16769527181SJed Brown ISColoringValue *colors; 16869527181SJed Brown 16969527181SJed Brown PetscFunctionBegin; 17069527181SJed Brown switch (ctype) { 1719371c9d4SSatish Balay case IS_COLORING_GLOBAL: nloc = red->n; break; 1729371c9d4SSatish Balay case IS_COLORING_LOCAL: nloc = red->N; break; 17398921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown ISColoringType %d", (int)ctype); 17469527181SJed Brown } 1759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &colors)); 17669527181SJed Brown for (i = 0; i < nloc; i++) colors[i] = i; 1779566063dSJacob Faibussowitsch PetscCall(ISColoringCreate(PetscObjectComm((PetscObject)dm), red->N, nloc, colors, PETSC_OWN_POINTER, coloring)); 1789566063dSJacob Faibussowitsch PetscCall(ISColoringSetType(*coloring, ctype)); 17969527181SJed Brown PetscFunctionReturn(0); 18069527181SJed Brown } 1818ac4e037SJed Brown 1829371c9d4SSatish Balay static PetscErrorCode DMRefine_Redundant(DM dmc, MPI_Comm comm, DM *dmf) { 1838ac4e037SJed Brown PetscMPIInt flag; 1848ac4e037SJed Brown DM_Redundant *redc = (DM_Redundant *)dmc->data; 1858ac4e037SJed Brown 1868ac4e037SJed Brown PetscFunctionBegin; 18748a46eb9SPierre Jolivet if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmc, &comm)); 1889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)dmc), comm, &flag)); 1891dca8a05SBarry Smith PetscCheck(flag == MPI_CONGRUENT || flag == MPI_IDENT, PetscObjectComm((PetscObject)dmc), PETSC_ERR_SUP, "cannot change communicators"); 1909566063dSJacob Faibussowitsch PetscCall(DMRedundantCreate(comm, redc->rank, redc->N, dmf)); 1918ac4e037SJed Brown PetscFunctionReturn(0); 1928ac4e037SJed Brown } 1938ac4e037SJed Brown 1949371c9d4SSatish Balay static PetscErrorCode DMCoarsen_Redundant(DM dmf, MPI_Comm comm, DM *dmc) { 1958ac4e037SJed Brown PetscMPIInt flag; 1968ac4e037SJed Brown DM_Redundant *redf = (DM_Redundant *)dmf->data; 1978ac4e037SJed Brown 1988ac4e037SJed Brown PetscFunctionBegin; 19948a46eb9SPierre Jolivet if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmf, &comm)); 2009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)dmf), comm, &flag)); 2011dca8a05SBarry Smith PetscCheck(flag == MPI_CONGRUENT || flag == MPI_IDENT, PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "cannot change communicators"); 2029566063dSJacob Faibussowitsch PetscCall(DMRedundantCreate(comm, redf->rank, redf->N, dmc)); 2038ac4e037SJed Brown PetscFunctionReturn(0); 2048ac4e037SJed Brown } 2058ac4e037SJed Brown 2069371c9d4SSatish Balay static PetscErrorCode DMCreateInterpolation_Redundant(DM dmc, DM dmf, Mat *P, Vec *scale) { 2078ac4e037SJed Brown DM_Redundant *redc = (DM_Redundant *)dmc->data; 2088ac4e037SJed Brown DM_Redundant *redf = (DM_Redundant *)dmf->data; 2098ac4e037SJed Brown PetscMPIInt flag; 2108ac4e037SJed Brown PetscInt i, rstart, rend; 2118ac4e037SJed Brown 2128ac4e037SJed Brown PetscFunctionBegin; 2139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)dmc), PetscObjectComm((PetscObject)dmf), &flag)); 2141dca8a05SBarry Smith PetscCheck(flag == MPI_CONGRUENT || flag == MPI_IDENT, PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "cannot change communicators"); 21508401ef6SPierre Jolivet PetscCheck(redc->rank == redf->rank, PetscObjectComm((PetscObject)dmf), PETSC_ERR_ARG_INCOMP, "Owning rank does not match"); 21608401ef6SPierre Jolivet PetscCheck(redc->N == redf->N, PetscObjectComm((PetscObject)dmf), PETSC_ERR_ARG_INCOMP, "Global size does not match"); 2179566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dmc), P)); 2189566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*P, redc->n, redc->n, redc->N, redc->N)); 2199566063dSJacob Faibussowitsch PetscCall(MatSetType(*P, MATAIJ)); 2209566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*P, 1, NULL)); 2219566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*P, 1, NULL, 0, NULL)); 2229566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(*P, &rstart, &rend)); 2239566063dSJacob Faibussowitsch for (i = rstart; i < rend; i++) PetscCall(MatSetValue(*P, i, i, 1.0, INSERT_VALUES)); 2249566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY)); 2259566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY)); 2269566063dSJacob Faibussowitsch if (scale) PetscCall(DMCreateInterpolationScale(dmc, dmf, *P, scale)); 2278ac4e037SJed Brown PetscFunctionReturn(0); 2288ac4e037SJed Brown } 2298ac4e037SJed Brown 2308ac4e037SJed Brown /*@ 2318ac4e037SJed Brown DMRedundantSetSize - Sets the size of a densely coupled redundant object 2328ac4e037SJed Brown 233d083f849SBarry Smith Collective on dm 2348ac4e037SJed Brown 235d8d19677SJose E. Roman Input Parameters: 2368ac4e037SJed Brown + dm - redundant DM 2378ac4e037SJed Brown . rank - rank of process to own redundant degrees of freedom 2388ac4e037SJed Brown - N - total number of redundant degrees of freedom 2398ac4e037SJed Brown 2408ac4e037SJed Brown Level: advanced 2418ac4e037SJed Brown 242db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCreateGlobalVector()`, `DMRedundantCreate()`, `DMRedundantGetSize()` 2438ac4e037SJed Brown @*/ 2449371c9d4SSatish Balay PetscErrorCode DMRedundantSetSize(DM dm, PetscMPIInt rank, PetscInt N) { 2458ac4e037SJed Brown PetscFunctionBegin; 2468ac4e037SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2478ac4e037SJed Brown PetscValidType(dm, 1); 24869fbec6eSBarry Smith PetscValidLogicalCollectiveMPIInt(dm, rank, 2); 2498ac4e037SJed Brown PetscValidLogicalCollectiveInt(dm, N, 3); 250cac4c232SBarry Smith PetscTryMethod(dm, "DMRedundantSetSize_C", (DM, PetscMPIInt, PetscInt), (dm, rank, N)); 2518ac4e037SJed Brown PetscFunctionReturn(0); 2528ac4e037SJed Brown } 2538ac4e037SJed Brown 2548ac4e037SJed Brown /*@ 2558ac4e037SJed Brown DMRedundantGetSize - Gets the size of a densely coupled redundant object 2568ac4e037SJed Brown 2578ac4e037SJed Brown Not Collective 2588ac4e037SJed Brown 2598ac4e037SJed Brown Input Parameter: 260a2b725a8SWilliam Gropp . dm - redundant DM 2618ac4e037SJed Brown 2628ac4e037SJed Brown Output Parameters: 2630298fd71SBarry Smith + rank - rank of process to own redundant degrees of freedom (or NULL) 2640298fd71SBarry Smith - N - total number of redundant degrees of freedom (or NULL) 2658ac4e037SJed Brown 2668ac4e037SJed Brown Level: advanced 2678ac4e037SJed Brown 268db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCreateGlobalVector()`, `DMRedundantCreate()`, `DMRedundantSetSize()` 2698ac4e037SJed Brown @*/ 2709371c9d4SSatish Balay PetscErrorCode DMRedundantGetSize(DM dm, PetscMPIInt *rank, PetscInt *N) { 2718ac4e037SJed Brown PetscFunctionBegin; 2728ac4e037SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2738ac4e037SJed Brown PetscValidType(dm, 1); 274cac4c232SBarry Smith PetscUseMethod(dm, "DMRedundantGetSize_C", (DM, PetscMPIInt *, PetscInt *), (dm, rank, N)); 2758ac4e037SJed Brown PetscFunctionReturn(0); 2768ac4e037SJed Brown } 2778ac4e037SJed Brown 2789371c9d4SSatish Balay static PetscErrorCode DMRedundantSetSize_Redundant(DM dm, PetscMPIInt rank, PetscInt N) { 2798ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data; 2808ac4e037SJed Brown PetscMPIInt myrank; 28173c9e547SStefano Zampini PetscInt i, *globals; 2828ac4e037SJed Brown 2838ac4e037SJed Brown PetscFunctionBegin; 2849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &myrank)); 2858ac4e037SJed Brown red->rank = rank; 2868ac4e037SJed Brown red->N = N; 2878ac4e037SJed Brown red->n = (myrank == rank) ? N : 0; 28873c9e547SStefano Zampini 28973c9e547SStefano Zampini /* mapping is setup here */ 2909566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(red->N, &globals)); 29173c9e547SStefano Zampini for (i = 0; i < red->N; i++) globals[i] = i; 2929566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&dm->ltogmap)); 2939566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), 1, red->N, globals, PETSC_OWN_POINTER, &dm->ltogmap)); 2948ac4e037SJed Brown PetscFunctionReturn(0); 2958ac4e037SJed Brown } 2968ac4e037SJed Brown 2979371c9d4SSatish Balay static PetscErrorCode DMRedundantGetSize_Redundant(DM dm, PetscInt *rank, PetscInt *N) { 2988ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data; 2998ac4e037SJed Brown 3008ac4e037SJed Brown PetscFunctionBegin; 3018ac4e037SJed Brown if (rank) *rank = red->rank; 3028ac4e037SJed Brown if (N) *N = red->N; 3038ac4e037SJed Brown PetscFunctionReturn(0); 3048ac4e037SJed Brown } 3058ac4e037SJed Brown 3069371c9d4SSatish Balay static PetscErrorCode DMSetUpGLVisViewer_Redundant(PetscObject odm, PetscViewer viewer) { 30736623e74SStefano Zampini PetscFunctionBegin; 30836623e74SStefano Zampini PetscFunctionReturn(0); 30936623e74SStefano Zampini } 31036623e74SStefano Zampini 3113efe6655SBarry Smith /*MC 3123efe6655SBarry Smith DMREDUNDANT = "redundant" - A DM object that is used to manage data for a small set of dense globally coupled variables. 3133efe6655SBarry Smith In the global representation of the vector the variables are all stored on a single MPI process (all the other MPI processes have 3143efe6655SBarry Smith no variables) in the local representation all the variables are stored on ALL the MPI processes (because they are all needed for each 3153efe6655SBarry Smith processes local computations). 3163efe6655SBarry Smith 3173efe6655SBarry Smith This DM is generally used inside a DMCOMPOSITE object. For example, it may be used to store continuation parameters for a bifurcation problem. 3183efe6655SBarry Smith 3193efe6655SBarry Smith Level: intermediate 3203efe6655SBarry Smith 321db781477SPatrick Sanan .seealso: `DMType`, `DMCOMPOSITE`, `DMCreate()`, `DMRedundantSetSize()`, `DMRedundantGetSize()` 3223efe6655SBarry Smith M*/ 3233efe6655SBarry Smith 3249371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMCreate_Redundant(DM dm) { 3258ac4e037SJed Brown DM_Redundant *red; 3268ac4e037SJed Brown 3278ac4e037SJed Brown PetscFunctionBegin; 328*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&red)); 3298ac4e037SJed Brown dm->data = red; 3308ac4e037SJed Brown 3318ac4e037SJed Brown dm->ops->setup = DMSetUp_Redundant; 3328ac4e037SJed Brown dm->ops->view = DMView_Redundant; 3338ac4e037SJed Brown dm->ops->createglobalvector = DMCreateGlobalVector_Redundant; 3348ac4e037SJed Brown dm->ops->createlocalvector = DMCreateLocalVector_Redundant; 33525296bd5SBarry Smith dm->ops->creatematrix = DMCreateMatrix_Redundant; 3368ac4e037SJed Brown dm->ops->destroy = DMDestroy_Redundant; 3378ac4e037SJed Brown dm->ops->globaltolocalbegin = DMGlobalToLocalBegin_Redundant; 3388ac4e037SJed Brown dm->ops->globaltolocalend = DMGlobalToLocalEnd_Redundant; 3398ac4e037SJed Brown dm->ops->localtoglobalbegin = DMLocalToGlobalBegin_Redundant; 3408ac4e037SJed Brown dm->ops->localtoglobalend = DMLocalToGlobalEnd_Redundant; 3418ac4e037SJed Brown dm->ops->refine = DMRefine_Redundant; 3428ac4e037SJed Brown dm->ops->coarsen = DMCoarsen_Redundant; 34325296bd5SBarry Smith dm->ops->createinterpolation = DMCreateInterpolation_Redundant; 344e727c939SJed Brown dm->ops->getcoloring = DMCreateColoring_Redundant; 3458865f1eaSKarl Rupp 3469566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMRedundantSetSize_C", DMRedundantSetSize_Redundant)); 3479566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMRedundantGetSize_C", DMRedundantGetSize_Redundant)); 3489566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_Redundant)); 3498ac4e037SJed Brown PetscFunctionReturn(0); 3508ac4e037SJed Brown } 3518ac4e037SJed Brown 3528ac4e037SJed Brown /*@C 3538ac4e037SJed Brown DMRedundantCreate - Creates a DM object, used to manage data for dense globally coupled variables 3548ac4e037SJed Brown 355d083f849SBarry Smith Collective 3568ac4e037SJed Brown 357d8d19677SJose E. Roman Input Parameters: 3588ac4e037SJed Brown + comm - the processors that will share the global vector 3598ac4e037SJed Brown . rank - rank to own the redundant values 3608ac4e037SJed Brown - N - total number of degrees of freedom 3618ac4e037SJed Brown 3628ac4e037SJed Brown Output Parameters: 363907376e6SBarry Smith . dm - the redundant DM 3648ac4e037SJed Brown 3658ac4e037SJed Brown Level: advanced 3668ac4e037SJed Brown 367db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCreateGlobalVector()`, `DMCreateMatrix()`, `DMCompositeAddDM()`, `DMREDUNDANT`, `DMSetType()`, `DMRedundantSetSize()`, `DMRedundantGetSize()` 3688ac4e037SJed Brown 3698ac4e037SJed Brown @*/ 3709371c9d4SSatish Balay PetscErrorCode DMRedundantCreate(MPI_Comm comm, PetscMPIInt rank, PetscInt N, DM *dm) { 3718ac4e037SJed Brown PetscFunctionBegin; 372064a246eSJacob Faibussowitsch PetscValidPointer(dm, 4); 3739566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 3749566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMREDUNDANT)); 3759566063dSJacob Faibussowitsch PetscCall(DMRedundantSetSize(*dm, rank, N)); 3769566063dSJacob Faibussowitsch PetscCall(DMSetUp(*dm)); 3778ac4e037SJed Brown PetscFunctionReturn(0); 3788ac4e037SJed Brown } 379