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 @*/ 16*9371c9d4SSatish Balay PetscErrorCode DMPlexSetMigrationSF(DM dm, PetscSF migrationSF) { 17f94b4a02SBlaise Bourdin PetscFunctionBegin; 18f94b4a02SBlaise Bourdin dm->sfMigration = migrationSF; 199566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)migrationSF)); 20f94b4a02SBlaise Bourdin PetscFunctionReturn(0); 21f94b4a02SBlaise Bourdin } 22f94b4a02SBlaise Bourdin 23f94b4a02SBlaise Bourdin /*@ 24f94b4a02SBlaise Bourdin DMPlexGetMigrationSF - Gets the SF for migrating from a parent DM into this DM 25f94b4a02SBlaise Bourdin 265d3b26e6SMatthew G. Knepley Input Parameter: 275d3b26e6SMatthew G. Knepley . dm - The DM 285d3b26e6SMatthew G. Knepley 295d3b26e6SMatthew G. Knepley Output Parameter: 305d3b26e6SMatthew G. Knepley . migrationSF - The PetscSF 315d3b26e6SMatthew G. Knepley 32f94b4a02SBlaise Bourdin Level: intermediate 33f94b4a02SBlaise Bourdin 34db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateMigrationSF()`, `DMPlexSetMigrationSF` 35f94b4a02SBlaise Bourdin @*/ 36*9371c9d4SSatish Balay PetscErrorCode DMPlexGetMigrationSF(DM dm, PetscSF *migrationSF) { 37f94b4a02SBlaise Bourdin PetscFunctionBegin; 38f94b4a02SBlaise Bourdin *migrationSF = dm->sfMigration; 39f94b4a02SBlaise Bourdin PetscFunctionReturn(0); 40f94b4a02SBlaise Bourdin } 41f94b4a02SBlaise Bourdin 42f94b4a02SBlaise Bourdin /*@ 43f94b4a02SBlaise Bourdin DMPlexSetGlobalToNaturalSF - Sets the SF for mapping Global Vec to the Natural Vec 44f94b4a02SBlaise Bourdin 455d3b26e6SMatthew G. Knepley Input Parameters: 46f94b4a02SBlaise Bourdin + dm - The DM 475d3b26e6SMatthew G. Knepley - naturalSF - The PetscSF 485d3b26e6SMatthew G. Knepley 49f94b4a02SBlaise Bourdin Level: intermediate 50f94b4a02SBlaise Bourdin 51db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateGlobalToNaturalSF()`, `DMPlexGetGlobaltoNaturalSF()` 52f94b4a02SBlaise Bourdin @*/ 53*9371c9d4SSatish Balay PetscErrorCode DMPlexSetGlobalToNaturalSF(DM dm, PetscSF naturalSF) { 54f94b4a02SBlaise Bourdin PetscFunctionBegin; 55f94b4a02SBlaise Bourdin dm->sfNatural = naturalSF; 569566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)naturalSF)); 57f94b4a02SBlaise Bourdin dm->useNatural = PETSC_TRUE; 58f94b4a02SBlaise Bourdin PetscFunctionReturn(0); 59f94b4a02SBlaise Bourdin } 60f94b4a02SBlaise Bourdin 61f94b4a02SBlaise Bourdin /*@ 62f94b4a02SBlaise Bourdin DMPlexGetGlobalToNaturalSF - Gets the SF for mapping Global Vec to the Natural Vec 63f94b4a02SBlaise Bourdin 645d3b26e6SMatthew G. Knepley Input Parameter: 655d3b26e6SMatthew G. Knepley . dm - The DM 665d3b26e6SMatthew G. Knepley 675d3b26e6SMatthew G. Knepley Output Parameter: 685d3b26e6SMatthew G. Knepley . naturalSF - The PetscSF 695d3b26e6SMatthew G. Knepley 70f94b4a02SBlaise Bourdin Level: intermediate 71f94b4a02SBlaise Bourdin 72db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexCreateGlobalToNaturalSF()`, `DMPlexSetGlobaltoNaturalSF` 73f94b4a02SBlaise Bourdin @*/ 74*9371c9d4SSatish Balay PetscErrorCode DMPlexGetGlobalToNaturalSF(DM dm, PetscSF *naturalSF) { 75f94b4a02SBlaise Bourdin PetscFunctionBegin; 76f94b4a02SBlaise Bourdin *naturalSF = dm->sfNatural; 77f94b4a02SBlaise Bourdin PetscFunctionReturn(0); 78f94b4a02SBlaise Bourdin } 79f94b4a02SBlaise Bourdin 80f94b4a02SBlaise Bourdin /*@ 81fa534816SMatthew G. Knepley DMPlexCreateGlobalToNaturalSF - Creates the SF for mapping Global Vec to the Natural Vec 82fa534816SMatthew G. Knepley 83fa534816SMatthew G. Knepley Input Parameters: 84e7858701SMatthew G. Knepley + dm - The redistributed DM 85e7858701SMatthew G. Knepley . section - The local PetscSection describing the Vec before the mesh was distributed, or NULL if not available 865c3f5608SAlexisMarb - sfMigration - The PetscSF used to distribute the mesh, or NULL if it cannot be computed 87fa534816SMatthew G. Knepley 885d3b26e6SMatthew G. Knepley Output Parameter: 89fa534816SMatthew G. Knepley . sfNatural - PetscSF for mapping the Vec in PETSc ordering to the canonical ordering 90fa534816SMatthew G. Knepley 915d3b26e6SMatthew G. Knepley Note: This is not typically called by the user. 925d3b26e6SMatthew G. Knepley 93fa534816SMatthew G. Knepley Level: intermediate 94fa534816SMatthew G. Knepley 95db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()` 96fa534816SMatthew G. Knepley @*/ 97*9371c9d4SSatish Balay PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM dm, PetscSection section, PetscSF sfMigration, PetscSF *sfNatural) { 98fa534816SMatthew G. Knepley MPI_Comm comm; 99e7858701SMatthew G. Knepley PetscSF sf, sfEmbed, sfField; 100e5b44f4fSMatthew G. Knepley PetscSection gSection, sectionDist, gLocSection; 101fa534816SMatthew G. Knepley PetscInt *spoints, *remoteOffsets; 1025c3f5608SAlexisMarb PetscInt ssize, pStart, pEnd, p, globalSize; 103e7858701SMatthew G. Knepley PetscBool destroyFlag = PETSC_FALSE, debug = PETSC_FALSE; 104fa534816SMatthew G. Knepley 105fa534816SMatthew G. Knepley PetscFunctionBegin; 1069566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 1075c3f5608SAlexisMarb if (!sfMigration) { 108e7858701SMatthew G. Knepley /* If sfMigration is missing, sfNatural cannot be computed and is set to NULL */ 1095c3f5608SAlexisMarb *sfNatural = NULL; 1105c3f5608SAlexisMarb PetscFunctionReturn(0); 1115c3f5608SAlexisMarb } else if (!section) { 112e7858701SMatthew G. Knepley /* If the sequential section is not provided (NULL), it is reconstructed from the parallel section */ 1135c3f5608SAlexisMarb PetscSF sfMigrationInv; 1145c3f5608SAlexisMarb PetscSection localSection; 1155c3f5608SAlexisMarb 1169566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, &localSection)); 1179566063dSJacob Faibussowitsch PetscCall(PetscSFCreateInverseSF(sfMigration, &sfMigrationInv)); 1189566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ion)); 1199566063dSJacob Faibussowitsch PetscCall(PetscSFDistributeSection(sfMigrationInv, localSection, NULL, section)); 1209566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfMigrationInv)); 1215c3f5608SAlexisMarb destroyFlag = PETSC_TRUE; 1225c3f5608SAlexisMarb } 123e7858701SMatthew G. Knepley if (debug) PetscCall(PetscSFView(sfMigration, NULL)); 124e5b44f4fSMatthew G. Knepley /* Create a new section from distributing the original section */ 1259566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, §ionDist)); 1269566063dSJacob Faibussowitsch PetscCall(PetscSFDistributeSection(sfMigration, section, &remoteOffsets, sectionDist)); 127e7858701SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)sectionDist, "Migrated Section")); 128e7858701SMatthew G. Knepley if (debug) PetscCall(PetscSectionView(sectionDist, NULL)); 1299566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(dm, sectionDist)); 130e7858701SMatthew G. Knepley /* If a sequential section is provided but no dof is affected, sfNatural cannot be computed and is set to NULL */ 13199a026b9SAlexis Marboeuf PetscCall(PetscSectionGetStorageSize(sectionDist, &globalSize)); 13214744f87SBlaise Bourdin if (globalSize) { 133e7858701SMatthew G. Knepley const PetscInt *leaves; 134e7858701SMatthew G. Knepley PetscInt *sortleaves, *indices; 135e7858701SMatthew G. Knepley PetscInt Nl; 136e7858701SMatthew G. Knepley 137fa534816SMatthew G. Knepley /* Get a pruned version of migration SF */ 1389566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(dm, &gSection)); 139e7858701SMatthew G. Knepley if (debug) PetscCall(PetscSectionView(gSection, NULL)); 140e7858701SMatthew G. Knepley PetscCall(PetscSFGetGraph(sfMigration, NULL, &Nl, &leaves, NULL)); 1419566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(gSection, &pStart, &pEnd)); 142fa534816SMatthew G. Knepley for (p = pStart, ssize = 0; p < pEnd; ++p) { 143fa534816SMatthew G. Knepley PetscInt dof, off; 144fa534816SMatthew G. Knepley 1459566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(gSection, p, &dof)); 1469566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(gSection, p, &off)); 147fa534816SMatthew G. Knepley if ((dof > 0) && (off >= 0)) ++ssize; 148fa534816SMatthew G. Knepley } 149e7858701SMatthew G. Knepley PetscCall(PetscMalloc3(ssize, &spoints, Nl, &sortleaves, Nl, &indices)); 150*9371c9d4SSatish Balay for (p = 0; p < Nl; ++p) { 151*9371c9d4SSatish Balay sortleaves[p] = leaves ? leaves[p] : p; 152*9371c9d4SSatish Balay indices[p] = p; 153*9371c9d4SSatish Balay } 154e7858701SMatthew G. Knepley PetscCall(PetscSortIntWithArray(Nl, sortleaves, indices)); 155fa534816SMatthew G. Knepley for (p = pStart, ssize = 0; p < pEnd; ++p) { 156e7858701SMatthew G. Knepley PetscInt dof, off, loc; 157fa534816SMatthew G. Knepley 1589566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(gSection, p, &dof)); 1599566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(gSection, p, &off)); 160e7858701SMatthew G. Knepley if ((dof > 0) && (off >= 0)) { 161e7858701SMatthew G. Knepley PetscCall(PetscFindInt(p, Nl, sortleaves, &loc)); 162e7858701SMatthew 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); 163e7858701SMatthew G. Knepley spoints[ssize++] = indices[loc]; 164e7858701SMatthew G. Knepley } 165fa534816SMatthew G. Knepley } 1669566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedLeafSF(sfMigration, ssize, spoints, &sfEmbed)); 167e7858701SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)sfEmbed, "Embedded SF")); 168e7858701SMatthew G. Knepley PetscCall(PetscFree3(spoints, sortleaves, indices)); 169e7858701SMatthew G. Knepley if (debug) PetscCall(PetscSFView(sfEmbed, NULL)); 170e7858701SMatthew G. Knepley /* Create the SF associated with this section 171e7858701SMatthew G. Knepley Roots are natural dofs, leaves are global dofs */ 1729566063dSJacob Faibussowitsch PetscCall(DMGetPointSF(dm, &sf)); 17399a026b9SAlexis Marboeuf PetscCall(PetscSectionCreateGlobalSection(sectionDist, sf, PETSC_TRUE, PETSC_TRUE, &gLocSection)); 1749566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, gLocSection, &sfField)); 1759566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfEmbed)); 1769566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&gLocSection)); 177e7858701SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)sfField, "Natural-to-Global SF")); 178e7858701SMatthew G. Knepley if (debug) PetscCall(PetscSFView(sfField, NULL)); 179e7858701SMatthew G. Knepley /* Invert the field SF 180e7858701SMatthew G. Knepley Roots are global dofs, leaves are natural dofs */ 181e7858701SMatthew G. Knepley PetscCall(PetscSFCreateInverseSF(sfField, sfNatural)); 182e7858701SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)*sfNatural, "Global-to-Natural SF")); 1839566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)*sfNatural, NULL, "-globaltonatural_sf_view")); 184fa534816SMatthew G. Knepley /* Clean up */ 185e7858701SMatthew G. Knepley PetscCall(PetscSFDestroy(&sfField)); 18614744f87SBlaise Bourdin } else { 18714744f87SBlaise Bourdin *sfNatural = NULL; 18814744f87SBlaise Bourdin } 1899566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ionDist)); 1909566063dSJacob Faibussowitsch PetscCall(PetscFree(remoteOffsets)); 1919566063dSJacob Faibussowitsch if (destroyFlag) PetscCall(PetscSectionDestroy(§ion)); 192fa534816SMatthew G. Knepley PetscFunctionReturn(0); 193fa534816SMatthew G. Knepley } 194fa534816SMatthew G. Knepley 195fa534816SMatthew G. Knepley /*@ 196fa534816SMatthew G. Knepley DMPlexGlobalToNaturalBegin - Rearranges a global Vector in the natural order. 197fa534816SMatthew G. Knepley 198fa534816SMatthew G. Knepley Collective on dm 199fa534816SMatthew G. Knepley 200fa534816SMatthew G. Knepley Input Parameters: 201fa534816SMatthew G. Knepley + dm - The distributed DMPlex 202fa534816SMatthew G. Knepley - gv - The global Vec 203fa534816SMatthew G. Knepley 204fa534816SMatthew G. Knepley Output Parameters: 205fa534816SMatthew G. Knepley . nv - Vec in the canonical ordering distributed over all processors associated with gv 206fa534816SMatthew G. Knepley 20742ea106eSTristan Konolige Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 208fa534816SMatthew G. Knepley 209fa534816SMatthew G. Knepley Level: intermediate 210fa534816SMatthew G. Knepley 211db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()` 212fa534816SMatthew G. Knepley @*/ 213*9371c9d4SSatish Balay PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv) { 214fa534816SMatthew G. Knepley const PetscScalar *inarray; 215fa534816SMatthew G. Knepley PetscScalar *outarray; 216a6a55facSMatthew G. Knepley PetscMPIInt size; 217fa534816SMatthew G. Knepley 218fa534816SMatthew G. Knepley PetscFunctionBegin; 2199566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0)); 2209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 221fa534816SMatthew G. Knepley if (dm->sfNatural) { 2229566063dSJacob Faibussowitsch PetscCall(VecGetArray(nv, &outarray)); 2239566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gv, &inarray)); 2249566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE)); 2259566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gv, &inarray)); 2269566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(nv, &outarray)); 227a6a55facSMatthew G. Knepley } else if (size == 1) { 2289566063dSJacob Faibussowitsch PetscCall(VecCopy(gv, nv)); 229f7d195e4SLawrence Mitchell } else { 230f7d195e4SLawrence 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."); 231f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 232f7d195e4SLawrence Mitchell } 2339566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin, dm, 0, 0, 0)); 234fa534816SMatthew G. Knepley PetscFunctionReturn(0); 235fa534816SMatthew G. Knepley } 236fa534816SMatthew G. Knepley 237fa534816SMatthew G. Knepley /*@ 238fa534816SMatthew G. Knepley DMPlexGlobalToNaturalEnd - Rearranges a global Vector in the natural order. 239fa534816SMatthew G. Knepley 240fa534816SMatthew G. Knepley Collective on dm 241fa534816SMatthew G. Knepley 242fa534816SMatthew G. Knepley Input Parameters: 243fa534816SMatthew G. Knepley + dm - The distributed DMPlex 244fa534816SMatthew G. Knepley - gv - The global Vec 245fa534816SMatthew G. Knepley 24697bb3fdcSJose E. Roman Output Parameter: 247fa534816SMatthew G. Knepley . nv - The natural Vec 248fa534816SMatthew G. Knepley 24942ea106eSTristan Konolige Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 250fa534816SMatthew G. Knepley 251fa534816SMatthew G. Knepley Level: intermediate 252fa534816SMatthew G. Knepley 253db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()` 254fa534816SMatthew G. Knepley @*/ 255*9371c9d4SSatish Balay PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv) { 256fa534816SMatthew G. Knepley const PetscScalar *inarray; 257fa534816SMatthew G. Knepley PetscScalar *outarray; 258a6a55facSMatthew G. Knepley PetscMPIInt size; 259fa534816SMatthew G. Knepley 260fa534816SMatthew G. Knepley PetscFunctionBegin; 2619566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0)); 2629566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 263fa534816SMatthew G. Knepley if (dm->sfNatural) { 2649566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gv, &inarray)); 2659566063dSJacob Faibussowitsch PetscCall(VecGetArray(nv, &outarray)); 2669566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_REPLACE)); 2679566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gv, &inarray)); 2689566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(nv, &outarray)); 269a6a55facSMatthew G. Knepley } else if (size == 1) { 270f7d195e4SLawrence Mitchell } else { 271f7d195e4SLawrence 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."); 272f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 273f7d195e4SLawrence Mitchell } 2749566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd, dm, 0, 0, 0)); 275fa534816SMatthew G. Knepley PetscFunctionReturn(0); 276fa534816SMatthew G. Knepley } 277fa534816SMatthew G. Knepley 278fa534816SMatthew G. Knepley /*@ 279fa534816SMatthew G. Knepley DMPlexNaturalToGlobalBegin - Rearranges a Vector in the natural order to the Global order. 280fa534816SMatthew G. Knepley 281fa534816SMatthew G. Knepley Collective on dm 282fa534816SMatthew G. Knepley 283fa534816SMatthew G. Knepley Input Parameters: 284fa534816SMatthew G. Knepley + dm - The distributed DMPlex 285fa534816SMatthew G. Knepley - nv - The natural Vec 286fa534816SMatthew G. Knepley 287fa534816SMatthew G. Knepley Output Parameters: 288fa534816SMatthew G. Knepley . gv - The global Vec 289fa534816SMatthew G. Knepley 29042ea106eSTristan Konolige Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 291fa534816SMatthew G. Knepley 292fa534816SMatthew G. Knepley Level: intermediate 293fa534816SMatthew G. Knepley 294c2e3fba1SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()` 295fa534816SMatthew G. Knepley @*/ 296*9371c9d4SSatish Balay PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv) { 297fa534816SMatthew G. Knepley const PetscScalar *inarray; 298fa534816SMatthew G. Knepley PetscScalar *outarray; 299a6a55facSMatthew G. Knepley PetscMPIInt size; 300fa534816SMatthew G. Knepley 301fa534816SMatthew G. Knepley PetscFunctionBegin; 3029566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0)); 3039566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 304fa534816SMatthew G. Knepley if (dm->sfNatural) { 305a5b23f4aSJose E. Roman /* We only have access to the SF that goes from Global to Natural. 306b64f75a9SBlaise Bourdin Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM. 307b64f75a9SBlaise Bourdin Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */ 3089566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(gv)); 3099566063dSJacob Faibussowitsch PetscCall(VecGetArray(gv, &outarray)); 3109566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(nv, &inarray)); 3119566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM)); 3129566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(nv, &inarray)); 3139566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gv, &outarray)); 314a6a55facSMatthew G. Knepley } else if (size == 1) { 3159566063dSJacob Faibussowitsch PetscCall(VecCopy(nv, gv)); 316f7d195e4SLawrence Mitchell } else { 317f7d195e4SLawrence 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."); 318f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 319f7d195e4SLawrence Mitchell } 3209566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin, dm, 0, 0, 0)); 321fa534816SMatthew G. Knepley PetscFunctionReturn(0); 322fa534816SMatthew G. Knepley } 323fa534816SMatthew G. Knepley 324fa534816SMatthew G. Knepley /*@ 325fa534816SMatthew G. Knepley DMPlexNaturalToGlobalEnd - Rearranges a Vector in the natural order to the Global order. 326fa534816SMatthew G. Knepley 327fa534816SMatthew G. Knepley Collective on dm 328fa534816SMatthew G. Knepley 329fa534816SMatthew G. Knepley Input Parameters: 330fa534816SMatthew G. Knepley + dm - The distributed DMPlex 331fa534816SMatthew G. Knepley - nv - The natural Vec 332fa534816SMatthew G. Knepley 333fa534816SMatthew G. Knepley Output Parameters: 334fa534816SMatthew G. Knepley . gv - The global Vec 335fa534816SMatthew G. Knepley 33642ea106eSTristan Konolige Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 337fa534816SMatthew G. Knepley 338fa534816SMatthew G. Knepley Level: intermediate 339fa534816SMatthew G. Knepley 340db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()` 341fa534816SMatthew G. Knepley @*/ 342*9371c9d4SSatish Balay PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv) { 343fa534816SMatthew G. Knepley const PetscScalar *inarray; 344fa534816SMatthew G. Knepley PetscScalar *outarray; 345a6a55facSMatthew G. Knepley PetscMPIInt size; 346fa534816SMatthew G. Knepley 347fa534816SMatthew G. Knepley PetscFunctionBegin; 3489566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0)); 3499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 350fa534816SMatthew G. Knepley if (dm->sfNatural) { 3519566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(nv, &inarray)); 3529566063dSJacob Faibussowitsch PetscCall(VecGetArray(gv, &outarray)); 3539566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *)inarray, outarray, MPI_SUM)); 3549566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(nv, &inarray)); 3559566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gv, &outarray)); 356a6a55facSMatthew G. Knepley } else if (size == 1) { 357f7d195e4SLawrence Mitchell } else { 358f7d195e4SLawrence 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."); 359f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 360f7d195e4SLawrence Mitchell } 3619566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd, dm, 0, 0, 0)); 362fa534816SMatthew G. Knepley PetscFunctionReturn(0); 363fa534816SMatthew G. Knepley } 36409cf685aSAlexis Marboeuf 36509cf685aSAlexis Marboeuf /*@ 3665cb80ecdSBarry Smith DMPlexCreateNaturalVector - Provide a `Vec` capable of holding the natural ordering and distribution. 36709cf685aSAlexis Marboeuf 36809cf685aSAlexis Marboeuf Collective on dm 36909cf685aSAlexis Marboeuf 3705cb80ecdSBarry Smith Input Parameter: 3715cb80ecdSBarry Smith . dm - The distributed `DMPlex` 37209cf685aSAlexis Marboeuf 3735cb80ecdSBarry Smith Output Parameter: 3745cb80ecdSBarry Smith . nv - The natural `Vec` 37509cf685aSAlexis Marboeuf 3765cb80ecdSBarry Smith Note: The user must call `DMSetUseNatural`(dm, PETSC_TRUE) before `DMPlexDistribute()`. 37709cf685aSAlexis Marboeuf 37809cf685aSAlexis Marboeuf Level: intermediate 37909cf685aSAlexis Marboeuf 38009cf685aSAlexis Marboeuf .seealso: `DMPlexDistribute()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()` 38109cf685aSAlexis Marboeuf @*/ 382*9371c9d4SSatish Balay PetscErrorCode DMPlexCreateNaturalVector(DM dm, Vec *nv) { 38309cf685aSAlexis Marboeuf PetscMPIInt size; 38409cf685aSAlexis Marboeuf 38509cf685aSAlexis Marboeuf PetscFunctionBegin; 38609cf685aSAlexis Marboeuf PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 38709cf685aSAlexis Marboeuf if (dm->sfNatural) { 3882e8d78feSBlaise Bourdin PetscInt nleaves, bs; 3892e8d78feSBlaise Bourdin Vec v; 3902e8d78feSBlaise Bourdin PetscCall(DMGetLocalVector(dm, &v)); 3912e8d78feSBlaise Bourdin PetscCall(VecGetBlockSize(v, &bs)); 3922e8d78feSBlaise Bourdin PetscCall(DMRestoreLocalVector(dm, &v)); 3932e8d78feSBlaise Bourdin 39409cf685aSAlexis Marboeuf PetscCall(PetscSFGetGraph(dm->sfNatural, NULL, &nleaves, NULL, NULL)); 39509cf685aSAlexis Marboeuf PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), nv)); 39609cf685aSAlexis Marboeuf PetscCall(VecSetSizes(*nv, nleaves, PETSC_DETERMINE)); 39709cf685aSAlexis Marboeuf PetscCall(VecSetBlockSize(*nv, bs)); 39809cf685aSAlexis Marboeuf PetscCall(VecSetType(*nv, dm->vectype)); 39909cf685aSAlexis Marboeuf PetscCall(VecSetDM(*nv, dm)); 40009cf685aSAlexis Marboeuf } else if (size == 1) { 4012e8d78feSBlaise Bourdin PetscCall(DMCreateLocalVector(dm, nv)); 40209cf685aSAlexis Marboeuf } else { 40309cf685aSAlexis 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."); 40409cf685aSAlexis Marboeuf SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 40509cf685aSAlexis Marboeuf } 40609cf685aSAlexis Marboeuf PetscFunctionReturn(0); 40709cf685aSAlexis Marboeuf } 408