1fa534816SMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2fa534816SMatthew G. Knepley 3fa534816SMatthew G. Knepley #undef __FUNCT__ 4fa534816SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateGlobalToNaturalSF" 5fa534816SMatthew G. Knepley /*@ 6fa534816SMatthew G. Knepley DMPlexCreateGlobalToNaturalSF - Creates the SF for mapping Global Vec to the Natural Vec 7fa534816SMatthew G. Knepley 8fa534816SMatthew G. Knepley Input Parameters: 9fa534816SMatthew G. Knepley + dm - The DM 10fa534816SMatthew G. Knepley . section - The PetscSection before the mesh was distributed 11fa534816SMatthew G. Knepley - sfMigration - The PetscSF used to distribute the mesh 12fa534816SMatthew G. Knepley 13fa534816SMatthew G. Knepley Output Parameters: 14fa534816SMatthew G. Knepley . sfNatural - PetscSF for mapping the Vec in PETSc ordering to the canonical ordering 15fa534816SMatthew G. Knepley 16fa534816SMatthew G. Knepley Level: intermediate 17fa534816SMatthew G. Knepley 18fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField() 19fa534816SMatthew G. Knepley @*/ 20fa534816SMatthew G. Knepley PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM dm, PetscSection section, PetscSF sfMigration, PetscSF *sfNatural) 21fa534816SMatthew G. Knepley { 22fa534816SMatthew G. Knepley MPI_Comm comm; 23fa534816SMatthew G. Knepley Vec gv; 24e5b44f4fSMatthew G. Knepley PetscSF sf, sfEmbed, sfSeqToNatural, sfField, sfFieldInv; 25e5b44f4fSMatthew G. Knepley PetscSection gSection, sectionDist, gLocSection; 26fa534816SMatthew G. Knepley PetscInt *spoints, *remoteOffsets; 27fa534816SMatthew G. Knepley PetscInt ssize, pStart, pEnd, p; 28fa534816SMatthew G. Knepley PetscErrorCode ierr; 29fa534816SMatthew G. Knepley 30fa534816SMatthew G. Knepley PetscFunctionBegin; 31fa534816SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 32e5b44f4fSMatthew G. Knepley /* ierr = PetscPrintf(comm, "Point migration SF\n");CHKERRQ(ierr); 33e5b44f4fSMatthew G. Knepley ierr = PetscSFView(sfMigration, 0);CHKERRQ(ierr); */ 34e5b44f4fSMatthew G. Knepley /* Create a new section from distributing the original section */ 35e5b44f4fSMatthew G. Knepley ierr = PetscSectionCreate(comm, §ionDist);CHKERRQ(ierr); 36e5b44f4fSMatthew G. Knepley ierr = PetscSFDistributeSection(sfMigration, section, &remoteOffsets, sectionDist);CHKERRQ(ierr); 37e5b44f4fSMatthew G. Knepley /* ierr = PetscPrintf(comm, "Distributed Section\n");CHKERRQ(ierr); 38e5b44f4fSMatthew G. Knepley ierr = PetscSectionView(sectionDist, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 39e5b44f4fSMatthew G. Knepley ierr = DMSetDefaultSection(dm, sectionDist);CHKERRQ(ierr); 40fa534816SMatthew G. Knepley /* Get a pruned version of migration SF */ 41fa534816SMatthew G. Knepley ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr); 42fa534816SMatthew G. Knepley ierr = PetscSectionGetChart(gSection, &pStart, &pEnd);CHKERRQ(ierr); 43fa534816SMatthew G. Knepley for (p = pStart, ssize = 0; p < pEnd; ++p) { 44fa534816SMatthew G. Knepley PetscInt dof, off; 45fa534816SMatthew G. Knepley 46fa534816SMatthew G. Knepley ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr); 47fa534816SMatthew G. Knepley ierr = PetscSectionGetOffset(gSection, p, &off);CHKERRQ(ierr); 48fa534816SMatthew G. Knepley if ((dof > 0) && (off >= 0)) ++ssize; 49fa534816SMatthew G. Knepley } 50fa534816SMatthew G. Knepley ierr = PetscMalloc1(ssize, &spoints);CHKERRQ(ierr); 51fa534816SMatthew G. Knepley for (p = pStart, ssize = 0; p < pEnd; ++p) { 52fa534816SMatthew G. Knepley PetscInt dof, off; 53fa534816SMatthew G. Knepley 54fa534816SMatthew G. Knepley ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr); 55fa534816SMatthew G. Knepley ierr = PetscSectionGetOffset(gSection, p, &off);CHKERRQ(ierr); 56fa534816SMatthew G. Knepley if ((dof > 0) && (off >= 0)) spoints[ssize++] = p; 57fa534816SMatthew G. Knepley } 58fa534816SMatthew G. Knepley ierr = PetscSFCreateEmbeddedLeafSF(sfMigration, ssize, spoints, &sfEmbed);CHKERRQ(ierr); 59e5b44f4fSMatthew G. Knepley ierr = PetscFree(spoints);CHKERRQ(ierr); 60e5b44f4fSMatthew G. Knepley /* ierr = PetscPrintf(comm, "Embedded SF\n");CHKERRQ(ierr); 61e5b44f4fSMatthew G. Knepley ierr = PetscSFView(sfEmbed, 0);CHKERRQ(ierr); */ 62fa534816SMatthew G. Knepley /* Create the SF for seq to natural */ 63fa534816SMatthew G. Knepley ierr = DMGetGlobalVector(dm, &gv);CHKERRQ(ierr); 64fa534816SMatthew G. Knepley ierr = PetscSFCreateFromZero(comm, gv, &sfSeqToNatural);CHKERRQ(ierr); 65fa534816SMatthew G. Knepley ierr = DMRestoreGlobalVector(dm, &gv);CHKERRQ(ierr); 66e5b44f4fSMatthew G. Knepley /* ierr = PetscPrintf(comm, "Seq-to-Natural SF\n");CHKERRQ(ierr); 67e5b44f4fSMatthew G. Knepley ierr = PetscSFView(sfSeqToNatural, 0);CHKERRQ(ierr); */ 68fa534816SMatthew G. Knepley /* Create the SF associated with this section */ 69e5b44f4fSMatthew G. Knepley ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 70e5b44f4fSMatthew G. Knepley ierr = PetscSectionCreateGlobalSection(sectionDist, sf, PETSC_FALSE, PETSC_TRUE, &gLocSection);CHKERRQ(ierr); 71e5b44f4fSMatthew G. Knepley ierr = PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, gLocSection, &sfField);CHKERRQ(ierr); 72*0c374c54SMatthew G. Knepley ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 73fa534816SMatthew G. Knepley ierr = PetscSFDestroy(&sfEmbed);CHKERRQ(ierr); 74e5b44f4fSMatthew G. Knepley ierr = PetscSectionDestroy(&gLocSection);CHKERRQ(ierr); 75e5b44f4fSMatthew G. Knepley ierr = PetscSectionDestroy(§ionDist);CHKERRQ(ierr); 76e5b44f4fSMatthew G. Knepley /* ierr = PetscPrintf(comm, "Field SF\n");CHKERRQ(ierr); 77e5b44f4fSMatthew G. Knepley ierr = PetscSFView(sfField, 0);CHKERRQ(ierr); */ 78fa534816SMatthew G. Knepley /* Invert the field SF so it's now from distributed to sequential */ 79fa534816SMatthew G. Knepley ierr = PetscSFCreateInverseSF(sfField, &sfFieldInv);CHKERRQ(ierr); 80fa534816SMatthew G. Knepley ierr = PetscSFDestroy(&sfField);CHKERRQ(ierr); 81e5b44f4fSMatthew G. Knepley /* ierr = PetscPrintf(comm, "Inverse Field SF\n");CHKERRQ(ierr); 82e5b44f4fSMatthew G. Knepley ierr = PetscSFView(sfFieldInv, 0);CHKERRQ(ierr); */ 83fa534816SMatthew G. Knepley /* Multiply the sfFieldInv with the */ 84fa534816SMatthew G. Knepley ierr = PetscSFCompose(sfFieldInv, sfSeqToNatural, sfNatural);CHKERRQ(ierr); 85e5b44f4fSMatthew G. Knepley ierr = PetscObjectViewFromOptions((PetscObject) *sfNatural, NULL, "-globaltonatural_sf_view");CHKERRQ(ierr); 86fa534816SMatthew G. Knepley /* Clean up */ 87fa534816SMatthew G. Knepley ierr = PetscSFDestroy(&sfFieldInv);CHKERRQ(ierr); 88fa534816SMatthew G. Knepley ierr = PetscSFDestroy(&sfSeqToNatural);CHKERRQ(ierr); 89fa534816SMatthew G. Knepley PetscFunctionReturn(0); 90fa534816SMatthew G. Knepley } 91fa534816SMatthew G. Knepley 92fa534816SMatthew G. Knepley #undef __FUNCT__ 93fa534816SMatthew G. Knepley #define __FUNCT__ "DMPlexGlobalToNaturalBegin" 94fa534816SMatthew G. Knepley /*@ 95fa534816SMatthew G. Knepley DMPlexGlobalToNaturalBegin - Rearranges a global Vector in the natural order. 96fa534816SMatthew G. Knepley 97fa534816SMatthew G. Knepley Collective on dm 98fa534816SMatthew G. Knepley 99fa534816SMatthew G. Knepley Input Parameters: 100fa534816SMatthew G. Knepley + dm - The distributed DMPlex 101fa534816SMatthew G. Knepley - gv - The global Vec 102fa534816SMatthew G. Knepley 103fa534816SMatthew G. Knepley Output Parameters: 104fa534816SMatthew G. Knepley . nv - Vec in the canonical ordering distributed over all processors associated with gv 105fa534816SMatthew G. Knepley 106fa534816SMatthew G. Knepley Note: The user must call DMPlexSetUseNaturalSF(dm, PETSC_TRUE) before DMPlexDistribute(). 107fa534816SMatthew G. Knepley 108fa534816SMatthew G. Knepley Level: intermediate 109fa534816SMatthew G. Knepley 110fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalEnd() 111fa534816SMatthew G. Knepley @*/ 112fa534816SMatthew G. Knepley PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv) 113fa534816SMatthew G. Knepley { 114fa534816SMatthew G. Knepley const PetscScalar *inarray; 115fa534816SMatthew G. Knepley PetscScalar *outarray; 116fa534816SMatthew G. Knepley PetscErrorCode ierr; 117fa534816SMatthew G. Knepley 118fa534816SMatthew G. Knepley PetscFunctionBegin; 119fa534816SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin,dm,0,0,0);CHKERRQ(ierr); 120fa534816SMatthew G. Knepley if (dm->sfNatural) { 121fa534816SMatthew G. Knepley ierr = VecGetArray(nv, &outarray);CHKERRQ(ierr); 122fa534816SMatthew G. Knepley ierr = VecGetArrayRead(gv, &inarray);CHKERRQ(ierr); 123fa534816SMatthew G. Knepley ierr = PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray);CHKERRQ(ierr); 124fa534816SMatthew G. Knepley ierr = VecRestoreArrayRead(gv, &inarray);CHKERRQ(ierr); 125fa534816SMatthew G. Knepley ierr = VecRestoreArray(nv, &outarray);CHKERRQ(ierr); 126e5b44f4fSMatthew G. Knepley } else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMSetUseNatural() before DMPlexDistribute().\n"); 127fa534816SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin,dm,0,0,0);CHKERRQ(ierr); 128fa534816SMatthew G. Knepley PetscFunctionReturn(0); 129fa534816SMatthew G. Knepley } 130fa534816SMatthew G. Knepley 131fa534816SMatthew G. Knepley #undef __FUNCT__ 132fa534816SMatthew G. Knepley #define __FUNCT__ "DMPlexGlobalToNaturalEnd" 133fa534816SMatthew G. Knepley /*@ 134fa534816SMatthew G. Knepley DMPlexGlobalToNaturalEnd - Rearranges a global Vector in the natural order. 135fa534816SMatthew G. Knepley 136fa534816SMatthew G. Knepley Collective on dm 137fa534816SMatthew G. Knepley 138fa534816SMatthew G. Knepley Input Parameters: 139fa534816SMatthew G. Knepley + dm - The distributed DMPlex 140fa534816SMatthew G. Knepley - gv - The global Vec 141fa534816SMatthew G. Knepley 142fa534816SMatthew G. Knepley Output Parameters: 143fa534816SMatthew G. Knepley . nv - The natural Vec 144fa534816SMatthew G. Knepley 145fa534816SMatthew G. Knepley Note: The user must call DMPlexSetUseNaturalSF(dm, PETSC_TRUE) before DMPlexDistribute(). 146fa534816SMatthew G. Knepley 147fa534816SMatthew G. Knepley Level: intermediate 148fa534816SMatthew G. Knepley 149fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalBegin() 150fa534816SMatthew G. Knepley @*/ 151fa534816SMatthew G. Knepley PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv) 152fa534816SMatthew G. Knepley { 153fa534816SMatthew G. Knepley const PetscScalar *inarray; 154fa534816SMatthew G. Knepley PetscScalar *outarray; 155fa534816SMatthew G. Knepley PetscErrorCode ierr; 156fa534816SMatthew G. Knepley 157fa534816SMatthew G. Knepley PetscFunctionBegin; 158fa534816SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd,dm,0,0,0);CHKERRQ(ierr); 159fa534816SMatthew G. Knepley if (dm->sfNatural) { 160fa534816SMatthew G. Knepley ierr = VecGetArrayRead(gv, &inarray);CHKERRQ(ierr); 161fa534816SMatthew G. Knepley ierr = VecGetArray(nv, &outarray);CHKERRQ(ierr);CHKERRQ(ierr); 162fa534816SMatthew G. Knepley ierr = PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray);CHKERRQ(ierr); 163fa534816SMatthew G. Knepley ierr = VecRestoreArrayRead(gv, &inarray);CHKERRQ(ierr); 164fa534816SMatthew G. Knepley ierr = VecRestoreArray(nv, &outarray);CHKERRQ(ierr); 165fa534816SMatthew G. Knepley } 166fa534816SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd,dm,0,0,0);CHKERRQ(ierr); 167fa534816SMatthew G. Knepley PetscFunctionReturn(0); 168fa534816SMatthew G. Knepley } 169fa534816SMatthew G. Knepley 170fa534816SMatthew G. Knepley #undef __FUNCT__ 171fa534816SMatthew G. Knepley #define __FUNCT__ "DMPlexNaturalToGlobalBegin" 172fa534816SMatthew G. Knepley /*@ 173fa534816SMatthew G. Knepley DMPlexNaturalToGlobalBegin - Rearranges a Vector in the natural order to the Global order. 174fa534816SMatthew G. Knepley 175fa534816SMatthew G. Knepley Collective on dm 176fa534816SMatthew G. Knepley 177fa534816SMatthew G. Knepley Input Parameters: 178fa534816SMatthew G. Knepley + dm - The distributed DMPlex 179fa534816SMatthew G. Knepley - nv - The natural Vec 180fa534816SMatthew G. Knepley 181fa534816SMatthew G. Knepley Output Parameters: 182fa534816SMatthew G. Knepley . gv - The global Vec 183fa534816SMatthew G. Knepley 184fa534816SMatthew G. Knepley Note: The user must call DMPlexSetUseNaturalSF(dm, PETSC_TRUE) before DMPlexDistribute(). 185fa534816SMatthew G. Knepley 186fa534816SMatthew G. Knepley Level: intermediate 187fa534816SMatthew G. Knepley 188fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(),DMPlexGlobalToNaturalEnd() 189fa534816SMatthew G. Knepley @*/ 190fa534816SMatthew G. Knepley PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv) 191fa534816SMatthew G. Knepley { 192fa534816SMatthew G. Knepley const PetscScalar *inarray; 193fa534816SMatthew G. Knepley PetscScalar *outarray; 194fa534816SMatthew G. Knepley PetscErrorCode ierr; 195fa534816SMatthew G. Knepley 196fa534816SMatthew G. Knepley PetscFunctionBegin; 197fa534816SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin,dm,0,0,0);CHKERRQ(ierr); 198fa534816SMatthew G. Knepley if (dm->sfNatural) { 199fa534816SMatthew G. Knepley ierr = VecGetArray(gv, &outarray);CHKERRQ(ierr); 200fa534816SMatthew G. Knepley ierr = VecGetArrayRead(nv, &inarray);CHKERRQ(ierr); 201fa534816SMatthew G. Knepley ierr = PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM);CHKERRQ(ierr); 202fa534816SMatthew G. Knepley ierr = VecRestoreArrayRead(nv, &inarray);CHKERRQ(ierr); 203fa534816SMatthew G. Knepley ierr = VecRestoreArray(gv, &outarray);CHKERRQ(ierr); 204fa534816SMatthew G. Knepley } else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created.\nYou must call DMPlexSetUseNaturalSF() before DMPlexDistribute().\n"); 205fa534816SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin,dm,0,0,0);CHKERRQ(ierr); 206fa534816SMatthew G. Knepley PetscFunctionReturn(0); 207fa534816SMatthew G. Knepley } 208fa534816SMatthew G. Knepley 209fa534816SMatthew G. Knepley #undef __FUNCT__ 210fa534816SMatthew G. Knepley #define __FUNCT__ "DMPlexNaturalToGlobalEnd" 211fa534816SMatthew G. Knepley /*@ 212fa534816SMatthew G. Knepley DMPlexNaturalToGlobalEnd - Rearranges a Vector in the natural order to the Global order. 213fa534816SMatthew G. Knepley 214fa534816SMatthew G. Knepley Collective on dm 215fa534816SMatthew G. Knepley 216fa534816SMatthew G. Knepley Input Parameters: 217fa534816SMatthew G. Knepley + dm - The distributed DMPlex 218fa534816SMatthew G. Knepley - nv - The natural Vec 219fa534816SMatthew G. Knepley 220fa534816SMatthew G. Knepley Output Parameters: 221fa534816SMatthew G. Knepley . gv - The global Vec 222fa534816SMatthew G. Knepley 223fa534816SMatthew G. Knepley Note: The user must call DMPlexSetUseNaturalSF(dm, PETSC_TRUE) before DMPlexDistribute(). 224fa534816SMatthew G. Knepley 225fa534816SMatthew G. Knepley Level: intermediate 226fa534816SMatthew G. Knepley 227fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalBegin() 228fa534816SMatthew G. Knepley @*/ 229fa534816SMatthew G. Knepley PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv) 230fa534816SMatthew G. Knepley { 231fa534816SMatthew G. Knepley const PetscScalar *inarray; 232fa534816SMatthew G. Knepley PetscScalar *outarray; 233fa534816SMatthew G. Knepley PetscErrorCode ierr; 234fa534816SMatthew G. Knepley 235fa534816SMatthew G. Knepley PetscFunctionBegin; 236fa534816SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd,dm,0,0,0);CHKERRQ(ierr); 237fa534816SMatthew G. Knepley if (dm->sfNatural) { 238fa534816SMatthew G. Knepley ierr = VecGetArrayRead(nv, &inarray);CHKERRQ(ierr); 239fa534816SMatthew G. Knepley ierr = VecGetArray(gv, &outarray);CHKERRQ(ierr);CHKERRQ(ierr); 240fa534816SMatthew G. Knepley ierr = PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM);CHKERRQ(ierr); 241fa534816SMatthew G. Knepley ierr = VecRestoreArrayRead(nv, &inarray);CHKERRQ(ierr); 242fa534816SMatthew G. Knepley ierr = VecRestoreArray(gv, &outarray);CHKERRQ(ierr); 243fa534816SMatthew G. Knepley } 244fa534816SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd,dm,0,0,0);CHKERRQ(ierr); 245fa534816SMatthew G. Knepley PetscFunctionReturn(0); 246fa534816SMatthew G. Knepley } 247