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: 88fa534816SMatthew G. Knepley + dm - The DM 895c3f5608SAlexisMarb . section - The PetscSection describing the Vec before the mesh was distributed, 905c3f5608SAlexisMarb or NULL if not available 915c3f5608SAlexisMarb - sfMigration - The PetscSF used to distribute the mesh, or NULL if it cannot be computed 92fa534816SMatthew G. Knepley 935d3b26e6SMatthew G. Knepley Output Parameter: 94fa534816SMatthew G. Knepley . sfNatural - PetscSF for mapping the Vec in PETSc ordering to the canonical ordering 95fa534816SMatthew G. Knepley 965d3b26e6SMatthew G. Knepley Note: This is not typically called by the user. 975d3b26e6SMatthew G. Knepley 98fa534816SMatthew G. Knepley Level: intermediate 99fa534816SMatthew G. Knepley 100db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()` 101fa534816SMatthew G. Knepley @*/ 102fa534816SMatthew G. Knepley PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM dm, PetscSection section, PetscSF sfMigration, PetscSF *sfNatural) 103fa534816SMatthew G. Knepley { 104fa534816SMatthew G. Knepley MPI_Comm comm; 1055c3f5608SAlexisMarb Vec gv, tmpVec; 106e5b44f4fSMatthew G. Knepley PetscSF sf, sfEmbed, sfSeqToNatural, sfField, sfFieldInv; 107e5b44f4fSMatthew G. Knepley PetscSection gSection, sectionDist, gLocSection; 108fa534816SMatthew G. Knepley PetscInt *spoints, *remoteOffsets; 1095c3f5608SAlexisMarb PetscInt ssize, pStart, pEnd, p, globalSize; 1103a350576SJunchao Zhang PetscLayout map; 1115c3f5608SAlexisMarb PetscBool destroyFlag = PETSC_FALSE; 112fa534816SMatthew G. Knepley 113fa534816SMatthew G. Knepley PetscFunctionBegin; 1149566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 1155c3f5608SAlexisMarb if (!sfMigration) { 1165c3f5608SAlexisMarb /* If sfMigration is missing, 1175c3f5608SAlexisMarb sfNatural cannot be computed and is set to NULL */ 1185c3f5608SAlexisMarb *sfNatural = NULL; 1195c3f5608SAlexisMarb PetscFunctionReturn(0); 1205c3f5608SAlexisMarb } else if (!section) { 1215c3f5608SAlexisMarb /* If the sequential section is not provided (NULL), 1225c3f5608SAlexisMarb it is reconstructed from the parallel section */ 1235c3f5608SAlexisMarb PetscSF sfMigrationInv; 1245c3f5608SAlexisMarb PetscSection localSection; 1255c3f5608SAlexisMarb 1269566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, &localSection)); 1279566063dSJacob Faibussowitsch PetscCall(PetscSFCreateInverseSF(sfMigration, &sfMigrationInv)); 1289566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ion)); 1299566063dSJacob Faibussowitsch PetscCall(PetscSFDistributeSection(sfMigrationInv, localSection, NULL, section)); 1309566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfMigrationInv)); 1315c3f5608SAlexisMarb destroyFlag = PETSC_TRUE; 1325c3f5608SAlexisMarb } 1339566063dSJacob Faibussowitsch /* PetscCall(PetscPrintf(comm, "Point migration SF\n")); 1349566063dSJacob Faibussowitsch PetscCall(PetscSFView(sfMigration, 0)); */ 135e5b44f4fSMatthew G. Knepley /* Create a new section from distributing the original section */ 1369566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, §ionDist)); 1379566063dSJacob Faibussowitsch PetscCall(PetscSFDistributeSection(sfMigration, section, &remoteOffsets, sectionDist)); 1389566063dSJacob Faibussowitsch /* PetscCall(PetscPrintf(comm, "Distributed Section\n")); 1399566063dSJacob Faibussowitsch PetscCall(PetscSectionView(sectionDist, PETSC_VIEWER_STDOUT_WORLD)); */ 1409566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(dm, sectionDist)); 1415c3f5608SAlexisMarb /* If a sequential section is provided but no dof is affected, 1425c3f5608SAlexisMarb sfNatural cannot be computed and is set to NULL */ 1439566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &tmpVec)); 1449566063dSJacob Faibussowitsch PetscCall(VecGetSize(tmpVec, &globalSize)); 1459566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &tmpVec)); 14614744f87SBlaise Bourdin if (globalSize) { 147fa534816SMatthew G. Knepley /* Get a pruned version of migration SF */ 1489566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(dm, &gSection)); 1499566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(gSection, &pStart, &pEnd)); 150fa534816SMatthew G. Knepley for (p = pStart, ssize = 0; p < pEnd; ++p) { 151fa534816SMatthew G. Knepley PetscInt dof, off; 152fa534816SMatthew G. Knepley 1539566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(gSection, p, &dof)); 1549566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(gSection, p, &off)); 155fa534816SMatthew G. Knepley if ((dof > 0) && (off >= 0)) ++ssize; 156fa534816SMatthew G. Knepley } 1579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ssize, &spoints)); 158fa534816SMatthew G. Knepley for (p = pStart, ssize = 0; p < pEnd; ++p) { 159fa534816SMatthew G. Knepley PetscInt dof, off; 160fa534816SMatthew G. Knepley 1619566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(gSection, p, &dof)); 1629566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(gSection, p, &off)); 163fa534816SMatthew G. Knepley if ((dof > 0) && (off >= 0)) spoints[ssize++] = p; 164fa534816SMatthew G. Knepley } 1659566063dSJacob Faibussowitsch PetscCall(PetscSFCreateEmbeddedLeafSF(sfMigration, ssize, spoints, &sfEmbed)); 1669566063dSJacob Faibussowitsch PetscCall(PetscFree(spoints)); 1679566063dSJacob Faibussowitsch /* PetscCall(PetscPrintf(comm, "Embedded SF\n")); 1689566063dSJacob Faibussowitsch PetscCall(PetscSFView(sfEmbed, 0)); */ 169fa534816SMatthew G. Knepley /* Create the SF for seq to natural */ 1709566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &gv)); 1719566063dSJacob Faibussowitsch PetscCall(VecGetLayout(gv,&map)); 1723a350576SJunchao Zhang /* Note that entries of gv are leaves in sfSeqToNatural, entries of the seq vec are roots */ 1739566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, &sfSeqToNatural)); 1749566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraphWithPattern(sfSeqToNatural, map, PETSCSF_PATTERN_GATHER)); 1759566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &gv)); 1769566063dSJacob Faibussowitsch /* PetscCall(PetscPrintf(comm, "Seq-to-Natural SF\n")); 1779566063dSJacob Faibussowitsch PetscCall(PetscSFView(sfSeqToNatural, 0)); */ 178fa534816SMatthew G. Knepley /* Create the SF associated with this section */ 1799566063dSJacob Faibussowitsch PetscCall(DMGetPointSF(dm, &sf)); 1809566063dSJacob Faibussowitsch PetscCall(PetscSectionCreateGlobalSection(sectionDist, sf, PETSC_FALSE, PETSC_TRUE, &gLocSection)); 1819566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, gLocSection, &sfField)); 1829566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfEmbed)); 1839566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&gLocSection)); 1849566063dSJacob Faibussowitsch /* PetscCall(PetscPrintf(comm, "Field SF\n")); 1859566063dSJacob Faibussowitsch PetscCall(PetscSFView(sfField, 0)); */ 186fa534816SMatthew G. Knepley /* Invert the field SF so it's now from distributed to sequential */ 1879566063dSJacob Faibussowitsch PetscCall(PetscSFCreateInverseSF(sfField, &sfFieldInv)); 1889566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfField)); 1899566063dSJacob Faibussowitsch /* PetscCall(PetscPrintf(comm, "Inverse Field SF\n")); 1909566063dSJacob Faibussowitsch PetscCall(PetscSFView(sfFieldInv, 0)); */ 191fa534816SMatthew G. Knepley /* Multiply the sfFieldInv with the */ 1929566063dSJacob Faibussowitsch PetscCall(PetscSFComposeInverse(sfFieldInv, sfSeqToNatural, sfNatural)); 1939566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject) *sfNatural, NULL, "-globaltonatural_sf_view")); 194fa534816SMatthew G. Knepley /* Clean up */ 1959566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfFieldInv)); 1969566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfSeqToNatural)); 19714744f87SBlaise Bourdin } else { 19814744f87SBlaise Bourdin *sfNatural = NULL; 19914744f87SBlaise Bourdin } 2009566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ionDist)); 2019566063dSJacob Faibussowitsch PetscCall(PetscFree(remoteOffsets)); 2029566063dSJacob Faibussowitsch if (destroyFlag) PetscCall(PetscSectionDestroy(§ion)); 203fa534816SMatthew G. Knepley PetscFunctionReturn(0); 204fa534816SMatthew G. Knepley } 205fa534816SMatthew G. Knepley 206fa534816SMatthew G. Knepley /*@ 207fa534816SMatthew G. Knepley DMPlexGlobalToNaturalBegin - Rearranges a global Vector in the natural order. 208fa534816SMatthew G. Knepley 209fa534816SMatthew G. Knepley Collective on dm 210fa534816SMatthew G. Knepley 211fa534816SMatthew G. Knepley Input Parameters: 212fa534816SMatthew G. Knepley + dm - The distributed DMPlex 213fa534816SMatthew G. Knepley - gv - The global Vec 214fa534816SMatthew G. Knepley 215fa534816SMatthew G. Knepley Output Parameters: 216fa534816SMatthew G. Knepley . nv - Vec in the canonical ordering distributed over all processors associated with gv 217fa534816SMatthew G. Knepley 21842ea106eSTristan Konolige Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 219fa534816SMatthew G. Knepley 220fa534816SMatthew G. Knepley Level: intermediate 221fa534816SMatthew G. Knepley 222db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()` 223fa534816SMatthew G. Knepley @*/ 224fa534816SMatthew G. Knepley PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv) 225fa534816SMatthew G. Knepley { 226fa534816SMatthew G. Knepley const PetscScalar *inarray; 227fa534816SMatthew G. Knepley PetscScalar *outarray; 228a6a55facSMatthew G. Knepley PetscMPIInt size; 229fa534816SMatthew G. Knepley 230fa534816SMatthew G. Knepley PetscFunctionBegin; 2319566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin,dm,0,0,0)); 2329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size)); 233fa534816SMatthew G. Knepley if (dm->sfNatural) { 2349566063dSJacob Faibussowitsch PetscCall(VecGetArray(nv, &outarray)); 2359566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gv, &inarray)); 2369566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray,MPI_REPLACE)); 2379566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gv, &inarray)); 2389566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(nv, &outarray)); 239a6a55facSMatthew G. Knepley } else if (size == 1) { 2409566063dSJacob Faibussowitsch PetscCall(VecCopy(gv, nv)); 241*f7d195e4SLawrence Mitchell } else { 242*f7d195e4SLawrence 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."); 243*f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 244*f7d195e4SLawrence Mitchell } 2459566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin,dm,0,0,0)); 246fa534816SMatthew G. Knepley PetscFunctionReturn(0); 247fa534816SMatthew G. Knepley } 248fa534816SMatthew G. Knepley 249fa534816SMatthew G. Knepley /*@ 250fa534816SMatthew G. Knepley DMPlexGlobalToNaturalEnd - Rearranges a global Vector in the natural order. 251fa534816SMatthew G. Knepley 252fa534816SMatthew G. Knepley Collective on dm 253fa534816SMatthew G. Knepley 254fa534816SMatthew G. Knepley Input Parameters: 255fa534816SMatthew G. Knepley + dm - The distributed DMPlex 256fa534816SMatthew G. Knepley - gv - The global Vec 257fa534816SMatthew G. Knepley 25897bb3fdcSJose E. Roman Output Parameter: 259fa534816SMatthew G. Knepley . nv - The natural Vec 260fa534816SMatthew G. Knepley 26142ea106eSTristan Konolige Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 262fa534816SMatthew G. Knepley 263fa534816SMatthew G. Knepley Level: intermediate 264fa534816SMatthew G. Knepley 265db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()` 266fa534816SMatthew G. Knepley @*/ 267fa534816SMatthew G. Knepley PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv) 268fa534816SMatthew G. Knepley { 269fa534816SMatthew G. Knepley const PetscScalar *inarray; 270fa534816SMatthew G. Knepley PetscScalar *outarray; 271a6a55facSMatthew G. Knepley PetscMPIInt size; 272fa534816SMatthew G. Knepley 273fa534816SMatthew G. Knepley PetscFunctionBegin; 2749566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd,dm,0,0,0)); 2759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size)); 276fa534816SMatthew G. Knepley if (dm->sfNatural) { 2779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(gv, &inarray)); 2789566063dSJacob Faibussowitsch PetscCall(VecGetArray(nv, &outarray)); 2799566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray,MPI_REPLACE)); 2809566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(gv, &inarray)); 2819566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(nv, &outarray)); 282a6a55facSMatthew G. Knepley } else if (size == 1) { 283*f7d195e4SLawrence Mitchell } else { 284*f7d195e4SLawrence 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."); 285*f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 286*f7d195e4SLawrence Mitchell } 2879566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd,dm,0,0,0)); 288fa534816SMatthew G. Knepley PetscFunctionReturn(0); 289fa534816SMatthew G. Knepley } 290fa534816SMatthew G. Knepley 291fa534816SMatthew G. Knepley /*@ 292fa534816SMatthew G. Knepley DMPlexNaturalToGlobalBegin - Rearranges a Vector in the natural order to the Global order. 293fa534816SMatthew G. Knepley 294fa534816SMatthew G. Knepley Collective on dm 295fa534816SMatthew G. Knepley 296fa534816SMatthew G. Knepley Input Parameters: 297fa534816SMatthew G. Knepley + dm - The distributed DMPlex 298fa534816SMatthew G. Knepley - nv - The natural Vec 299fa534816SMatthew G. Knepley 300fa534816SMatthew G. Knepley Output Parameters: 301fa534816SMatthew G. Knepley . gv - The global Vec 302fa534816SMatthew G. Knepley 30342ea106eSTristan Konolige Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 304fa534816SMatthew G. Knepley 305fa534816SMatthew G. Knepley Level: intermediate 306fa534816SMatthew G. Knepley 307c2e3fba1SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalEnd()` 308fa534816SMatthew G. Knepley @*/ 309fa534816SMatthew G. Knepley PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv) 310fa534816SMatthew G. Knepley { 311fa534816SMatthew G. Knepley const PetscScalar *inarray; 312fa534816SMatthew G. Knepley PetscScalar *outarray; 313a6a55facSMatthew G. Knepley PetscMPIInt size; 314fa534816SMatthew G. Knepley 315fa534816SMatthew G. Knepley PetscFunctionBegin; 3169566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin,dm,0,0,0)); 3179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size)); 318fa534816SMatthew G. Knepley if (dm->sfNatural) { 319a5b23f4aSJose E. Roman /* We only have access to the SF that goes from Global to Natural. 320b64f75a9SBlaise Bourdin Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM. 321b64f75a9SBlaise Bourdin Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */ 3229566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(gv)); 3239566063dSJacob Faibussowitsch PetscCall(VecGetArray(gv, &outarray)); 3249566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(nv, &inarray)); 3259566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM)); 3269566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(nv, &inarray)); 3279566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gv, &outarray)); 328a6a55facSMatthew G. Knepley } else if (size == 1) { 3299566063dSJacob Faibussowitsch PetscCall(VecCopy(nv, gv)); 330*f7d195e4SLawrence Mitchell } else { 331*f7d195e4SLawrence 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."); 332*f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 333*f7d195e4SLawrence Mitchell } 3349566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin,dm,0,0,0)); 335fa534816SMatthew G. Knepley PetscFunctionReturn(0); 336fa534816SMatthew G. Knepley } 337fa534816SMatthew G. Knepley 338fa534816SMatthew G. Knepley /*@ 339fa534816SMatthew G. Knepley DMPlexNaturalToGlobalEnd - Rearranges a Vector in the natural order to the Global order. 340fa534816SMatthew G. Knepley 341fa534816SMatthew G. Knepley Collective on dm 342fa534816SMatthew G. Knepley 343fa534816SMatthew G. Knepley Input Parameters: 344fa534816SMatthew G. Knepley + dm - The distributed DMPlex 345fa534816SMatthew G. Knepley - nv - The natural Vec 346fa534816SMatthew G. Knepley 347fa534816SMatthew G. Knepley Output Parameters: 348fa534816SMatthew G. Knepley . gv - The global Vec 349fa534816SMatthew G. Knepley 35042ea106eSTristan Konolige Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 351fa534816SMatthew G. Knepley 352fa534816SMatthew G. Knepley Level: intermediate 353fa534816SMatthew G. Knepley 354db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexNaturalToGlobalBegin()`, `DMPlexGlobalToNaturalBegin()` 355fa534816SMatthew G. Knepley @*/ 356fa534816SMatthew G. Knepley PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv) 357fa534816SMatthew G. Knepley { 358fa534816SMatthew G. Knepley const PetscScalar *inarray; 359fa534816SMatthew G. Knepley PetscScalar *outarray; 360a6a55facSMatthew G. Knepley PetscMPIInt size; 361fa534816SMatthew G. Knepley 362fa534816SMatthew G. Knepley PetscFunctionBegin; 3639566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd,dm,0,0,0)); 3649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size)); 365fa534816SMatthew G. Knepley if (dm->sfNatural) { 3669566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(nv, &inarray)); 3679566063dSJacob Faibussowitsch PetscCall(VecGetArray(gv, &outarray)); 3689566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM)); 3699566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(nv, &inarray)); 3709566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(gv, &outarray)); 371a6a55facSMatthew G. Knepley } else if (size == 1) { 372*f7d195e4SLawrence Mitchell } else { 373*f7d195e4SLawrence 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."); 374*f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 375*f7d195e4SLawrence Mitchell } 3769566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd,dm,0,0,0)); 377fa534816SMatthew G. Knepley PetscFunctionReturn(0); 378fa534816SMatthew G. Knepley } 379