1fa534816SMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2fa534816SMatthew G. Knepley 3fa534816SMatthew G. Knepley /*@ 4f94b4a02SBlaise Bourdin DMPlexSetMigrationSF - Sets the SF for migrating from a parent DM into this DM 5f94b4a02SBlaise Bourdin 65d3b26e6SMatthew G. Knepley Input Parameters: 7f94b4a02SBlaise Bourdin + dm - The DM 85d3b26e6SMatthew G. Knepley - naturalSF - The PetscSF 95d3b26e6SMatthew G. Knepley 105d3b26e6SMatthew G. Knepley Note: It is necessary to call this in order to have DMCreateSubDM() or DMCreateSuperDM() build the Global-To-Natural map 115d3b26e6SMatthew G. Knepley 12f94b4a02SBlaise Bourdin Level: intermediate 13f94b4a02SBlaise Bourdin 14db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateMigrationSF()`, `DMPlexGetMigrationSF()` 15f94b4a02SBlaise Bourdin @*/ 16d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetMigrationSF(DM dm, PetscSF migrationSF) 17d71ae5a4SJacob Faibussowitsch { 18f94b4a02SBlaise Bourdin PetscFunctionBegin; 19f94b4a02SBlaise Bourdin dm->sfMigration = migrationSF; 209566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)migrationSF)); 21f94b4a02SBlaise Bourdin PetscFunctionReturn(0); 22f94b4a02SBlaise Bourdin } 23f94b4a02SBlaise Bourdin 24f94b4a02SBlaise Bourdin /*@ 25f94b4a02SBlaise Bourdin DMPlexGetMigrationSF - Gets the SF for migrating from a parent DM into this DM 26f94b4a02SBlaise Bourdin 275d3b26e6SMatthew G. Knepley Input Parameter: 285d3b26e6SMatthew G. Knepley . dm - The DM 295d3b26e6SMatthew G. Knepley 305d3b26e6SMatthew G. Knepley Output Parameter: 315d3b26e6SMatthew G. Knepley . migrationSF - The PetscSF 325d3b26e6SMatthew G. Knepley 33f94b4a02SBlaise Bourdin Level: intermediate 34f94b4a02SBlaise Bourdin 35db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateMigrationSF()`, `DMPlexSetMigrationSF` 36f94b4a02SBlaise Bourdin @*/ 37d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetMigrationSF(DM dm, PetscSF *migrationSF) 38d71ae5a4SJacob Faibussowitsch { 39f94b4a02SBlaise Bourdin PetscFunctionBegin; 40f94b4a02SBlaise Bourdin *migrationSF = dm->sfMigration; 41f94b4a02SBlaise Bourdin PetscFunctionReturn(0); 42f94b4a02SBlaise Bourdin } 43f94b4a02SBlaise Bourdin 44f94b4a02SBlaise Bourdin /*@ 45f94b4a02SBlaise Bourdin DMPlexSetGlobalToNaturalSF - Sets the SF for mapping Global Vec to the Natural Vec 46f94b4a02SBlaise Bourdin 475d3b26e6SMatthew G. Knepley Input Parameters: 48f94b4a02SBlaise Bourdin + dm - The DM 495d3b26e6SMatthew G. Knepley - naturalSF - The PetscSF 505d3b26e6SMatthew G. Knepley 51f94b4a02SBlaise Bourdin Level: intermediate 52f94b4a02SBlaise Bourdin 53db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateGlobalToNaturalSF()`, `DMPlexGetGlobaltoNaturalSF()` 54f94b4a02SBlaise Bourdin @*/ 55d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetGlobalToNaturalSF(DM dm, PetscSF naturalSF) 56d71ae5a4SJacob Faibussowitsch { 57f94b4a02SBlaise Bourdin PetscFunctionBegin; 58f94b4a02SBlaise Bourdin dm->sfNatural = naturalSF; 599566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)naturalSF)); 60f94b4a02SBlaise Bourdin dm->useNatural = PETSC_TRUE; 61f94b4a02SBlaise Bourdin PetscFunctionReturn(0); 62f94b4a02SBlaise Bourdin } 63f94b4a02SBlaise Bourdin 64f94b4a02SBlaise Bourdin /*@ 65f94b4a02SBlaise Bourdin DMPlexGetGlobalToNaturalSF - Gets the SF for mapping Global Vec to the Natural Vec 66f94b4a02SBlaise Bourdin 675d3b26e6SMatthew G. Knepley Input Parameter: 685d3b26e6SMatthew G. Knepley . dm - The DM 695d3b26e6SMatthew G. Knepley 705d3b26e6SMatthew G. Knepley Output Parameter: 715d3b26e6SMatthew G. Knepley . naturalSF - The PetscSF 725d3b26e6SMatthew G. Knepley 73f94b4a02SBlaise Bourdin Level: intermediate 74f94b4a02SBlaise Bourdin 75db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateGlobalToNaturalSF()`, `DMPlexSetGlobaltoNaturalSF` 76f94b4a02SBlaise Bourdin @*/ 77d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetGlobalToNaturalSF(DM dm, PetscSF *naturalSF) 78d71ae5a4SJacob Faibussowitsch { 79f94b4a02SBlaise Bourdin PetscFunctionBegin; 80f94b4a02SBlaise Bourdin *naturalSF = dm->sfNatural; 81f94b4a02SBlaise Bourdin PetscFunctionReturn(0); 82f94b4a02SBlaise Bourdin } 83f94b4a02SBlaise Bourdin 84f94b4a02SBlaise Bourdin /*@ 85fa534816SMatthew G. Knepley DMPlexCreateGlobalToNaturalSF - Creates the SF for mapping Global Vec to the Natural Vec 86fa534816SMatthew G. Knepley 87fa534816SMatthew G. Knepley Input Parameters: 88e7858701SMatthew G. Knepley + dm - The redistributed DM 89e7858701SMatthew G. Knepley . section - The local PetscSection describing the Vec before the mesh was distributed, or NULL if not available 905c3f5608SAlexisMarb - sfMigration - The PetscSF used to distribute the mesh, or NULL if it cannot be computed 91fa534816SMatthew G. Knepley 925d3b26e6SMatthew G. Knepley Output Parameter: 93fa534816SMatthew G. Knepley . sfNatural - PetscSF for mapping the Vec in PETSc ordering to the canonical ordering 94fa534816SMatthew G. Knepley 955d3b26e6SMatthew G. Knepley Note: This is not typically called by the user. 965d3b26e6SMatthew G. Knepley 97fa534816SMatthew G. Knepley Level: intermediate 98fa534816SMatthew G. Knepley 99db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()` 100fa534816SMatthew G. Knepley @*/ 101d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM dm, PetscSection section, PetscSF sfMigration, PetscSF *sfNatural) 102d71ae5a4SJacob Faibussowitsch { 103fa534816SMatthew G. Knepley MPI_Comm comm; 104e7858701SMatthew G. Knepley PetscSF sf, sfEmbed, sfField; 105e5b44f4fSMatthew G. Knepley PetscSection gSection, sectionDist, gLocSection; 106fa534816SMatthew G. Knepley PetscInt *spoints, *remoteOffsets; 107*57a2b1f8SBlaise Bourdin PetscInt ssize, pStart, pEnd, p, localSize, maxStorageSize; 108e7858701SMatthew G. Knepley PetscBool destroyFlag = PETSC_FALSE, debug = PETSC_FALSE; 109fa534816SMatthew G. Knepley 110fa534816SMatthew G. Knepley PetscFunctionBegin; 1119566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 1125c3f5608SAlexisMarb if (!sfMigration) { 113e7858701SMatthew G. Knepley /* If sfMigration is missing, sfNatural cannot be computed and is set to NULL */ 1145c3f5608SAlexisMarb *sfNatural = NULL; 1155c3f5608SAlexisMarb PetscFunctionReturn(0); 1165c3f5608SAlexisMarb } else if (!section) { 117e7858701SMatthew G. Knepley /* If the sequential section is not provided (NULL), it is reconstructed from the parallel section */ 1185c3f5608SAlexisMarb PetscSF sfMigrationInv; 1195c3f5608SAlexisMarb PetscSection localSection; 1205c3f5608SAlexisMarb 1219566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, &localSection)); 1229566063dSJacob Faibussowitsch PetscCall(PetscSFCreateInverseSF(sfMigration, &sfMigrationInv)); 1239566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ion)); 1249566063dSJacob Faibussowitsch PetscCall(PetscSFDistributeSection(sfMigrationInv, localSection, NULL, section)); 1259566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfMigrationInv)); 1265c3f5608SAlexisMarb destroyFlag = PETSC_TRUE; 1275c3f5608SAlexisMarb } 128e7858701SMatthew G. Knepley if (debug) PetscCall(PetscSFView(sfMigration, NULL)); 129e5b44f4fSMatthew G. Knepley /* Create a new section from distributing the original section */ 1309566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, §ionDist)); 1319566063dSJacob Faibussowitsch PetscCall(PetscSFDistributeSection(sfMigration, section, &remoteOffsets, sectionDist)); 132e7858701SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)sectionDist, "Migrated Section")); 133e7858701SMatthew G. Knepley if (debug) PetscCall(PetscSectionView(sectionDist, NULL)); 1349566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(dm, sectionDist)); 135e7858701SMatthew G. Knepley /* If a sequential section is provided but no dof is affected, sfNatural cannot be computed and is set to NULL */ 136*57a2b1f8SBlaise Bourdin PetscCall(PetscSectionGetStorageSize(sectionDist, &localSize)); 137*57a2b1f8SBlaise Bourdin PetscCallMPI(MPI_Allreduce(&localSize, &maxStorageSize, 1, MPI_INT, MPI_MAX, PetscObjectComm((PetscObject)dm))); 138*57a2b1f8SBlaise Bourdin if (maxStorageSize) { 139e7858701SMatthew G. Knepley const PetscInt *leaves; 140e7858701SMatthew G. Knepley PetscInt *sortleaves, *indices; 141e7858701SMatthew G. Knepley PetscInt Nl; 142e7858701SMatthew G. Knepley 143fa534816SMatthew G. Knepley /* Get a pruned version of migration SF */ 1449566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(dm, &gSection)); 145e7858701SMatthew G. Knepley if (debug) PetscCall(PetscSectionView(gSection, NULL)); 146e7858701SMatthew G. Knepley PetscCall(PetscSFGetGraph(sfMigration, NULL, &Nl, &leaves, NULL)); 1479566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(gSection, &pStart, &pEnd)); 148fa534816SMatthew G. Knepley for (p = pStart, ssize = 0; p < pEnd; ++p) { 149fa534816SMatthew G. Knepley PetscInt dof, off; 150fa534816SMatthew G. Knepley 1519566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(gSection, p, &dof)); 1529566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(gSection, p, &off)); 153fa534816SMatthew G. Knepley if ((dof > 0) && (off >= 0)) ++ssize; 154fa534816SMatthew G. Knepley } 155e7858701SMatthew G. Knepley PetscCall(PetscMalloc3(ssize, &spoints, Nl, &sortleaves, Nl, &indices)); 1569371c9d4SSatish Balay for (p = 0; p < Nl; ++p) { 1579371c9d4SSatish Balay sortleaves[p] = leaves ? leaves[p] : p; 1589371c9d4SSatish Balay indices[p] = p; 1599371c9d4SSatish Balay } 160e7858701SMatthew G. Knepley PetscCall(PetscSortIntWithArray(Nl, sortleaves, indices)); 161fa534816SMatthew G. Knepley for (p = pStart, ssize = 0; p < pEnd; ++p) { 162e7858701SMatthew G. Knepley PetscInt dof, off, loc; 163fa534816SMatthew G. Knepley 1649566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(gSection, p, &dof)); 1659566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(gSection, p, &off)); 166e7858701SMatthew G. Knepley if ((dof > 0) && (off >= 0)) { 167e7858701SMatthew G. Knepley PetscCall(PetscFindInt(p, Nl, sortleaves, &loc)); 168e7858701SMatthew G. Knepley PetscCheck(loc >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %" PetscInt_FMT " with nonzero dof is not a leaf of the migration SF", p); 169e7858701SMatthew G. Knepley spoints[ssize++] = indices[loc]; 170e7858701SMatthew G. Knepley } 171fa534816SMatthew G. Knepley } 1729566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedLeafSF(sfMigration, ssize, spoints, &sfEmbed)); 173e7858701SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)sfEmbed, "Embedded SF")); 174e7858701SMatthew G. Knepley PetscCall(PetscFree3(spoints, sortleaves, indices)); 175e7858701SMatthew G. Knepley if (debug) PetscCall(PetscSFView(sfEmbed, NULL)); 176e7858701SMatthew G. Knepley /* Create the SF associated with this section 177e7858701SMatthew G. Knepley Roots are natural dofs, leaves are global dofs */ 1789566063dSJacob Faibussowitsch PetscCall(DMGetPointSF(dm, &sf)); 17999a026b9SAlexis Marboeuf PetscCall(PetscSectionCreateGlobalSection(sectionDist, sf, PETSC_TRUE, PETSC_TRUE, &gLocSection)); 1809566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, gLocSection, &sfField)); 1819566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfEmbed)); 1829566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&gLocSection)); 183e7858701SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)sfField, "Natural-to-Global SF")); 184e7858701SMatthew G. Knepley if (debug) PetscCall(PetscSFView(sfField, NULL)); 185e7858701SMatthew G. Knepley /* Invert the field SF 186e7858701SMatthew G. Knepley Roots are global dofs, leaves are natural dofs */ 187e7858701SMatthew G. Knepley PetscCall(PetscSFCreateInverseSF(sfField, sfNatural)); 188e7858701SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)*sfNatural, "Global-to-Natural SF")); 1899566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)*sfNatural, NULL, "-globaltonatural_sf_view")); 190fa534816SMatthew G. Knepley /* Clean up */ 191e7858701SMatthew G. Knepley PetscCall(PetscSFDestroy(&sfField)); 19214744f87SBlaise Bourdin } else { 19314744f87SBlaise Bourdin *sfNatural = NULL; 19414744f87SBlaise Bourdin } 1959566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ionDist)); 1969566063dSJacob Faibussowitsch PetscCall(PetscFree(remoteOffsets)); 1979566063dSJacob Faibussowitsch if (destroyFlag) PetscCall(PetscSectionDestroy(§ion)); 198fa534816SMatthew G. Knepley PetscFunctionReturn(0); 199fa534816SMatthew G. Knepley } 200fa534816SMatthew G. Knepley 201fa534816SMatthew G. Knepley /*@ 202fa534816SMatthew G. Knepley DMPlexGlobalToNaturalBegin - Rearranges a global Vector in the natural order. 203fa534816SMatthew G. Knepley 204fa534816SMatthew G. Knepley Collective on dm 205fa534816SMatthew G. Knepley 206fa534816SMatthew G. Knepley Input Parameters: 207fa534816SMatthew G. Knepley + dm - The distributed DMPlex 208fa534816SMatthew G. Knepley - gv - The global Vec 209fa534816SMatthew G. Knepley 210fa534816SMatthew G. Knepley Output Parameters: 211fa534816SMatthew G. Knepley . nv - Vec in the canonical ordering distributed over all processors associated with gv 212fa534816SMatthew G. Knepley 21342ea106eSTristan Konolige Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 214fa534816SMatthew G. Knepley 215fa534816SMatthew G. Knepley Level: intermediate 216fa534816SMatthew G. Knepley 217db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()` 218fa534816SMatthew G. Knepley @*/ 219d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv) 220d71ae5a4SJacob Faibussowitsch { 221fa534816SMatthew G. Knepley const PetscScalar *inarray; 222fa534816SMatthew G. Knepley PetscScalar *outarray; 223a6a55facSMatthew G. Knepley PetscMPIInt size; 224fa534816SMatthew G. Knepley 225fa534816SMatthew G. Knepley PetscFunctionBegin; 2269566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0)); 2279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 228fa534816SMatthew G. Knepley if (dm->sfNatural) { 2299566063dSJacob Faibussowitsch PetscCall(VecGetArray(nv, &outarray)); 2309566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gv, &inarray)); 2319566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE)); 2329566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gv, &inarray)); 2339566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(nv, &outarray)); 234a6a55facSMatthew G. Knepley } else if (size == 1) { 2359566063dSJacob Faibussowitsch PetscCall(VecCopy(gv, nv)); 236f7d195e4SLawrence Mitchell } else { 237f7d195e4SLawrence Mitchell PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present.\nIf DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov."); 238f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 239f7d195e4SLawrence Mitchell } 2409566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0)); 241fa534816SMatthew G. Knepley PetscFunctionReturn(0); 242fa534816SMatthew G. Knepley } 243fa534816SMatthew G. Knepley 244fa534816SMatthew G. Knepley /*@ 245fa534816SMatthew G. Knepley DMPlexGlobalToNaturalEnd - Rearranges a global Vector in the natural order. 246fa534816SMatthew G. Knepley 247fa534816SMatthew G. Knepley Collective on dm 248fa534816SMatthew G. Knepley 249fa534816SMatthew G. Knepley Input Parameters: 250fa534816SMatthew G. Knepley + dm - The distributed DMPlex 251fa534816SMatthew G. Knepley - gv - The global Vec 252fa534816SMatthew G. Knepley 25397bb3fdcSJose E. Roman Output Parameter: 254fa534816SMatthew G. Knepley . nv - The natural Vec 255fa534816SMatthew G. Knepley 25642ea106eSTristan Konolige Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 257fa534816SMatthew G. Knepley 258fa534816SMatthew G. Knepley Level: intermediate 259fa534816SMatthew G. Knepley 260db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()` 261fa534816SMatthew G. Knepley @*/ 262d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv) 263d71ae5a4SJacob Faibussowitsch { 264fa534816SMatthew G. Knepley const PetscScalar *inarray; 265fa534816SMatthew G. Knepley PetscScalar *outarray; 266a6a55facSMatthew G. Knepley PetscMPIInt size; 267fa534816SMatthew G. Knepley 268fa534816SMatthew G. Knepley PetscFunctionBegin; 2699566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0)); 2709566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 271fa534816SMatthew G. Knepley if (dm->sfNatural) { 2729566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gv, &inarray)); 2739566063dSJacob Faibussowitsch PetscCall(VecGetArray(nv, &outarray)); 2749566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE)); 2759566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gv, &inarray)); 2769566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(nv, &outarray)); 277a6a55facSMatthew G. Knepley } else if (size == 1) { 278f7d195e4SLawrence Mitchell } else { 279f7d195e4SLawrence Mitchell PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present.\nIf DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov."); 280f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 281f7d195e4SLawrence Mitchell } 2829566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0)); 283fa534816SMatthew G. Knepley PetscFunctionReturn(0); 284fa534816SMatthew G. Knepley } 285fa534816SMatthew G. Knepley 286fa534816SMatthew G. Knepley /*@ 287fa534816SMatthew G. Knepley DMPlexNaturalToGlobalBegin - Rearranges a Vector in the natural order to the Global order. 288fa534816SMatthew G. Knepley 289fa534816SMatthew G. Knepley Collective on dm 290fa534816SMatthew G. Knepley 291fa534816SMatthew G. Knepley Input Parameters: 292fa534816SMatthew G. Knepley + dm - The distributed DMPlex 293fa534816SMatthew G. Knepley - nv - The natural Vec 294fa534816SMatthew G. Knepley 295fa534816SMatthew G. Knepley Output Parameters: 296fa534816SMatthew G. Knepley . gv - The global Vec 297fa534816SMatthew G. Knepley 29842ea106eSTristan Konolige Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 299fa534816SMatthew G. Knepley 300fa534816SMatthew G. Knepley Level: intermediate 301fa534816SMatthew G. Knepley 302c2e3fba1SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()` 303fa534816SMatthew G. Knepley @*/ 304d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv) 305d71ae5a4SJacob Faibussowitsch { 306fa534816SMatthew G. Knepley const PetscScalar *inarray; 307fa534816SMatthew G. Knepley PetscScalar *outarray; 308a6a55facSMatthew G. Knepley PetscMPIInt size; 309fa534816SMatthew G. Knepley 310fa534816SMatthew G. Knepley PetscFunctionBegin; 3119566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0)); 3129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 313fa534816SMatthew G. Knepley if (dm->sfNatural) { 314a5b23f4aSJose E. Roman /* We only have access to the SF that goes from Global to Natural. 315b64f75a9SBlaise Bourdin Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM. 316b64f75a9SBlaise Bourdin Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */ 3179566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(gv)); 3189566063dSJacob Faibussowitsch PetscCall(VecGetArray(gv, &outarray)); 3199566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(nv, &inarray)); 3209566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM)); 3219566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(nv, &inarray)); 3229566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gv, &outarray)); 323a6a55facSMatthew G. Knepley } else if (size == 1) { 3249566063dSJacob Faibussowitsch PetscCall(VecCopy(nv, gv)); 325f7d195e4SLawrence Mitchell } else { 326f7d195e4SLawrence Mitchell PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present.\nIf DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov."); 327f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 328f7d195e4SLawrence Mitchell } 3299566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0)); 330fa534816SMatthew G. Knepley PetscFunctionReturn(0); 331fa534816SMatthew G. Knepley } 332fa534816SMatthew G. Knepley 333fa534816SMatthew G. Knepley /*@ 334fa534816SMatthew G. Knepley DMPlexNaturalToGlobalEnd - Rearranges a Vector in the natural order to the Global order. 335fa534816SMatthew G. Knepley 336fa534816SMatthew G. Knepley Collective on dm 337fa534816SMatthew G. Knepley 338fa534816SMatthew G. Knepley Input Parameters: 339fa534816SMatthew G. Knepley + dm - The distributed DMPlex 340fa534816SMatthew G. Knepley - nv - The natural Vec 341fa534816SMatthew G. Knepley 342fa534816SMatthew G. Knepley Output Parameters: 343fa534816SMatthew G. Knepley . gv - The global Vec 344fa534816SMatthew G. Knepley 34542ea106eSTristan Konolige Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 346fa534816SMatthew G. Knepley 347fa534816SMatthew G. Knepley Level: intermediate 348fa534816SMatthew G. Knepley 349db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()` 350fa534816SMatthew G. Knepley @*/ 351d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv) 352d71ae5a4SJacob Faibussowitsch { 353fa534816SMatthew G. Knepley const PetscScalar *inarray; 354fa534816SMatthew G. Knepley PetscScalar *outarray; 355a6a55facSMatthew G. Knepley PetscMPIInt size; 356fa534816SMatthew G. Knepley 357fa534816SMatthew G. Knepley PetscFunctionBegin; 3589566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0)); 3599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 360fa534816SMatthew G. Knepley if (dm->sfNatural) { 3619566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(nv, &inarray)); 3629566063dSJacob Faibussowitsch PetscCall(VecGetArray(gv, &outarray)); 3639566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM)); 3649566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(nv, &inarray)); 3659566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gv, &outarray)); 366a6a55facSMatthew G. Knepley } else if (size == 1) { 367f7d195e4SLawrence Mitchell } else { 368f7d195e4SLawrence Mitchell PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present.\nIf DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov."); 369f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 370f7d195e4SLawrence Mitchell } 3719566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0)); 372fa534816SMatthew G. Knepley PetscFunctionReturn(0); 373fa534816SMatthew G. Knepley } 37409cf685aSAlexis Marboeuf 37509cf685aSAlexis Marboeuf /*@ 3765cb80ecdSBarry Smith DMPlexCreateNaturalVector - Provide a `Vec` capable of holding the natural ordering and distribution. 37709cf685aSAlexis Marboeuf 37809cf685aSAlexis Marboeuf Collective on dm 37909cf685aSAlexis Marboeuf 3805cb80ecdSBarry Smith Input Parameter: 3815cb80ecdSBarry Smith . dm - The distributed `DMPlex` 38209cf685aSAlexis Marboeuf 3835cb80ecdSBarry Smith Output Parameter: 3845cb80ecdSBarry Smith . nv - The natural `Vec` 38509cf685aSAlexis Marboeuf 3865cb80ecdSBarry Smith Note: The user must call `DMSetUseNatural`(dm, PETSC_TRUE) before `DMPlexDistribute()`. 38709cf685aSAlexis Marboeuf 38809cf685aSAlexis Marboeuf Level: intermediate 38909cf685aSAlexis Marboeuf 39009cf685aSAlexis Marboeuf .seealso: `DMPlexDistribute()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()` 39109cf685aSAlexis Marboeuf @*/ 392d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateNaturalVector(DM dm, Vec *nv) 393d71ae5a4SJacob Faibussowitsch { 39409cf685aSAlexis Marboeuf PetscMPIInt size; 39509cf685aSAlexis Marboeuf 39609cf685aSAlexis Marboeuf PetscFunctionBegin; 39709cf685aSAlexis Marboeuf PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 39809cf685aSAlexis Marboeuf if (dm->sfNatural) { 3992e8d78feSBlaise Bourdin PetscInt nleaves, bs; 4002e8d78feSBlaise Bourdin Vec v; 4012e8d78feSBlaise Bourdin PetscCall(DMGetLocalVector(dm, &v)); 4022e8d78feSBlaise Bourdin PetscCall(VecGetBlockSize(v, &bs)); 4032e8d78feSBlaise Bourdin PetscCall(DMRestoreLocalVector(dm, &v)); 4042e8d78feSBlaise Bourdin 40509cf685aSAlexis Marboeuf PetscCall(PetscSFGetGraph(dm->sfNatural, NULL, &nleaves, NULL, NULL)); 40609cf685aSAlexis Marboeuf PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), nv)); 40709cf685aSAlexis Marboeuf PetscCall(VecSetSizes(*nv, nleaves, PETSC_DETERMINE)); 40809cf685aSAlexis Marboeuf PetscCall(VecSetBlockSize(*nv, bs)); 40909cf685aSAlexis Marboeuf PetscCall(VecSetType(*nv, dm->vectype)); 41009cf685aSAlexis Marboeuf PetscCall(VecSetDM(*nv, dm)); 41109cf685aSAlexis Marboeuf } else if (size == 1) { 4122e8d78feSBlaise Bourdin PetscCall(DMCreateLocalVector(dm, nv)); 41309cf685aSAlexis Marboeuf } else { 41409cf685aSAlexis Marboeuf PetscCheck(!dm->useNatural, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "DM global to natural SF not present.\nIf DMPlexDistribute() was called and a section was defined, report to petsc-maint@mcs.anl.gov."); 41509cf685aSAlexis Marboeuf SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 41609cf685aSAlexis Marboeuf } 41709cf685aSAlexis Marboeuf PetscFunctionReturn(0); 41809cf685aSAlexis Marboeuf } 419