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