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 10d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateMatrix_Redundant(DM dm, Mat *J) 11d71ae5a4SJacob Faibussowitsch { 128ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data; 1345b6f7e9SBarry Smith ISLocalToGlobalMapping ltog; 14d2e2b695SJed Brown PetscInt i, rstart, rend, *cols; 15d2e2b695SJed Brown PetscScalar *vals; 168ac4e037SJed Brown 178ac4e037SJed Brown PetscFunctionBegin; 189566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), J)); 199566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*J, red->n, red->n, red->N, red->N)); 209566063dSJacob Faibussowitsch PetscCall(MatSetType(*J, dm->mattype)); 219566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*J, red->n, NULL)); 229566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(*J, 1, red->n, NULL)); 239566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*J, red->n, NULL, red->N - red->n, NULL)); 249566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(*J, 1, red->n, NULL, red->N - red->n, NULL)); 258ac4e037SJed Brown 269566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dm, <og)); 279566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(*J, ltog, ltog)); 289566063dSJacob Faibussowitsch PetscCall(MatSetDM(*J, dm)); 29d2e2b695SJed Brown 309566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(red->N, &cols, red->N, &vals)); 31d2e2b695SJed Brown for (i = 0; i < red->N; i++) { 32d2e2b695SJed Brown cols[i] = i; 33d2e2b695SJed Brown vals[i] = 0.0; 34d2e2b695SJed Brown } 359566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(*J, &rstart, &rend)); 3648a46eb9SPierre Jolivet for (i = rstart; i < rend; i++) PetscCall(MatSetValues(*J, 1, &i, red->N, cols, vals, INSERT_VALUES)); 379566063dSJacob Faibussowitsch PetscCall(PetscFree2(cols, vals)); 389566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY)); 399566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY)); 403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 418ac4e037SJed Brown } 428ac4e037SJed Brown 43d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDestroy_Redundant(DM dm) 44d71ae5a4SJacob Faibussowitsch { 458ac4e037SJed Brown PetscFunctionBegin; 469566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMRedundantSetSize_C", NULL)); 479566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMRedundantGetSize_C", NULL)); 489566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", NULL)); 49435a35e8SMatthew G Knepley /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */ 509566063dSJacob Faibussowitsch PetscCall(PetscFree(dm->data)); 513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 528ac4e037SJed Brown } 538ac4e037SJed Brown 54d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateGlobalVector_Redundant(DM dm, Vec *gvec) 55d71ae5a4SJacob Faibussowitsch { 568ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data; 5745b6f7e9SBarry Smith ISLocalToGlobalMapping ltog; 588ac4e037SJed Brown 598ac4e037SJed Brown PetscFunctionBegin; 608ac4e037SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 614f572ea9SToby Isaac PetscAssertPointer(gvec, 2); 62ea78f98cSLisandro Dalcin *gvec = NULL; 639566063dSJacob Faibussowitsch PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), gvec)); 649566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*gvec, red->n, red->N)); 659566063dSJacob Faibussowitsch PetscCall(VecSetType(*gvec, dm->vectype)); 669566063dSJacob Faibussowitsch PetscCall(DMGetLocalToGlobalMapping(dm, <og)); 679566063dSJacob Faibussowitsch PetscCall(VecSetLocalToGlobalMapping(*gvec, ltog)); 689566063dSJacob Faibussowitsch PetscCall(VecSetDM(*gvec, dm)); 693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 708ac4e037SJed Brown } 718ac4e037SJed Brown 72d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateLocalVector_Redundant(DM dm, Vec *lvec) 73d71ae5a4SJacob Faibussowitsch { 748ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data; 758ac4e037SJed Brown 768ac4e037SJed Brown PetscFunctionBegin; 778ac4e037SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 784f572ea9SToby Isaac PetscAssertPointer(lvec, 2); 79ea78f98cSLisandro Dalcin *lvec = NULL; 809566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, lvec)); 819566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*lvec, red->N, red->N)); 829566063dSJacob Faibussowitsch PetscCall(VecSetType(*lvec, dm->vectype)); 839566063dSJacob Faibussowitsch PetscCall(VecSetDM(*lvec, dm)); 843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 858ac4e037SJed Brown } 868ac4e037SJed Brown 87d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMLocalToGlobalBegin_Redundant(DM dm, Vec l, InsertMode imode, Vec g) 88d71ae5a4SJacob Faibussowitsch { 898ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data; 908ac4e037SJed Brown const PetscScalar *lv; 918ac4e037SJed Brown PetscScalar *gv; 928ac4e037SJed Brown PetscMPIInt rank; 938ac4e037SJed Brown 948ac4e037SJed Brown PetscFunctionBegin; 959566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 969566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(l, &lv)); 979566063dSJacob Faibussowitsch PetscCall(VecGetArray(g, &gv)); 988ac4e037SJed Brown switch (imode) { 998ac4e037SJed Brown case ADD_VALUES: 1009371c9d4SSatish Balay case MAX_VALUES: { 1018ac4e037SJed Brown void *source; 1028ac4e037SJed Brown PetscScalar *buffer; 1038ac4e037SJed Brown PetscInt i; 1048ac4e037SJed Brown if (rank == red->rank) { 1058ac4e037SJed Brown buffer = gv; 1068ac4e037SJed Brown source = MPI_IN_PLACE; 1079371c9d4SSatish Balay if (imode == ADD_VALUES) 1089371c9d4SSatish Balay for (i = 0; i < red->N; i++) buffer[i] = gv[i] + lv[i]; 1098ac4e037SJed Brown #if !defined(PETSC_USE_COMPLEX) 1109371c9d4SSatish Balay if (imode == MAX_VALUES) 1119371c9d4SSatish Balay for (i = 0; i < red->N; i++) buffer[i] = PetscMax(gv[i], lv[i]); 1128ac4e037SJed Brown #endif 1138865f1eaSKarl Rupp } else source = (void *)lv; 1149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Reduce(source, gv, red->N, MPIU_SCALAR, (imode == ADD_VALUES) ? MPIU_SUM : MPIU_MAX, red->rank, PetscObjectComm((PetscObject)dm))); 1158ac4e037SJed Brown } break; 116d71ae5a4SJacob Faibussowitsch case INSERT_VALUES: 117d71ae5a4SJacob Faibussowitsch PetscCall(PetscArraycpy(gv, lv, red->n)); 118d71ae5a4SJacob Faibussowitsch break; 119d71ae5a4SJacob Faibussowitsch default: 120d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "InsertMode not supported"); 1218ac4e037SJed Brown } 1229566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(l, &lv)); 1239566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(g, &gv)); 1243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1258ac4e037SJed Brown } 1268ac4e037SJed Brown 127d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMLocalToGlobalEnd_Redundant(DM dm, Vec l, InsertMode imode, Vec g) 128d71ae5a4SJacob Faibussowitsch { 1298ac4e037SJed Brown PetscFunctionBegin; 1303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1318ac4e037SJed Brown } 1328ac4e037SJed Brown 133d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalBegin_Redundant(DM dm, Vec g, InsertMode imode, Vec l) 134d71ae5a4SJacob Faibussowitsch { 1358ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data; 1368ac4e037SJed Brown const PetscScalar *gv; 1378ac4e037SJed Brown PetscScalar *lv; 1388ac4e037SJed Brown 1398ac4e037SJed Brown PetscFunctionBegin; 1409566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(g, &gv)); 1419566063dSJacob Faibussowitsch PetscCall(VecGetArray(l, &lv)); 1428ac4e037SJed Brown switch (imode) { 1438ac4e037SJed Brown case INSERT_VALUES: 1449566063dSJacob Faibussowitsch if (red->n) PetscCall(PetscArraycpy(lv, gv, red->n)); 1459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(lv, red->N, MPIU_SCALAR, red->rank, PetscObjectComm((PetscObject)dm))); 1468ac4e037SJed Brown break; 147d71ae5a4SJacob Faibussowitsch default: 148d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "InsertMode not supported"); 1498ac4e037SJed Brown } 1509566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(g, &gv)); 1519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(l, &lv)); 1523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1538ac4e037SJed Brown } 1548ac4e037SJed Brown 155d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalEnd_Redundant(DM dm, Vec g, InsertMode imode, Vec l) 156d71ae5a4SJacob Faibussowitsch { 1578ac4e037SJed Brown PetscFunctionBegin; 1583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1598ac4e037SJed Brown } 1608ac4e037SJed Brown 161d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSetUp_Redundant(DM dm) 162d71ae5a4SJacob Faibussowitsch { 1638ac4e037SJed Brown PetscFunctionBegin; 1643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1658ac4e037SJed Brown } 1668ac4e037SJed Brown 167d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMView_Redundant(DM dm, PetscViewer viewer) 168d71ae5a4SJacob Faibussowitsch { 1698ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data; 1708ac4e037SJed Brown PetscBool iascii; 1718ac4e037SJed Brown 1728ac4e037SJed Brown PetscFunctionBegin; 1739566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 17448a46eb9SPierre Jolivet if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, "redundant: rank=%d N=%" PetscInt_FMT "\n", red->rank, red->N)); 1753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1768ac4e037SJed Brown } 1778ac4e037SJed Brown 178d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateColoring_Redundant(DM dm, ISColoringType ctype, ISColoring *coloring) 179d71ae5a4SJacob Faibussowitsch { 18069527181SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data; 18169527181SJed Brown PetscInt i, nloc; 18269527181SJed Brown ISColoringValue *colors; 18369527181SJed Brown 18469527181SJed Brown PetscFunctionBegin; 18569527181SJed Brown switch (ctype) { 186d71ae5a4SJacob Faibussowitsch case IS_COLORING_GLOBAL: 187d71ae5a4SJacob Faibussowitsch nloc = red->n; 188d71ae5a4SJacob Faibussowitsch break; 189d71ae5a4SJacob Faibussowitsch case IS_COLORING_LOCAL: 190d71ae5a4SJacob Faibussowitsch nloc = red->N; 191d71ae5a4SJacob Faibussowitsch break; 192d71ae5a4SJacob Faibussowitsch default: 193d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown ISColoringType %d", (int)ctype); 19469527181SJed Brown } 1959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &colors)); 19669527181SJed Brown for (i = 0; i < nloc; i++) colors[i] = i; 1979566063dSJacob Faibussowitsch PetscCall(ISColoringCreate(PetscObjectComm((PetscObject)dm), red->N, nloc, colors, PETSC_OWN_POINTER, coloring)); 1989566063dSJacob Faibussowitsch PetscCall(ISColoringSetType(*coloring, ctype)); 1993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 20069527181SJed Brown } 2018ac4e037SJed Brown 202d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRefine_Redundant(DM dmc, MPI_Comm comm, DM *dmf) 203d71ae5a4SJacob Faibussowitsch { 2048ac4e037SJed Brown PetscMPIInt flag; 2058ac4e037SJed Brown DM_Redundant *redc = (DM_Redundant *)dmc->data; 2068ac4e037SJed Brown 2078ac4e037SJed Brown PetscFunctionBegin; 20848a46eb9SPierre Jolivet if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmc, &comm)); 2099566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)dmc), comm, &flag)); 2101dca8a05SBarry Smith PetscCheck(flag == MPI_CONGRUENT || flag == MPI_IDENT, PetscObjectComm((PetscObject)dmc), PETSC_ERR_SUP, "cannot change communicators"); 2119566063dSJacob Faibussowitsch PetscCall(DMRedundantCreate(comm, redc->rank, redc->N, dmf)); 2123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2138ac4e037SJed Brown } 2148ac4e037SJed Brown 215d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsen_Redundant(DM dmf, MPI_Comm comm, DM *dmc) 216d71ae5a4SJacob Faibussowitsch { 2178ac4e037SJed Brown PetscMPIInt flag; 2188ac4e037SJed Brown DM_Redundant *redf = (DM_Redundant *)dmf->data; 2198ac4e037SJed Brown 2208ac4e037SJed Brown PetscFunctionBegin; 22148a46eb9SPierre Jolivet if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmf, &comm)); 2229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)dmf), comm, &flag)); 2231dca8a05SBarry Smith PetscCheck(flag == MPI_CONGRUENT || flag == MPI_IDENT, PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "cannot change communicators"); 2249566063dSJacob Faibussowitsch PetscCall(DMRedundantCreate(comm, redf->rank, redf->N, dmc)); 2253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2268ac4e037SJed Brown } 2278ac4e037SJed Brown 228d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateInterpolation_Redundant(DM dmc, DM dmf, Mat *P, Vec *scale) 229d71ae5a4SJacob Faibussowitsch { 2308ac4e037SJed Brown DM_Redundant *redc = (DM_Redundant *)dmc->data; 2318ac4e037SJed Brown DM_Redundant *redf = (DM_Redundant *)dmf->data; 2328ac4e037SJed Brown PetscMPIInt flag; 2338ac4e037SJed Brown PetscInt i, rstart, rend; 2348ac4e037SJed Brown 2358ac4e037SJed Brown PetscFunctionBegin; 2369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)dmc), PetscObjectComm((PetscObject)dmf), &flag)); 2371dca8a05SBarry Smith PetscCheck(flag == MPI_CONGRUENT || flag == MPI_IDENT, PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "cannot change communicators"); 23808401ef6SPierre Jolivet PetscCheck(redc->rank == redf->rank, PetscObjectComm((PetscObject)dmf), PETSC_ERR_ARG_INCOMP, "Owning rank does not match"); 23908401ef6SPierre Jolivet PetscCheck(redc->N == redf->N, PetscObjectComm((PetscObject)dmf), PETSC_ERR_ARG_INCOMP, "Global size does not match"); 2409566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dmc), P)); 2419566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*P, redc->n, redc->n, redc->N, redc->N)); 2429566063dSJacob Faibussowitsch PetscCall(MatSetType(*P, MATAIJ)); 2439566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*P, 1, NULL)); 2449566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*P, 1, NULL, 0, NULL)); 2459566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(*P, &rstart, &rend)); 2469566063dSJacob Faibussowitsch for (i = rstart; i < rend; i++) PetscCall(MatSetValue(*P, i, i, 1.0, INSERT_VALUES)); 2479566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY)); 2489566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY)); 2499566063dSJacob Faibussowitsch if (scale) PetscCall(DMCreateInterpolationScale(dmc, dmf, *P, scale)); 2503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2518ac4e037SJed Brown } 2528ac4e037SJed Brown 2538ac4e037SJed Brown /*@ 2548ac4e037SJed Brown DMRedundantSetSize - Sets the size of a densely coupled redundant object 2558ac4e037SJed Brown 25627430b45SBarry Smith Collective 2578ac4e037SJed Brown 258d8d19677SJose E. Roman Input Parameters: 25927430b45SBarry Smith + dm - `DM` object of type `DMREDUNDANT` 26027430b45SBarry Smith . rank - rank of process to own the redundant degrees of freedom 2618ac4e037SJed Brown - N - total number of redundant degrees of freedom 2628ac4e037SJed Brown 2638ac4e037SJed Brown Level: advanced 2648ac4e037SJed Brown 265*a4e35b19SJacob Faibussowitsch .seealso: `DM`, `DMREDUNDANT`, `DMDestroy()`, `DMCreateGlobalVector()`, `DMRedundantCreate()`, `DMRedundantGetSize()` 2668ac4e037SJed Brown @*/ 267d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRedundantSetSize(DM dm, PetscMPIInt rank, PetscInt N) 268d71ae5a4SJacob Faibussowitsch { 2698ac4e037SJed Brown PetscFunctionBegin; 2708ac4e037SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2718ac4e037SJed Brown PetscValidType(dm, 1); 27269fbec6eSBarry Smith PetscValidLogicalCollectiveMPIInt(dm, rank, 2); 2738ac4e037SJed Brown PetscValidLogicalCollectiveInt(dm, N, 3); 274cac4c232SBarry Smith PetscTryMethod(dm, "DMRedundantSetSize_C", (DM, PetscMPIInt, PetscInt), (dm, rank, N)); 2753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2768ac4e037SJed Brown } 2778ac4e037SJed Brown 2788ac4e037SJed Brown /*@ 2798ac4e037SJed Brown DMRedundantGetSize - Gets the size of a densely coupled redundant object 2808ac4e037SJed Brown 2818ac4e037SJed Brown Not Collective 2828ac4e037SJed Brown 2838ac4e037SJed Brown Input Parameter: 28427430b45SBarry Smith . dm - `DM` object of type `DMREDUNDANT` 2858ac4e037SJed Brown 2868ac4e037SJed Brown Output Parameters: 28727430b45SBarry Smith + rank - rank of process to own the redundant degrees of freedom (or `NULL`) 28827430b45SBarry Smith - N - total number of redundant degrees of freedom (or `NULL`) 2898ac4e037SJed Brown 2908ac4e037SJed Brown Level: advanced 2918ac4e037SJed Brown 292*a4e35b19SJacob Faibussowitsch .seealso: `DM`, `DMREDUNDANT`, `DMDestroy()`, `DMCreateGlobalVector()`, `DMRedundantCreate()`, `DMRedundantSetSize()` 2938ac4e037SJed Brown @*/ 294d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRedundantGetSize(DM dm, PetscMPIInt *rank, PetscInt *N) 295d71ae5a4SJacob Faibussowitsch { 2968ac4e037SJed Brown PetscFunctionBegin; 2978ac4e037SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2988ac4e037SJed Brown PetscValidType(dm, 1); 299cac4c232SBarry Smith PetscUseMethod(dm, "DMRedundantGetSize_C", (DM, PetscMPIInt *, PetscInt *), (dm, rank, N)); 3003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3018ac4e037SJed Brown } 3028ac4e037SJed Brown 303d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRedundantSetSize_Redundant(DM dm, PetscMPIInt rank, PetscInt N) 304d71ae5a4SJacob Faibussowitsch { 3058ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data; 3068ac4e037SJed Brown PetscMPIInt myrank; 30773c9e547SStefano Zampini PetscInt i, *globals; 3088ac4e037SJed Brown 3098ac4e037SJed Brown PetscFunctionBegin; 3109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &myrank)); 3118ac4e037SJed Brown red->rank = rank; 3128ac4e037SJed Brown red->N = N; 3138ac4e037SJed Brown red->n = (myrank == rank) ? N : 0; 31473c9e547SStefano Zampini 31573c9e547SStefano Zampini /* mapping is setup here */ 3169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(red->N, &globals)); 31773c9e547SStefano Zampini for (i = 0; i < red->N; i++) globals[i] = i; 3189566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&dm->ltogmap)); 3199566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), 1, red->N, globals, PETSC_OWN_POINTER, &dm->ltogmap)); 3203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3218ac4e037SJed Brown } 3228ac4e037SJed Brown 323d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRedundantGetSize_Redundant(DM dm, PetscInt *rank, PetscInt *N) 324d71ae5a4SJacob Faibussowitsch { 3258ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data; 3268ac4e037SJed Brown 3278ac4e037SJed Brown PetscFunctionBegin; 3288ac4e037SJed Brown if (rank) *rank = red->rank; 3298ac4e037SJed Brown if (N) *N = red->N; 3303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3318ac4e037SJed Brown } 3328ac4e037SJed Brown 333d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSetUpGLVisViewer_Redundant(PetscObject odm, PetscViewer viewer) 334d71ae5a4SJacob Faibussowitsch { 33536623e74SStefano Zampini PetscFunctionBegin; 3363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 33736623e74SStefano Zampini } 33836623e74SStefano Zampini 3393efe6655SBarry Smith /*MC 3403efe6655SBarry Smith DMREDUNDANT = "redundant" - A DM object that is used to manage data for a small set of dense globally coupled variables. 3413efe6655SBarry Smith In the global representation of the vector the variables are all stored on a single MPI process (all the other MPI processes have 3423efe6655SBarry Smith no variables) in the local representation all the variables are stored on ALL the MPI processes (because they are all needed for each 3433efe6655SBarry Smith processes local computations). 3443efe6655SBarry Smith 3453efe6655SBarry Smith This DM is generally used inside a DMCOMPOSITE object. For example, it may be used to store continuation parameters for a bifurcation problem. 3463efe6655SBarry Smith 3473efe6655SBarry Smith Level: intermediate 3483efe6655SBarry Smith 349db781477SPatrick Sanan .seealso: `DMType`, `DMCOMPOSITE`, `DMCreate()`, `DMRedundantSetSize()`, `DMRedundantGetSize()` 3503efe6655SBarry Smith M*/ 3513efe6655SBarry Smith 352d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Redundant(DM dm) 353d71ae5a4SJacob Faibussowitsch { 3548ac4e037SJed Brown DM_Redundant *red; 3558ac4e037SJed Brown 3568ac4e037SJed Brown PetscFunctionBegin; 3574dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&red)); 3588ac4e037SJed Brown dm->data = red; 3598ac4e037SJed Brown 3608ac4e037SJed Brown dm->ops->setup = DMSetUp_Redundant; 3618ac4e037SJed Brown dm->ops->view = DMView_Redundant; 3628ac4e037SJed Brown dm->ops->createglobalvector = DMCreateGlobalVector_Redundant; 3638ac4e037SJed Brown dm->ops->createlocalvector = DMCreateLocalVector_Redundant; 36425296bd5SBarry Smith dm->ops->creatematrix = DMCreateMatrix_Redundant; 3658ac4e037SJed Brown dm->ops->destroy = DMDestroy_Redundant; 3668ac4e037SJed Brown dm->ops->globaltolocalbegin = DMGlobalToLocalBegin_Redundant; 3678ac4e037SJed Brown dm->ops->globaltolocalend = DMGlobalToLocalEnd_Redundant; 3688ac4e037SJed Brown dm->ops->localtoglobalbegin = DMLocalToGlobalBegin_Redundant; 3698ac4e037SJed Brown dm->ops->localtoglobalend = DMLocalToGlobalEnd_Redundant; 3708ac4e037SJed Brown dm->ops->refine = DMRefine_Redundant; 3718ac4e037SJed Brown dm->ops->coarsen = DMCoarsen_Redundant; 37225296bd5SBarry Smith dm->ops->createinterpolation = DMCreateInterpolation_Redundant; 373e727c939SJed Brown dm->ops->getcoloring = DMCreateColoring_Redundant; 3748865f1eaSKarl Rupp 3759566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMRedundantSetSize_C", DMRedundantSetSize_Redundant)); 3769566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMRedundantGetSize_C", DMRedundantGetSize_Redundant)); 3779566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_Redundant)); 3783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3798ac4e037SJed Brown } 3808ac4e037SJed Brown 3818ac4e037SJed Brown /*@C 38227430b45SBarry Smith DMRedundantCreate - Creates a `DM` object, used to manage data for dense globally coupled variables 3838ac4e037SJed Brown 384d083f849SBarry Smith Collective 3858ac4e037SJed Brown 386d8d19677SJose E. Roman Input Parameters: 3878ac4e037SJed Brown + comm - the processors that will share the global vector 38827430b45SBarry Smith . rank - the MPI rank to own the redundant values 3898ac4e037SJed Brown - N - total number of degrees of freedom 3908ac4e037SJed Brown 3912fe279fdSBarry Smith Output Parameter: 39227430b45SBarry Smith . dm - the `DM` object of type `DMREDUNDANT` 3938ac4e037SJed Brown 3948ac4e037SJed Brown Level: advanced 3958ac4e037SJed Brown 396*a4e35b19SJacob Faibussowitsch .seealso: `DM`, `DMREDUNDANT`, `DMDestroy()`, `DMCreateGlobalVector()`, `DMCreateMatrix()`, `DMCompositeAddDM()`, `DMSetType()`, `DMRedundantSetSize()`, `DMRedundantGetSize()` 3978ac4e037SJed Brown @*/ 398d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRedundantCreate(MPI_Comm comm, PetscMPIInt rank, PetscInt N, DM *dm) 399d71ae5a4SJacob Faibussowitsch { 4008ac4e037SJed Brown PetscFunctionBegin; 4014f572ea9SToby Isaac PetscAssertPointer(dm, 4); 4029566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 4039566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMREDUNDANT)); 4049566063dSJacob Faibussowitsch PetscCall(DMRedundantSetSize(*dm, rank, N)); 4059566063dSJacob Faibussowitsch PetscCall(DMSetUp(*dm)); 4063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4078ac4e037SJed Brown } 408