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 14f94b4a02SBlaise Bourdin .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexCreateMigrationSF(), DMPlexGetMigrationSF() 15f94b4a02SBlaise Bourdin @*/ 16f94b4a02SBlaise Bourdin PetscErrorCode DMPlexSetMigrationSF(DM dm, PetscSF migrationSF) 17f94b4a02SBlaise Bourdin { 18736995cdSBlaise Bourdin PetscErrorCode ierr; 19f94b4a02SBlaise Bourdin PetscFunctionBegin; 20f94b4a02SBlaise Bourdin dm->sfMigration = migrationSF; 21736995cdSBlaise Bourdin ierr = PetscObjectReference((PetscObject) migrationSF);CHKERRQ(ierr); 22f94b4a02SBlaise Bourdin PetscFunctionReturn(0); 23f94b4a02SBlaise Bourdin } 24f94b4a02SBlaise Bourdin 25f94b4a02SBlaise Bourdin /*@ 26f94b4a02SBlaise Bourdin DMPlexGetMigrationSF - Gets the SF for migrating from a parent DM into this DM 27f94b4a02SBlaise Bourdin 285d3b26e6SMatthew G. Knepley Input Parameter: 295d3b26e6SMatthew G. Knepley . dm - The DM 305d3b26e6SMatthew G. Knepley 315d3b26e6SMatthew G. Knepley Output Parameter: 325d3b26e6SMatthew G. Knepley . migrationSF - The PetscSF 335d3b26e6SMatthew G. Knepley 34f94b4a02SBlaise Bourdin Level: intermediate 35f94b4a02SBlaise Bourdin 36f94b4a02SBlaise Bourdin .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexCreateMigrationSF(), DMPlexSetMigrationSF 37f94b4a02SBlaise Bourdin @*/ 38f94b4a02SBlaise Bourdin PetscErrorCode DMPlexGetMigrationSF(DM dm, PetscSF *migrationSF) 39f94b4a02SBlaise Bourdin { 40f94b4a02SBlaise Bourdin PetscFunctionBegin; 41f94b4a02SBlaise Bourdin *migrationSF = dm->sfMigration; 42f94b4a02SBlaise Bourdin PetscFunctionReturn(0); 43f94b4a02SBlaise Bourdin } 44f94b4a02SBlaise Bourdin 45f94b4a02SBlaise Bourdin /*@ 46f94b4a02SBlaise Bourdin DMPlexSetGlobalToNaturalSF - Sets the SF for mapping Global Vec to the Natural Vec 47f94b4a02SBlaise Bourdin 485d3b26e6SMatthew G. Knepley Input Parameters: 49f94b4a02SBlaise Bourdin + dm - The DM 505d3b26e6SMatthew G. Knepley - naturalSF - The PetscSF 515d3b26e6SMatthew G. Knepley 52f94b4a02SBlaise Bourdin Level: intermediate 53f94b4a02SBlaise Bourdin 54f94b4a02SBlaise Bourdin .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexCreateGlobalToNaturalSF(), DMPlexGetGlobaltoNaturalSF() 55f94b4a02SBlaise Bourdin @*/ 56f94b4a02SBlaise Bourdin PetscErrorCode DMPlexSetGlobalToNaturalSF(DM dm, PetscSF naturalSF) 57f94b4a02SBlaise Bourdin { 58736995cdSBlaise Bourdin PetscErrorCode ierr; 59f94b4a02SBlaise Bourdin PetscFunctionBegin; 60f94b4a02SBlaise Bourdin dm->sfNatural = naturalSF; 61736995cdSBlaise Bourdin ierr = PetscObjectReference((PetscObject) naturalSF);CHKERRQ(ierr); 62f94b4a02SBlaise Bourdin dm->useNatural = PETSC_TRUE; 63f94b4a02SBlaise Bourdin PetscFunctionReturn(0); 64f94b4a02SBlaise Bourdin } 65f94b4a02SBlaise Bourdin 66f94b4a02SBlaise Bourdin /*@ 67f94b4a02SBlaise Bourdin DMPlexGetGlobalToNaturalSF - Gets the SF for mapping Global Vec to the Natural Vec 68f94b4a02SBlaise Bourdin 695d3b26e6SMatthew G. Knepley Input Parameter: 705d3b26e6SMatthew G. Knepley . dm - The DM 715d3b26e6SMatthew G. Knepley 725d3b26e6SMatthew G. Knepley Output Parameter: 735d3b26e6SMatthew G. Knepley . naturalSF - The PetscSF 745d3b26e6SMatthew G. Knepley 75f94b4a02SBlaise Bourdin Level: intermediate 76f94b4a02SBlaise Bourdin 77f94b4a02SBlaise Bourdin .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexCreateGlobalToNaturalSF(), DMPlexSetGlobaltoNaturalSF 78f94b4a02SBlaise Bourdin @*/ 79f94b4a02SBlaise Bourdin PetscErrorCode DMPlexGetGlobalToNaturalSF(DM dm, PetscSF *naturalSF) 80f94b4a02SBlaise Bourdin { 81f94b4a02SBlaise Bourdin PetscFunctionBegin; 82f94b4a02SBlaise Bourdin *naturalSF = dm->sfNatural; 83f94b4a02SBlaise Bourdin PetscFunctionReturn(0); 84f94b4a02SBlaise Bourdin } 85f94b4a02SBlaise Bourdin 86f94b4a02SBlaise Bourdin /*@ 87fa534816SMatthew G. Knepley DMPlexCreateGlobalToNaturalSF - Creates the SF for mapping Global Vec to the Natural Vec 88fa534816SMatthew G. Knepley 89fa534816SMatthew G. Knepley Input Parameters: 90fa534816SMatthew G. Knepley + dm - The DM 915c3f5608SAlexisMarb . section - The PetscSection describing the Vec before the mesh was distributed, 925c3f5608SAlexisMarb or NULL if not available 935c3f5608SAlexisMarb - sfMigration - The PetscSF used to distribute the mesh, or NULL if it cannot be computed 94fa534816SMatthew G. Knepley 955d3b26e6SMatthew G. Knepley Output Parameter: 96fa534816SMatthew G. Knepley . sfNatural - PetscSF for mapping the Vec in PETSc ordering to the canonical ordering 97fa534816SMatthew G. Knepley 985d3b26e6SMatthew G. Knepley Note: This is not typically called by the user. 995d3b26e6SMatthew G. Knepley 100fa534816SMatthew G. Knepley Level: intermediate 101fa534816SMatthew G. Knepley 102fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField() 103fa534816SMatthew G. Knepley @*/ 104fa534816SMatthew G. Knepley PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM dm, PetscSection section, PetscSF sfMigration, PetscSF *sfNatural) 105fa534816SMatthew G. Knepley { 106fa534816SMatthew G. Knepley MPI_Comm comm; 1075c3f5608SAlexisMarb Vec gv, tmpVec; 108e5b44f4fSMatthew G. Knepley PetscSF sf, sfEmbed, sfSeqToNatural, sfField, sfFieldInv; 109e5b44f4fSMatthew G. Knepley PetscSection gSection, sectionDist, gLocSection; 110fa534816SMatthew G. Knepley PetscInt *spoints, *remoteOffsets; 1115c3f5608SAlexisMarb PetscInt ssize, pStart, pEnd, p, globalSize; 1123a350576SJunchao Zhang PetscLayout map; 1135c3f5608SAlexisMarb PetscBool destroyFlag = PETSC_FALSE; 114fa534816SMatthew G. Knepley PetscErrorCode ierr; 115fa534816SMatthew G. Knepley 116fa534816SMatthew G. Knepley PetscFunctionBegin; 117fa534816SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 1185c3f5608SAlexisMarb if (!sfMigration) { 1195c3f5608SAlexisMarb /* If sfMigration is missing, 1205c3f5608SAlexisMarb sfNatural cannot be computed and is set to NULL */ 1215c3f5608SAlexisMarb *sfNatural = NULL; 1225c3f5608SAlexisMarb PetscFunctionReturn(0); 1235c3f5608SAlexisMarb } else if (!section) { 1245c3f5608SAlexisMarb /* If the sequential section is not provided (NULL), 1255c3f5608SAlexisMarb it is reconstructed from the parallel section */ 1265c3f5608SAlexisMarb PetscSF sfMigrationInv; 1275c3f5608SAlexisMarb PetscSection localSection; 1285c3f5608SAlexisMarb 1295c3f5608SAlexisMarb ierr = DMGetLocalSection(dm, &localSection);CHKERRQ(ierr); 1305c3f5608SAlexisMarb ierr = PetscSFCreateInverseSF(sfMigration, &sfMigrationInv);CHKERRQ(ierr); 1315c3f5608SAlexisMarb ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ion);CHKERRQ(ierr); 1325c3f5608SAlexisMarb ierr = PetscSFDistributeSection(sfMigrationInv, localSection, NULL, section);CHKERRQ(ierr); 1335c3f5608SAlexisMarb ierr = PetscSFDestroy(&sfMigrationInv);CHKERRQ(ierr); 1345c3f5608SAlexisMarb destroyFlag = PETSC_TRUE; 1355c3f5608SAlexisMarb } 136e5b44f4fSMatthew G. Knepley /* ierr = PetscPrintf(comm, "Point migration SF\n");CHKERRQ(ierr); 137e5b44f4fSMatthew G. Knepley ierr = PetscSFView(sfMigration, 0);CHKERRQ(ierr); */ 138e5b44f4fSMatthew G. Knepley /* Create a new section from distributing the original section */ 139e5b44f4fSMatthew G. Knepley ierr = PetscSectionCreate(comm, §ionDist);CHKERRQ(ierr); 140e5b44f4fSMatthew G. Knepley ierr = PetscSFDistributeSection(sfMigration, section, &remoteOffsets, sectionDist);CHKERRQ(ierr); 141e5b44f4fSMatthew G. Knepley /* ierr = PetscPrintf(comm, "Distributed Section\n");CHKERRQ(ierr); 142e5b44f4fSMatthew G. Knepley ierr = PetscSectionView(sectionDist, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 14392fd8e1eSJed Brown ierr = DMSetLocalSection(dm, sectionDist);CHKERRQ(ierr); 1445c3f5608SAlexisMarb /* If a sequential section is provided but no dof is affected, 1455c3f5608SAlexisMarb sfNatural cannot be computed and is set to NULL */ 1465c3f5608SAlexisMarb ierr = DMCreateGlobalVector(dm, &tmpVec);CHKERRQ(ierr); 1475c3f5608SAlexisMarb ierr = VecGetSize(tmpVec, &globalSize);CHKERRQ(ierr); 1485c3f5608SAlexisMarb ierr = DMRestoreGlobalVector(dm, &tmpVec);CHKERRQ(ierr); 14914744f87SBlaise Bourdin if (globalSize) { 150fa534816SMatthew G. Knepley /* Get a pruned version of migration SF */ 151e87a4003SBarry Smith ierr = DMGetGlobalSection(dm, &gSection);CHKERRQ(ierr); 152fa534816SMatthew G. Knepley ierr = PetscSectionGetChart(gSection, &pStart, &pEnd);CHKERRQ(ierr); 153fa534816SMatthew G. Knepley for (p = pStart, ssize = 0; p < pEnd; ++p) { 154fa534816SMatthew G. Knepley PetscInt dof, off; 155fa534816SMatthew G. Knepley 156fa534816SMatthew G. Knepley ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr); 157fa534816SMatthew G. Knepley ierr = PetscSectionGetOffset(gSection, p, &off);CHKERRQ(ierr); 158fa534816SMatthew G. Knepley if ((dof > 0) && (off >= 0)) ++ssize; 159fa534816SMatthew G. Knepley } 160fa534816SMatthew G. Knepley ierr = PetscMalloc1(ssize, &spoints);CHKERRQ(ierr); 161fa534816SMatthew G. Knepley for (p = pStart, ssize = 0; p < pEnd; ++p) { 162fa534816SMatthew G. Knepley PetscInt dof, off; 163fa534816SMatthew G. Knepley 164fa534816SMatthew G. Knepley ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr); 165fa534816SMatthew G. Knepley ierr = PetscSectionGetOffset(gSection, p, &off);CHKERRQ(ierr); 166fa534816SMatthew G. Knepley if ((dof > 0) && (off >= 0)) spoints[ssize++] = p; 167fa534816SMatthew G. Knepley } 168fa534816SMatthew G. Knepley ierr = PetscSFCreateEmbeddedLeafSF(sfMigration, ssize, spoints, &sfEmbed);CHKERRQ(ierr); 169e5b44f4fSMatthew G. Knepley ierr = PetscFree(spoints);CHKERRQ(ierr); 170e5b44f4fSMatthew G. Knepley /* ierr = PetscPrintf(comm, "Embedded SF\n");CHKERRQ(ierr); 171e5b44f4fSMatthew G. Knepley ierr = PetscSFView(sfEmbed, 0);CHKERRQ(ierr); */ 172fa534816SMatthew G. Knepley /* Create the SF for seq to natural */ 173fa534816SMatthew G. Knepley ierr = DMGetGlobalVector(dm, &gv);CHKERRQ(ierr); 1743a350576SJunchao Zhang ierr = VecGetLayout(gv,&map);CHKERRQ(ierr); 1753a350576SJunchao Zhang /* Note that entries of gv are leaves in sfSeqToNatural, entries of the seq vec are roots */ 1763a350576SJunchao Zhang ierr = PetscSFCreate(comm, &sfSeqToNatural);CHKERRQ(ierr); 1773a350576SJunchao Zhang ierr = PetscSFSetGraphWithPattern(sfSeqToNatural, map, PETSCSF_PATTERN_GATHER);CHKERRQ(ierr); 178fa534816SMatthew G. Knepley ierr = DMRestoreGlobalVector(dm, &gv);CHKERRQ(ierr); 179e5b44f4fSMatthew G. Knepley /* ierr = PetscPrintf(comm, "Seq-to-Natural SF\n");CHKERRQ(ierr); 180e5b44f4fSMatthew G. Knepley ierr = PetscSFView(sfSeqToNatural, 0);CHKERRQ(ierr); */ 181fa534816SMatthew G. Knepley /* Create the SF associated with this section */ 182e5b44f4fSMatthew G. Knepley ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 183e5b44f4fSMatthew G. Knepley ierr = PetscSectionCreateGlobalSection(sectionDist, sf, PETSC_FALSE, PETSC_TRUE, &gLocSection);CHKERRQ(ierr); 184e5b44f4fSMatthew G. Knepley ierr = PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, gLocSection, &sfField);CHKERRQ(ierr); 185fa534816SMatthew G. Knepley ierr = PetscSFDestroy(&sfEmbed);CHKERRQ(ierr); 186e5b44f4fSMatthew G. Knepley ierr = PetscSectionDestroy(&gLocSection);CHKERRQ(ierr); 187e5b44f4fSMatthew G. Knepley /* ierr = PetscPrintf(comm, "Field SF\n");CHKERRQ(ierr); 188e5b44f4fSMatthew G. Knepley ierr = PetscSFView(sfField, 0);CHKERRQ(ierr); */ 189fa534816SMatthew G. Knepley /* Invert the field SF so it's now from distributed to sequential */ 190fa534816SMatthew G. Knepley ierr = PetscSFCreateInverseSF(sfField, &sfFieldInv);CHKERRQ(ierr); 191fa534816SMatthew G. Knepley ierr = PetscSFDestroy(&sfField);CHKERRQ(ierr); 192e5b44f4fSMatthew G. Knepley /* ierr = PetscPrintf(comm, "Inverse Field SF\n");CHKERRQ(ierr); 193e5b44f4fSMatthew G. Knepley ierr = PetscSFView(sfFieldInv, 0);CHKERRQ(ierr); */ 194fa534816SMatthew G. Knepley /* Multiply the sfFieldInv with the */ 1953a350576SJunchao Zhang ierr = PetscSFComposeInverse(sfFieldInv, sfSeqToNatural, sfNatural);CHKERRQ(ierr); 196e5b44f4fSMatthew G. Knepley ierr = PetscObjectViewFromOptions((PetscObject) *sfNatural, NULL, "-globaltonatural_sf_view");CHKERRQ(ierr); 197fa534816SMatthew G. Knepley /* Clean up */ 198fa534816SMatthew G. Knepley ierr = PetscSFDestroy(&sfFieldInv);CHKERRQ(ierr); 199fa534816SMatthew G. Knepley ierr = PetscSFDestroy(&sfSeqToNatural);CHKERRQ(ierr); 20014744f87SBlaise Bourdin } else { 20114744f87SBlaise Bourdin *sfNatural = NULL; 20214744f87SBlaise Bourdin } 20314744f87SBlaise Bourdin ierr = PetscSectionDestroy(§ionDist);CHKERRQ(ierr); 20414744f87SBlaise Bourdin ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 2055c3f5608SAlexisMarb if (destroyFlag) {ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr);} 206fa534816SMatthew G. Knepley PetscFunctionReturn(0); 207fa534816SMatthew G. Knepley } 208fa534816SMatthew G. Knepley 209fa534816SMatthew G. Knepley /*@ 210fa534816SMatthew G. Knepley DMPlexGlobalToNaturalBegin - Rearranges a global Vector in the natural order. 211fa534816SMatthew G. Knepley 212fa534816SMatthew G. Knepley Collective on dm 213fa534816SMatthew G. Knepley 214fa534816SMatthew G. Knepley Input Parameters: 215fa534816SMatthew G. Knepley + dm - The distributed DMPlex 216fa534816SMatthew G. Knepley - gv - The global Vec 217fa534816SMatthew G. Knepley 218fa534816SMatthew G. Knepley Output Parameters: 219fa534816SMatthew G. Knepley . nv - Vec in the canonical ordering distributed over all processors associated with gv 220fa534816SMatthew G. Knepley 22142ea106eSTristan Konolige Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 222fa534816SMatthew G. Knepley 223fa534816SMatthew G. Knepley Level: intermediate 224fa534816SMatthew G. Knepley 225fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalEnd() 226fa534816SMatthew G. Knepley @*/ 227fa534816SMatthew G. Knepley PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv) 228fa534816SMatthew G. Knepley { 229fa534816SMatthew G. Knepley const PetscScalar *inarray; 230fa534816SMatthew G. Knepley PetscScalar *outarray; 231a6a55facSMatthew G. Knepley PetscMPIInt size; 232fa534816SMatthew G. Knepley PetscErrorCode ierr; 233fa534816SMatthew G. Knepley 234fa534816SMatthew G. Knepley PetscFunctionBegin; 235fa534816SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin,dm,0,0,0);CHKERRQ(ierr); 236ffc4695bSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRMPI(ierr); 237fa534816SMatthew G. Knepley if (dm->sfNatural) { 238fa534816SMatthew G. Knepley ierr = VecGetArray(nv, &outarray);CHKERRQ(ierr); 239fa534816SMatthew G. Knepley ierr = VecGetArrayRead(gv, &inarray);CHKERRQ(ierr); 240ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray,MPI_REPLACE);CHKERRQ(ierr); 241fa534816SMatthew G. Knepley ierr = VecRestoreArrayRead(gv, &inarray);CHKERRQ(ierr); 242fa534816SMatthew G. Knepley ierr = VecRestoreArray(nv, &outarray);CHKERRQ(ierr); 243a6a55facSMatthew G. Knepley } else if (size == 1) { 244cafebca8SMatthew G. Knepley ierr = VecCopy(gv, nv);CHKERRQ(ierr); 245*2c71b3e2SJacob Faibussowitsch } else PetscCheckFalse(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."); 246546078acSJacob Faibussowitsch else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 247fa534816SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin,dm,0,0,0);CHKERRQ(ierr); 248fa534816SMatthew G. Knepley PetscFunctionReturn(0); 249fa534816SMatthew G. Knepley } 250fa534816SMatthew G. Knepley 251fa534816SMatthew G. Knepley /*@ 252fa534816SMatthew G. Knepley DMPlexGlobalToNaturalEnd - Rearranges a global Vector in the natural order. 253fa534816SMatthew G. Knepley 254fa534816SMatthew G. Knepley Collective on dm 255fa534816SMatthew G. Knepley 256fa534816SMatthew G. Knepley Input Parameters: 257fa534816SMatthew G. Knepley + dm - The distributed DMPlex 258fa534816SMatthew G. Knepley - gv - The global Vec 259fa534816SMatthew G. Knepley 26097bb3fdcSJose E. Roman Output Parameter: 261fa534816SMatthew G. Knepley . nv - The natural Vec 262fa534816SMatthew G. Knepley 26342ea106eSTristan Konolige Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 264fa534816SMatthew G. Knepley 265fa534816SMatthew G. Knepley Level: intermediate 266fa534816SMatthew G. Knepley 267fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalBegin() 268fa534816SMatthew G. Knepley @*/ 269fa534816SMatthew G. Knepley PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv) 270fa534816SMatthew G. Knepley { 271fa534816SMatthew G. Knepley const PetscScalar *inarray; 272fa534816SMatthew G. Knepley PetscScalar *outarray; 273a6a55facSMatthew G. Knepley PetscMPIInt size; 274fa534816SMatthew G. Knepley PetscErrorCode ierr; 275fa534816SMatthew G. Knepley 276fa534816SMatthew G. Knepley PetscFunctionBegin; 277fa534816SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd,dm,0,0,0);CHKERRQ(ierr); 278ffc4695bSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRMPI(ierr); 279fa534816SMatthew G. Knepley if (dm->sfNatural) { 280fa534816SMatthew G. Knepley ierr = VecGetArrayRead(gv, &inarray);CHKERRQ(ierr); 281c3b366b1Sprj- ierr = VecGetArray(nv, &outarray);CHKERRQ(ierr); 282ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray,MPI_REPLACE);CHKERRQ(ierr); 283fa534816SMatthew G. Knepley ierr = VecRestoreArrayRead(gv, &inarray);CHKERRQ(ierr); 284fa534816SMatthew G. Knepley ierr = VecRestoreArray(nv, &outarray);CHKERRQ(ierr); 285a6a55facSMatthew G. Knepley } else if (size == 1) { 286*2c71b3e2SJacob Faibussowitsch } else PetscCheckFalse(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."); 287546078acSJacob Faibussowitsch else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 288fa534816SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd,dm,0,0,0);CHKERRQ(ierr); 289fa534816SMatthew G. Knepley PetscFunctionReturn(0); 290fa534816SMatthew G. Knepley } 291fa534816SMatthew G. Knepley 292fa534816SMatthew G. Knepley /*@ 293fa534816SMatthew G. Knepley DMPlexNaturalToGlobalBegin - Rearranges a Vector in the natural order to the Global order. 294fa534816SMatthew G. Knepley 295fa534816SMatthew G. Knepley Collective on dm 296fa534816SMatthew G. Knepley 297fa534816SMatthew G. Knepley Input Parameters: 298fa534816SMatthew G. Knepley + dm - The distributed DMPlex 299fa534816SMatthew G. Knepley - nv - The natural Vec 300fa534816SMatthew G. Knepley 301fa534816SMatthew G. Knepley Output Parameters: 302fa534816SMatthew G. Knepley . gv - The global Vec 303fa534816SMatthew G. Knepley 30442ea106eSTristan Konolige Note: The user must call DMSetUseNatural(dm, PETSC_TRUE) before DMPlexDistribute(). 305fa534816SMatthew G. Knepley 306fa534816SMatthew G. Knepley Level: intermediate 307fa534816SMatthew G. Knepley 308fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(),DMPlexGlobalToNaturalEnd() 309fa534816SMatthew G. Knepley @*/ 310fa534816SMatthew G. Knepley PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv) 311fa534816SMatthew G. Knepley { 312fa534816SMatthew G. Knepley const PetscScalar *inarray; 313fa534816SMatthew G. Knepley PetscScalar *outarray; 314a6a55facSMatthew G. Knepley PetscMPIInt size; 315fa534816SMatthew G. Knepley PetscErrorCode ierr; 316fa534816SMatthew G. Knepley 317fa534816SMatthew G. Knepley PetscFunctionBegin; 318fa534816SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin,dm,0,0,0);CHKERRQ(ierr); 319ffc4695bSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRMPI(ierr); 320fa534816SMatthew G. Knepley if (dm->sfNatural) { 321a5b23f4aSJose E. Roman /* We only have access to the SF that goes from Global to Natural. 322b64f75a9SBlaise Bourdin Instead of inverting dm->sfNatural, we can call PetscSFReduceBegin/End with MPI_Op MPI_SUM. 323b64f75a9SBlaise Bourdin Here the SUM really does nothing since sfNatural is one to one, as long as gV is set to zero first. */ 324b64f75a9SBlaise Bourdin ierr = VecZeroEntries(gv);CHKERRQ(ierr); 325fa534816SMatthew G. Knepley ierr = VecGetArray(gv, &outarray);CHKERRQ(ierr); 326fa534816SMatthew G. Knepley ierr = VecGetArrayRead(nv, &inarray);CHKERRQ(ierr); 327fa534816SMatthew G. Knepley ierr = PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM);CHKERRQ(ierr); 328fa534816SMatthew G. Knepley ierr = VecRestoreArrayRead(nv, &inarray);CHKERRQ(ierr); 329fa534816SMatthew G. Knepley ierr = VecRestoreArray(gv, &outarray);CHKERRQ(ierr); 330a6a55facSMatthew G. Knepley } else if (size == 1) { 331a6a55facSMatthew G. Knepley ierr = VecCopy(nv, gv);CHKERRQ(ierr); 332*2c71b3e2SJacob Faibussowitsch } else PetscCheckFalse(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."); 333546078acSJacob Faibussowitsch else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 334fa534816SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin,dm,0,0,0);CHKERRQ(ierr); 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 354fa534816SMatthew G. Knepley .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 PetscErrorCode ierr; 362fa534816SMatthew G. Knepley 363fa534816SMatthew G. Knepley PetscFunctionBegin; 364fa534816SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd,dm,0,0,0);CHKERRQ(ierr); 365ffc4695bSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRMPI(ierr); 366fa534816SMatthew G. Knepley if (dm->sfNatural) { 367fa534816SMatthew G. Knepley ierr = VecGetArrayRead(nv, &inarray);CHKERRQ(ierr); 368c3b366b1Sprj- ierr = VecGetArray(gv, &outarray);CHKERRQ(ierr); 369fa534816SMatthew G. Knepley ierr = PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM);CHKERRQ(ierr); 370fa534816SMatthew G. Knepley ierr = VecRestoreArrayRead(nv, &inarray);CHKERRQ(ierr); 371fa534816SMatthew G. Knepley ierr = VecRestoreArray(gv, &outarray);CHKERRQ(ierr); 372a6a55facSMatthew G. Knepley } else if (size == 1) { 373*2c71b3e2SJacob Faibussowitsch } else PetscCheckFalse(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."); 374546078acSJacob Faibussowitsch else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute()."); 375fa534816SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd,dm,0,0,0);CHKERRQ(ierr); 376fa534816SMatthew G. Knepley PetscFunctionReturn(0); 377fa534816SMatthew G. Knepley } 378