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; 24*e5b44f4fSMatthew G. Knepley PetscSF sf, sfEmbed, sfSeqToNatural, sfField, sfFieldInv; 25*e5b44f4fSMatthew 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); 32*e5b44f4fSMatthew G. Knepley /* ierr = PetscPrintf(comm, "Point migration SF\n");CHKERRQ(ierr); 33*e5b44f4fSMatthew G. Knepley ierr = PetscSFView(sfMigration, 0);CHKERRQ(ierr); */ 34*e5b44f4fSMatthew G. Knepley /* Create a new section from distributing the original section */ 35*e5b44f4fSMatthew G. Knepley ierr = PetscSectionCreate(comm, §ionDist);CHKERRQ(ierr); 36*e5b44f4fSMatthew G. Knepley ierr = PetscSFDistributeSection(sfMigration, section, &remoteOffsets, sectionDist);CHKERRQ(ierr); 37*e5b44f4fSMatthew G. Knepley /* ierr = PetscPrintf(comm, "Distributed Section\n");CHKERRQ(ierr); 38*e5b44f4fSMatthew G. Knepley ierr = PetscSectionView(sectionDist, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 39*e5b44f4fSMatthew 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); 59*e5b44f4fSMatthew G. Knepley ierr = PetscFree(spoints);CHKERRQ(ierr); 60*e5b44f4fSMatthew G. Knepley /* ierr = PetscPrintf(comm, "Embedded SF\n");CHKERRQ(ierr); 61*e5b44f4fSMatthew 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); 66*e5b44f4fSMatthew G. Knepley /* ierr = PetscPrintf(comm, "Seq-to-Natural SF\n");CHKERRQ(ierr); 67*e5b44f4fSMatthew G. Knepley ierr = PetscSFView(sfSeqToNatural, 0);CHKERRQ(ierr); */ 68fa534816SMatthew G. Knepley /* Create the SF associated with this section */ 69*e5b44f4fSMatthew G. Knepley ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 70*e5b44f4fSMatthew G. Knepley ierr = PetscSectionCreateGlobalSection(sectionDist, sf, PETSC_FALSE, PETSC_TRUE, &gLocSection);CHKERRQ(ierr); 71*e5b44f4fSMatthew G. Knepley ierr = PetscSFCreateSectionSF(sfEmbed, section, remoteOffsets, gLocSection, &sfField);CHKERRQ(ierr); 72fa534816SMatthew G. Knepley ierr = PetscSFDestroy(&sfEmbed);CHKERRQ(ierr); 73*e5b44f4fSMatthew G. Knepley ierr = PetscSectionDestroy(&gLocSection);CHKERRQ(ierr); 74*e5b44f4fSMatthew G. Knepley ierr = PetscSectionDestroy(§ionDist);CHKERRQ(ierr); 75*e5b44f4fSMatthew G. Knepley /* ierr = PetscPrintf(comm, "Field SF\n");CHKERRQ(ierr); 76*e5b44f4fSMatthew G. Knepley ierr = PetscSFView(sfField, 0);CHKERRQ(ierr); */ 77fa534816SMatthew G. Knepley /* Invert the field SF so it's now from distributed to sequential */ 78fa534816SMatthew G. Knepley ierr = PetscSFCreateInverseSF(sfField, &sfFieldInv);CHKERRQ(ierr); 79fa534816SMatthew G. Knepley ierr = PetscSFDestroy(&sfField);CHKERRQ(ierr); 80*e5b44f4fSMatthew G. Knepley /* ierr = PetscPrintf(comm, "Inverse Field SF\n");CHKERRQ(ierr); 81*e5b44f4fSMatthew G. Knepley ierr = PetscSFView(sfFieldInv, 0);CHKERRQ(ierr); */ 82fa534816SMatthew G. Knepley /* Multiply the sfFieldInv with the */ 83fa534816SMatthew G. Knepley ierr = PetscSFCompose(sfFieldInv, sfSeqToNatural, sfNatural);CHKERRQ(ierr); 84*e5b44f4fSMatthew G. Knepley ierr = PetscObjectViewFromOptions((PetscObject) *sfNatural, NULL, "-globaltonatural_sf_view");CHKERRQ(ierr); 85fa534816SMatthew G. Knepley /* Clean up */ 86fa534816SMatthew G. Knepley ierr = PetscSFDestroy(&sfFieldInv);CHKERRQ(ierr); 87fa534816SMatthew G. Knepley ierr = PetscSFDestroy(&sfSeqToNatural);CHKERRQ(ierr); 88fa534816SMatthew G. Knepley PetscFunctionReturn(0); 89fa534816SMatthew G. Knepley } 90fa534816SMatthew G. Knepley 91fa534816SMatthew G. Knepley #undef __FUNCT__ 92fa534816SMatthew G. Knepley #define __FUNCT__ "DMPlexGlobalToNaturalBegin" 93fa534816SMatthew G. Knepley /*@ 94fa534816SMatthew G. Knepley DMPlexGlobalToNaturalBegin - Rearranges a global Vector in the natural order. 95fa534816SMatthew G. Knepley 96fa534816SMatthew G. Knepley Collective on dm 97fa534816SMatthew G. Knepley 98fa534816SMatthew G. Knepley Input Parameters: 99fa534816SMatthew G. Knepley + dm - The distributed DMPlex 100fa534816SMatthew G. Knepley - gv - The global Vec 101fa534816SMatthew G. Knepley 102fa534816SMatthew G. Knepley Output Parameters: 103fa534816SMatthew G. Knepley . nv - Vec in the canonical ordering distributed over all processors associated with gv 104fa534816SMatthew G. Knepley 105fa534816SMatthew G. Knepley Note: The user must call DMPlexSetUseNaturalSF(dm, PETSC_TRUE) before DMPlexDistribute(). 106fa534816SMatthew G. Knepley 107fa534816SMatthew G. Knepley Level: intermediate 108fa534816SMatthew G. Knepley 109fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalEnd() 110fa534816SMatthew G. Knepley @*/ 111fa534816SMatthew G. Knepley PetscErrorCode DMPlexGlobalToNaturalBegin(DM dm, Vec gv, Vec nv) 112fa534816SMatthew G. Knepley { 113fa534816SMatthew G. Knepley const PetscScalar *inarray; 114fa534816SMatthew G. Knepley PetscScalar *outarray; 115fa534816SMatthew G. Knepley PetscErrorCode ierr; 116fa534816SMatthew G. Knepley 117fa534816SMatthew G. Knepley PetscFunctionBegin; 118fa534816SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_GlobalToNaturalBegin,dm,0,0,0);CHKERRQ(ierr); 119fa534816SMatthew G. Knepley if (dm->sfNatural) { 120fa534816SMatthew G. Knepley ierr = VecGetArray(nv, &outarray);CHKERRQ(ierr); 121fa534816SMatthew G. Knepley ierr = VecGetArrayRead(gv, &inarray);CHKERRQ(ierr); 122fa534816SMatthew G. Knepley ierr = PetscSFBcastBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray);CHKERRQ(ierr); 123fa534816SMatthew G. Knepley ierr = VecRestoreArrayRead(gv, &inarray);CHKERRQ(ierr); 124fa534816SMatthew G. Knepley ierr = VecRestoreArray(nv, &outarray);CHKERRQ(ierr); 125*e5b44f4fSMatthew 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"); 126fa534816SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_GlobalToNaturalBegin,dm,0,0,0);CHKERRQ(ierr); 127fa534816SMatthew G. Knepley PetscFunctionReturn(0); 128fa534816SMatthew G. Knepley } 129fa534816SMatthew G. Knepley 130fa534816SMatthew G. Knepley #undef __FUNCT__ 131fa534816SMatthew G. Knepley #define __FUNCT__ "DMPlexGlobalToNaturalEnd" 132fa534816SMatthew G. Knepley /*@ 133fa534816SMatthew G. Knepley DMPlexGlobalToNaturalEnd - Rearranges a global Vector in the natural order. 134fa534816SMatthew G. Knepley 135fa534816SMatthew G. Knepley Collective on dm 136fa534816SMatthew G. Knepley 137fa534816SMatthew G. Knepley Input Parameters: 138fa534816SMatthew G. Knepley + dm - The distributed DMPlex 139fa534816SMatthew G. Knepley - gv - The global Vec 140fa534816SMatthew G. Knepley 141fa534816SMatthew G. Knepley Output Parameters: 142fa534816SMatthew G. Knepley . nv - The natural Vec 143fa534816SMatthew G. Knepley 144fa534816SMatthew G. Knepley Note: The user must call DMPlexSetUseNaturalSF(dm, PETSC_TRUE) before DMPlexDistribute(). 145fa534816SMatthew G. Knepley 146fa534816SMatthew G. Knepley Level: intermediate 147fa534816SMatthew G. Knepley 148fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalBegin() 149fa534816SMatthew G. Knepley @*/ 150fa534816SMatthew G. Knepley PetscErrorCode DMPlexGlobalToNaturalEnd(DM dm, Vec gv, Vec nv) 151fa534816SMatthew G. Knepley { 152fa534816SMatthew G. Knepley const PetscScalar *inarray; 153fa534816SMatthew G. Knepley PetscScalar *outarray; 154fa534816SMatthew G. Knepley PetscErrorCode ierr; 155fa534816SMatthew G. Knepley 156fa534816SMatthew G. Knepley PetscFunctionBegin; 157fa534816SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_GlobalToNaturalEnd,dm,0,0,0);CHKERRQ(ierr); 158fa534816SMatthew G. Knepley if (dm->sfNatural) { 159fa534816SMatthew G. Knepley ierr = VecGetArrayRead(gv, &inarray);CHKERRQ(ierr); 160fa534816SMatthew G. Knepley ierr = VecGetArray(nv, &outarray);CHKERRQ(ierr);CHKERRQ(ierr); 161fa534816SMatthew G. Knepley ierr = PetscSFBcastEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray);CHKERRQ(ierr); 162fa534816SMatthew G. Knepley ierr = VecRestoreArrayRead(gv, &inarray);CHKERRQ(ierr); 163fa534816SMatthew G. Knepley ierr = VecRestoreArray(nv, &outarray);CHKERRQ(ierr); 164fa534816SMatthew G. Knepley } 165fa534816SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_GlobalToNaturalEnd,dm,0,0,0);CHKERRQ(ierr); 166fa534816SMatthew G. Knepley PetscFunctionReturn(0); 167fa534816SMatthew G. Knepley } 168fa534816SMatthew G. Knepley 169fa534816SMatthew G. Knepley #undef __FUNCT__ 170fa534816SMatthew G. Knepley #define __FUNCT__ "DMPlexNaturalToGlobalBegin" 171fa534816SMatthew G. Knepley /*@ 172fa534816SMatthew G. Knepley DMPlexNaturalToGlobalBegin - Rearranges a Vector in the natural order to the Global order. 173fa534816SMatthew G. Knepley 174fa534816SMatthew G. Knepley Collective on dm 175fa534816SMatthew G. Knepley 176fa534816SMatthew G. Knepley Input Parameters: 177fa534816SMatthew G. Knepley + dm - The distributed DMPlex 178fa534816SMatthew G. Knepley - nv - The natural Vec 179fa534816SMatthew G. Knepley 180fa534816SMatthew G. Knepley Output Parameters: 181fa534816SMatthew G. Knepley . gv - The global Vec 182fa534816SMatthew G. Knepley 183fa534816SMatthew G. Knepley Note: The user must call DMPlexSetUseNaturalSF(dm, PETSC_TRUE) before DMPlexDistribute(). 184fa534816SMatthew G. Knepley 185fa534816SMatthew G. Knepley Level: intermediate 186fa534816SMatthew G. Knepley 187fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(),DMPlexGlobalToNaturalEnd() 188fa534816SMatthew G. Knepley @*/ 189fa534816SMatthew G. Knepley PetscErrorCode DMPlexNaturalToGlobalBegin(DM dm, Vec nv, Vec gv) 190fa534816SMatthew G. Knepley { 191fa534816SMatthew G. Knepley const PetscScalar *inarray; 192fa534816SMatthew G. Knepley PetscScalar *outarray; 193fa534816SMatthew G. Knepley PetscErrorCode ierr; 194fa534816SMatthew G. Knepley 195fa534816SMatthew G. Knepley PetscFunctionBegin; 196fa534816SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_NaturalToGlobalBegin,dm,0,0,0);CHKERRQ(ierr); 197fa534816SMatthew G. Knepley if (dm->sfNatural) { 198fa534816SMatthew G. Knepley ierr = VecGetArray(gv, &outarray);CHKERRQ(ierr); 199fa534816SMatthew G. Knepley ierr = VecGetArrayRead(nv, &inarray);CHKERRQ(ierr); 200fa534816SMatthew G. Knepley ierr = PetscSFReduceBegin(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM);CHKERRQ(ierr); 201fa534816SMatthew G. Knepley ierr = VecRestoreArrayRead(nv, &inarray);CHKERRQ(ierr); 202fa534816SMatthew G. Knepley ierr = VecRestoreArray(gv, &outarray);CHKERRQ(ierr); 203fa534816SMatthew 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"); 204fa534816SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_NaturalToGlobalBegin,dm,0,0,0);CHKERRQ(ierr); 205fa534816SMatthew G. Knepley PetscFunctionReturn(0); 206fa534816SMatthew G. Knepley } 207fa534816SMatthew G. Knepley 208fa534816SMatthew G. Knepley #undef __FUNCT__ 209fa534816SMatthew G. Knepley #define __FUNCT__ "DMPlexNaturalToGlobalEnd" 210fa534816SMatthew G. Knepley /*@ 211fa534816SMatthew G. Knepley DMPlexNaturalToGlobalEnd - Rearranges a Vector in the natural order to the Global order. 212fa534816SMatthew G. Knepley 213fa534816SMatthew G. Knepley Collective on dm 214fa534816SMatthew G. Knepley 215fa534816SMatthew G. Knepley Input Parameters: 216fa534816SMatthew G. Knepley + dm - The distributed DMPlex 217fa534816SMatthew G. Knepley - nv - The natural Vec 218fa534816SMatthew G. Knepley 219fa534816SMatthew G. Knepley Output Parameters: 220fa534816SMatthew G. Knepley . gv - The global Vec 221fa534816SMatthew G. Knepley 222fa534816SMatthew G. Knepley Note: The user must call DMPlexSetUseNaturalSF(dm, PETSC_TRUE) before DMPlexDistribute(). 223fa534816SMatthew G. Knepley 224fa534816SMatthew G. Knepley Level: intermediate 225fa534816SMatthew G. Knepley 226fa534816SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexNaturalToGlobalBegin(), DMPlexGlobalToNaturalBegin() 227fa534816SMatthew G. Knepley @*/ 228fa534816SMatthew G. Knepley PetscErrorCode DMPlexNaturalToGlobalEnd(DM dm, Vec nv, Vec gv) 229fa534816SMatthew G. Knepley { 230fa534816SMatthew G. Knepley const PetscScalar *inarray; 231fa534816SMatthew G. Knepley PetscScalar *outarray; 232fa534816SMatthew G. Knepley PetscErrorCode ierr; 233fa534816SMatthew G. Knepley 234fa534816SMatthew G. Knepley PetscFunctionBegin; 235fa534816SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_NaturalToGlobalEnd,dm,0,0,0);CHKERRQ(ierr); 236fa534816SMatthew G. Knepley if (dm->sfNatural) { 237fa534816SMatthew G. Knepley ierr = VecGetArrayRead(nv, &inarray);CHKERRQ(ierr); 238fa534816SMatthew G. Knepley ierr = VecGetArray(gv, &outarray);CHKERRQ(ierr);CHKERRQ(ierr); 239fa534816SMatthew G. Knepley ierr = PetscSFReduceEnd(dm->sfNatural, MPIU_SCALAR, (PetscScalar *) inarray, outarray, MPI_SUM);CHKERRQ(ierr); 240fa534816SMatthew G. Knepley ierr = VecRestoreArrayRead(nv, &inarray);CHKERRQ(ierr); 241fa534816SMatthew G. Knepley ierr = VecRestoreArray(gv, &outarray);CHKERRQ(ierr); 242fa534816SMatthew G. Knepley } 243fa534816SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_NaturalToGlobalEnd,dm,0,0,0);CHKERRQ(ierr); 244fa534816SMatthew G. Knepley PetscFunctionReturn(0); 245fa534816SMatthew G. Knepley } 246