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 @*/ 16f94b4a02SBlaise Bourdin PetscErrorCode DMPlexSetMigrationSF(DM dm, PetscSF migrationSF) 17f94b4a02SBlaise Bourdin { 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 @*/ 37f94b4a02SBlaise Bourdin PetscErrorCode DMPlexGetMigrationSF(DM dm, PetscSF *migrationSF) 38f94b4a02SBlaise Bourdin { 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 @*/ 55f94b4a02SBlaise Bourdin PetscErrorCode DMPlexSetGlobalToNaturalSF(DM dm, PetscSF naturalSF) 56f94b4a02SBlaise Bourdin { 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 @*/ 77f94b4a02SBlaise Bourdin PetscErrorCode DMPlexGetGlobalToNaturalSF(DM dm, PetscSF *naturalSF) 78f94b4a02SBlaise Bourdin { 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 @*/ 101fa534816SMatthew G. Knepley PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM dm, PetscSection section, PetscSF sfMigration, PetscSF *sfNatural) 102fa534816SMatthew G. Knepley { 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; 1075c3f5608SAlexisMarb PetscInt ssize, pStart, pEnd, p, globalSize; 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 */ 13699a026b9SAlexis Marboeuf PetscCall(PetscSectionGetStorageSize(sectionDist, &globalSize)); 13714744f87SBlaise Bourdin if (globalSize) { 138e7858701SMatthew G. Knepley const PetscInt *leaves; 139e7858701SMatthew G. Knepley PetscInt *sortleaves, *indices; 140e7858701SMatthew G. Knepley PetscInt Nl; 141e7858701SMatthew G. Knepley 142fa534816SMatthew G. Knepley /* Get a pruned version of migration SF */ 1439566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(dm, &gSection)); 144e7858701SMatthew G. Knepley if (debug) PetscCall(PetscSectionView(gSection, NULL)); 145e7858701SMatthew G. Knepley PetscCall(PetscSFGetGraph(sfMigration, NULL, &Nl, &leaves, NULL)); 1469566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(gSection, &pStart, &pEnd)); 147fa534816SMatthew G. Knepley for (p = pStart, ssize = 0; p < pEnd; ++p) { 148fa534816SMatthew G. Knepley PetscInt dof, off; 149fa534816SMatthew G. Knepley 1509566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(gSection, p, &dof)); 1519566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(gSection, p, &off)); 152fa534816SMatthew G. Knepley if ((dof > 0) && (off >= 0)) ++ssize; 153fa534816SMatthew G. Knepley } 154e7858701SMatthew G. Knepley PetscCall(PetscMalloc3(ssize, &spoints, Nl, &sortleaves, Nl, &indices)); 155e7858701SMatthew G. Knepley for (p = 0; p < Nl; ++p) {sortleaves[p] = leaves ? leaves[p] : p; indices[p] = p;} 156e7858701SMatthew G. Knepley PetscCall(PetscSortIntWithArray(Nl, sortleaves, indices)); 157fa534816SMatthew G. Knepley for (p = pStart, ssize = 0; p < pEnd; ++p) { 158e7858701SMatthew G. Knepley PetscInt dof, off, loc; 159fa534816SMatthew G. Knepley 1609566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(gSection, p, &dof)); 1619566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(gSection, p, &off)); 162e7858701SMatthew G. Knepley if ((dof > 0) && (off >= 0)) { 163e7858701SMatthew G. Knepley PetscCall(PetscFindInt(p, Nl, sortleaves, &loc)); 164e7858701SMatthew 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); 165e7858701SMatthew G. Knepley spoints[ssize++] = indices[loc]; 166e7858701SMatthew G. Knepley } 167fa534816SMatthew G. Knepley } 1689566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedLeafSF(sfMigration, ssize, spoints, &sfEmbed)); 169e7858701SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject) sfEmbed, "Embedded SF")); 170e7858701SMatthew G. Knepley PetscCall(PetscFree3(spoints, sortleaves, indices)); 171e7858701SMatthew G. Knepley if (debug) PetscCall(PetscSFView(sfEmbed, NULL)); 172e7858701SMatthew G. Knepley /* Create the SF associated with this section 173e7858701SMatthew G. Knepley Roots are natural dofs, leaves are global dofs */ 1749566063dSJacob Faibussowitsch PetscCall(DMGetPointSF(dm, &sf)); 17599a026b9SAlexis Marboeuf PetscCall(PetscSectionCreateGlobalSection(sectionDist, sf, PETSC_TRUE, PETSC_TRUE, &gLocSection)); 1769566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, gLocSection, &sfField)); 1779566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfEmbed)); 1789566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&gLocSection)); 179e7858701SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject) sfField, "Natural-to-Global SF")); 180e7858701SMatthew G. Knepley if (debug) PetscCall(PetscSFView(sfField, NULL)); 181e7858701SMatthew G. Knepley /* Invert the field SF 182e7858701SMatthew G. Knepley Roots are global dofs, leaves are natural dofs */ 183e7858701SMatthew G. Knepley PetscCall(PetscSFCreateInverseSF(sfField, sfNatural)); 184e7858701SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject) *sfNatural, "Global-to-Natural SF")); 1859566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject) *sfNatural, NULL, "-globaltonatural_sf_view")); 186fa534816SMatthew G. Knepley /* Clean up */ 187e7858701SMatthew G. Knepley PetscCall(PetscSFDestroy(&sfField)); 18814744f87SBlaise Bourdin } else { 18914744f87SBlaise Bourdin *sfNatural = NULL; 19014744f87SBlaise Bourdin } 1919566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ionDist)); 1929566063dSJacob Faibussowitsch PetscCall(PetscFree(remoteOffsets)); 1939566063dSJacob Faibussowitsch if (destroyFlag) PetscCall(PetscSectionDestroy(§ion)); 194fa534816SMatthew G. Knepley PetscFunctionReturn(0); 195fa534816SMatthew G. Knepley } 196fa534816SMatthew G. Knepley 197fa534816SMatthew G. Knepley /*@ 198fa534816SMatthew G. Knepley DMPlexGlobalToNaturalBegin - Rearranges a global Vector in the natural order. 199fa534816SMatthew G. Knepley 200fa534816SMatthew G. Knepley Collective on dm 201fa534816SMatthew G. Knepley 202fa534816SMatthew G. Knepley Input Parameters: 203fa534816SMatthew G. Knepley + dm - The distributed DMPlex 204fa534816SMatthew G. Knepley - gv - The global Vec 205fa534816SMatthew G. Knepley 206fa534816SMatthew G. Knepley Output Parameters: 207fa534816SMatthew G. Knepley . nv - Vec in the canonical ordering distributed over all processors associated with gv 208fa534816SMatthew G. Knepley 20942ea106eSTristan Konolige Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 210fa534816SMatthew G. Knepley 211fa534816SMatthew G. Knepley Level: intermediate 212fa534816SMatthew G. Knepley 213db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()` 214fa534816SMatthew G. Knepley @*/ 215fa534816SMatthew G. Knepley PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv) 216fa534816SMatthew G. Knepley { 217fa534816SMatthew G. Knepley const PetscScalar *inarray; 218fa534816SMatthew G. Knepley PetscScalar *outarray; 219a6a55facSMatthew G. Knepley PetscMPIInt size; 220fa534816SMatthew G. Knepley 221fa534816SMatthew G. Knepley PetscFunctionBegin; 2229566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin,dm,0,0,0)); 2239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size)); 224fa534816SMatthew G. Knepley if (dm->sfNatural) { 2259566063dSJacob Faibussowitsch PetscCall(VecGetArray(nv, &outarray)); 2269566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gv, &inarray)); 2279566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray,MPI_REPLACE)); 2289566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gv, &inarray)); 2299566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(nv, &outarray)); 230a6a55facSMatthew G. Knepley } else if (size == 1) { 2319566063dSJacob Faibussowitsch PetscCall(VecCopy(gv, nv)); 232f7d195e4SLawrence Mitchell } else { 233f7d195e4SLawrence 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."); 234f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 235f7d195e4SLawrence Mitchell } 2369566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin,dm,0,0,0)); 237fa534816SMatthew G. Knepley PetscFunctionReturn(0); 238fa534816SMatthew G. Knepley } 239fa534816SMatthew G. Knepley 240fa534816SMatthew G. Knepley /*@ 241fa534816SMatthew G. Knepley DMPlexGlobalToNaturalEnd - Rearranges a global Vector in the natural order. 242fa534816SMatthew G. Knepley 243fa534816SMatthew G. Knepley Collective on dm 244fa534816SMatthew G. Knepley 245fa534816SMatthew G. Knepley Input Parameters: 246fa534816SMatthew G. Knepley + dm - The distributed DMPlex 247fa534816SMatthew G. Knepley - gv - The global Vec 248fa534816SMatthew G. Knepley 24997bb3fdcSJose E. Roman Output Parameter: 250fa534816SMatthew G. Knepley . nv - The natural Vec 251fa534816SMatthew G. Knepley 25242ea106eSTristan Konolige Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 253fa534816SMatthew G. Knepley 254fa534816SMatthew G. Knepley Level: intermediate 255fa534816SMatthew G. Knepley 256db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()` 257fa534816SMatthew G. Knepley @*/ 258fa534816SMatthew G. Knepley PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv) 259fa534816SMatthew G. Knepley { 260fa534816SMatthew G. Knepley const PetscScalar *inarray; 261fa534816SMatthew G. Knepley PetscScalar *outarray; 262a6a55facSMatthew G. Knepley PetscMPIInt size; 263fa534816SMatthew G. Knepley 264fa534816SMatthew G. Knepley PetscFunctionBegin; 2659566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd,dm,0,0,0)); 2669566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size)); 267fa534816SMatthew G. Knepley if (dm->sfNatural) { 2689566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gv, &inarray)); 2699566063dSJacob Faibussowitsch PetscCall(VecGetArray(nv, &outarray)); 2709566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray,MPI_REPLACE)); 2719566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gv, &inarray)); 2729566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(nv, &outarray)); 273a6a55facSMatthew G. Knepley } else if (size == 1) { 274f7d195e4SLawrence Mitchell } else { 275f7d195e4SLawrence 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."); 276f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 277f7d195e4SLawrence Mitchell } 2789566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd,dm,0,0,0)); 279fa534816SMatthew G. Knepley PetscFunctionReturn(0); 280fa534816SMatthew G. Knepley } 281fa534816SMatthew G. Knepley 282fa534816SMatthew G. Knepley /*@ 283fa534816SMatthew G. Knepley DMPlexNaturalToGlobalBegin - Rearranges a Vector in the natural order to the Global order. 284fa534816SMatthew G. Knepley 285fa534816SMatthew G. Knepley Collective on dm 286fa534816SMatthew G. Knepley 287fa534816SMatthew G. Knepley Input Parameters: 288fa534816SMatthew G. Knepley + dm - The distributed DMPlex 289fa534816SMatthew G. Knepley - nv - The natural Vec 290fa534816SMatthew G. Knepley 291fa534816SMatthew G. Knepley Output Parameters: 292fa534816SMatthew G. Knepley . gv - The global Vec 293fa534816SMatthew G. Knepley 29442ea106eSTristan Konolige Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 295fa534816SMatthew G. Knepley 296fa534816SMatthew G. Knepley Level: intermediate 297fa534816SMatthew G. Knepley 298c2e3fba1SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()` 299fa534816SMatthew G. Knepley @*/ 300fa534816SMatthew G. Knepley PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv) 301fa534816SMatthew G. Knepley { 302fa534816SMatthew G. Knepley const PetscScalar *inarray; 303fa534816SMatthew G. Knepley PetscScalar *outarray; 304a6a55facSMatthew G. Knepley PetscMPIInt size; 305fa534816SMatthew G. Knepley 306fa534816SMatthew G. Knepley PetscFunctionBegin; 3079566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin,dm,0,0,0)); 3089566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size)); 309fa534816SMatthew G. Knepley if (dm->sfNatural) { 310a5b23f4aSJose E. Roman /* We only have access to the SF that goes from Global to Natural. 311b64f75a9SBlaise Bourdin Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM. 312b64f75a9SBlaise Bourdin Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */ 3139566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(gv)); 3149566063dSJacob Faibussowitsch PetscCall(VecGetArray(gv, &outarray)); 3159566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(nv, &inarray)); 3169566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM)); 3179566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(nv, &inarray)); 3189566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gv, &outarray)); 319a6a55facSMatthew G. Knepley } else if (size == 1) { 3209566063dSJacob Faibussowitsch PetscCall(VecCopy(nv, gv)); 321f7d195e4SLawrence Mitchell } else { 322f7d195e4SLawrence 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."); 323f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 324f7d195e4SLawrence Mitchell } 3259566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin,dm,0,0,0)); 326fa534816SMatthew G. Knepley PetscFunctionReturn(0); 327fa534816SMatthew G. Knepley } 328fa534816SMatthew G. Knepley 329fa534816SMatthew G. Knepley /*@ 330fa534816SMatthew G. Knepley DMPlexNaturalToGlobalEnd - Rearranges a Vector in the natural order to the Global order. 331fa534816SMatthew G. Knepley 332fa534816SMatthew G. Knepley Collective on dm 333fa534816SMatthew G. Knepley 334fa534816SMatthew G. Knepley Input Parameters: 335fa534816SMatthew G. Knepley + dm - The distributed DMPlex 336fa534816SMatthew G. Knepley - nv - The natural Vec 337fa534816SMatthew G. Knepley 338fa534816SMatthew G. Knepley Output Parameters: 339fa534816SMatthew G. Knepley . gv - The global Vec 340fa534816SMatthew G. Knepley 34142ea106eSTristan Konolige Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 342fa534816SMatthew G. Knepley 343fa534816SMatthew G. Knepley Level: intermediate 344fa534816SMatthew G. Knepley 345db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()` 346fa534816SMatthew G. Knepley @*/ 347fa534816SMatthew G. Knepley PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv) 348fa534816SMatthew G. Knepley { 349fa534816SMatthew G. Knepley const PetscScalar *inarray; 350fa534816SMatthew G. Knepley PetscScalar *outarray; 351a6a55facSMatthew G. Knepley PetscMPIInt size; 352fa534816SMatthew G. Knepley 353fa534816SMatthew G. Knepley PetscFunctionBegin; 3549566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd,dm,0,0,0)); 3559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size)); 356fa534816SMatthew G. Knepley if (dm->sfNatural) { 3579566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(nv, &inarray)); 3589566063dSJacob Faibussowitsch PetscCall(VecGetArray(gv, &outarray)); 3599566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM)); 3609566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(nv, &inarray)); 3619566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gv, &outarray)); 362a6a55facSMatthew G. Knepley } else if (size == 1) { 363f7d195e4SLawrence Mitchell } else { 364f7d195e4SLawrence 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."); 365f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 366f7d195e4SLawrence Mitchell } 3679566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd,dm,0,0,0)); 368fa534816SMatthew G. Knepley PetscFunctionReturn(0); 369fa534816SMatthew G. Knepley } 37009cf685aSAlexis Marboeuf 37109cf685aSAlexis Marboeuf /*@ 372*5cb80ecdSBarry Smith DMPlexCreateNaturalVector - Provide a `Vec` capable of holding the natural ordering and distribution. 37309cf685aSAlexis Marboeuf 37409cf685aSAlexis Marboeuf Collective on dm 37509cf685aSAlexis Marboeuf 376*5cb80ecdSBarry Smith Input Parameter: 377*5cb80ecdSBarry Smith . dm - The distributed `DMPlex` 37809cf685aSAlexis Marboeuf 379*5cb80ecdSBarry Smith Output Parameter: 380*5cb80ecdSBarry Smith . nv - The natural `Vec` 38109cf685aSAlexis Marboeuf 382*5cb80ecdSBarry Smith Note: The user must call `DMSetUseNatural`(dm, PETSC_TRUE) before `DMPlexDistribute()`. 38309cf685aSAlexis Marboeuf 38409cf685aSAlexis Marboeuf Level: intermediate 38509cf685aSAlexis Marboeuf 38609cf685aSAlexis Marboeuf .seealso: `DMPlexDistribute()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()` 38709cf685aSAlexis Marboeuf @*/ 38809cf685aSAlexis Marboeuf PetscErrorCode DMPlexCreateNaturalVector(DM dm, Vec *nv) 38909cf685aSAlexis Marboeuf { 39009cf685aSAlexis Marboeuf PetscMPIInt size; 39109cf685aSAlexis Marboeuf 39209cf685aSAlexis Marboeuf PetscFunctionBegin; 39309cf685aSAlexis Marboeuf PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size)); 39409cf685aSAlexis Marboeuf if (dm->sfNatural) { 3952e8d78feSBlaise Bourdin PetscInt nleaves, bs; 3962e8d78feSBlaise Bourdin Vec v; 3972e8d78feSBlaise Bourdin PetscCall(DMGetLocalVector(dm,&v)); 3982e8d78feSBlaise Bourdin PetscCall(VecGetBlockSize(v,&bs)); 3992e8d78feSBlaise Bourdin PetscCall(DMRestoreLocalVector(dm,&v)); 4002e8d78feSBlaise Bourdin 40109cf685aSAlexis Marboeuf PetscCall(PetscSFGetGraph(dm->sfNatural,NULL,&nleaves,NULL,NULL)); 40209cf685aSAlexis Marboeuf PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), nv)); 40309cf685aSAlexis Marboeuf PetscCall(VecSetSizes(*nv, nleaves, PETSC_DETERMINE)); 40409cf685aSAlexis Marboeuf PetscCall(VecSetBlockSize(*nv, bs)); 40509cf685aSAlexis Marboeuf PetscCall(VecSetType(*nv,dm->vectype)); 40609cf685aSAlexis Marboeuf PetscCall(VecSetDM(*nv, dm)); 40709cf685aSAlexis Marboeuf } else if (size == 1) { 4082e8d78feSBlaise Bourdin PetscCall(DMCreateLocalVector(dm,nv)); 40909cf685aSAlexis Marboeuf } else { 41009cf685aSAlexis 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."); 41109cf685aSAlexis Marboeuf SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 41209cf685aSAlexis Marboeuf } 41309cf685aSAlexis Marboeuf PetscFunctionReturn(0); 41409cf685aSAlexis Marboeuf } 415