16858538eSMatthew G. Knepley #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/ 26858538eSMatthew G. Knepley 36858538eSMatthew G. Knepley #include <petscdmplex.h> /* For DMProjectCoordinates() */ 46858538eSMatthew G. Knepley #include <petscsf.h> /* For DMLocatePoints() */ 56858538eSMatthew G. Knepley 69371c9d4SSatish Balay PetscErrorCode DMRestrictHook_Coordinates(DM dm, DM dmc, void *ctx) { 76858538eSMatthew G. Knepley DM dm_coord, dmc_coord; 86858538eSMatthew G. Knepley Vec coords, ccoords; 96858538eSMatthew G. Knepley Mat inject; 106858538eSMatthew G. Knepley PetscFunctionBegin; 116858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 126858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dmc, &dmc_coord)); 136858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(dm, &coords)); 146858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(dmc, &ccoords)); 156858538eSMatthew G. Knepley if (coords && !ccoords) { 166858538eSMatthew G. Knepley PetscCall(DMCreateGlobalVector(dmc_coord, &ccoords)); 176858538eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)ccoords, "coordinates")); 186858538eSMatthew G. Knepley PetscCall(DMCreateInjection(dmc_coord, dm_coord, &inject)); 196858538eSMatthew G. Knepley PetscCall(MatRestrict(inject, coords, ccoords)); 206858538eSMatthew G. Knepley PetscCall(MatDestroy(&inject)); 216858538eSMatthew G. Knepley PetscCall(DMSetCoordinates(dmc, ccoords)); 226858538eSMatthew G. Knepley PetscCall(VecDestroy(&ccoords)); 236858538eSMatthew G. Knepley } 246858538eSMatthew G. Knepley PetscFunctionReturn(0); 256858538eSMatthew G. Knepley } 266858538eSMatthew G. Knepley 279371c9d4SSatish Balay static PetscErrorCode DMSubDomainHook_Coordinates(DM dm, DM subdm, void *ctx) { 286858538eSMatthew G. Knepley DM dm_coord, subdm_coord; 296858538eSMatthew G. Knepley Vec coords, ccoords, clcoords; 306858538eSMatthew G. Knepley VecScatter *scat_i, *scat_g; 316858538eSMatthew G. Knepley PetscFunctionBegin; 326858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 336858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(subdm, &subdm_coord)); 346858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(dm, &coords)); 356858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(subdm, &ccoords)); 366858538eSMatthew G. Knepley if (coords && !ccoords) { 376858538eSMatthew G. Knepley PetscCall(DMCreateGlobalVector(subdm_coord, &ccoords)); 386858538eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)ccoords, "coordinates")); 396858538eSMatthew G. Knepley PetscCall(DMCreateLocalVector(subdm_coord, &clcoords)); 406858538eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)clcoords, "coordinates")); 416858538eSMatthew G. Knepley PetscCall(DMCreateDomainDecompositionScatters(dm_coord, 1, &subdm_coord, NULL, &scat_i, &scat_g)); 426858538eSMatthew G. Knepley PetscCall(VecScatterBegin(scat_i[0], coords, ccoords, INSERT_VALUES, SCATTER_FORWARD)); 436858538eSMatthew G. Knepley PetscCall(VecScatterEnd(scat_i[0], coords, ccoords, INSERT_VALUES, SCATTER_FORWARD)); 446858538eSMatthew G. Knepley PetscCall(VecScatterBegin(scat_g[0], coords, clcoords, INSERT_VALUES, SCATTER_FORWARD)); 456858538eSMatthew G. Knepley PetscCall(VecScatterEnd(scat_g[0], coords, clcoords, INSERT_VALUES, SCATTER_FORWARD)); 466858538eSMatthew G. Knepley PetscCall(DMSetCoordinates(subdm, ccoords)); 476858538eSMatthew G. Knepley PetscCall(DMSetCoordinatesLocal(subdm, clcoords)); 486858538eSMatthew G. Knepley PetscCall(VecScatterDestroy(&scat_i[0])); 496858538eSMatthew G. Knepley PetscCall(VecScatterDestroy(&scat_g[0])); 506858538eSMatthew G. Knepley PetscCall(VecDestroy(&ccoords)); 516858538eSMatthew G. Knepley PetscCall(VecDestroy(&clcoords)); 526858538eSMatthew G. Knepley PetscCall(PetscFree(scat_i)); 536858538eSMatthew G. Knepley PetscCall(PetscFree(scat_g)); 546858538eSMatthew G. Knepley } 556858538eSMatthew G. Knepley PetscFunctionReturn(0); 566858538eSMatthew G. Knepley } 576858538eSMatthew G. Knepley 586858538eSMatthew G. Knepley /*@ 596858538eSMatthew G. Knepley DMGetCoordinateDM - Gets the DM that prescribes coordinate layout and scatters between global and local coordinates 606858538eSMatthew G. Knepley 616858538eSMatthew G. Knepley Collective on dm 626858538eSMatthew G. Knepley 636858538eSMatthew G. Knepley Input Parameter: 646858538eSMatthew G. Knepley . dm - the DM 656858538eSMatthew G. Knepley 666858538eSMatthew G. Knepley Output Parameter: 676858538eSMatthew G. Knepley . cdm - coordinate DM 686858538eSMatthew G. Knepley 696858538eSMatthew G. Knepley Level: intermediate 706858538eSMatthew G. Knepley 716858538eSMatthew G. Knepley .seealso: `DMSetCoordinateDM()`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()` 726858538eSMatthew G. Knepley @*/ 739371c9d4SSatish Balay PetscErrorCode DMGetCoordinateDM(DM dm, DM *cdm) { 746858538eSMatthew G. Knepley PetscFunctionBegin; 756858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 766858538eSMatthew G. Knepley PetscValidPointer(cdm, 2); 776858538eSMatthew G. Knepley if (!dm->coordinates[0].dm) { 786858538eSMatthew G. Knepley DM cdm; 796858538eSMatthew G. Knepley 80dbbe0bcdSBarry Smith PetscUseTypeMethod(dm, createcoordinatedm, &cdm); 816858538eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)cdm, "coordinateDM")); 826858538eSMatthew G. Knepley /* Just in case the DM sets the coordinate DM when creating it (DMP4est can do this, because it may not setup 836858538eSMatthew G. Knepley * until the call to CreateCoordinateDM) */ 846858538eSMatthew G. Knepley PetscCall(DMDestroy(&dm->coordinates[0].dm)); 856858538eSMatthew G. Knepley dm->coordinates[0].dm = cdm; 866858538eSMatthew G. Knepley } 876858538eSMatthew G. Knepley *cdm = dm->coordinates[0].dm; 886858538eSMatthew G. Knepley PetscFunctionReturn(0); 896858538eSMatthew G. Knepley } 906858538eSMatthew G. Knepley 916858538eSMatthew G. Knepley /*@ 926858538eSMatthew G. Knepley DMSetCoordinateDM - Sets the DM that prescribes coordinate layout and scatters between global and local coordinates 936858538eSMatthew G. Knepley 946858538eSMatthew G. Knepley Logically Collective on dm 956858538eSMatthew G. Knepley 966858538eSMatthew G. Knepley Input Parameters: 976858538eSMatthew G. Knepley + dm - the DM 986858538eSMatthew G. Knepley - cdm - coordinate DM 996858538eSMatthew G. Knepley 1006858538eSMatthew G. Knepley Level: intermediate 1016858538eSMatthew G. Knepley 1026858538eSMatthew G. Knepley .seealso: `DMGetCoordinateDM()`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()` 1036858538eSMatthew G. Knepley @*/ 1049371c9d4SSatish Balay PetscErrorCode DMSetCoordinateDM(DM dm, DM cdm) { 1056858538eSMatthew G. Knepley PetscFunctionBegin; 1066858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1076858538eSMatthew G. Knepley PetscValidHeaderSpecific(cdm, DM_CLASSID, 2); 1086858538eSMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)cdm)); 1096858538eSMatthew G. Knepley PetscCall(DMDestroy(&dm->coordinates[0].dm)); 1106858538eSMatthew G. Knepley dm->coordinates[0].dm = cdm; 1116858538eSMatthew G. Knepley PetscFunctionReturn(0); 1126858538eSMatthew G. Knepley } 1136858538eSMatthew G. Knepley 1146858538eSMatthew G. Knepley /*@ 1156858538eSMatthew G. Knepley DMGetCellCoordinateDM - Gets the DM that prescribes cellwise coordinate layout and scatters between global and local cellwise coordinates 1166858538eSMatthew G. Knepley 1176858538eSMatthew G. Knepley Collective on dm 1186858538eSMatthew G. Knepley 1196858538eSMatthew G. Knepley Input Parameter: 1206858538eSMatthew G. Knepley . dm - the DM 1216858538eSMatthew G. Knepley 1226858538eSMatthew G. Knepley Output Parameter: 1236858538eSMatthew G. Knepley . cdm - cellwise coordinate DM, or NULL if they are not defined 1246858538eSMatthew G. Knepley 1256858538eSMatthew G. Knepley Note: 1266858538eSMatthew G. Knepley Call DMLocalizeCoordinates() to automatically create cellwise coordinates for periodic geometries. 1276858538eSMatthew G. Knepley 1286858538eSMatthew G. Knepley Level: intermediate 1296858538eSMatthew G. Knepley 1306858538eSMatthew G. Knepley .seealso: `DMSetCellCoordinateDM()`, `DMSetCellCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMLocalizeCoordinates()` 1316858538eSMatthew G. Knepley @*/ 1329371c9d4SSatish Balay PetscErrorCode DMGetCellCoordinateDM(DM dm, DM *cdm) { 1336858538eSMatthew G. Knepley PetscFunctionBegin; 1346858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1356858538eSMatthew G. Knepley PetscValidPointer(cdm, 2); 1366858538eSMatthew G. Knepley *cdm = dm->coordinates[1].dm; 1376858538eSMatthew G. Knepley PetscFunctionReturn(0); 1386858538eSMatthew G. Knepley } 1396858538eSMatthew G. Knepley 1406858538eSMatthew G. Knepley /*@ 1416858538eSMatthew G. Knepley DMSetCellCoordinateDM - Sets the DM that prescribes cellwise coordinate layout and scatters between global and local cellwise coordinates 1426858538eSMatthew G. Knepley 1436858538eSMatthew G. Knepley Logically Collective on dm 1446858538eSMatthew G. Knepley 1456858538eSMatthew G. Knepley Input Parameters: 1466858538eSMatthew G. Knepley + dm - the DM 1476858538eSMatthew G. Knepley - cdm - cellwise coordinate DM 1486858538eSMatthew G. Knepley 1496858538eSMatthew G. Knepley Level: intermediate 1506858538eSMatthew G. Knepley 1516858538eSMatthew G. Knepley .seealso: `DMGetCellCoordinateDM()`, `DMSetCellCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()` 1526858538eSMatthew G. Knepley @*/ 1539371c9d4SSatish Balay PetscErrorCode DMSetCellCoordinateDM(DM dm, DM cdm) { 1546858538eSMatthew G. Knepley PetscInt dim; 1556858538eSMatthew G. Knepley 1566858538eSMatthew G. Knepley PetscFunctionBegin; 1576858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1586858538eSMatthew G. Knepley if (cdm) { 1596858538eSMatthew G. Knepley PetscValidHeaderSpecific(cdm, DM_CLASSID, 2); 1606858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &dim)); 1616858538eSMatthew G. Knepley dm->coordinates[1].dim = dim; 1626858538eSMatthew G. Knepley } 1636858538eSMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)cdm)); 1646858538eSMatthew G. Knepley PetscCall(DMDestroy(&dm->coordinates[1].dm)); 1656858538eSMatthew G. Knepley dm->coordinates[1].dm = cdm; 1666858538eSMatthew G. Knepley PetscFunctionReturn(0); 1676858538eSMatthew G. Knepley } 1686858538eSMatthew G. Knepley 1696858538eSMatthew G. Knepley /*@ 1706858538eSMatthew G. Knepley DMGetCoordinateDim - Retrieve the dimension of embedding space for coordinate values. 1716858538eSMatthew G. Knepley 1726858538eSMatthew G. Knepley Not Collective 1736858538eSMatthew G. Knepley 1746858538eSMatthew G. Knepley Input Parameter: 1756858538eSMatthew G. Knepley . dm - The DM object 1766858538eSMatthew G. Knepley 1776858538eSMatthew G. Knepley Output Parameter: 1786858538eSMatthew G. Knepley . dim - The embedding dimension 1796858538eSMatthew G. Knepley 1806858538eSMatthew G. Knepley Level: intermediate 1816858538eSMatthew G. Knepley 1826858538eSMatthew G. Knepley .seealso: `DMSetCoordinateDim()`, `DMGetCoordinateSection()`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()` 1836858538eSMatthew G. Knepley @*/ 1849371c9d4SSatish Balay PetscErrorCode DMGetCoordinateDim(DM dm, PetscInt *dim) { 1856858538eSMatthew G. Knepley PetscFunctionBegin; 1866858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1876858538eSMatthew G. Knepley PetscValidIntPointer(dim, 2); 1886858538eSMatthew G. Knepley if (dm->coordinates[0].dim == PETSC_DEFAULT) dm->coordinates[0].dim = dm->dim; 1896858538eSMatthew G. Knepley *dim = dm->coordinates[0].dim; 1906858538eSMatthew G. Knepley PetscFunctionReturn(0); 1916858538eSMatthew G. Knepley } 1926858538eSMatthew G. Knepley 1936858538eSMatthew G. Knepley /*@ 1946858538eSMatthew G. Knepley DMSetCoordinateDim - Set the dimension of the embedding space for coordinate values. 1956858538eSMatthew G. Knepley 1966858538eSMatthew G. Knepley Not Collective 1976858538eSMatthew G. Knepley 1986858538eSMatthew G. Knepley Input Parameters: 1996858538eSMatthew G. Knepley + dm - The DM object 2006858538eSMatthew G. Knepley - dim - The embedding dimension 2016858538eSMatthew G. Knepley 2026858538eSMatthew G. Knepley Level: intermediate 2036858538eSMatthew G. Knepley 2046858538eSMatthew G. Knepley .seealso: `DMGetCoordinateDim()`, `DMSetCoordinateSection()`, `DMGetCoordinateSection()`, `DMGetLocalSection()`, `DMSetLocalSection()` 2056858538eSMatthew G. Knepley @*/ 2069371c9d4SSatish Balay PetscErrorCode DMSetCoordinateDim(DM dm, PetscInt dim) { 2076858538eSMatthew G. Knepley PetscDS ds; 2086858538eSMatthew G. Knepley PetscInt Nds, n; 2096858538eSMatthew G. Knepley 2106858538eSMatthew G. Knepley PetscFunctionBegin; 2116858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2126858538eSMatthew G. Knepley dm->coordinates[0].dim = dim; 2136858538eSMatthew G. Knepley if (dm->dim >= 0) { 2146858538eSMatthew G. Knepley PetscCall(DMGetNumDS(dm, &Nds)); 2156858538eSMatthew G. Knepley for (n = 0; n < Nds; ++n) { 2166858538eSMatthew G. Knepley PetscCall(DMGetRegionNumDS(dm, n, NULL, NULL, &ds)); 2176858538eSMatthew G. Knepley PetscCall(PetscDSSetCoordinateDimension(ds, dim)); 2186858538eSMatthew G. Knepley } 2196858538eSMatthew G. Knepley } 2206858538eSMatthew G. Knepley PetscFunctionReturn(0); 2216858538eSMatthew G. Knepley } 2226858538eSMatthew G. Knepley 2236858538eSMatthew G. Knepley /*@ 2246858538eSMatthew G. Knepley DMGetCoordinateSection - Retrieve the layout of coordinate values over the mesh. 2256858538eSMatthew G. Knepley 2266858538eSMatthew G. Knepley Collective on dm 2276858538eSMatthew G. Knepley 2286858538eSMatthew G. Knepley Input Parameter: 2296858538eSMatthew G. Knepley . dm - The DM object 2306858538eSMatthew G. Knepley 2316858538eSMatthew G. Knepley Output Parameter: 2326858538eSMatthew G. Knepley . section - The PetscSection object 2336858538eSMatthew G. Knepley 2346858538eSMatthew G. Knepley Level: intermediate 2356858538eSMatthew G. Knepley 2366858538eSMatthew G. Knepley .seealso: `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()` 2376858538eSMatthew G. Knepley @*/ 2389371c9d4SSatish Balay PetscErrorCode DMGetCoordinateSection(DM dm, PetscSection *section) { 2396858538eSMatthew G. Knepley DM cdm; 2406858538eSMatthew G. Knepley 2416858538eSMatthew G. Knepley PetscFunctionBegin; 2426858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2436858538eSMatthew G. Knepley PetscValidPointer(section, 2); 2446858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &cdm)); 2456858538eSMatthew G. Knepley PetscCall(DMGetLocalSection(cdm, section)); 2466858538eSMatthew G. Knepley PetscFunctionReturn(0); 2476858538eSMatthew G. Knepley } 2486858538eSMatthew G. Knepley 2496858538eSMatthew G. Knepley /*@ 2506858538eSMatthew G. Knepley DMSetCoordinateSection - Set the layout of coordinate values over the mesh. 2516858538eSMatthew G. Knepley 2526858538eSMatthew G. Knepley Not Collective 2536858538eSMatthew G. Knepley 2546858538eSMatthew G. Knepley Input Parameters: 2556858538eSMatthew G. Knepley + dm - The DM object 2566858538eSMatthew G. Knepley . dim - The embedding dimension, or PETSC_DETERMINE 2576858538eSMatthew G. Knepley - section - The PetscSection object 2586858538eSMatthew G. Knepley 2596858538eSMatthew G. Knepley Level: intermediate 2606858538eSMatthew G. Knepley 2616858538eSMatthew G. Knepley .seealso: `DMGetCoordinateSection()`, `DMGetLocalSection()`, `DMSetLocalSection()` 2626858538eSMatthew G. Knepley @*/ 2639371c9d4SSatish Balay PetscErrorCode DMSetCoordinateSection(DM dm, PetscInt dim, PetscSection section) { 2646858538eSMatthew G. Knepley DM cdm; 2656858538eSMatthew G. Knepley 2666858538eSMatthew G. Knepley PetscFunctionBegin; 2676858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2686858538eSMatthew G. Knepley PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 3); 2696858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &cdm)); 2706858538eSMatthew G. Knepley PetscCall(DMSetLocalSection(cdm, section)); 2716858538eSMatthew G. Knepley if (dim == PETSC_DETERMINE) { 2726858538eSMatthew G. Knepley PetscInt d = PETSC_DEFAULT; 2736858538eSMatthew G. Knepley PetscInt pStart, pEnd, vStart, vEnd, v, dd; 2746858538eSMatthew G. Knepley 2756858538eSMatthew G. Knepley PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 2766858538eSMatthew G. Knepley PetscCall(DMGetDimPoints(dm, 0, &vStart, &vEnd)); 2776858538eSMatthew G. Knepley pStart = PetscMax(vStart, pStart); 2786858538eSMatthew G. Knepley pEnd = PetscMin(vEnd, pEnd); 2796858538eSMatthew G. Knepley for (v = pStart; v < pEnd; ++v) { 2806858538eSMatthew G. Knepley PetscCall(PetscSectionGetDof(section, v, &dd)); 2819371c9d4SSatish Balay if (dd) { 2829371c9d4SSatish Balay d = dd; 2839371c9d4SSatish Balay break; 2849371c9d4SSatish Balay } 2856858538eSMatthew G. Knepley } 2866858538eSMatthew G. Knepley if (d >= 0) PetscCall(DMSetCoordinateDim(dm, d)); 2876858538eSMatthew G. Knepley } 2886858538eSMatthew G. Knepley PetscFunctionReturn(0); 2896858538eSMatthew G. Knepley } 2906858538eSMatthew G. Knepley 2916858538eSMatthew G. Knepley /*@ 2926858538eSMatthew G. Knepley DMGetCellCoordinateSection - Retrieve the layout of cellwise coordinate values over the mesh. 2936858538eSMatthew G. Knepley 2946858538eSMatthew G. Knepley Collective on dm 2956858538eSMatthew G. Knepley 2966858538eSMatthew G. Knepley Input Parameter: 2976858538eSMatthew G. Knepley . dm - The DM object 2986858538eSMatthew G. Knepley 2996858538eSMatthew G. Knepley Output Parameter: 3006858538eSMatthew G. Knepley . section - The PetscSection object, or NULL if no cellwise coordinates are defined 3016858538eSMatthew G. Knepley 3026858538eSMatthew G. Knepley Level: intermediate 3036858538eSMatthew G. Knepley 3046858538eSMatthew G. Knepley .seealso: `DMGetCoordinateSection()`, `DMSetCellCoordinateSection()`, `DMGetCellCoordinateDM()`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()` 3056858538eSMatthew G. Knepley @*/ 3069371c9d4SSatish Balay PetscErrorCode DMGetCellCoordinateSection(DM dm, PetscSection *section) { 3076858538eSMatthew G. Knepley DM cdm; 3086858538eSMatthew G. Knepley 3096858538eSMatthew G. Knepley PetscFunctionBegin; 3106858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3116858538eSMatthew G. Knepley PetscValidPointer(section, 2); 3126858538eSMatthew G. Knepley *section = NULL; 3136858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 3146858538eSMatthew G. Knepley if (cdm) PetscCall(DMGetLocalSection(cdm, section)); 3156858538eSMatthew G. Knepley PetscFunctionReturn(0); 3166858538eSMatthew G. Knepley } 3176858538eSMatthew G. Knepley 3186858538eSMatthew G. Knepley /*@ 3196858538eSMatthew G. Knepley DMSetCellCoordinateSection - Set the layout of cellwise coordinate values over the mesh. 3206858538eSMatthew G. Knepley 3216858538eSMatthew G. Knepley Not Collective 3226858538eSMatthew G. Knepley 3236858538eSMatthew G. Knepley Input Parameters: 3246858538eSMatthew G. Knepley + dm - The DM object 3256858538eSMatthew G. Knepley . dim - The embedding dimension, or PETSC_DETERMINE 3266858538eSMatthew G. Knepley - section - The PetscSection object for a cellwise layout 3276858538eSMatthew G. Knepley 3286858538eSMatthew G. Knepley Level: intermediate 3296858538eSMatthew G. Knepley 3306858538eSMatthew G. Knepley .seealso: `DMSetCoordinateSection()`, `DMGetCellCoordinateSection()`, `DMGetCoordinateSection()`, `DMGetCellCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()` 3316858538eSMatthew G. Knepley @*/ 3329371c9d4SSatish Balay PetscErrorCode DMSetCellCoordinateSection(DM dm, PetscInt dim, PetscSection section) { 3336858538eSMatthew G. Knepley DM cdm; 3346858538eSMatthew G. Knepley 3356858538eSMatthew G. Knepley PetscFunctionBegin; 3366858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3376858538eSMatthew G. Knepley PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 3); 3386858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 3396858538eSMatthew G. Knepley PetscCheck(cdm, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "No DM defined for cellwise coordinates"); 3406858538eSMatthew G. Knepley PetscCall(DMSetLocalSection(cdm, section)); 3416858538eSMatthew G. Knepley if (dim == PETSC_DETERMINE) { 3426858538eSMatthew G. Knepley PetscInt d = PETSC_DEFAULT; 3436858538eSMatthew G. Knepley PetscInt pStart, pEnd, vStart, vEnd, v, dd; 3446858538eSMatthew G. Knepley 3456858538eSMatthew G. Knepley PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 3466858538eSMatthew G. Knepley PetscCall(DMGetDimPoints(dm, 0, &vStart, &vEnd)); 3476858538eSMatthew G. Knepley pStart = PetscMax(vStart, pStart); 3486858538eSMatthew G. Knepley pEnd = PetscMin(vEnd, pEnd); 3496858538eSMatthew G. Knepley for (v = pStart; v < pEnd; ++v) { 3506858538eSMatthew G. Knepley PetscCall(PetscSectionGetDof(section, v, &dd)); 3519371c9d4SSatish Balay if (dd) { 3529371c9d4SSatish Balay d = dd; 3539371c9d4SSatish Balay break; 3549371c9d4SSatish Balay } 3556858538eSMatthew G. Knepley } 3566858538eSMatthew G. Knepley if (d >= 0) PetscCall(DMSetCoordinateDim(dm, d)); 3576858538eSMatthew G. Knepley } 3586858538eSMatthew G. Knepley PetscFunctionReturn(0); 3596858538eSMatthew G. Knepley } 3606858538eSMatthew G. Knepley 3616858538eSMatthew G. Knepley /*@ 3626858538eSMatthew G. Knepley DMGetCoordinates - Gets a global vector with the coordinates associated with the DM. 3636858538eSMatthew G. Knepley 3646858538eSMatthew G. Knepley Collective on dm 3656858538eSMatthew G. Knepley 3666858538eSMatthew G. Knepley Input Parameter: 3676858538eSMatthew G. Knepley . dm - the DM 3686858538eSMatthew G. Knepley 3696858538eSMatthew G. Knepley Output Parameter: 3706858538eSMatthew G. Knepley . c - global coordinate vector 3716858538eSMatthew G. Knepley 3726858538eSMatthew G. Knepley Note: 3736858538eSMatthew G. Knepley This is a borrowed reference, so the user should NOT destroy this vector. When the DM is 3746858538eSMatthew G. Knepley destroyed the array will no longer be valid. 3756858538eSMatthew G. Knepley 3766858538eSMatthew G. Knepley Each process has only the locally-owned portion of the global coordinates (does NOT have the ghost coordinates). 3776858538eSMatthew G. Knepley 3786858538eSMatthew G. Knepley For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...) 3796858538eSMatthew G. Knepley and (x_0,y_0,z_0,x_1,y_1,z_1...) 3806858538eSMatthew G. Knepley 3816858538eSMatthew G. Knepley Level: intermediate 3826858538eSMatthew G. Knepley 3836858538eSMatthew G. Knepley .seealso: `DMSetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCoordinateDM()`, `DMDASetUniformCoordinates()` 3846858538eSMatthew G. Knepley @*/ 3859371c9d4SSatish Balay PetscErrorCode DMGetCoordinates(DM dm, Vec *c) { 3866858538eSMatthew G. Knepley PetscFunctionBegin; 3876858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3886858538eSMatthew G. Knepley PetscValidPointer(c, 2); 3896858538eSMatthew G. Knepley if (!dm->coordinates[0].x && dm->coordinates[0].xl) { 3906858538eSMatthew G. Knepley DM cdm = NULL; 3916858538eSMatthew G. Knepley 3926858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &cdm)); 3936858538eSMatthew G. Knepley PetscCall(DMCreateGlobalVector(cdm, &dm->coordinates[0].x)); 3946858538eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].x, "coordinates")); 3956858538eSMatthew G. Knepley PetscCall(DMLocalToGlobalBegin(cdm, dm->coordinates[0].xl, INSERT_VALUES, dm->coordinates[0].x)); 3966858538eSMatthew G. Knepley PetscCall(DMLocalToGlobalEnd(cdm, dm->coordinates[0].xl, INSERT_VALUES, dm->coordinates[0].x)); 3976858538eSMatthew G. Knepley } 3986858538eSMatthew G. Knepley *c = dm->coordinates[0].x; 3996858538eSMatthew G. Knepley PetscFunctionReturn(0); 4006858538eSMatthew G. Knepley } 4016858538eSMatthew G. Knepley 4026858538eSMatthew G. Knepley /*@ 4036858538eSMatthew G. Knepley DMSetCoordinates - Sets into the DM a global vector that holds the coordinates 4046858538eSMatthew G. Knepley 4056858538eSMatthew G. Knepley Collective on dm 4066858538eSMatthew G. Knepley 4076858538eSMatthew G. Knepley Input Parameters: 4086858538eSMatthew G. Knepley + dm - the DM 4096858538eSMatthew G. Knepley - c - coordinate vector 4106858538eSMatthew G. Knepley 4116858538eSMatthew G. Knepley Notes: 4126858538eSMatthew G. Knepley The coordinates do include those for ghost points, which are in the local vector. 4136858538eSMatthew G. Knepley 4146858538eSMatthew G. Knepley The vector c should be destroyed by the caller. 4156858538eSMatthew G. Knepley 4166858538eSMatthew G. Knepley Level: intermediate 4176858538eSMatthew G. Knepley 4186858538eSMatthew G. Knepley .seealso: `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCoordinateDM()`, `DMDASetUniformCoordinates()` 4196858538eSMatthew G. Knepley @*/ 4209371c9d4SSatish Balay PetscErrorCode DMSetCoordinates(DM dm, Vec c) { 4216858538eSMatthew G. Knepley PetscFunctionBegin; 4226858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4236858538eSMatthew G. Knepley if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2); 4246858538eSMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)c)); 4256858538eSMatthew G. Knepley PetscCall(VecDestroy(&dm->coordinates[0].x)); 4266858538eSMatthew G. Knepley dm->coordinates[0].x = c; 4276858538eSMatthew G. Knepley PetscCall(VecDestroy(&dm->coordinates[0].xl)); 4286858538eSMatthew G. Knepley PetscCall(DMCoarsenHookAdd(dm, DMRestrictHook_Coordinates, NULL, NULL)); 4296858538eSMatthew G. Knepley PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_Coordinates, NULL, NULL)); 4306858538eSMatthew G. Knepley PetscFunctionReturn(0); 4316858538eSMatthew G. Knepley } 4326858538eSMatthew G. Knepley 4336858538eSMatthew G. Knepley /*@ 4346858538eSMatthew G. Knepley DMGetCellCoordinates - Gets a global vector with the cellwise coordinates associated with the DM. 4356858538eSMatthew G. Knepley 4366858538eSMatthew G. Knepley Collective on dm 4376858538eSMatthew G. Knepley 4386858538eSMatthew G. Knepley Input Parameter: 4396858538eSMatthew G. Knepley . dm - the DM 4406858538eSMatthew G. Knepley 4416858538eSMatthew G. Knepley Output Parameter: 4426858538eSMatthew G. Knepley . c - global coordinate vector 4436858538eSMatthew G. Knepley 4446858538eSMatthew G. Knepley Note: 4456858538eSMatthew G. Knepley This is a borrowed reference, so the user should NOT destroy this vector. When the DM is 4466858538eSMatthew G. Knepley destroyed the array will no longer be valid. 4476858538eSMatthew G. Knepley 4486858538eSMatthew G. Knepley Each process has only the locally-owned portion of the global coordinates (does NOT have the ghost coordinates). 4496858538eSMatthew G. Knepley 4506858538eSMatthew G. Knepley Level: intermediate 4516858538eSMatthew G. Knepley 4526858538eSMatthew G. Knepley .seealso: `DMSetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMGetCellCoordinateDM()` 4536858538eSMatthew G. Knepley @*/ 4549371c9d4SSatish Balay PetscErrorCode DMGetCellCoordinates(DM dm, Vec *c) { 4556858538eSMatthew G. Knepley PetscFunctionBegin; 4566858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4576858538eSMatthew G. Knepley PetscValidPointer(c, 2); 4586858538eSMatthew G. Knepley if (!dm->coordinates[1].x && dm->coordinates[1].xl) { 4596858538eSMatthew G. Knepley DM cdm = NULL; 4606858538eSMatthew G. Knepley 46111ea91a7SMatthew G. Knepley PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 4626858538eSMatthew G. Knepley PetscCall(DMCreateGlobalVector(cdm, &dm->coordinates[1].x)); 4636858538eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].x, "DG coordinates")); 4646858538eSMatthew G. Knepley PetscCall(DMLocalToGlobalBegin(cdm, dm->coordinates[1].xl, INSERT_VALUES, dm->coordinates[1].x)); 4656858538eSMatthew G. Knepley PetscCall(DMLocalToGlobalEnd(cdm, dm->coordinates[1].xl, INSERT_VALUES, dm->coordinates[1].x)); 4666858538eSMatthew G. Knepley } 4676858538eSMatthew G. Knepley *c = dm->coordinates[1].x; 4686858538eSMatthew G. Knepley PetscFunctionReturn(0); 4696858538eSMatthew G. Knepley } 4706858538eSMatthew G. Knepley 4716858538eSMatthew G. Knepley /*@ 4726858538eSMatthew G. Knepley DMSetCellCoordinates - Sets into the DM a global vector that holds the cellwise coordinates 4736858538eSMatthew G. Knepley 4746858538eSMatthew G. Knepley Collective on dm 4756858538eSMatthew G. Knepley 4766858538eSMatthew G. Knepley Input Parameters: 4776858538eSMatthew G. Knepley + dm - the DM 4786858538eSMatthew G. Knepley - c - cellwise coordinate vector 4796858538eSMatthew G. Knepley 4806858538eSMatthew G. Knepley Notes: 4816858538eSMatthew G. Knepley The coordinates do include those for ghost points, which are in the local vector. 4826858538eSMatthew G. Knepley 4836858538eSMatthew G. Knepley The vector c should be destroyed by the caller. 4846858538eSMatthew G. Knepley 4856858538eSMatthew G. Knepley Level: intermediate 4866858538eSMatthew G. Knepley 4876858538eSMatthew G. Knepley .seealso: `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMGetCellCoordinateDM()` 4886858538eSMatthew G. Knepley @*/ 4899371c9d4SSatish Balay PetscErrorCode DMSetCellCoordinates(DM dm, Vec c) { 4906858538eSMatthew G. Knepley PetscFunctionBegin; 4916858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4926858538eSMatthew G. Knepley if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2); 4936858538eSMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)c)); 4946858538eSMatthew G. Knepley PetscCall(VecDestroy(&dm->coordinates[1].x)); 4956858538eSMatthew G. Knepley dm->coordinates[1].x = c; 4966858538eSMatthew G. Knepley PetscCall(VecDestroy(&dm->coordinates[1].xl)); 4976858538eSMatthew G. Knepley PetscFunctionReturn(0); 4986858538eSMatthew G. Knepley } 4996858538eSMatthew G. Knepley 5006858538eSMatthew G. Knepley /*@ 5016858538eSMatthew G. Knepley DMGetCoordinatesLocalSetUp - Prepares a local vector of coordinates, so that DMGetCoordinatesLocalNoncollective() can be used as non-collective afterwards. 5026858538eSMatthew G. Knepley 5036858538eSMatthew G. Knepley Collective on dm 5046858538eSMatthew G. Knepley 5056858538eSMatthew G. Knepley Input Parameter: 5066858538eSMatthew G. Knepley . dm - the DM 5076858538eSMatthew G. Knepley 5086858538eSMatthew G. Knepley Level: advanced 5096858538eSMatthew G. Knepley 5106858538eSMatthew G. Knepley .seealso: `DMGetCoordinatesLocalNoncollective()` 5116858538eSMatthew G. Knepley @*/ 5129371c9d4SSatish Balay PetscErrorCode DMGetCoordinatesLocalSetUp(DM dm) { 5136858538eSMatthew G. Knepley PetscFunctionBegin; 5146858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5156858538eSMatthew G. Knepley if (!dm->coordinates[0].xl && dm->coordinates[0].x) { 5166858538eSMatthew G. Knepley DM cdm = NULL; 5176858538eSMatthew G. Knepley 5186858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &cdm)); 5196858538eSMatthew G. Knepley PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[0].xl)); 5206858538eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].xl, "coordinates")); 5216858538eSMatthew G. Knepley PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl)); 5226858538eSMatthew G. Knepley PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl)); 5236858538eSMatthew G. Knepley } 5246858538eSMatthew G. Knepley PetscFunctionReturn(0); 5256858538eSMatthew G. Knepley } 5266858538eSMatthew G. Knepley 5276858538eSMatthew G. Knepley /*@ 5286858538eSMatthew G. Knepley DMGetCoordinatesLocal - Gets a local vector with the coordinates associated with the DM. 5296858538eSMatthew G. Knepley 5306858538eSMatthew G. Knepley Collective on dm 5316858538eSMatthew G. Knepley 5326858538eSMatthew G. Knepley Input Parameter: 5336858538eSMatthew G. Knepley . dm - the DM 5346858538eSMatthew G. Knepley 5356858538eSMatthew G. Knepley Output Parameter: 5366858538eSMatthew G. Knepley . c - coordinate vector 5376858538eSMatthew G. Knepley 5386858538eSMatthew G. Knepley Note: 5396858538eSMatthew G. Knepley This is a borrowed reference, so the user should NOT destroy this vector 5406858538eSMatthew G. Knepley 5416858538eSMatthew G. Knepley Each process has the local and ghost coordinates 5426858538eSMatthew G. Knepley 5436858538eSMatthew G. Knepley For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...) 5446858538eSMatthew G. Knepley and (x_0,y_0,z_0,x_1,y_1,z_1...) 5456858538eSMatthew G. Knepley 5466858538eSMatthew G. Knepley Level: intermediate 5476858538eSMatthew G. Knepley 5486858538eSMatthew G. Knepley .seealso: `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`, `DMGetCoordinatesLocalNoncollective()` 5496858538eSMatthew G. Knepley @*/ 5509371c9d4SSatish Balay PetscErrorCode DMGetCoordinatesLocal(DM dm, Vec *c) { 5516858538eSMatthew G. Knepley PetscFunctionBegin; 5526858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5536858538eSMatthew G. Knepley PetscValidPointer(c, 2); 5546858538eSMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(dm)); 5556858538eSMatthew G. Knepley *c = dm->coordinates[0].xl; 5566858538eSMatthew G. Knepley PetscFunctionReturn(0); 5576858538eSMatthew G. Knepley } 5586858538eSMatthew G. Knepley 5596858538eSMatthew G. Knepley /*@ 5606858538eSMatthew G. Knepley DMGetCoordinatesLocalNoncollective - Non-collective version of DMGetCoordinatesLocal(). Fails if global coordinates have been set and DMGetCoordinatesLocalSetUp() not called. 5616858538eSMatthew G. Knepley 5626858538eSMatthew G. Knepley Not collective 5636858538eSMatthew G. Knepley 5646858538eSMatthew G. Knepley Input Parameter: 5656858538eSMatthew G. Knepley . dm - the DM 5666858538eSMatthew G. Knepley 5676858538eSMatthew G. Knepley Output Parameter: 5686858538eSMatthew G. Knepley . c - coordinate vector 5696858538eSMatthew G. Knepley 5706858538eSMatthew G. Knepley Level: advanced 5716858538eSMatthew G. Knepley 5726858538eSMatthew G. Knepley .seealso: `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinatesLocal()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()` 5736858538eSMatthew G. Knepley @*/ 5749371c9d4SSatish Balay PetscErrorCode DMGetCoordinatesLocalNoncollective(DM dm, Vec *c) { 5756858538eSMatthew G. Knepley PetscFunctionBegin; 5766858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5776858538eSMatthew G. Knepley PetscValidPointer(c, 2); 5786858538eSMatthew G. Knepley PetscCheck(dm->coordinates[0].xl || !dm->coordinates[0].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called"); 5796858538eSMatthew G. Knepley *c = dm->coordinates[0].xl; 5806858538eSMatthew G. Knepley PetscFunctionReturn(0); 5816858538eSMatthew G. Knepley } 5826858538eSMatthew G. Knepley 5836858538eSMatthew G. Knepley /*@ 5846858538eSMatthew G. Knepley DMGetCoordinatesLocalTuple - Gets a local vector with the coordinates of specified points and section describing its layout. 5856858538eSMatthew G. Knepley 5866858538eSMatthew G. Knepley Not collective 5876858538eSMatthew G. Knepley 5886858538eSMatthew G. Knepley Input Parameters: 5896858538eSMatthew G. Knepley + dm - the DM 5906858538eSMatthew G. Knepley - p - the IS of points whose coordinates will be returned 5916858538eSMatthew G. Knepley 5926858538eSMatthew G. Knepley Output Parameters: 5936858538eSMatthew G. Knepley + pCoordSection - the PetscSection describing the layout of pCoord, i.e. each point corresponds to one point in p, and DOFs correspond to coordinates 5946858538eSMatthew G. Knepley - pCoord - the Vec with coordinates of points in p 5956858538eSMatthew G. Knepley 5966858538eSMatthew G. Knepley Note: 5976858538eSMatthew G. Knepley DMGetCoordinatesLocalSetUp() must be called first. This function employs DMGetCoordinatesLocalNoncollective() so it is not collective. 5986858538eSMatthew G. Knepley 5996858538eSMatthew G. Knepley This creates a new vector, so the user SHOULD destroy this vector 6006858538eSMatthew G. Knepley 6016858538eSMatthew G. Knepley Each process has the local and ghost coordinates 6026858538eSMatthew G. Knepley 6036858538eSMatthew G. Knepley For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...) 6046858538eSMatthew G. Knepley and (x_0,y_0,z_0,x_1,y_1,z_1...) 6056858538eSMatthew G. Knepley 6066858538eSMatthew G. Knepley Level: advanced 6076858538eSMatthew G. Knepley 6086858538eSMatthew G. Knepley .seealso: `DMSetCoordinatesLocal()`, `DMGetCoordinatesLocal()`, `DMGetCoordinatesLocalNoncollective()`, `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()` 6096858538eSMatthew G. Knepley @*/ 6109371c9d4SSatish Balay PetscErrorCode DMGetCoordinatesLocalTuple(DM dm, IS p, PetscSection *pCoordSection, Vec *pCoord) { 6116858538eSMatthew G. Knepley DM cdm; 6126858538eSMatthew G. Knepley PetscSection cs, newcs; 6136858538eSMatthew G. Knepley Vec coords; 6146858538eSMatthew G. Knepley const PetscScalar *arr; 6156858538eSMatthew G. Knepley PetscScalar *newarr = NULL; 6166858538eSMatthew G. Knepley PetscInt n; 6176858538eSMatthew G. Knepley 6186858538eSMatthew G. Knepley PetscFunctionBegin; 6196858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 6206858538eSMatthew G. Knepley PetscValidHeaderSpecific(p, IS_CLASSID, 2); 6216858538eSMatthew G. Knepley if (pCoordSection) PetscValidPointer(pCoordSection, 3); 6226858538eSMatthew G. Knepley if (pCoord) PetscValidPointer(pCoord, 4); 6236858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &cdm)); 6246858538eSMatthew G. Knepley PetscCall(DMGetLocalSection(cdm, &cs)); 6256858538eSMatthew G. Knepley PetscCall(DMGetCoordinatesLocal(dm, &coords)); 6266858538eSMatthew G. Knepley PetscCheck(coords, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called or coordinates not set"); 6276858538eSMatthew G. Knepley PetscCheck(cdm && cs, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM not supported"); 6286858538eSMatthew G. Knepley PetscCall(VecGetArrayRead(coords, &arr)); 6296858538eSMatthew G. Knepley PetscCall(PetscSectionExtractDofsFromArray(cs, MPIU_SCALAR, arr, p, &newcs, pCoord ? ((void **)&newarr) : NULL)); 6306858538eSMatthew G. Knepley PetscCall(VecRestoreArrayRead(coords, &arr)); 6316858538eSMatthew G. Knepley if (pCoord) { 6326858538eSMatthew G. Knepley PetscCall(PetscSectionGetStorageSize(newcs, &n)); 6336858538eSMatthew G. Knepley /* set array in two steps to mimic PETSC_OWN_POINTER */ 6346858538eSMatthew G. Knepley PetscCall(VecCreateSeqWithArray(PetscObjectComm((PetscObject)p), 1, n, NULL, pCoord)); 6356858538eSMatthew G. Knepley PetscCall(VecReplaceArray(*pCoord, newarr)); 6366858538eSMatthew G. Knepley } else { 6376858538eSMatthew G. Knepley PetscCall(PetscFree(newarr)); 6386858538eSMatthew G. Knepley } 6399371c9d4SSatish Balay if (pCoordSection) { 6409371c9d4SSatish Balay *pCoordSection = newcs; 6419371c9d4SSatish Balay } else PetscCall(PetscSectionDestroy(&newcs)); 6426858538eSMatthew G. Knepley PetscFunctionReturn(0); 6436858538eSMatthew G. Knepley } 6446858538eSMatthew G. Knepley 6456858538eSMatthew G. Knepley /*@ 6466858538eSMatthew G. Knepley DMSetCoordinatesLocal - Sets into the DM a local vector that holds the coordinates 6476858538eSMatthew G. Knepley 6486858538eSMatthew G. Knepley Not collective 6496858538eSMatthew G. Knepley 6506858538eSMatthew G. Knepley Input Parameters: 6516858538eSMatthew G. Knepley + dm - the DM 6526858538eSMatthew G. Knepley - c - coordinate vector 6536858538eSMatthew G. Knepley 6546858538eSMatthew G. Knepley Notes: 6556858538eSMatthew G. Knepley The coordinates of ghost points can be set using DMSetCoordinates() 6566858538eSMatthew G. Knepley followed by DMGetCoordinatesLocal(). This is intended to enable the 6576858538eSMatthew G. Knepley setting of ghost coordinates outside of the domain. 6586858538eSMatthew G. Knepley 6596858538eSMatthew G. Knepley The vector c should be destroyed by the caller. 6606858538eSMatthew G. Knepley 6616858538eSMatthew G. Knepley Level: intermediate 6626858538eSMatthew G. Knepley 6636858538eSMatthew G. Knepley .seealso: `DMGetCoordinatesLocal()`, `DMSetCoordinates()`, `DMGetCoordinates()`, `DMGetCoordinateDM()` 6646858538eSMatthew G. Knepley @*/ 6659371c9d4SSatish Balay PetscErrorCode DMSetCoordinatesLocal(DM dm, Vec c) { 6666858538eSMatthew G. Knepley PetscFunctionBegin; 6676858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 6686858538eSMatthew G. Knepley if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2); 6696858538eSMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)c)); 6706858538eSMatthew G. Knepley PetscCall(VecDestroy(&dm->coordinates[0].xl)); 6716858538eSMatthew G. Knepley dm->coordinates[0].xl = c; 6726858538eSMatthew G. Knepley PetscCall(VecDestroy(&dm->coordinates[0].x)); 6736858538eSMatthew G. Knepley PetscFunctionReturn(0); 6746858538eSMatthew G. Knepley } 6756858538eSMatthew G. Knepley 6766858538eSMatthew G. Knepley /*@ 6776858538eSMatthew G. Knepley DMGetCellCoordinatesLocalSetUp - Prepares a local vector of cellwise coordinates, so that DMGetCellCoordinatesLocalNoncollective() can be used as non-collective afterwards. 6786858538eSMatthew G. Knepley 6796858538eSMatthew G. Knepley Collective on dm 6806858538eSMatthew G. Knepley 6816858538eSMatthew G. Knepley Input Parameter: 6826858538eSMatthew G. Knepley . dm - the DM 6836858538eSMatthew G. Knepley 6846858538eSMatthew G. Knepley Level: advanced 6856858538eSMatthew G. Knepley 6866858538eSMatthew G. Knepley .seealso: `DMGetCellCoordinatesLocalNoncollective()` 6876858538eSMatthew G. Knepley @*/ 6889371c9d4SSatish Balay PetscErrorCode DMGetCellCoordinatesLocalSetUp(DM dm) { 6896858538eSMatthew G. Knepley PetscFunctionBegin; 6906858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 6916858538eSMatthew G. Knepley if (!dm->coordinates[1].xl && dm->coordinates[1].x) { 6926858538eSMatthew G. Knepley DM cdm = NULL; 6936858538eSMatthew G. Knepley 69411ea91a7SMatthew G. Knepley PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 6956858538eSMatthew G. Knepley PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[1].xl)); 6966858538eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].xl, "DG coordinates")); 6976858538eSMatthew G. Knepley PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl)); 6986858538eSMatthew G. Knepley PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl)); 6996858538eSMatthew G. Knepley } 7006858538eSMatthew G. Knepley PetscFunctionReturn(0); 7016858538eSMatthew G. Knepley } 7026858538eSMatthew G. Knepley 7036858538eSMatthew G. Knepley /*@ 7046858538eSMatthew G. Knepley DMGetCellCoordinatesLocal - Gets a local vector with the cellwise coordinates associated with the DM. 7056858538eSMatthew G. Knepley 7066858538eSMatthew G. Knepley Collective on dm 7076858538eSMatthew G. Knepley 7086858538eSMatthew G. Knepley Input Parameter: 7096858538eSMatthew G. Knepley . dm - the DM 7106858538eSMatthew G. Knepley 7116858538eSMatthew G. Knepley Output Parameter: 7126858538eSMatthew G. Knepley . c - coordinate vector 7136858538eSMatthew G. Knepley 7146858538eSMatthew G. Knepley Note: 7156858538eSMatthew G. Knepley This is a borrowed reference, so the user should NOT destroy this vector 7166858538eSMatthew G. Knepley 7176858538eSMatthew G. Knepley Each process has the local and ghost coordinates 7186858538eSMatthew G. Knepley 7196858538eSMatthew G. Knepley Level: intermediate 7206858538eSMatthew G. Knepley 7216858538eSMatthew G. Knepley .seealso: `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()`, `DMGetCellCoordinatesLocalNoncollective()` 7226858538eSMatthew G. Knepley @*/ 7239371c9d4SSatish Balay PetscErrorCode DMGetCellCoordinatesLocal(DM dm, Vec *c) { 7246858538eSMatthew G. Knepley PetscFunctionBegin; 7256858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 7266858538eSMatthew G. Knepley PetscValidPointer(c, 2); 7276858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinatesLocalSetUp(dm)); 7286858538eSMatthew G. Knepley *c = dm->coordinates[1].xl; 7296858538eSMatthew G. Knepley PetscFunctionReturn(0); 7306858538eSMatthew G. Knepley } 7316858538eSMatthew G. Knepley 7326858538eSMatthew G. Knepley /*@ 7336858538eSMatthew G. Knepley DMGetCellCoordinatesLocalNoncollective - Non-collective version of DMGetCellCoordinatesLocal(). Fails if global cellwise coordinates have been set and DMGetCellCoordinatesLocalSetUp() not called. 7346858538eSMatthew G. Knepley 7356858538eSMatthew G. Knepley Not collective 7366858538eSMatthew G. Knepley 7376858538eSMatthew G. Knepley Input Parameter: 7386858538eSMatthew G. Knepley . dm - the DM 7396858538eSMatthew G. Knepley 7406858538eSMatthew G. Knepley Output Parameter: 7416858538eSMatthew G. Knepley . c - cellwise coordinate vector 7426858538eSMatthew G. Knepley 7436858538eSMatthew G. Knepley Level: advanced 7446858538eSMatthew G. Knepley 7456858538eSMatthew G. Knepley .seealso: `DMGetCellCoordinatesLocalSetUp()`, `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()` 7466858538eSMatthew G. Knepley @*/ 7479371c9d4SSatish Balay PetscErrorCode DMGetCellCoordinatesLocalNoncollective(DM dm, Vec *c) { 7486858538eSMatthew G. Knepley PetscFunctionBegin; 7496858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 7506858538eSMatthew G. Knepley PetscValidPointer(c, 2); 7516858538eSMatthew G. Knepley PetscCheck(dm->coordinates[1].xl || !dm->coordinates[1].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCellCoordinatesLocalSetUp() has not been called"); 7526858538eSMatthew G. Knepley *c = dm->coordinates[1].xl; 7536858538eSMatthew G. Knepley PetscFunctionReturn(0); 7546858538eSMatthew G. Knepley } 7556858538eSMatthew G. Knepley 7566858538eSMatthew G. Knepley /*@ 7576858538eSMatthew G. Knepley DMSetCellCoordinatesLocal - Sets into the DM a local vector that holds the cellwise coordinates 7586858538eSMatthew G. Knepley 7596858538eSMatthew G. Knepley Not collective 7606858538eSMatthew G. Knepley 7616858538eSMatthew G. Knepley Input Parameters: 7626858538eSMatthew G. Knepley + dm - the DM 7636858538eSMatthew G. Knepley - c - cellwise coordinate vector 7646858538eSMatthew G. Knepley 7656858538eSMatthew G. Knepley Notes: 7666858538eSMatthew G. Knepley The coordinates of ghost points can be set using DMSetCoordinates() 7676858538eSMatthew G. Knepley followed by DMGetCoordinatesLocal(). This is intended to enable the 7686858538eSMatthew G. Knepley setting of ghost coordinates outside of the domain. 7696858538eSMatthew G. Knepley 7706858538eSMatthew G. Knepley The vector c should be destroyed by the caller. 7716858538eSMatthew G. Knepley 7726858538eSMatthew G. Knepley Level: intermediate 7736858538eSMatthew G. Knepley 7746858538eSMatthew G. Knepley .seealso: `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinates()`, `DMGetCellCoordinates()`, `DMGetCellCoordinateDM()` 7756858538eSMatthew G. Knepley @*/ 7769371c9d4SSatish Balay PetscErrorCode DMSetCellCoordinatesLocal(DM dm, Vec c) { 7776858538eSMatthew G. Knepley PetscFunctionBegin; 7786858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 7796858538eSMatthew G. Knepley if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2); 7806858538eSMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)c)); 7816858538eSMatthew G. Knepley PetscCall(VecDestroy(&dm->coordinates[1].xl)); 7826858538eSMatthew G. Knepley dm->coordinates[1].xl = c; 7836858538eSMatthew G. Knepley PetscCall(VecDestroy(&dm->coordinates[1].x)); 7846858538eSMatthew G. Knepley PetscFunctionReturn(0); 7856858538eSMatthew G. Knepley } 7866858538eSMatthew G. Knepley 7879371c9d4SSatish Balay PetscErrorCode DMGetCoordinateField(DM dm, DMField *field) { 7886858538eSMatthew G. Knepley PetscFunctionBegin; 7896858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 7906858538eSMatthew G. Knepley PetscValidPointer(field, 2); 7916858538eSMatthew G. Knepley if (!dm->coordinates[0].field) { 7926858538eSMatthew G. Knepley if (dm->ops->createcoordinatefield) PetscCall((*dm->ops->createcoordinatefield)(dm, &dm->coordinates[0].field)); 7936858538eSMatthew G. Knepley } 7946858538eSMatthew G. Knepley *field = dm->coordinates[0].field; 7956858538eSMatthew G. Knepley PetscFunctionReturn(0); 7966858538eSMatthew G. Knepley } 7976858538eSMatthew G. Knepley 7989371c9d4SSatish Balay PetscErrorCode DMSetCoordinateField(DM dm, DMField field) { 7996858538eSMatthew G. Knepley PetscFunctionBegin; 8006858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8016858538eSMatthew G. Knepley if (field) PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 2); 8026858538eSMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)field)); 8036858538eSMatthew G. Knepley PetscCall(DMFieldDestroy(&dm->coordinates[0].field)); 8046858538eSMatthew G. Knepley dm->coordinates[0].field = field; 8056858538eSMatthew G. Knepley PetscFunctionReturn(0); 8066858538eSMatthew G. Knepley } 8076858538eSMatthew G. Knepley 8086858538eSMatthew G. Knepley /*@ 8096858538eSMatthew G. Knepley DMGetLocalBoundingBox - Returns the bounding box for the piece of the DM on this process. 8106858538eSMatthew G. Knepley 8116858538eSMatthew G. Knepley Not collective 8126858538eSMatthew G. Knepley 8136858538eSMatthew G. Knepley Input Parameter: 8146858538eSMatthew G. Knepley . dm - the DM 8156858538eSMatthew G. Knepley 8166858538eSMatthew G. Knepley Output Parameters: 8176858538eSMatthew G. Knepley + lmin - local minimum coordinates (length coord dim, optional) 8186858538eSMatthew G. Knepley - lmax - local maximim coordinates (length coord dim, optional) 8196858538eSMatthew G. Knepley 8206858538eSMatthew G. Knepley Level: beginner 8216858538eSMatthew G. Knepley 8226858538eSMatthew G. Knepley Note: If the DM is a DMDA and has no coordinates, the index bounds are returned instead. 8236858538eSMatthew G. Knepley 8246858538eSMatthew G. Knepley .seealso: `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetBoundingBox()` 8256858538eSMatthew G. Knepley @*/ 8269371c9d4SSatish Balay PetscErrorCode DMGetLocalBoundingBox(DM dm, PetscReal lmin[], PetscReal lmax[]) { 8276858538eSMatthew G. Knepley Vec coords = NULL; 8286858538eSMatthew G. Knepley PetscReal min[3] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL}; 8296858538eSMatthew G. Knepley PetscReal max[3] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 8306858538eSMatthew G. Knepley const PetscScalar *local_coords; 8316858538eSMatthew G. Knepley PetscInt N, Ni; 8326858538eSMatthew G. Knepley PetscInt cdim, i, j; 8336858538eSMatthew G. Knepley 8346858538eSMatthew G. Knepley PetscFunctionBegin; 8356858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8366858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &cdim)); 8376858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(dm, &coords)); 8386858538eSMatthew G. Knepley if (coords) { 8396858538eSMatthew G. Knepley PetscCall(VecGetArrayRead(coords, &local_coords)); 8406858538eSMatthew G. Knepley PetscCall(VecGetLocalSize(coords, &N)); 8416858538eSMatthew G. Knepley Ni = N / cdim; 8426858538eSMatthew G. Knepley for (i = 0; i < Ni; ++i) { 8436858538eSMatthew G. Knepley for (j = 0; j < 3; ++j) { 8446858538eSMatthew G. Knepley min[j] = j < cdim ? PetscMin(min[j], PetscRealPart(local_coords[i * cdim + j])) : 0; 8456858538eSMatthew G. Knepley max[j] = j < cdim ? PetscMax(max[j], PetscRealPart(local_coords[i * cdim + j])) : 0; 8466858538eSMatthew G. Knepley } 8476858538eSMatthew G. Knepley } 8486858538eSMatthew G. Knepley PetscCall(VecRestoreArrayRead(coords, &local_coords)); 8496858538eSMatthew G. Knepley } else { 8506858538eSMatthew G. Knepley PetscBool isda; 8516858538eSMatthew G. Knepley 8526858538eSMatthew G. Knepley PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &isda)); 8536858538eSMatthew G. Knepley if (isda) PetscCall(DMGetLocalBoundingIndices_DMDA(dm, min, max)); 8546858538eSMatthew G. Knepley } 8556858538eSMatthew G. Knepley if (lmin) PetscCall(PetscArraycpy(lmin, min, cdim)); 8566858538eSMatthew G. Knepley if (lmax) PetscCall(PetscArraycpy(lmax, max, cdim)); 8576858538eSMatthew G. Knepley PetscFunctionReturn(0); 8586858538eSMatthew G. Knepley } 8596858538eSMatthew G. Knepley 8606858538eSMatthew G. Knepley /*@ 8616858538eSMatthew G. Knepley DMGetBoundingBox - Returns the global bounding box for the DM. 8626858538eSMatthew G. Knepley 8636858538eSMatthew G. Knepley Collective 8646858538eSMatthew G. Knepley 8656858538eSMatthew G. Knepley Input Parameter: 8666858538eSMatthew G. Knepley . dm - the DM 8676858538eSMatthew G. Knepley 8686858538eSMatthew G. Knepley Output Parameters: 8696858538eSMatthew G. Knepley + gmin - global minimum coordinates (length coord dim, optional) 8706858538eSMatthew G. Knepley - gmax - global maximim coordinates (length coord dim, optional) 8716858538eSMatthew G. Knepley 8726858538eSMatthew G. Knepley Level: beginner 8736858538eSMatthew G. Knepley 8746858538eSMatthew G. Knepley .seealso: `DMGetLocalBoundingBox()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()` 8756858538eSMatthew G. Knepley @*/ 8769371c9d4SSatish Balay PetscErrorCode DMGetBoundingBox(DM dm, PetscReal gmin[], PetscReal gmax[]) { 8776858538eSMatthew G. Knepley PetscReal lmin[3], lmax[3]; 8786858538eSMatthew G. Knepley PetscInt cdim; 8796858538eSMatthew G. Knepley PetscMPIInt count; 8806858538eSMatthew G. Knepley 8816858538eSMatthew G. Knepley PetscFunctionBegin; 8826858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8836858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &cdim)); 8846858538eSMatthew G. Knepley PetscCall(PetscMPIIntCast(cdim, &count)); 8856858538eSMatthew G. Knepley PetscCall(DMGetLocalBoundingBox(dm, lmin, lmax)); 8866858538eSMatthew G. Knepley if (gmin) PetscCall(MPIU_Allreduce(lmin, gmin, count, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm))); 8876858538eSMatthew G. Knepley if (gmax) PetscCall(MPIU_Allreduce(lmax, gmax, count, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)dm))); 8886858538eSMatthew G. Knepley PetscFunctionReturn(0); 8896858538eSMatthew G. Knepley } 8906858538eSMatthew G. Knepley 8916858538eSMatthew G. Knepley /*@ 8926858538eSMatthew G. Knepley DMProjectCoordinates - Project coordinates to a different space 8936858538eSMatthew G. Knepley 8946858538eSMatthew G. Knepley Input Parameters: 8956858538eSMatthew G. Knepley + dm - The DM object 8966858538eSMatthew G. Knepley - disc - The new coordinate discretization or NULL to ensure a coordinate discretization exists 8976858538eSMatthew G. Knepley 8986858538eSMatthew G. Knepley Level: intermediate 8996858538eSMatthew G. Knepley 9006858538eSMatthew G. Knepley .seealso: `DMGetCoordinateField()` 9016858538eSMatthew G. Knepley @*/ 9029371c9d4SSatish Balay PetscErrorCode DMProjectCoordinates(DM dm, PetscFE disc) { 9036858538eSMatthew G. Knepley PetscFE discOld; 9046858538eSMatthew G. Knepley PetscClassId classid; 9056858538eSMatthew G. Knepley DM cdmOld, cdmNew; 9066858538eSMatthew G. Knepley Vec coordsOld, coordsNew; 9076858538eSMatthew G. Knepley Mat matInterp; 9086858538eSMatthew G. Knepley PetscBool same_space = PETSC_TRUE; 9096858538eSMatthew G. Knepley 9106858538eSMatthew G. Knepley PetscFunctionBegin; 9116858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 9126858538eSMatthew G. Knepley if (disc) PetscValidHeaderSpecific(disc, PETSCFE_CLASSID, 2); 9136858538eSMatthew G. Knepley 9146858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &cdmOld)); 9156858538eSMatthew G. Knepley /* Check current discretization is compatible */ 9166858538eSMatthew G. Knepley PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld)); 9176858538eSMatthew G. Knepley PetscCall(PetscObjectGetClassId((PetscObject)discOld, &classid)); 9186858538eSMatthew G. Knepley if (classid != PETSCFE_CLASSID) { 9196858538eSMatthew G. Knepley if (classid == PETSC_CONTAINER_CLASSID) { 9206858538eSMatthew G. Knepley PetscFE feLinear; 9216858538eSMatthew G. Knepley DMPolytopeType ct; 9226858538eSMatthew G. Knepley PetscInt dim, dE, cStart, cEnd; 9236858538eSMatthew G. Knepley PetscBool simplex; 9246858538eSMatthew G. Knepley 9256858538eSMatthew G. Knepley /* Assume linear vertex coordinates */ 9266858538eSMatthew G. Knepley PetscCall(DMGetDimension(dm, &dim)); 9276858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &dE)); 9286858538eSMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 9296858538eSMatthew G. Knepley if (cEnd > cStart) { 9306858538eSMatthew G. Knepley PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 9316858538eSMatthew G. Knepley switch (ct) { 9326858538eSMatthew G. Knepley case DM_POLYTOPE_TRI_PRISM: 9339371c9d4SSatish Balay case DM_POLYTOPE_TRI_PRISM_TENSOR: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot autoamtically create coordinate space for prisms"); 9346858538eSMatthew G. Knepley default: break; 9356858538eSMatthew G. Knepley } 9366858538eSMatthew G. Knepley } 9376858538eSMatthew G. Knepley PetscCall(DMPlexIsSimplex(dm, &simplex)); 9386858538eSMatthew G. Knepley PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, dE, simplex, 1, -1, &feLinear)); 9396858538eSMatthew G. Knepley PetscCall(DMSetField(cdmOld, 0, NULL, (PetscObject)feLinear)); 9406858538eSMatthew G. Knepley PetscCall(PetscFEDestroy(&feLinear)); 9416858538eSMatthew G. Knepley PetscCall(DMCreateDS(cdmOld)); 9426858538eSMatthew G. Knepley PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld)); 9436858538eSMatthew G. Knepley } else { 9446858538eSMatthew G. Knepley const char *discname; 9456858538eSMatthew G. Knepley 9466858538eSMatthew G. Knepley PetscCall(PetscObjectGetType((PetscObject)discOld, &discname)); 9476858538eSMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject)discOld), PETSC_ERR_SUP, "Discretization type %s not supported", discname); 9486858538eSMatthew G. Knepley } 9496858538eSMatthew G. Knepley } 9506858538eSMatthew G. Knepley if (!disc) PetscFunctionReturn(0); 9516858538eSMatthew G. Knepley { // Check if the new space is the same as the old modulo quadrature 9526858538eSMatthew G. Knepley PetscDualSpace dsOld, ds; 9536858538eSMatthew G. Knepley PetscCall(PetscFEGetDualSpace(discOld, &dsOld)); 9546858538eSMatthew G. Knepley PetscCall(PetscFEGetDualSpace(disc, &ds)); 9556858538eSMatthew G. Knepley PetscCall(PetscDualSpaceEqual(dsOld, ds, &same_space)); 9566858538eSMatthew G. Knepley } 9576858538eSMatthew G. Knepley /* Make a fresh clone of the coordinate DM */ 9586858538eSMatthew G. Knepley PetscCall(DMClone(cdmOld, &cdmNew)); 9596858538eSMatthew G. Knepley PetscCall(DMSetField(cdmNew, 0, NULL, (PetscObject)disc)); 9606858538eSMatthew G. Knepley PetscCall(DMCreateDS(cdmNew)); 9616858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(dm, &coordsOld)); 9626858538eSMatthew G. Knepley if (same_space) { 9636858538eSMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)coordsOld)); 9646858538eSMatthew G. Knepley coordsNew = coordsOld; 9656858538eSMatthew G. Knepley } else { // Project the coordinate vector from old to new space 9666858538eSMatthew G. Knepley PetscCall(DMCreateGlobalVector(cdmNew, &coordsNew)); 9676858538eSMatthew G. Knepley PetscCall(DMCreateInterpolation(cdmOld, cdmNew, &matInterp, NULL)); 9686858538eSMatthew G. Knepley PetscCall(MatInterpolate(matInterp, coordsOld, coordsNew)); 9696858538eSMatthew G. Knepley PetscCall(MatDestroy(&matInterp)); 9706858538eSMatthew G. Knepley } 9716858538eSMatthew G. Knepley /* Set new coordinate structures */ 9726858538eSMatthew G. Knepley PetscCall(DMSetCoordinateField(dm, NULL)); 9736858538eSMatthew G. Knepley PetscCall(DMSetCoordinateDM(dm, cdmNew)); 9746858538eSMatthew G. Knepley PetscCall(DMSetCoordinates(dm, coordsNew)); 9756858538eSMatthew G. Knepley PetscCall(VecDestroy(&coordsNew)); 9766858538eSMatthew G. Knepley PetscCall(DMDestroy(&cdmNew)); 9776858538eSMatthew G. Knepley PetscFunctionReturn(0); 9786858538eSMatthew G. Knepley } 9796858538eSMatthew G. Knepley 9806858538eSMatthew G. Knepley /*@ 9816858538eSMatthew G. Knepley DMLocatePoints - Locate the points in v in the mesh and return a PetscSF of the containing cells 9826858538eSMatthew G. Knepley 9836858538eSMatthew G. Knepley Collective on v (see explanation below) 9846858538eSMatthew G. Knepley 9856858538eSMatthew G. Knepley Input Parameters: 9866858538eSMatthew G. Knepley + dm - The DM 9876858538eSMatthew G. Knepley - ltype - The type of point location, e.g. DM_POINTLOCATION_NONE or DM_POINTLOCATION_NEAREST 9886858538eSMatthew G. Knepley 9896858538eSMatthew G. Knepley Input/Output Parameters: 9906858538eSMatthew G. Knepley + v - The Vec of points, on output contains the nearest mesh points to the given points if DM_POINTLOCATION_NEAREST is used 9916858538eSMatthew G. Knepley - cellSF - Points to either NULL, or a PetscSF with guesses for which cells contain each point; 9926858538eSMatthew G. Knepley on output, the PetscSF containing the ranks and local indices of the containing points 9936858538eSMatthew G. Knepley 9946858538eSMatthew G. Knepley Level: developer 9956858538eSMatthew G. Knepley 9966858538eSMatthew G. Knepley Notes: 9976858538eSMatthew G. Knepley To do a search of the local cells of the mesh, v should have PETSC_COMM_SELF as its communicator. 9986858538eSMatthew G. Knepley To do a search of all the cells in the distributed mesh, v should have the same communicator as dm. 9996858538eSMatthew G. Knepley 1000*d8206211SMatthew G. Knepley Points will only be located in owned cells, not overlap cells arising from `DMPlexDistribute()` or other overlapping distributions. 1001*d8206211SMatthew G. Knepley 10026858538eSMatthew G. Knepley If *cellSF is NULL on input, a PetscSF will be created. 10036858538eSMatthew G. Knepley If *cellSF is not NULL on input, it should point to an existing PetscSF, whose graph will be used as initial guesses. 10046858538eSMatthew G. Knepley 10056858538eSMatthew G. Knepley An array that maps each point to its containing cell can be obtained with 10066858538eSMatthew G. Knepley 10076858538eSMatthew G. Knepley $ const PetscSFNode *cells; 10086858538eSMatthew G. Knepley $ PetscInt nFound; 10096858538eSMatthew G. Knepley $ const PetscInt *found; 10106858538eSMatthew G. Knepley $ 10116858538eSMatthew G. Knepley $ PetscSFGetGraph(cellSF,NULL,&nFound,&found,&cells); 10126858538eSMatthew G. Knepley 10136858538eSMatthew G. Knepley Where cells[i].rank is the rank of the cell containing point found[i] (or i if found == NULL), and cells[i].index is 10146858538eSMatthew G. Knepley the index of the cell in its rank's local numbering. 10156858538eSMatthew G. Knepley 10166858538eSMatthew G. Knepley .seealso: `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMPointLocationType` 10176858538eSMatthew G. Knepley @*/ 10189371c9d4SSatish Balay PetscErrorCode DMLocatePoints(DM dm, Vec v, DMPointLocationType ltype, PetscSF *cellSF) { 10196858538eSMatthew G. Knepley PetscFunctionBegin; 10206858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 10216858538eSMatthew G. Knepley PetscValidHeaderSpecific(v, VEC_CLASSID, 2); 10226858538eSMatthew G. Knepley PetscValidPointer(cellSF, 4); 10236858538eSMatthew G. Knepley if (*cellSF) { 10246858538eSMatthew G. Knepley PetscMPIInt result; 10256858538eSMatthew G. Knepley 10266858538eSMatthew G. Knepley PetscValidHeaderSpecific(*cellSF, PETSCSF_CLASSID, 4); 10276858538eSMatthew G. Knepley PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)v), PetscObjectComm((PetscObject)*cellSF), &result)); 10286858538eSMatthew G. Knepley PetscCheck(result == MPI_IDENT || result == MPI_CONGRUENT, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "cellSF must have a communicator congruent to v's"); 10296858538eSMatthew G. Knepley } else { 10306858538eSMatthew G. Knepley PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)v), cellSF)); 10316858538eSMatthew G. Knepley } 10326858538eSMatthew G. Knepley PetscCall(PetscLogEventBegin(DM_LocatePoints, dm, 0, 0, 0)); 1033dbbe0bcdSBarry Smith PetscUseTypeMethod(dm, locatepoints, v, ltype, *cellSF); 10346858538eSMatthew G. Knepley PetscCall(PetscLogEventEnd(DM_LocatePoints, dm, 0, 0, 0)); 10356858538eSMatthew G. Knepley PetscFunctionReturn(0); 10366858538eSMatthew G. Knepley } 1037