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
DMCreateMatrix_Redundant(DM dm,Mat * J)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
DMDestroy_Redundant(DM dm)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
DMCreateGlobalVector_Redundant(DM dm,Vec * gvec)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
DMCreateLocalVector_Redundant(DM dm,Vec * lvec)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
DMLocalToGlobalBegin_Redundant(DM dm,Vec l,InsertMode imode,Vec g)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;
926497c311SBarry Smith PetscMPIInt rank, iN;
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;
1146497c311SBarry Smith PetscCall(PetscMPIIntCast(red->N, &iN));
1156497c311SBarry Smith PetscCallMPI(MPI_Reduce(source, gv, iN, MPIU_SCALAR, (imode == ADD_VALUES) ? MPIU_SUM : MPIU_MAX, red->rank, PetscObjectComm((PetscObject)dm)));
1168ac4e037SJed Brown } break;
117d71ae5a4SJacob Faibussowitsch case INSERT_VALUES:
118d71ae5a4SJacob Faibussowitsch PetscCall(PetscArraycpy(gv, lv, red->n));
119d71ae5a4SJacob Faibussowitsch break;
120d71ae5a4SJacob Faibussowitsch default:
121d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "InsertMode not supported");
1228ac4e037SJed Brown }
1239566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(l, &lv));
1249566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(g, &gv));
1253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1268ac4e037SJed Brown }
1278ac4e037SJed Brown
DMLocalToGlobalEnd_Redundant(DM dm,Vec l,InsertMode imode,Vec g)128d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMLocalToGlobalEnd_Redundant(DM dm, Vec l, InsertMode imode, Vec g)
129d71ae5a4SJacob Faibussowitsch {
1308ac4e037SJed Brown PetscFunctionBegin;
1313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1328ac4e037SJed Brown }
1338ac4e037SJed Brown
DMGlobalToLocalBegin_Redundant(DM dm,Vec g,InsertMode imode,Vec l)134d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalBegin_Redundant(DM dm, Vec g, InsertMode imode, Vec l)
135d71ae5a4SJacob Faibussowitsch {
1368ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data;
1378ac4e037SJed Brown const PetscScalar *gv;
1388ac4e037SJed Brown PetscScalar *lv;
1396497c311SBarry Smith PetscMPIInt iN;
1408ac4e037SJed Brown
1418ac4e037SJed Brown PetscFunctionBegin;
1429566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(g, &gv));
1439566063dSJacob Faibussowitsch PetscCall(VecGetArray(l, &lv));
1448ac4e037SJed Brown switch (imode) {
1458ac4e037SJed Brown case INSERT_VALUES:
1469566063dSJacob Faibussowitsch if (red->n) PetscCall(PetscArraycpy(lv, gv, red->n));
1476497c311SBarry Smith PetscCall(PetscMPIIntCast(red->N, &iN));
1486497c311SBarry Smith PetscCallMPI(MPI_Bcast(lv, iN, MPIU_SCALAR, red->rank, PetscObjectComm((PetscObject)dm)));
1498ac4e037SJed Brown break;
150d71ae5a4SJacob Faibussowitsch default:
151d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "InsertMode not supported");
1528ac4e037SJed Brown }
1539566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(g, &gv));
1549566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(l, &lv));
1553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1568ac4e037SJed Brown }
1578ac4e037SJed Brown
DMGlobalToLocalEnd_Redundant(DM dm,Vec g,InsertMode imode,Vec l)158d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalEnd_Redundant(DM dm, Vec g, InsertMode imode, Vec l)
159d71ae5a4SJacob Faibussowitsch {
1608ac4e037SJed Brown PetscFunctionBegin;
1613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1628ac4e037SJed Brown }
1638ac4e037SJed Brown
DMView_Redundant(DM dm,PetscViewer viewer)164d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMView_Redundant(DM dm, PetscViewer viewer)
165d71ae5a4SJacob Faibussowitsch {
1668ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data;
167*9f196a02SMartin Diehl PetscBool isascii;
1688ac4e037SJed Brown
1698ac4e037SJed Brown PetscFunctionBegin;
170*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
171*9f196a02SMartin Diehl if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "redundant: rank=%d N=%" PetscInt_FMT "\n", red->rank, red->N));
1723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1738ac4e037SJed Brown }
1748ac4e037SJed Brown
DMCreateColoring_Redundant(DM dm,ISColoringType ctype,ISColoring * coloring)175d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateColoring_Redundant(DM dm, ISColoringType ctype, ISColoring *coloring)
176d71ae5a4SJacob Faibussowitsch {
17769527181SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data;
17869527181SJed Brown PetscInt i, nloc;
17969527181SJed Brown ISColoringValue *colors;
18069527181SJed Brown
18169527181SJed Brown PetscFunctionBegin;
18269527181SJed Brown switch (ctype) {
183d71ae5a4SJacob Faibussowitsch case IS_COLORING_GLOBAL:
184d71ae5a4SJacob Faibussowitsch nloc = red->n;
185d71ae5a4SJacob Faibussowitsch break;
186d71ae5a4SJacob Faibussowitsch case IS_COLORING_LOCAL:
187d71ae5a4SJacob Faibussowitsch nloc = red->N;
188d71ae5a4SJacob Faibussowitsch break;
189d71ae5a4SJacob Faibussowitsch default:
190d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown ISColoringType %d", (int)ctype);
19169527181SJed Brown }
1929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc, &colors));
1936497c311SBarry Smith for (i = 0; i < nloc; i++) PetscCall(ISColoringValueCast(i, colors + i));
1949566063dSJacob Faibussowitsch PetscCall(ISColoringCreate(PetscObjectComm((PetscObject)dm), red->N, nloc, colors, PETSC_OWN_POINTER, coloring));
1959566063dSJacob Faibussowitsch PetscCall(ISColoringSetType(*coloring, ctype));
1963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
19769527181SJed Brown }
1988ac4e037SJed Brown
DMRefine_Redundant(DM dmc,MPI_Comm comm,DM * dmf)199d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRefine_Redundant(DM dmc, MPI_Comm comm, DM *dmf)
200d71ae5a4SJacob Faibussowitsch {
2018ac4e037SJed Brown PetscMPIInt flag;
2028ac4e037SJed Brown DM_Redundant *redc = (DM_Redundant *)dmc->data;
2038ac4e037SJed Brown
2048ac4e037SJed Brown PetscFunctionBegin;
20548a46eb9SPierre Jolivet if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmc, &comm));
2069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)dmc), comm, &flag));
2071dca8a05SBarry Smith PetscCheck(flag == MPI_CONGRUENT || flag == MPI_IDENT, PetscObjectComm((PetscObject)dmc), PETSC_ERR_SUP, "cannot change communicators");
2089566063dSJacob Faibussowitsch PetscCall(DMRedundantCreate(comm, redc->rank, redc->N, dmf));
2093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2108ac4e037SJed Brown }
2118ac4e037SJed Brown
DMCoarsen_Redundant(DM dmf,MPI_Comm comm,DM * dmc)212d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsen_Redundant(DM dmf, MPI_Comm comm, DM *dmc)
213d71ae5a4SJacob Faibussowitsch {
2148ac4e037SJed Brown PetscMPIInt flag;
2158ac4e037SJed Brown DM_Redundant *redf = (DM_Redundant *)dmf->data;
2168ac4e037SJed Brown
2178ac4e037SJed Brown PetscFunctionBegin;
21848a46eb9SPierre Jolivet if (comm == MPI_COMM_NULL) PetscCall(PetscObjectGetComm((PetscObject)dmf, &comm));
2199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)dmf), comm, &flag));
2201dca8a05SBarry Smith PetscCheck(flag == MPI_CONGRUENT || flag == MPI_IDENT, PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "cannot change communicators");
2219566063dSJacob Faibussowitsch PetscCall(DMRedundantCreate(comm, redf->rank, redf->N, dmc));
2223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2238ac4e037SJed Brown }
2248ac4e037SJed Brown
DMCreateInterpolation_Redundant(DM dmc,DM dmf,Mat * P,Vec * scale)225d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateInterpolation_Redundant(DM dmc, DM dmf, Mat *P, Vec *scale)
226d71ae5a4SJacob Faibussowitsch {
2278ac4e037SJed Brown DM_Redundant *redc = (DM_Redundant *)dmc->data;
2288ac4e037SJed Brown DM_Redundant *redf = (DM_Redundant *)dmf->data;
2298ac4e037SJed Brown PetscMPIInt flag;
2308ac4e037SJed Brown PetscInt i, rstart, rend;
2318ac4e037SJed Brown
2328ac4e037SJed Brown PetscFunctionBegin;
2339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)dmc), PetscObjectComm((PetscObject)dmf), &flag));
2341dca8a05SBarry Smith PetscCheck(flag == MPI_CONGRUENT || flag == MPI_IDENT, PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "cannot change communicators");
23508401ef6SPierre Jolivet PetscCheck(redc->rank == redf->rank, PetscObjectComm((PetscObject)dmf), PETSC_ERR_ARG_INCOMP, "Owning rank does not match");
23608401ef6SPierre Jolivet PetscCheck(redc->N == redf->N, PetscObjectComm((PetscObject)dmf), PETSC_ERR_ARG_INCOMP, "Global size does not match");
2379566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dmc), P));
2389566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*P, redc->n, redc->n, redc->N, redc->N));
2399566063dSJacob Faibussowitsch PetscCall(MatSetType(*P, MATAIJ));
2409566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*P, 1, NULL));
2419566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*P, 1, NULL, 0, NULL));
2429566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(*P, &rstart, &rend));
2439566063dSJacob Faibussowitsch for (i = rstart; i < rend; i++) PetscCall(MatSetValue(*P, i, i, 1.0, INSERT_VALUES));
2449566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*P, MAT_FINAL_ASSEMBLY));
2459566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*P, MAT_FINAL_ASSEMBLY));
2469566063dSJacob Faibussowitsch if (scale) PetscCall(DMCreateInterpolationScale(dmc, dmf, *P, scale));
2473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2488ac4e037SJed Brown }
2498ac4e037SJed Brown
2508ac4e037SJed Brown /*@
2518ac4e037SJed Brown DMRedundantSetSize - Sets the size of a densely coupled redundant object
2528ac4e037SJed Brown
25327430b45SBarry Smith Collective
2548ac4e037SJed Brown
255d8d19677SJose E. Roman Input Parameters:
25627430b45SBarry Smith + dm - `DM` object of type `DMREDUNDANT`
25727430b45SBarry Smith . rank - rank of process to own the redundant degrees of freedom
2588ac4e037SJed Brown - N - total number of redundant degrees of freedom
2598ac4e037SJed Brown
2608ac4e037SJed Brown Level: advanced
2618ac4e037SJed Brown
262a4e35b19SJacob Faibussowitsch .seealso: `DM`, `DMREDUNDANT`, `DMDestroy()`, `DMCreateGlobalVector()`, `DMRedundantCreate()`, `DMRedundantGetSize()`
2638ac4e037SJed Brown @*/
DMRedundantSetSize(DM dm,PetscMPIInt rank,PetscInt N)264d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRedundantSetSize(DM dm, PetscMPIInt rank, PetscInt N)
265d71ae5a4SJacob Faibussowitsch {
2668ac4e037SJed Brown PetscFunctionBegin;
2678ac4e037SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2688ac4e037SJed Brown PetscValidType(dm, 1);
26969fbec6eSBarry Smith PetscValidLogicalCollectiveMPIInt(dm, rank, 2);
2708ac4e037SJed Brown PetscValidLogicalCollectiveInt(dm, N, 3);
271cac4c232SBarry Smith PetscTryMethod(dm, "DMRedundantSetSize_C", (DM, PetscMPIInt, PetscInt), (dm, rank, N));
2723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2738ac4e037SJed Brown }
2748ac4e037SJed Brown
2758ac4e037SJed Brown /*@
2768ac4e037SJed Brown DMRedundantGetSize - Gets the size of a densely coupled redundant object
2778ac4e037SJed Brown
2788ac4e037SJed Brown Not Collective
2798ac4e037SJed Brown
2808ac4e037SJed Brown Input Parameter:
28127430b45SBarry Smith . dm - `DM` object of type `DMREDUNDANT`
2828ac4e037SJed Brown
2838ac4e037SJed Brown Output Parameters:
28427430b45SBarry Smith + rank - rank of process to own the redundant degrees of freedom (or `NULL`)
28527430b45SBarry Smith - N - total number of redundant degrees of freedom (or `NULL`)
2868ac4e037SJed Brown
2878ac4e037SJed Brown Level: advanced
2888ac4e037SJed Brown
289a4e35b19SJacob Faibussowitsch .seealso: `DM`, `DMREDUNDANT`, `DMDestroy()`, `DMCreateGlobalVector()`, `DMRedundantCreate()`, `DMRedundantSetSize()`
2908ac4e037SJed Brown @*/
DMRedundantGetSize(DM dm,PetscMPIInt * rank,PetscInt * N)291d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRedundantGetSize(DM dm, PetscMPIInt *rank, PetscInt *N)
292d71ae5a4SJacob Faibussowitsch {
2938ac4e037SJed Brown PetscFunctionBegin;
2948ac4e037SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2958ac4e037SJed Brown PetscValidType(dm, 1);
296cac4c232SBarry Smith PetscUseMethod(dm, "DMRedundantGetSize_C", (DM, PetscMPIInt *, PetscInt *), (dm, rank, N));
2973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2988ac4e037SJed Brown }
2998ac4e037SJed Brown
DMRedundantSetSize_Redundant(DM dm,PetscMPIInt rank,PetscInt N)300d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRedundantSetSize_Redundant(DM dm, PetscMPIInt rank, PetscInt N)
301d71ae5a4SJacob Faibussowitsch {
3028ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data;
3038ac4e037SJed Brown PetscMPIInt myrank;
30473c9e547SStefano Zampini PetscInt i, *globals;
3058ac4e037SJed Brown
3068ac4e037SJed Brown PetscFunctionBegin;
3079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &myrank));
3088ac4e037SJed Brown red->rank = rank;
3098ac4e037SJed Brown red->N = N;
3108ac4e037SJed Brown red->n = (myrank == rank) ? N : 0;
31173c9e547SStefano Zampini
31273c9e547SStefano Zampini /* mapping is setup here */
3139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(red->N, &globals));
31473c9e547SStefano Zampini for (i = 0; i < red->N; i++) globals[i] = i;
3159566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&dm->ltogmap));
3169566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), 1, red->N, globals, PETSC_OWN_POINTER, &dm->ltogmap));
3173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3188ac4e037SJed Brown }
3198ac4e037SJed Brown
DMRedundantGetSize_Redundant(DM dm,PetscInt * rank,PetscInt * N)320d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRedundantGetSize_Redundant(DM dm, PetscInt *rank, PetscInt *N)
321d71ae5a4SJacob Faibussowitsch {
3228ac4e037SJed Brown DM_Redundant *red = (DM_Redundant *)dm->data;
3238ac4e037SJed Brown
3248ac4e037SJed Brown PetscFunctionBegin;
3258ac4e037SJed Brown if (rank) *rank = red->rank;
3268ac4e037SJed Brown if (N) *N = red->N;
3273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3288ac4e037SJed Brown }
3298ac4e037SJed Brown
DMSetUpGLVisViewer_Redundant(PetscObject odm,PetscViewer viewer)330d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSetUpGLVisViewer_Redundant(PetscObject odm, PetscViewer viewer)
331d71ae5a4SJacob Faibussowitsch {
33236623e74SStefano Zampini PetscFunctionBegin;
3333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
33436623e74SStefano Zampini }
33536623e74SStefano Zampini
3363efe6655SBarry Smith /*MC
3373efe6655SBarry Smith DMREDUNDANT = "redundant" - A DM object that is used to manage data for a small set of dense globally coupled variables.
3383efe6655SBarry Smith In the global representation of the vector the variables are all stored on a single MPI process (all the other MPI processes have
3393efe6655SBarry Smith no variables) in the local representation all the variables are stored on ALL the MPI processes (because they are all needed for each
3403efe6655SBarry Smith processes local computations).
3413efe6655SBarry Smith
3423efe6655SBarry Smith This DM is generally used inside a DMCOMPOSITE object. For example, it may be used to store continuation parameters for a bifurcation problem.
3433efe6655SBarry Smith
3443efe6655SBarry Smith Level: intermediate
3453efe6655SBarry Smith
346db781477SPatrick Sanan .seealso: `DMType`, `DMCOMPOSITE`, `DMCreate()`, `DMRedundantSetSize()`, `DMRedundantGetSize()`
3473efe6655SBarry Smith M*/
3483efe6655SBarry Smith
DMCreate_Redundant(DM dm)349d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Redundant(DM dm)
350d71ae5a4SJacob Faibussowitsch {
3518ac4e037SJed Brown DM_Redundant *red;
3528ac4e037SJed Brown
3538ac4e037SJed Brown PetscFunctionBegin;
3544dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&red));
3558ac4e037SJed Brown dm->data = red;
3568ac4e037SJed Brown
3578ac4e037SJed Brown dm->ops->view = DMView_Redundant;
3588ac4e037SJed Brown dm->ops->createglobalvector = DMCreateGlobalVector_Redundant;
3598ac4e037SJed Brown dm->ops->createlocalvector = DMCreateLocalVector_Redundant;
36025296bd5SBarry Smith dm->ops->creatematrix = DMCreateMatrix_Redundant;
3618ac4e037SJed Brown dm->ops->destroy = DMDestroy_Redundant;
3628ac4e037SJed Brown dm->ops->globaltolocalbegin = DMGlobalToLocalBegin_Redundant;
3638ac4e037SJed Brown dm->ops->globaltolocalend = DMGlobalToLocalEnd_Redundant;
3648ac4e037SJed Brown dm->ops->localtoglobalbegin = DMLocalToGlobalBegin_Redundant;
3658ac4e037SJed Brown dm->ops->localtoglobalend = DMLocalToGlobalEnd_Redundant;
3668ac4e037SJed Brown dm->ops->refine = DMRefine_Redundant;
3678ac4e037SJed Brown dm->ops->coarsen = DMCoarsen_Redundant;
36825296bd5SBarry Smith dm->ops->createinterpolation = DMCreateInterpolation_Redundant;
369e727c939SJed Brown dm->ops->getcoloring = DMCreateColoring_Redundant;
3708865f1eaSKarl Rupp
3719566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMRedundantSetSize_C", DMRedundantSetSize_Redundant));
3729566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMRedundantGetSize_C", DMRedundantGetSize_Redundant));
3739566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_Redundant));
3743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3758ac4e037SJed Brown }
3768ac4e037SJed Brown
377cc4c1da9SBarry Smith /*@
37827430b45SBarry Smith DMRedundantCreate - Creates a `DM` object, used to manage data for dense globally coupled variables
3798ac4e037SJed Brown
380d083f849SBarry Smith Collective
3818ac4e037SJed Brown
382d8d19677SJose E. Roman Input Parameters:
3838ac4e037SJed Brown + comm - the processors that will share the global vector
38427430b45SBarry Smith . rank - the MPI rank to own the redundant values
3858ac4e037SJed Brown - N - total number of degrees of freedom
3868ac4e037SJed Brown
3872fe279fdSBarry Smith Output Parameter:
38827430b45SBarry Smith . dm - the `DM` object of type `DMREDUNDANT`
3898ac4e037SJed Brown
3908ac4e037SJed Brown Level: advanced
3918ac4e037SJed Brown
392a4e35b19SJacob Faibussowitsch .seealso: `DM`, `DMREDUNDANT`, `DMDestroy()`, `DMCreateGlobalVector()`, `DMCreateMatrix()`, `DMCompositeAddDM()`, `DMSetType()`, `DMRedundantSetSize()`, `DMRedundantGetSize()`
3938ac4e037SJed Brown @*/
DMRedundantCreate(MPI_Comm comm,PetscMPIInt rank,PetscInt N,DM * dm)394d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRedundantCreate(MPI_Comm comm, PetscMPIInt rank, PetscInt N, DM *dm)
395d71ae5a4SJacob Faibussowitsch {
3968ac4e037SJed Brown PetscFunctionBegin;
3974f572ea9SToby Isaac PetscAssertPointer(dm, 4);
3989566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm));
3999566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMREDUNDANT));
4009566063dSJacob Faibussowitsch PetscCall(DMRedundantSetSize(*dm, rank, N));
4019566063dSJacob Faibussowitsch PetscCall(DMSetUp(*dm));
4023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4038ac4e037SJed Brown }
404