16858538eSMatthew G. Knepley #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/ 26858538eSMatthew G. Knepley 3e44f6aebSMatthew G. Knepley #include <petscdmplex.h> /* For DMCreateAffineCoordinates_Internal() */ 46858538eSMatthew G. Knepley #include <petscsf.h> /* For DMLocatePoints() */ 56858538eSMatthew G. Knepley 6d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRestrictHook_Coordinates(DM dm, DM dmc, void *ctx) 7d71ae5a4SJacob Faibussowitsch { 86858538eSMatthew G. Knepley DM dm_coord, dmc_coord; 96858538eSMatthew G. Knepley Vec coords, ccoords; 106858538eSMatthew G. Knepley Mat inject; 11*4d86920dSPierre Jolivet 126858538eSMatthew G. Knepley PetscFunctionBegin; 136858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 146858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dmc, &dmc_coord)); 156858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(dm, &coords)); 166858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(dmc, &ccoords)); 176858538eSMatthew G. Knepley if (coords && !ccoords) { 186858538eSMatthew G. Knepley PetscCall(DMCreateGlobalVector(dmc_coord, &ccoords)); 196858538eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)ccoords, "coordinates")); 206858538eSMatthew G. Knepley PetscCall(DMCreateInjection(dmc_coord, dm_coord, &inject)); 216858538eSMatthew G. Knepley PetscCall(MatRestrict(inject, coords, ccoords)); 226858538eSMatthew G. Knepley PetscCall(MatDestroy(&inject)); 236858538eSMatthew G. Knepley PetscCall(DMSetCoordinates(dmc, ccoords)); 246858538eSMatthew G. Knepley PetscCall(VecDestroy(&ccoords)); 256858538eSMatthew G. Knepley } 263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 276858538eSMatthew G. Knepley } 286858538eSMatthew G. Knepley 29d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_Coordinates(DM dm, DM subdm, void *ctx) 30d71ae5a4SJacob Faibussowitsch { 316858538eSMatthew G. Knepley DM dm_coord, subdm_coord; 326858538eSMatthew G. Knepley Vec coords, ccoords, clcoords; 336858538eSMatthew G. Knepley VecScatter *scat_i, *scat_g; 34*4d86920dSPierre Jolivet 356858538eSMatthew G. Knepley PetscFunctionBegin; 366858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 376858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(subdm, &subdm_coord)); 386858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(dm, &coords)); 396858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(subdm, &ccoords)); 406858538eSMatthew G. Knepley if (coords && !ccoords) { 416858538eSMatthew G. Knepley PetscCall(DMCreateGlobalVector(subdm_coord, &ccoords)); 426858538eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)ccoords, "coordinates")); 436858538eSMatthew G. Knepley PetscCall(DMCreateLocalVector(subdm_coord, &clcoords)); 446858538eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)clcoords, "coordinates")); 456858538eSMatthew G. Knepley PetscCall(DMCreateDomainDecompositionScatters(dm_coord, 1, &subdm_coord, NULL, &scat_i, &scat_g)); 466858538eSMatthew G. Knepley PetscCall(VecScatterBegin(scat_i[0], coords, ccoords, INSERT_VALUES, SCATTER_FORWARD)); 476858538eSMatthew G. Knepley PetscCall(VecScatterEnd(scat_i[0], coords, ccoords, INSERT_VALUES, SCATTER_FORWARD)); 486858538eSMatthew G. Knepley PetscCall(VecScatterBegin(scat_g[0], coords, clcoords, INSERT_VALUES, SCATTER_FORWARD)); 496858538eSMatthew G. Knepley PetscCall(VecScatterEnd(scat_g[0], coords, clcoords, INSERT_VALUES, SCATTER_FORWARD)); 506858538eSMatthew G. Knepley PetscCall(DMSetCoordinates(subdm, ccoords)); 516858538eSMatthew G. Knepley PetscCall(DMSetCoordinatesLocal(subdm, clcoords)); 526858538eSMatthew G. Knepley PetscCall(VecScatterDestroy(&scat_i[0])); 536858538eSMatthew G. Knepley PetscCall(VecScatterDestroy(&scat_g[0])); 546858538eSMatthew G. Knepley PetscCall(VecDestroy(&ccoords)); 556858538eSMatthew G. Knepley PetscCall(VecDestroy(&clcoords)); 566858538eSMatthew G. Knepley PetscCall(PetscFree(scat_i)); 576858538eSMatthew G. Knepley PetscCall(PetscFree(scat_g)); 586858538eSMatthew G. Knepley } 593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 606858538eSMatthew G. Knepley } 616858538eSMatthew G. Knepley 626858538eSMatthew G. Knepley /*@ 63dce8aebaSBarry Smith DMGetCoordinateDM - Gets the `DM` that prescribes coordinate layout and scatters between global and local coordinates 646858538eSMatthew G. Knepley 6520f4b53cSBarry Smith Collective 666858538eSMatthew G. Knepley 676858538eSMatthew G. Knepley Input Parameter: 68dce8aebaSBarry Smith . dm - the `DM` 696858538eSMatthew G. Knepley 706858538eSMatthew G. Knepley Output Parameter: 71dce8aebaSBarry Smith . cdm - coordinate `DM` 726858538eSMatthew G. Knepley 736858538eSMatthew G. Knepley Level: intermediate 746858538eSMatthew G. Knepley 75dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinateDM()`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGSetCellCoordinateDM()`, 7660225df5SJacob Faibussowitsch 776858538eSMatthew G. Knepley @*/ 78d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinateDM(DM dm, DM *cdm) 79d71ae5a4SJacob Faibussowitsch { 806858538eSMatthew G. Knepley PetscFunctionBegin; 816858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 824f572ea9SToby Isaac PetscAssertPointer(cdm, 2); 836858538eSMatthew G. Knepley if (!dm->coordinates[0].dm) { 846858538eSMatthew G. Knepley DM cdm; 856858538eSMatthew G. Knepley 86dbbe0bcdSBarry Smith PetscUseTypeMethod(dm, createcoordinatedm, &cdm); 876858538eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)cdm, "coordinateDM")); 886858538eSMatthew G. Knepley /* Just in case the DM sets the coordinate DM when creating it (DMP4est can do this, because it may not setup 896858538eSMatthew G. Knepley * until the call to CreateCoordinateDM) */ 906858538eSMatthew G. Knepley PetscCall(DMDestroy(&dm->coordinates[0].dm)); 916858538eSMatthew G. Knepley dm->coordinates[0].dm = cdm; 926858538eSMatthew G. Knepley } 936858538eSMatthew G. Knepley *cdm = dm->coordinates[0].dm; 943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 956858538eSMatthew G. Knepley } 966858538eSMatthew G. Knepley 976858538eSMatthew G. Knepley /*@ 98dce8aebaSBarry Smith DMSetCoordinateDM - Sets the `DM` that prescribes coordinate layout and scatters between global and local coordinates 996858538eSMatthew G. Knepley 10020f4b53cSBarry Smith Logically Collective 1016858538eSMatthew G. Knepley 1026858538eSMatthew G. Knepley Input Parameters: 103dce8aebaSBarry Smith + dm - the `DM` 104dce8aebaSBarry Smith - cdm - coordinate `DM` 1056858538eSMatthew G. Knepley 1066858538eSMatthew G. Knepley Level: intermediate 1076858538eSMatthew G. Knepley 108dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateDM()`, `DMSetCoordinates()`, `DMGetCellCoordinateDM()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, 109dce8aebaSBarry Smith `DMGSetCellCoordinateDM()` 1106858538eSMatthew G. Knepley @*/ 111d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinateDM(DM dm, DM cdm) 112d71ae5a4SJacob Faibussowitsch { 1136858538eSMatthew G. Knepley PetscFunctionBegin; 1146858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 11590612307SJed Brown if (cdm) PetscValidHeaderSpecific(cdm, DM_CLASSID, 2); 1166858538eSMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)cdm)); 1176858538eSMatthew G. Knepley PetscCall(DMDestroy(&dm->coordinates[0].dm)); 1186858538eSMatthew G. Knepley dm->coordinates[0].dm = cdm; 1193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1206858538eSMatthew G. Knepley } 1216858538eSMatthew G. Knepley 1226858538eSMatthew G. Knepley /*@ 123dce8aebaSBarry Smith DMGetCellCoordinateDM - Gets the `DM` that prescribes cellwise coordinate layout and scatters between global and local cellwise coordinates 1246858538eSMatthew G. Knepley 12520f4b53cSBarry Smith Collective 1266858538eSMatthew G. Knepley 1276858538eSMatthew G. Knepley Input Parameter: 128dce8aebaSBarry Smith . dm - the `DM` 1296858538eSMatthew G. Knepley 1306858538eSMatthew G. Knepley Output Parameter: 1310b3275a6SBarry Smith . cdm - cellwise coordinate `DM`, or `NULL` if they are not defined 1326858538eSMatthew G. Knepley 1336858538eSMatthew G. Knepley Level: intermediate 1346858538eSMatthew G. Knepley 135dce8aebaSBarry Smith Note: 136dce8aebaSBarry Smith Call `DMLocalizeCoordinates()` to automatically create cellwise coordinates for periodic geometries. 137dce8aebaSBarry Smith 138dce8aebaSBarry Smith .seealso: `DM`, `DMSetCellCoordinateDM()`, `DMSetCellCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, 139dce8aebaSBarry Smith `DMLocalizeCoordinates()`, `DMSetCoordinateDM()`, `DMGetCoordinateDM()` 1406858538eSMatthew G. Knepley @*/ 141d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinateDM(DM dm, DM *cdm) 142d71ae5a4SJacob Faibussowitsch { 1436858538eSMatthew G. Knepley PetscFunctionBegin; 1446858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1454f572ea9SToby Isaac PetscAssertPointer(cdm, 2); 1466858538eSMatthew G. Knepley *cdm = dm->coordinates[1].dm; 1473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1486858538eSMatthew G. Knepley } 1496858538eSMatthew G. Knepley 1506858538eSMatthew G. Knepley /*@ 151dce8aebaSBarry Smith DMSetCellCoordinateDM - Sets the `DM` that prescribes cellwise coordinate layout and scatters between global and local cellwise coordinates 1526858538eSMatthew G. Knepley 15320f4b53cSBarry Smith Logically Collective 1546858538eSMatthew G. Knepley 1556858538eSMatthew G. Knepley Input Parameters: 156dce8aebaSBarry Smith + dm - the `DM` 157dce8aebaSBarry Smith - cdm - cellwise coordinate `DM` 1586858538eSMatthew G. Knepley 1596858538eSMatthew G. Knepley Level: intermediate 1606858538eSMatthew G. Knepley 161dce8aebaSBarry Smith Note: 16235cb6cd3SPierre Jolivet As opposed to `DMSetCoordinateDM()` these coordinates are useful for discontinuous Galerkin methods since they support coordinate fields that are discontinuous at cell boundaries. 163dce8aebaSBarry Smith 164dce8aebaSBarry Smith .seealso: `DMGetCellCoordinateDM()`, `DMSetCellCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, 165dce8aebaSBarry Smith `DMSetCoordinateDM()`, `DMGetCoordinateDM()` 1666858538eSMatthew G. Knepley @*/ 167d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCellCoordinateDM(DM dm, DM cdm) 168d71ae5a4SJacob Faibussowitsch { 1696858538eSMatthew G. Knepley PetscInt dim; 1706858538eSMatthew G. Knepley 1716858538eSMatthew G. Knepley PetscFunctionBegin; 1726858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1736858538eSMatthew G. Knepley if (cdm) { 1746858538eSMatthew G. Knepley PetscValidHeaderSpecific(cdm, DM_CLASSID, 2); 1756858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &dim)); 1766858538eSMatthew G. Knepley dm->coordinates[1].dim = dim; 1776858538eSMatthew G. Knepley } 1786858538eSMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)cdm)); 1796858538eSMatthew G. Knepley PetscCall(DMDestroy(&dm->coordinates[1].dm)); 1806858538eSMatthew G. Knepley dm->coordinates[1].dm = cdm; 1813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1826858538eSMatthew G. Knepley } 1836858538eSMatthew G. Knepley 1846858538eSMatthew G. Knepley /*@ 1850b3275a6SBarry Smith DMGetCoordinateDim - Retrieve the dimension of the embedding space for coordinate values. For example a mesh on the surface of a sphere would have a 3 dimensional embedding space 1866858538eSMatthew G. Knepley 1876858538eSMatthew G. Knepley Not Collective 1886858538eSMatthew G. Knepley 1896858538eSMatthew G. Knepley Input Parameter: 190dce8aebaSBarry Smith . dm - The `DM` object 1916858538eSMatthew G. Knepley 1926858538eSMatthew G. Knepley Output Parameter: 1936858538eSMatthew G. Knepley . dim - The embedding dimension 1946858538eSMatthew G. Knepley 1956858538eSMatthew G. Knepley Level: intermediate 1966858538eSMatthew G. Knepley 197dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinateDim()`, `DMGetCoordinateSection()`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()` 1986858538eSMatthew G. Knepley @*/ 199d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinateDim(DM dm, PetscInt *dim) 200d71ae5a4SJacob Faibussowitsch { 2016858538eSMatthew G. Knepley PetscFunctionBegin; 2026858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2034f572ea9SToby Isaac PetscAssertPointer(dim, 2); 2046858538eSMatthew G. Knepley if (dm->coordinates[0].dim == PETSC_DEFAULT) dm->coordinates[0].dim = dm->dim; 2056858538eSMatthew G. Knepley *dim = dm->coordinates[0].dim; 2063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2076858538eSMatthew G. Knepley } 2086858538eSMatthew G. Knepley 2096858538eSMatthew G. Knepley /*@ 2106858538eSMatthew G. Knepley DMSetCoordinateDim - Set the dimension of the embedding space for coordinate values. 2116858538eSMatthew G. Knepley 2126858538eSMatthew G. Knepley Not Collective 2136858538eSMatthew G. Knepley 2146858538eSMatthew G. Knepley Input Parameters: 215dce8aebaSBarry Smith + dm - The `DM` object 2166858538eSMatthew G. Knepley - dim - The embedding dimension 2176858538eSMatthew G. Knepley 2186858538eSMatthew G. Knepley Level: intermediate 2196858538eSMatthew G. Knepley 220dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateDim()`, `DMSetCoordinateSection()`, `DMGetCoordinateSection()`, `DMGetLocalSection()`, `DMSetLocalSection()` 2216858538eSMatthew G. Knepley @*/ 222d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinateDim(DM dm, PetscInt dim) 223d71ae5a4SJacob Faibussowitsch { 2246858538eSMatthew G. Knepley PetscDS ds; 2256858538eSMatthew G. Knepley PetscInt Nds, n; 2266858538eSMatthew G. Knepley 2276858538eSMatthew G. Knepley PetscFunctionBegin; 2286858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2296858538eSMatthew G. Knepley dm->coordinates[0].dim = dim; 2306858538eSMatthew G. Knepley if (dm->dim >= 0) { 2316858538eSMatthew G. Knepley PetscCall(DMGetNumDS(dm, &Nds)); 2326858538eSMatthew G. Knepley for (n = 0; n < Nds; ++n) { 23307218a29SMatthew G. Knepley PetscCall(DMGetRegionNumDS(dm, n, NULL, NULL, &ds, NULL)); 2346858538eSMatthew G. Knepley PetscCall(PetscDSSetCoordinateDimension(ds, dim)); 2356858538eSMatthew G. Knepley } 2366858538eSMatthew G. Knepley } 2373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2386858538eSMatthew G. Knepley } 2396858538eSMatthew G. Knepley 2406858538eSMatthew G. Knepley /*@ 2410b3275a6SBarry Smith DMGetCoordinateSection - Retrieve the `PetscSection` of coordinate values over the mesh. 2426858538eSMatthew G. Knepley 24320f4b53cSBarry Smith Collective 2446858538eSMatthew G. Knepley 2456858538eSMatthew G. Knepley Input Parameter: 246dce8aebaSBarry Smith . dm - The `DM` object 2476858538eSMatthew G. Knepley 2486858538eSMatthew G. Knepley Output Parameter: 249dce8aebaSBarry Smith . section - The `PetscSection` object 2506858538eSMatthew G. Knepley 2516858538eSMatthew G. Knepley Level: intermediate 2526858538eSMatthew G. Knepley 253dce8aebaSBarry Smith Note: 254dce8aebaSBarry Smith This just retrieves the local section from the coordinate `DM`. In other words, 255dce8aebaSBarry Smith .vb 256dce8aebaSBarry Smith DMGetCoordinateDM(dm, &cdm); 257dce8aebaSBarry Smith DMGetLocalSection(cdm, §ion); 258dce8aebaSBarry Smith .ve 259dce8aebaSBarry Smith 2600b3275a6SBarry Smith .seealso: `DM`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()` 2616858538eSMatthew G. Knepley @*/ 262d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinateSection(DM dm, PetscSection *section) 263d71ae5a4SJacob Faibussowitsch { 2646858538eSMatthew G. Knepley DM cdm; 2656858538eSMatthew G. Knepley 2666858538eSMatthew G. Knepley PetscFunctionBegin; 2676858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2684f572ea9SToby Isaac PetscAssertPointer(section, 2); 2696858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &cdm)); 2706858538eSMatthew G. Knepley PetscCall(DMGetLocalSection(cdm, section)); 2713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2726858538eSMatthew G. Knepley } 2736858538eSMatthew G. Knepley 2746858538eSMatthew G. Knepley /*@ 2750b3275a6SBarry Smith DMSetCoordinateSection - Set the `PetscSection` of coordinate values over the mesh. 2766858538eSMatthew G. Knepley 2776858538eSMatthew G. Knepley Not Collective 2786858538eSMatthew G. Knepley 2796858538eSMatthew G. Knepley Input Parameters: 280dce8aebaSBarry Smith + dm - The `DM` object 281dce8aebaSBarry Smith . dim - The embedding dimension, or `PETSC_DETERMINE` 282dce8aebaSBarry Smith - section - The `PetscSection` object 2836858538eSMatthew G. Knepley 2846858538eSMatthew G. Knepley Level: intermediate 2856858538eSMatthew G. Knepley 286dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateDim()`, `DMGetCoordinateSection()`, `DMGetLocalSection()`, `DMSetLocalSection()` 2876858538eSMatthew G. Knepley @*/ 288d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinateSection(DM dm, PetscInt dim, PetscSection section) 289d71ae5a4SJacob Faibussowitsch { 2906858538eSMatthew G. Knepley DM cdm; 2916858538eSMatthew G. Knepley 2926858538eSMatthew G. Knepley PetscFunctionBegin; 2936858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2946858538eSMatthew G. Knepley PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 3); 2956858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &cdm)); 2966858538eSMatthew G. Knepley PetscCall(DMSetLocalSection(cdm, section)); 2976858538eSMatthew G. Knepley if (dim == PETSC_DETERMINE) { 2986858538eSMatthew G. Knepley PetscInt d = PETSC_DEFAULT; 2996858538eSMatthew G. Knepley PetscInt pStart, pEnd, vStart, vEnd, v, dd; 3006858538eSMatthew G. Knepley 3016858538eSMatthew G. Knepley PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 3026858538eSMatthew G. Knepley PetscCall(DMGetDimPoints(dm, 0, &vStart, &vEnd)); 3036858538eSMatthew G. Knepley pStart = PetscMax(vStart, pStart); 3046858538eSMatthew G. Knepley pEnd = PetscMin(vEnd, pEnd); 3056858538eSMatthew G. Knepley for (v = pStart; v < pEnd; ++v) { 3066858538eSMatthew G. Knepley PetscCall(PetscSectionGetDof(section, v, &dd)); 3079371c9d4SSatish Balay if (dd) { 3089371c9d4SSatish Balay d = dd; 3099371c9d4SSatish Balay break; 3109371c9d4SSatish Balay } 3116858538eSMatthew G. Knepley } 3126858538eSMatthew G. Knepley if (d >= 0) PetscCall(DMSetCoordinateDim(dm, d)); 3136858538eSMatthew G. Knepley } 3143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3156858538eSMatthew G. Knepley } 3166858538eSMatthew G. Knepley 3176858538eSMatthew G. Knepley /*@ 3180b3275a6SBarry Smith DMGetCellCoordinateSection - Retrieve the `PetscSection` of cellwise coordinate values over the mesh. 3196858538eSMatthew G. Knepley 32020f4b53cSBarry Smith Collective 3216858538eSMatthew G. Knepley 3226858538eSMatthew G. Knepley Input Parameter: 323dce8aebaSBarry Smith . dm - The `DM` object 3246858538eSMatthew G. Knepley 3256858538eSMatthew G. Knepley Output Parameter: 32620f4b53cSBarry Smith . section - The `PetscSection` object, or `NULL` if no cellwise coordinates are defined 3276858538eSMatthew G. Knepley 3286858538eSMatthew G. Knepley Level: intermediate 3296858538eSMatthew G. Knepley 330dce8aebaSBarry Smith Note: 331dce8aebaSBarry Smith This just retrieves the local section from the cell coordinate `DM`. In other words, 332dce8aebaSBarry Smith .vb 333dce8aebaSBarry Smith DMGetCellCoordinateDM(dm, &cdm); 334dce8aebaSBarry Smith DMGetLocalSection(cdm, §ion); 335dce8aebaSBarry Smith .ve 336dce8aebaSBarry Smith 337dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateSection()`, `DMSetCellCoordinateSection()`, `DMGetCellCoordinateDM()`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()` 3386858538eSMatthew G. Knepley @*/ 339d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinateSection(DM dm, PetscSection *section) 340d71ae5a4SJacob Faibussowitsch { 3416858538eSMatthew G. Knepley DM cdm; 3426858538eSMatthew G. Knepley 3436858538eSMatthew G. Knepley PetscFunctionBegin; 3446858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3454f572ea9SToby Isaac PetscAssertPointer(section, 2); 3466858538eSMatthew G. Knepley *section = NULL; 3476858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 3486858538eSMatthew G. Knepley if (cdm) PetscCall(DMGetLocalSection(cdm, section)); 3493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3506858538eSMatthew G. Knepley } 3516858538eSMatthew G. Knepley 3526858538eSMatthew G. Knepley /*@ 3530b3275a6SBarry Smith DMSetCellCoordinateSection - Set the `PetscSection` of cellwise coordinate values over the mesh. 3546858538eSMatthew G. Knepley 3556858538eSMatthew G. Knepley Not Collective 3566858538eSMatthew G. Knepley 3576858538eSMatthew G. Knepley Input Parameters: 358dce8aebaSBarry Smith + dm - The `DM` object 359dce8aebaSBarry Smith . dim - The embedding dimension, or `PETSC_DETERMINE` 360dce8aebaSBarry Smith - section - The `PetscSection` object for a cellwise layout 3616858538eSMatthew G. Knepley 3626858538eSMatthew G. Knepley Level: intermediate 3636858538eSMatthew G. Knepley 364dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateDim()`, `DMSetCoordinateSection()`, `DMGetCellCoordinateSection()`, `DMGetCoordinateSection()`, `DMGetCellCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()` 3656858538eSMatthew G. Knepley @*/ 366d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCellCoordinateSection(DM dm, PetscInt dim, PetscSection section) 367d71ae5a4SJacob Faibussowitsch { 3686858538eSMatthew G. Knepley DM cdm; 3696858538eSMatthew G. Knepley 3706858538eSMatthew G. Knepley PetscFunctionBegin; 3716858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3726858538eSMatthew G. Knepley PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 3); 3736858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 3746858538eSMatthew G. Knepley PetscCheck(cdm, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "No DM defined for cellwise coordinates"); 3756858538eSMatthew G. Knepley PetscCall(DMSetLocalSection(cdm, section)); 3766858538eSMatthew G. Knepley if (dim == PETSC_DETERMINE) { 3776858538eSMatthew G. Knepley PetscInt d = PETSC_DEFAULT; 3786858538eSMatthew G. Knepley PetscInt pStart, pEnd, vStart, vEnd, v, dd; 3796858538eSMatthew G. Knepley 3806858538eSMatthew G. Knepley PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 3816858538eSMatthew G. Knepley PetscCall(DMGetDimPoints(dm, 0, &vStart, &vEnd)); 3826858538eSMatthew G. Knepley pStart = PetscMax(vStart, pStart); 3836858538eSMatthew G. Knepley pEnd = PetscMin(vEnd, pEnd); 3846858538eSMatthew G. Knepley for (v = pStart; v < pEnd; ++v) { 3856858538eSMatthew G. Knepley PetscCall(PetscSectionGetDof(section, v, &dd)); 3869371c9d4SSatish Balay if (dd) { 3879371c9d4SSatish Balay d = dd; 3889371c9d4SSatish Balay break; 3899371c9d4SSatish Balay } 3906858538eSMatthew G. Knepley } 3916858538eSMatthew G. Knepley if (d >= 0) PetscCall(DMSetCoordinateDim(dm, d)); 3926858538eSMatthew G. Knepley } 3933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3946858538eSMatthew G. Knepley } 3956858538eSMatthew G. Knepley 3966858538eSMatthew G. Knepley /*@ 397dce8aebaSBarry Smith DMGetCoordinates - Gets a global vector with the coordinates associated with the `DM`. 3986858538eSMatthew G. Knepley 39920f4b53cSBarry Smith Collective 4006858538eSMatthew G. Knepley 4016858538eSMatthew G. Knepley Input Parameter: 402dce8aebaSBarry Smith . dm - the `DM` 4036858538eSMatthew G. Knepley 4046858538eSMatthew G. Knepley Output Parameter: 4056858538eSMatthew G. Knepley . c - global coordinate vector 4066858538eSMatthew G. Knepley 407dce8aebaSBarry Smith Level: intermediate 408dce8aebaSBarry Smith 409dce8aebaSBarry Smith Notes: 410dce8aebaSBarry Smith This is a borrowed reference, so the user should NOT destroy this vector. When the `DM` is 4110b3275a6SBarry Smith destroyed `c` will no longer be valid. 4126858538eSMatthew G. Knepley 4136858538eSMatthew G. Knepley Each process has only the locally-owned portion of the global coordinates (does NOT have the ghost coordinates). 4146858538eSMatthew G. Knepley 415dce8aebaSBarry Smith For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...) 4166858538eSMatthew G. Knepley and (x_0,y_0,z_0,x_1,y_1,z_1...) 4176858538eSMatthew G. Knepley 418dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMSetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCoordinateDM()`, `DMDASetUniformCoordinates()` 4196858538eSMatthew G. Knepley @*/ 420d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinates(DM dm, Vec *c) 421d71ae5a4SJacob Faibussowitsch { 4226858538eSMatthew G. Knepley PetscFunctionBegin; 4236858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4244f572ea9SToby Isaac PetscAssertPointer(c, 2); 4256858538eSMatthew G. Knepley if (!dm->coordinates[0].x && dm->coordinates[0].xl) { 4266858538eSMatthew G. Knepley DM cdm = NULL; 4276858538eSMatthew G. Knepley 4286858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &cdm)); 4296858538eSMatthew G. Knepley PetscCall(DMCreateGlobalVector(cdm, &dm->coordinates[0].x)); 4306858538eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].x, "coordinates")); 4316858538eSMatthew G. Knepley PetscCall(DMLocalToGlobalBegin(cdm, dm->coordinates[0].xl, INSERT_VALUES, dm->coordinates[0].x)); 4326858538eSMatthew G. Knepley PetscCall(DMLocalToGlobalEnd(cdm, dm->coordinates[0].xl, INSERT_VALUES, dm->coordinates[0].x)); 4336858538eSMatthew G. Knepley } 4346858538eSMatthew G. Knepley *c = dm->coordinates[0].x; 4353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4366858538eSMatthew G. Knepley } 4376858538eSMatthew G. Knepley 4386858538eSMatthew G. Knepley /*@ 439dce8aebaSBarry Smith DMSetCoordinates - Sets into the `DM` a global vector that holds the coordinates 4406858538eSMatthew G. Knepley 44120f4b53cSBarry Smith Collective 4426858538eSMatthew G. Knepley 4436858538eSMatthew G. Knepley Input Parameters: 444dce8aebaSBarry Smith + dm - the `DM` 4456858538eSMatthew G. Knepley - c - coordinate vector 4466858538eSMatthew G. Knepley 4476858538eSMatthew G. Knepley Level: intermediate 4486858538eSMatthew G. Knepley 449dce8aebaSBarry Smith Notes: 450dce8aebaSBarry Smith The coordinates do not include those for ghost points, which are in the local vector. 451dce8aebaSBarry Smith 45220f4b53cSBarry Smith The vector `c` can be destroyed after the call 453dce8aebaSBarry Smith 454dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCoordinateDM()`, `DMDASetUniformCoordinates()` 4556858538eSMatthew G. Knepley @*/ 456d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinates(DM dm, Vec c) 457d71ae5a4SJacob Faibussowitsch { 4586858538eSMatthew G. Knepley PetscFunctionBegin; 4596858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4606858538eSMatthew G. Knepley if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2); 4616858538eSMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)c)); 4626858538eSMatthew G. Knepley PetscCall(VecDestroy(&dm->coordinates[0].x)); 4636858538eSMatthew G. Knepley dm->coordinates[0].x = c; 4646858538eSMatthew G. Knepley PetscCall(VecDestroy(&dm->coordinates[0].xl)); 4656858538eSMatthew G. Knepley PetscCall(DMCoarsenHookAdd(dm, DMRestrictHook_Coordinates, NULL, NULL)); 4666858538eSMatthew G. Knepley PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_Coordinates, NULL, NULL)); 4673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4686858538eSMatthew G. Knepley } 4696858538eSMatthew G. Knepley 4706858538eSMatthew G. Knepley /*@ 471dce8aebaSBarry Smith DMGetCellCoordinates - Gets a global vector with the cellwise coordinates associated with the `DM`. 4726858538eSMatthew G. Knepley 47320f4b53cSBarry Smith Collective 4746858538eSMatthew G. Knepley 4756858538eSMatthew G. Knepley Input Parameter: 476dce8aebaSBarry Smith . dm - the `DM` 4776858538eSMatthew G. Knepley 4786858538eSMatthew G. Knepley Output Parameter: 4796858538eSMatthew G. Knepley . c - global coordinate vector 4806858538eSMatthew G. Knepley 481dce8aebaSBarry Smith Level: intermediate 482dce8aebaSBarry Smith 483dce8aebaSBarry Smith Notes: 484dce8aebaSBarry Smith This is a borrowed reference, so the user should NOT destroy this vector. When the `DM` is 4850b3275a6SBarry Smith destroyed `c` will no longer be valid. 4866858538eSMatthew G. Knepley 4876858538eSMatthew G. Knepley Each process has only the locally-owned portion of the global coordinates (does NOT have the ghost coordinates). 4886858538eSMatthew G. Knepley 489dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMGetCellCoordinateDM()` 4906858538eSMatthew G. Knepley @*/ 491d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinates(DM dm, Vec *c) 492d71ae5a4SJacob Faibussowitsch { 4936858538eSMatthew G. Knepley PetscFunctionBegin; 4946858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4954f572ea9SToby Isaac PetscAssertPointer(c, 2); 4966858538eSMatthew G. Knepley if (!dm->coordinates[1].x && dm->coordinates[1].xl) { 4976858538eSMatthew G. Knepley DM cdm = NULL; 4986858538eSMatthew G. Knepley 49911ea91a7SMatthew G. Knepley PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 5006858538eSMatthew G. Knepley PetscCall(DMCreateGlobalVector(cdm, &dm->coordinates[1].x)); 5016858538eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].x, "DG coordinates")); 5026858538eSMatthew G. Knepley PetscCall(DMLocalToGlobalBegin(cdm, dm->coordinates[1].xl, INSERT_VALUES, dm->coordinates[1].x)); 5036858538eSMatthew G. Knepley PetscCall(DMLocalToGlobalEnd(cdm, dm->coordinates[1].xl, INSERT_VALUES, dm->coordinates[1].x)); 5046858538eSMatthew G. Knepley } 5056858538eSMatthew G. Knepley *c = dm->coordinates[1].x; 5063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5076858538eSMatthew G. Knepley } 5086858538eSMatthew G. Knepley 5096858538eSMatthew G. Knepley /*@ 510dce8aebaSBarry Smith DMSetCellCoordinates - Sets into the `DM` a global vector that holds the cellwise coordinates 5116858538eSMatthew G. Knepley 51220f4b53cSBarry Smith Collective 5136858538eSMatthew G. Knepley 5146858538eSMatthew G. Knepley Input Parameters: 515dce8aebaSBarry Smith + dm - the `DM` 5166858538eSMatthew G. Knepley - c - cellwise coordinate vector 5176858538eSMatthew G. Knepley 5186858538eSMatthew G. Knepley Level: intermediate 5196858538eSMatthew G. Knepley 520dce8aebaSBarry Smith Notes: 521dce8aebaSBarry Smith The coordinates do not include those for ghost points, which are in the local vector. 522dce8aebaSBarry Smith 52320f4b53cSBarry Smith The vector `c` should be destroyed by the caller. 524dce8aebaSBarry Smith 525dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMGetCellCoordinateDM()` 5266858538eSMatthew G. Knepley @*/ 527d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCellCoordinates(DM dm, Vec c) 528d71ae5a4SJacob Faibussowitsch { 5296858538eSMatthew G. Knepley PetscFunctionBegin; 5306858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5316858538eSMatthew G. Knepley if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2); 5326858538eSMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)c)); 5336858538eSMatthew G. Knepley PetscCall(VecDestroy(&dm->coordinates[1].x)); 5346858538eSMatthew G. Knepley dm->coordinates[1].x = c; 5356858538eSMatthew G. Knepley PetscCall(VecDestroy(&dm->coordinates[1].xl)); 5363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5376858538eSMatthew G. Knepley } 5386858538eSMatthew G. Knepley 5396858538eSMatthew G. Knepley /*@ 540dce8aebaSBarry Smith DMGetCoordinatesLocalSetUp - Prepares a local vector of coordinates, so that `DMGetCoordinatesLocalNoncollective()` can be used as non-collective afterwards. 5416858538eSMatthew G. Knepley 54220f4b53cSBarry Smith Collective 5436858538eSMatthew G. Knepley 5446858538eSMatthew G. Knepley Input Parameter: 545dce8aebaSBarry Smith . dm - the `DM` 5466858538eSMatthew G. Knepley 5476858538eSMatthew G. Knepley Level: advanced 5486858538eSMatthew G. Knepley 549dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinates()`, `DMGetCoordinatesLocalNoncollective()` 5506858538eSMatthew G. Knepley @*/ 551d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocalSetUp(DM dm) 552d71ae5a4SJacob Faibussowitsch { 5536858538eSMatthew G. Knepley PetscFunctionBegin; 5546858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5556858538eSMatthew G. Knepley if (!dm->coordinates[0].xl && dm->coordinates[0].x) { 5566858538eSMatthew G. Knepley DM cdm = NULL; 5579d92f5ffSMatthew G. Knepley PetscInt bs; 5586858538eSMatthew G. Knepley 5596858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &cdm)); 5606858538eSMatthew G. Knepley PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[0].xl)); 5619d92f5ffSMatthew G. Knepley // If the size of the vector is 0, it will not get the right block size 5629d92f5ffSMatthew G. Knepley PetscCall(VecGetBlockSize(dm->coordinates[0].x, &bs)); 5639d92f5ffSMatthew G. Knepley PetscCall(VecSetBlockSize(dm->coordinates[0].xl, bs)); 5646858538eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].xl, "coordinates")); 5656858538eSMatthew G. Knepley PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl)); 5666858538eSMatthew G. Knepley PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl)); 5676858538eSMatthew G. Knepley } 5683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5696858538eSMatthew G. Knepley } 5706858538eSMatthew G. Knepley 5716858538eSMatthew G. Knepley /*@ 572dce8aebaSBarry Smith DMGetCoordinatesLocal - Gets a local vector with the coordinates associated with the `DM`. 5736858538eSMatthew G. Knepley 57420f4b53cSBarry Smith Collective the first time it is called 5756858538eSMatthew G. Knepley 5766858538eSMatthew G. Knepley Input Parameter: 577dce8aebaSBarry Smith . dm - the `DM` 5786858538eSMatthew G. Knepley 5796858538eSMatthew G. Knepley Output Parameter: 5806858538eSMatthew G. Knepley . c - coordinate vector 5816858538eSMatthew G. Knepley 582dce8aebaSBarry Smith Level: intermediate 583dce8aebaSBarry Smith 584dce8aebaSBarry Smith Notes: 5850b3275a6SBarry Smith This is a borrowed reference, so the user should NOT destroy `c` 5866858538eSMatthew G. Knepley 5876858538eSMatthew G. Knepley Each process has the local and ghost coordinates 5886858538eSMatthew G. Knepley 589dce8aebaSBarry Smith For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...) 5906858538eSMatthew G. Knepley and (x_0,y_0,z_0,x_1,y_1,z_1...) 5916858538eSMatthew G. Knepley 592dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`, `DMGetCoordinatesLocalNoncollective()` 5936858538eSMatthew G. Knepley @*/ 594d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocal(DM dm, Vec *c) 595d71ae5a4SJacob Faibussowitsch { 5966858538eSMatthew G. Knepley PetscFunctionBegin; 5976858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5984f572ea9SToby Isaac PetscAssertPointer(c, 2); 5996858538eSMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(dm)); 6006858538eSMatthew G. Knepley *c = dm->coordinates[0].xl; 6013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6026858538eSMatthew G. Knepley } 6036858538eSMatthew G. Knepley 6046858538eSMatthew G. Knepley /*@ 605dce8aebaSBarry Smith DMGetCoordinatesLocalNoncollective - Non-collective version of `DMGetCoordinatesLocal()`. Fails if global coordinates have been set and `DMGetCoordinatesLocalSetUp()` not called. 6066858538eSMatthew G. Knepley 60720f4b53cSBarry Smith Not Collective 6086858538eSMatthew G. Knepley 6096858538eSMatthew G. Knepley Input Parameter: 610dce8aebaSBarry Smith . dm - the `DM` 6116858538eSMatthew G. Knepley 6126858538eSMatthew G. Knepley Output Parameter: 6136858538eSMatthew G. Knepley . c - coordinate vector 6146858538eSMatthew G. Knepley 6156858538eSMatthew G. Knepley Level: advanced 6166858538eSMatthew G. Knepley 617dce8aebaSBarry Smith Note: 618dce8aebaSBarry Smith A previous call to `DMGetCoordinatesLocal()` or `DMGetCoordinatesLocalSetUp()` ensures that a call to this function will not error. 619dce8aebaSBarry Smith 620dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinatesLocal()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()` 6216858538eSMatthew G. Knepley @*/ 622d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocalNoncollective(DM dm, Vec *c) 623d71ae5a4SJacob Faibussowitsch { 6246858538eSMatthew G. Knepley PetscFunctionBegin; 6256858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 6264f572ea9SToby Isaac PetscAssertPointer(c, 2); 6276858538eSMatthew G. Knepley PetscCheck(dm->coordinates[0].xl || !dm->coordinates[0].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called"); 6286858538eSMatthew G. Knepley *c = dm->coordinates[0].xl; 6293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6306858538eSMatthew G. Knepley } 6316858538eSMatthew G. Knepley 6326858538eSMatthew G. Knepley /*@ 633dce8aebaSBarry Smith DMGetCoordinatesLocalTuple - Gets a local vector with the coordinates of specified points and the section describing its layout. 6346858538eSMatthew G. Knepley 63520f4b53cSBarry Smith Not Collective 6366858538eSMatthew G. Knepley 6376858538eSMatthew G. Knepley Input Parameters: 638dce8aebaSBarry Smith + dm - the `DM` 639dce8aebaSBarry Smith - p - the `IS` of points whose coordinates will be returned 6406858538eSMatthew G. Knepley 6416858538eSMatthew G. Knepley Output Parameters: 6420b3275a6SBarry Smith + pCoordSection - the `PetscSection` describing the layout of pCoord, i.e. each point corresponds to one point in `p`, and DOFs correspond to coordinates 6430b3275a6SBarry Smith - pCoord - the `Vec` with coordinates of points in `p` 6446858538eSMatthew G. Knepley 645dce8aebaSBarry Smith Level: advanced 646dce8aebaSBarry Smith 647dce8aebaSBarry Smith Notes: 648dce8aebaSBarry Smith `DMGetCoordinatesLocalSetUp()` must be called first. This function employs `DMGetCoordinatesLocalNoncollective()` so it is not collective. 6496858538eSMatthew G. Knepley 6506858538eSMatthew G. Knepley This creates a new vector, so the user SHOULD destroy this vector 6516858538eSMatthew G. Knepley 6526858538eSMatthew G. Knepley Each process has the local and ghost coordinates 6536858538eSMatthew G. Knepley 654dce8aebaSBarry Smith For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...) 6556858538eSMatthew G. Knepley and (x_0,y_0,z_0,x_1,y_1,z_1...) 6566858538eSMatthew G. Knepley 657dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMSetCoordinatesLocal()`, `DMGetCoordinatesLocal()`, `DMGetCoordinatesLocalNoncollective()`, `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()` 6586858538eSMatthew G. Knepley @*/ 659d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocalTuple(DM dm, IS p, PetscSection *pCoordSection, Vec *pCoord) 660d71ae5a4SJacob Faibussowitsch { 6616858538eSMatthew G. Knepley DM cdm; 6626858538eSMatthew G. Knepley PetscSection cs, newcs; 6636858538eSMatthew G. Knepley Vec coords; 6646858538eSMatthew G. Knepley const PetscScalar *arr; 6656858538eSMatthew G. Knepley PetscScalar *newarr = NULL; 6666858538eSMatthew G. Knepley PetscInt n; 6676858538eSMatthew G. Knepley 6686858538eSMatthew G. Knepley PetscFunctionBegin; 6696858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 6706858538eSMatthew G. Knepley PetscValidHeaderSpecific(p, IS_CLASSID, 2); 6714f572ea9SToby Isaac if (pCoordSection) PetscAssertPointer(pCoordSection, 3); 6724f572ea9SToby Isaac if (pCoord) PetscAssertPointer(pCoord, 4); 6736858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &cdm)); 6746858538eSMatthew G. Knepley PetscCall(DMGetLocalSection(cdm, &cs)); 6756858538eSMatthew G. Knepley PetscCall(DMGetCoordinatesLocal(dm, &coords)); 6766858538eSMatthew G. Knepley PetscCheck(coords, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called or coordinates not set"); 6776858538eSMatthew G. Knepley PetscCheck(cdm && cs, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM not supported"); 6786858538eSMatthew G. Knepley PetscCall(VecGetArrayRead(coords, &arr)); 6796858538eSMatthew G. Knepley PetscCall(PetscSectionExtractDofsFromArray(cs, MPIU_SCALAR, arr, p, &newcs, pCoord ? ((void **)&newarr) : NULL)); 6806858538eSMatthew G. Knepley PetscCall(VecRestoreArrayRead(coords, &arr)); 6816858538eSMatthew G. Knepley if (pCoord) { 6826858538eSMatthew G. Knepley PetscCall(PetscSectionGetStorageSize(newcs, &n)); 6836858538eSMatthew G. Knepley /* set array in two steps to mimic PETSC_OWN_POINTER */ 6846858538eSMatthew G. Knepley PetscCall(VecCreateSeqWithArray(PetscObjectComm((PetscObject)p), 1, n, NULL, pCoord)); 6856858538eSMatthew G. Knepley PetscCall(VecReplaceArray(*pCoord, newarr)); 6866858538eSMatthew G. Knepley } else { 6876858538eSMatthew G. Knepley PetscCall(PetscFree(newarr)); 6886858538eSMatthew G. Knepley } 6899371c9d4SSatish Balay if (pCoordSection) { 6909371c9d4SSatish Balay *pCoordSection = newcs; 6919371c9d4SSatish Balay } else PetscCall(PetscSectionDestroy(&newcs)); 6923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6936858538eSMatthew G. Knepley } 6946858538eSMatthew G. Knepley 6956858538eSMatthew G. Knepley /*@ 696dce8aebaSBarry Smith DMSetCoordinatesLocal - Sets into the `DM` a local vector, including ghost points, that holds the coordinates 6976858538eSMatthew G. Knepley 69820f4b53cSBarry Smith Not Collective 6996858538eSMatthew G. Knepley 7006858538eSMatthew G. Knepley Input Parameters: 701dce8aebaSBarry Smith + dm - the `DM` 7026858538eSMatthew G. Knepley - c - coordinate vector 7036858538eSMatthew G. Knepley 704dce8aebaSBarry Smith Level: intermediate 705dce8aebaSBarry Smith 7066858538eSMatthew G. Knepley Notes: 707dce8aebaSBarry Smith The coordinates of ghost points can be set using `DMSetCoordinates()` 708dce8aebaSBarry Smith followed by `DMGetCoordinatesLocal()`. This is intended to enable the 7096858538eSMatthew G. Knepley setting of ghost coordinates outside of the domain. 7106858538eSMatthew G. Knepley 7110b3275a6SBarry Smith The vector `c` should be destroyed by the caller. 7126858538eSMatthew G. Knepley 713dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinatesLocal()`, `DMSetCoordinates()`, `DMGetCoordinates()`, `DMGetCoordinateDM()` 7146858538eSMatthew G. Knepley @*/ 715d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinatesLocal(DM dm, Vec c) 716d71ae5a4SJacob Faibussowitsch { 7176858538eSMatthew G. Knepley PetscFunctionBegin; 7186858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 7196858538eSMatthew G. Knepley if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2); 7206858538eSMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)c)); 7216858538eSMatthew G. Knepley PetscCall(VecDestroy(&dm->coordinates[0].xl)); 7226858538eSMatthew G. Knepley dm->coordinates[0].xl = c; 7236858538eSMatthew G. Knepley PetscCall(VecDestroy(&dm->coordinates[0].x)); 7243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7256858538eSMatthew G. Knepley } 7266858538eSMatthew G. Knepley 7276858538eSMatthew G. Knepley /*@ 728dce8aebaSBarry Smith DMGetCellCoordinatesLocalSetUp - Prepares a local vector of cellwise coordinates, so that `DMGetCellCoordinatesLocalNoncollective()` can be used as non-collective afterwards. 7296858538eSMatthew G. Knepley 73020f4b53cSBarry Smith Collective 7316858538eSMatthew G. Knepley 7326858538eSMatthew G. Knepley Input Parameter: 733dce8aebaSBarry Smith . dm - the `DM` 7346858538eSMatthew G. Knepley 7356858538eSMatthew G. Knepley Level: advanced 7366858538eSMatthew G. Knepley 737dce8aebaSBarry Smith .seealso: `DM`, `DMGetCellCoordinatesLocalNoncollective()` 7386858538eSMatthew G. Knepley @*/ 739d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinatesLocalSetUp(DM dm) 740d71ae5a4SJacob Faibussowitsch { 7416858538eSMatthew G. Knepley PetscFunctionBegin; 7426858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 7436858538eSMatthew G. Knepley if (!dm->coordinates[1].xl && dm->coordinates[1].x) { 7446858538eSMatthew G. Knepley DM cdm = NULL; 7456858538eSMatthew G. Knepley 74611ea91a7SMatthew G. Knepley PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 7476858538eSMatthew G. Knepley PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[1].xl)); 7486858538eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].xl, "DG coordinates")); 7496858538eSMatthew G. Knepley PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl)); 7506858538eSMatthew G. Knepley PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl)); 7516858538eSMatthew G. Knepley } 7523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7536858538eSMatthew G. Knepley } 7546858538eSMatthew G. Knepley 7556858538eSMatthew G. Knepley /*@ 756dce8aebaSBarry Smith DMGetCellCoordinatesLocal - Gets a local vector with the cellwise coordinates associated with the `DM`. 7576858538eSMatthew G. Knepley 75820f4b53cSBarry Smith Collective 7596858538eSMatthew G. Knepley 7606858538eSMatthew G. Knepley Input Parameter: 761dce8aebaSBarry Smith . dm - the `DM` 7626858538eSMatthew G. Knepley 7636858538eSMatthew G. Knepley Output Parameter: 7646858538eSMatthew G. Knepley . c - coordinate vector 7656858538eSMatthew G. Knepley 766dce8aebaSBarry Smith Level: intermediate 767dce8aebaSBarry Smith 768dce8aebaSBarry Smith Notes: 7696858538eSMatthew G. Knepley This is a borrowed reference, so the user should NOT destroy this vector 7706858538eSMatthew G. Knepley 7716858538eSMatthew G. Knepley Each process has the local and ghost coordinates 7726858538eSMatthew G. Knepley 773dce8aebaSBarry Smith .seealso: `DM`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()`, `DMGetCellCoordinatesLocalNoncollective()` 7746858538eSMatthew G. Knepley @*/ 775d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinatesLocal(DM dm, Vec *c) 776d71ae5a4SJacob Faibussowitsch { 7776858538eSMatthew G. Knepley PetscFunctionBegin; 7786858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 7794f572ea9SToby Isaac PetscAssertPointer(c, 2); 7806858538eSMatthew G. Knepley PetscCall(DMGetCellCoordinatesLocalSetUp(dm)); 7816858538eSMatthew G. Knepley *c = dm->coordinates[1].xl; 7823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7836858538eSMatthew G. Knepley } 7846858538eSMatthew G. Knepley 7856858538eSMatthew G. Knepley /*@ 786dce8aebaSBarry Smith DMGetCellCoordinatesLocalNoncollective - Non-collective version of `DMGetCellCoordinatesLocal()`. Fails if global cellwise coordinates have been set and `DMGetCellCoordinatesLocalSetUp()` not called. 7876858538eSMatthew G. Knepley 78820f4b53cSBarry Smith Not Collective 7896858538eSMatthew G. Knepley 7906858538eSMatthew G. Knepley Input Parameter: 791dce8aebaSBarry Smith . dm - the `DM` 7926858538eSMatthew G. Knepley 7936858538eSMatthew G. Knepley Output Parameter: 7946858538eSMatthew G. Knepley . c - cellwise coordinate vector 7956858538eSMatthew G. Knepley 7966858538eSMatthew G. Knepley Level: advanced 7976858538eSMatthew G. Knepley 798dce8aebaSBarry Smith .seealso: `DM`, `DMGetCellCoordinatesLocalSetUp()`, `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()` 7996858538eSMatthew G. Knepley @*/ 800d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinatesLocalNoncollective(DM dm, Vec *c) 801d71ae5a4SJacob Faibussowitsch { 8026858538eSMatthew G. Knepley PetscFunctionBegin; 8036858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8044f572ea9SToby Isaac PetscAssertPointer(c, 2); 8056858538eSMatthew G. Knepley PetscCheck(dm->coordinates[1].xl || !dm->coordinates[1].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCellCoordinatesLocalSetUp() has not been called"); 8066858538eSMatthew G. Knepley *c = dm->coordinates[1].xl; 8073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8086858538eSMatthew G. Knepley } 8096858538eSMatthew G. Knepley 8106858538eSMatthew G. Knepley /*@ 811dce8aebaSBarry Smith DMSetCellCoordinatesLocal - Sets into the `DM` a local vector including ghost points that holds the cellwise coordinates 8126858538eSMatthew G. Knepley 81320f4b53cSBarry Smith Not Collective 8146858538eSMatthew G. Knepley 8156858538eSMatthew G. Knepley Input Parameters: 816dce8aebaSBarry Smith + dm - the `DM` 8176858538eSMatthew G. Knepley - c - cellwise coordinate vector 8186858538eSMatthew G. Knepley 819dce8aebaSBarry Smith Level: intermediate 820dce8aebaSBarry Smith 8216858538eSMatthew G. Knepley Notes: 822dce8aebaSBarry Smith The coordinates of ghost points can be set using `DMSetCoordinates()` 823dce8aebaSBarry Smith followed by `DMGetCoordinatesLocal()`. This is intended to enable the 8246858538eSMatthew G. Knepley setting of ghost coordinates outside of the domain. 8256858538eSMatthew G. Knepley 8260b3275a6SBarry Smith The vector `c` should be destroyed by the caller. 8276858538eSMatthew G. Knepley 828dce8aebaSBarry Smith .seealso: `DM`, `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinates()`, `DMGetCellCoordinates()`, `DMGetCellCoordinateDM()` 8296858538eSMatthew G. Knepley @*/ 830d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCellCoordinatesLocal(DM dm, Vec c) 831d71ae5a4SJacob Faibussowitsch { 8326858538eSMatthew G. Knepley PetscFunctionBegin; 8336858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8346858538eSMatthew G. Knepley if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2); 8356858538eSMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)c)); 8366858538eSMatthew G. Knepley PetscCall(VecDestroy(&dm->coordinates[1].xl)); 8376858538eSMatthew G. Knepley dm->coordinates[1].xl = c; 8386858538eSMatthew G. Knepley PetscCall(VecDestroy(&dm->coordinates[1].x)); 8393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8406858538eSMatthew G. Knepley } 8416858538eSMatthew G. Knepley 842d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinateField(DM dm, DMField *field) 843d71ae5a4SJacob Faibussowitsch { 8446858538eSMatthew G. Knepley PetscFunctionBegin; 8456858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8464f572ea9SToby Isaac PetscAssertPointer(field, 2); 847213acdd3SPierre Jolivet if (!dm->coordinates[0].field) PetscTryTypeMethod(dm, createcoordinatefield, &dm->coordinates[0].field); 8486858538eSMatthew G. Knepley *field = dm->coordinates[0].field; 8493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8506858538eSMatthew G. Knepley } 8516858538eSMatthew G. Knepley 852d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinateField(DM dm, DMField field) 853d71ae5a4SJacob Faibussowitsch { 8546858538eSMatthew G. Knepley PetscFunctionBegin; 8556858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8566858538eSMatthew G. Knepley if (field) PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 2); 8576858538eSMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)field)); 8586858538eSMatthew G. Knepley PetscCall(DMFieldDestroy(&dm->coordinates[0].field)); 8596858538eSMatthew G. Knepley dm->coordinates[0].field = field; 8603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8616858538eSMatthew G. Knepley } 8626858538eSMatthew G. Knepley 8636858538eSMatthew G. Knepley /*@ 864dce8aebaSBarry Smith DMGetLocalBoundingBox - Returns the bounding box for the piece of the `DM` on this process. 8656858538eSMatthew G. Knepley 86620f4b53cSBarry Smith Not Collective 8676858538eSMatthew G. Knepley 8686858538eSMatthew G. Knepley Input Parameter: 869dce8aebaSBarry Smith . dm - the `DM` 8706858538eSMatthew G. Knepley 8716858538eSMatthew G. Knepley Output Parameters: 8726858538eSMatthew G. Knepley + lmin - local minimum coordinates (length coord dim, optional) 87335cb6cd3SPierre Jolivet - lmax - local maximum coordinates (length coord dim, optional) 8746858538eSMatthew G. Knepley 8756858538eSMatthew G. Knepley Level: beginner 8766858538eSMatthew G. Knepley 877dce8aebaSBarry Smith Note: 878dce8aebaSBarry Smith If the `DM` is a `DMDA` and has no coordinates, the index bounds are returned instead. 8796858538eSMatthew G. Knepley 880dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetBoundingBox()` 8816858538eSMatthew G. Knepley @*/ 882d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetLocalBoundingBox(DM dm, PetscReal lmin[], PetscReal lmax[]) 883d71ae5a4SJacob Faibussowitsch { 8846858538eSMatthew G. Knepley Vec coords = NULL; 8856858538eSMatthew G. Knepley PetscReal min[3] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL}; 8866858538eSMatthew G. Knepley PetscReal max[3] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 8876858538eSMatthew G. Knepley PetscInt cdim, i, j; 8886858538eSMatthew G. Knepley 8896858538eSMatthew G. Knepley PetscFunctionBegin; 8906858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8916858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &cdim)); 892bdf15c88SMatthew G. Knepley PetscCall(DMGetCoordinatesLocal(dm, &coords)); 893bdf15c88SMatthew G. Knepley if (coords) { 894bdf15c88SMatthew G. Knepley const PetscScalar *local_coords; 895bdf15c88SMatthew G. Knepley PetscInt N, Ni; 896bdf15c88SMatthew G. Knepley 897bdf15c88SMatthew G. Knepley for (j = cdim; j < 3; ++j) { 898bdf15c88SMatthew G. Knepley min[j] = 0; 899bdf15c88SMatthew G. Knepley max[j] = 0; 900bdf15c88SMatthew G. Knepley } 901bdf15c88SMatthew G. Knepley PetscCall(VecGetArrayRead(coords, &local_coords)); 902bdf15c88SMatthew G. Knepley PetscCall(VecGetLocalSize(coords, &N)); 903bdf15c88SMatthew G. Knepley Ni = N / cdim; 904bdf15c88SMatthew G. Knepley for (i = 0; i < Ni; ++i) { 905bdf15c88SMatthew G. Knepley for (j = 0; j < cdim; ++j) { 906bdf15c88SMatthew G. Knepley min[j] = PetscMin(min[j], PetscRealPart(local_coords[i * cdim + j])); 907bdf15c88SMatthew G. Knepley max[j] = PetscMax(max[j], PetscRealPart(local_coords[i * cdim + j])); 908bdf15c88SMatthew G. Knepley } 909bdf15c88SMatthew G. Knepley } 910bdf15c88SMatthew G. Knepley PetscCall(VecRestoreArrayRead(coords, &local_coords)); 911bdf15c88SMatthew G. Knepley PetscCall(DMGetCellCoordinatesLocal(dm, &coords)); 9126858538eSMatthew G. Knepley if (coords) { 9136858538eSMatthew G. Knepley PetscCall(VecGetArrayRead(coords, &local_coords)); 9146858538eSMatthew G. Knepley PetscCall(VecGetLocalSize(coords, &N)); 9156858538eSMatthew G. Knepley Ni = N / cdim; 9166858538eSMatthew G. Knepley for (i = 0; i < Ni; ++i) { 917bdf15c88SMatthew G. Knepley for (j = 0; j < cdim; ++j) { 918bdf15c88SMatthew G. Knepley min[j] = PetscMin(min[j], PetscRealPart(local_coords[i * cdim + j])); 919bdf15c88SMatthew G. Knepley max[j] = PetscMax(max[j], PetscRealPart(local_coords[i * cdim + j])); 9206858538eSMatthew G. Knepley } 9216858538eSMatthew G. Knepley } 9226858538eSMatthew G. Knepley PetscCall(VecRestoreArrayRead(coords, &local_coords)); 923bdf15c88SMatthew G. Knepley } 9246858538eSMatthew G. Knepley } else { 9256858538eSMatthew G. Knepley PetscBool isda; 9266858538eSMatthew G. Knepley 9276858538eSMatthew G. Knepley PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &isda)); 9286858538eSMatthew G. Knepley if (isda) PetscCall(DMGetLocalBoundingIndices_DMDA(dm, min, max)); 9296858538eSMatthew G. Knepley } 9306858538eSMatthew G. Knepley if (lmin) PetscCall(PetscArraycpy(lmin, min, cdim)); 9316858538eSMatthew G. Knepley if (lmax) PetscCall(PetscArraycpy(lmax, max, cdim)); 9323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9336858538eSMatthew G. Knepley } 9346858538eSMatthew G. Knepley 9356858538eSMatthew G. Knepley /*@ 936dce8aebaSBarry Smith DMGetBoundingBox - Returns the global bounding box for the `DM`. 9376858538eSMatthew G. Knepley 9386858538eSMatthew G. Knepley Collective 9396858538eSMatthew G. Knepley 9406858538eSMatthew G. Knepley Input Parameter: 941dce8aebaSBarry Smith . dm - the `DM` 9426858538eSMatthew G. Knepley 9436858538eSMatthew G. Knepley Output Parameters: 9446858538eSMatthew G. Knepley + gmin - global minimum coordinates (length coord dim, optional) 94535cb6cd3SPierre Jolivet - gmax - global maximum coordinates (length coord dim, optional) 9466858538eSMatthew G. Knepley 9476858538eSMatthew G. Knepley Level: beginner 9486858538eSMatthew G. Knepley 949dce8aebaSBarry Smith .seealso: `DM`, `DMGetLocalBoundingBox()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()` 9506858538eSMatthew G. Knepley @*/ 951d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetBoundingBox(DM dm, PetscReal gmin[], PetscReal gmax[]) 952d71ae5a4SJacob Faibussowitsch { 9536858538eSMatthew G. Knepley PetscReal lmin[3], lmax[3]; 9546858538eSMatthew G. Knepley PetscInt cdim; 9556858538eSMatthew G. Knepley PetscMPIInt count; 9566858538eSMatthew G. Knepley 9576858538eSMatthew G. Knepley PetscFunctionBegin; 9586858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 9596858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &cdim)); 9606858538eSMatthew G. Knepley PetscCall(PetscMPIIntCast(cdim, &count)); 9616858538eSMatthew G. Knepley PetscCall(DMGetLocalBoundingBox(dm, lmin, lmax)); 9626858538eSMatthew G. Knepley if (gmin) PetscCall(MPIU_Allreduce(lmin, gmin, count, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm))); 9636858538eSMatthew G. Knepley if (gmax) PetscCall(MPIU_Allreduce(lmax, gmax, count, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)dm))); 9643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9656858538eSMatthew G. Knepley } 9666858538eSMatthew G. Knepley 967e44f6aebSMatthew G. Knepley static PetscErrorCode DMCreateAffineCoordinates_Internal(DM dm) 96890612307SJed Brown { 969e44f6aebSMatthew G. Knepley DM cdm; 970e44f6aebSMatthew G. Knepley PetscFE feLinear; 971e44f6aebSMatthew G. Knepley DMPolytopeType ct; 9721df12153SMatthew G. Knepley PetscInt dim, dE, height, cStart, cEnd, gct; 973e44f6aebSMatthew G. Knepley 974e44f6aebSMatthew G. Knepley PetscFunctionBegin; 975e44f6aebSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &cdm)); 976e44f6aebSMatthew G. Knepley PetscCall(DMGetDimension(dm, &dim)); 977e44f6aebSMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &dE)); 9781df12153SMatthew G. Knepley PetscCall(DMPlexGetVTKCellHeight(dm, &height)); 9791df12153SMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dm, height, &cStart, &cEnd)); 980e44f6aebSMatthew G. Knepley if (cEnd > cStart) PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 981e44f6aebSMatthew G. Knepley else ct = DM_POLYTOPE_UNKNOWN; 982e44f6aebSMatthew G. Knepley gct = (PetscInt)ct; 983e44f6aebSMatthew G. Knepley PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &gct, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm))); 984e44f6aebSMatthew G. Knepley ct = (DMPolytopeType)gct; 985e44f6aebSMatthew G. Knepley // Work around current bug in PetscDualSpaceSetUp_Lagrange() 986e44f6aebSMatthew G. Knepley // Can be seen in plex_tutorials-ex10_1 987e44f6aebSMatthew G. Knepley if (ct != DM_POLYTOPE_SEG_PRISM_TENSOR && ct != DM_POLYTOPE_TRI_PRISM_TENSOR && ct != DM_POLYTOPE_QUAD_PRISM_TENSOR) { 988e44f6aebSMatthew G. Knepley PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, dE, ct, 1, -1, &feLinear)); 989e44f6aebSMatthew G. Knepley PetscCall(DMSetField(cdm, 0, NULL, (PetscObject)feLinear)); 990e44f6aebSMatthew G. Knepley PetscCall(PetscFEDestroy(&feLinear)); 991e44f6aebSMatthew G. Knepley PetscCall(DMCreateDS(cdm)); 992e44f6aebSMatthew G. Knepley } 993e44f6aebSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 994e44f6aebSMatthew G. Knepley } 995e44f6aebSMatthew G. Knepley 996e44f6aebSMatthew G. Knepley PetscErrorCode DMGetCoordinateDegree_Internal(DM dm, PetscInt *degree) 997e44f6aebSMatthew G. Knepley { 998e44f6aebSMatthew G. Knepley DM cdm; 999e44f6aebSMatthew G. Knepley PetscFE fe; 1000e44f6aebSMatthew G. Knepley PetscSpace sp; 1001e44f6aebSMatthew G. Knepley PetscClassId id; 1002e44f6aebSMatthew G. Knepley 1003e44f6aebSMatthew G. Knepley PetscFunctionBegin; 1004e44f6aebSMatthew G. Knepley *degree = 1; 1005e44f6aebSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &cdm)); 1006e44f6aebSMatthew G. Knepley PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe)); 1007e44f6aebSMatthew G. Knepley PetscCall(PetscObjectGetClassId((PetscObject)fe, &id)); 1008e44f6aebSMatthew G. Knepley if (id != PETSCFE_CLASSID) PetscFunctionReturn(PETSC_SUCCESS); 1009e44f6aebSMatthew G. Knepley PetscCall(PetscFEGetBasisSpace(fe, &sp)); 1010e44f6aebSMatthew G. Knepley PetscCall(PetscSpaceGetDegree(sp, degree, NULL)); 1011e44f6aebSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 101290612307SJed Brown } 101390612307SJed Brown 10146858538eSMatthew G. Knepley /*@ 1015e44f6aebSMatthew G. Knepley DMSetCoordinateDisc - Set a coordinate space 10166858538eSMatthew G. Knepley 10176858538eSMatthew G. Knepley Input Parameters: 1018dce8aebaSBarry Smith + dm - The `DM` object 10194688473bSSatish Balay . disc - The new coordinate discretization or `NULL` to ensure a coordinate discretization exists 1020e44f6aebSMatthew G. Knepley - project - Project coordinates to new discretization 10216858538eSMatthew G. Knepley 10226858538eSMatthew G. Knepley Level: intermediate 10236858538eSMatthew G. Knepley 1024dce8aebaSBarry Smith Notes: 1025e44f6aebSMatthew G. Knepley A `PetscFE` defines an approximation space using a `PetscSpace`, which represents the basis functions, and a `PetscDualSpace`, which defines the interpolation operation in the space. 1026dce8aebaSBarry Smith 1027dce8aebaSBarry Smith This function takes the current mesh coordinates, which are discretized using some `PetscFE` space, and projects this function into a new `PetscFE` space. 1028e44f6aebSMatthew G. Knepley The coordinate projection is done on the continuous coordinates, but the discontinuous coordinates are not updated. 1029dce8aebaSBarry Smith 10300b3275a6SBarry Smith Developer Note: 1031dce8aebaSBarry Smith With more effort, we could directly project the discontinuous coordinates also. 1032dce8aebaSBarry Smith 1033dce8aebaSBarry Smith .seealso: `DM`, `PetscFE`, `DMGetCoordinateField()` 10346858538eSMatthew G. Knepley @*/ 1035e44f6aebSMatthew G. Knepley PetscErrorCode DMSetCoordinateDisc(DM dm, PetscFE disc, PetscBool project) 1036d71ae5a4SJacob Faibussowitsch { 1037e44f6aebSMatthew G. Knepley DM cdmOld, cdmNew; 10386858538eSMatthew G. Knepley PetscFE discOld; 10396858538eSMatthew G. Knepley PetscClassId classid; 10406858538eSMatthew G. Knepley Vec coordsOld, coordsNew; 10416858538eSMatthew G. Knepley PetscBool same_space = PETSC_TRUE; 1042dd4c3f67SMatthew G. Knepley const char *prefix; 10436858538eSMatthew G. Knepley 10446858538eSMatthew G. Knepley PetscFunctionBegin; 10456858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 10466858538eSMatthew G. Knepley if (disc) PetscValidHeaderSpecific(disc, PETSCFE_CLASSID, 2); 10476858538eSMatthew G. Knepley 10486858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &cdmOld)); 10496858538eSMatthew G. Knepley /* Check current discretization is compatible */ 10506858538eSMatthew G. Knepley PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld)); 10516858538eSMatthew G. Knepley PetscCall(PetscObjectGetClassId((PetscObject)discOld, &classid)); 10526858538eSMatthew G. Knepley if (classid != PETSCFE_CLASSID) { 10536858538eSMatthew G. Knepley if (classid == PETSC_CONTAINER_CLASSID) { 1054e44f6aebSMatthew G. Knepley PetscCall(DMCreateAffineCoordinates_Internal(dm)); 10556858538eSMatthew G. Knepley PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld)); 10566858538eSMatthew G. Knepley } else { 10576858538eSMatthew G. Knepley const char *discname; 10586858538eSMatthew G. Knepley 10596858538eSMatthew G. Knepley PetscCall(PetscObjectGetType((PetscObject)discOld, &discname)); 10606858538eSMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject)discOld), PETSC_ERR_SUP, "Discretization type %s not supported", discname); 10616858538eSMatthew G. Knepley } 10626858538eSMatthew G. Knepley } 1063e44f6aebSMatthew G. Knepley // Linear space has been created by now 10643ba16761SJacob Faibussowitsch if (!disc) PetscFunctionReturn(PETSC_SUCCESS); 1065e44f6aebSMatthew G. Knepley // Check if the new space is the same as the old modulo quadrature 1066e44f6aebSMatthew G. Knepley { 10676858538eSMatthew G. Knepley PetscDualSpace dsOld, ds; 10686858538eSMatthew G. Knepley PetscCall(PetscFEGetDualSpace(discOld, &dsOld)); 10696858538eSMatthew G. Knepley PetscCall(PetscFEGetDualSpace(disc, &ds)); 10706858538eSMatthew G. Knepley PetscCall(PetscDualSpaceEqual(dsOld, ds, &same_space)); 10716858538eSMatthew G. Knepley } 1072e44f6aebSMatthew G. Knepley // Make a fresh clone of the coordinate DM 10736858538eSMatthew G. Knepley PetscCall(DMClone(cdmOld, &cdmNew)); 1074dd4c3f67SMatthew G. Knepley cdmNew->cloneOpts = PETSC_TRUE; 1075dd4c3f67SMatthew G. Knepley PetscCall(PetscObjectGetOptionsPrefix((PetscObject)cdmOld, &prefix)); 1076dd4c3f67SMatthew G. Knepley PetscCall(PetscObjectSetOptionsPrefix((PetscObject)cdmNew, prefix)); 10776858538eSMatthew G. Knepley PetscCall(DMSetField(cdmNew, 0, NULL, (PetscObject)disc)); 10786858538eSMatthew G. Knepley PetscCall(DMCreateDS(cdmNew)); 1079e44f6aebSMatthew G. Knepley { 1080e44f6aebSMatthew G. Knepley PetscDS ds, nds; 1081e44f6aebSMatthew G. Knepley 1082e44f6aebSMatthew G. Knepley PetscCall(DMGetDS(cdmOld, &ds)); 1083e44f6aebSMatthew G. Knepley PetscCall(DMGetDS(cdmNew, &nds)); 1084e44f6aebSMatthew G. Knepley PetscCall(PetscDSCopyConstants(ds, nds)); 1085e44f6aebSMatthew G. Knepley } 10866725e60dSJed Brown if (cdmOld->periodic.setup) { 10876725e60dSJed Brown cdmNew->periodic.setup = cdmOld->periodic.setup; 10886725e60dSJed Brown PetscCall(cdmNew->periodic.setup(cdmNew)); 10896725e60dSJed Brown } 1090dd4c3f67SMatthew G. Knepley if (dm->setfromoptionscalled) PetscCall(DMSetFromOptions(cdmNew)); 1091e44f6aebSMatthew G. Knepley if (project) { 10926858538eSMatthew G. Knepley PetscCall(DMGetCoordinates(dm, &coordsOld)); 10936858538eSMatthew G. Knepley PetscCall(DMCreateGlobalVector(cdmNew, &coordsNew)); 1094f0ac6e35SMatthew G. Knepley if (same_space) { 1095f0ac6e35SMatthew G. Knepley // Need to copy so that the new vector has the right dm 1096f0ac6e35SMatthew G. Knepley PetscCall(VecCopy(coordsOld, coordsNew)); 1097e44f6aebSMatthew G. Knepley } else { 1098e44f6aebSMatthew G. Knepley Mat In; 1099e44f6aebSMatthew G. Knepley 1100e44f6aebSMatthew G. Knepley PetscCall(DMCreateInterpolation(cdmOld, cdmNew, &In, NULL)); 1101e44f6aebSMatthew G. Knepley PetscCall(MatMult(In, coordsOld, coordsNew)); 1102e44f6aebSMatthew G. Knepley PetscCall(MatDestroy(&In)); 1103e44f6aebSMatthew G. Knepley } 1104e44f6aebSMatthew G. Knepley PetscCall(DMSetCoordinates(dm, coordsNew)); 1105e44f6aebSMatthew G. Knepley PetscCall(VecDestroy(&coordsNew)); 11066858538eSMatthew G. Knepley } 11076858538eSMatthew G. Knepley /* Set new coordinate structures */ 11086858538eSMatthew G. Knepley PetscCall(DMSetCoordinateField(dm, NULL)); 11096858538eSMatthew G. Knepley PetscCall(DMSetCoordinateDM(dm, cdmNew)); 11106858538eSMatthew G. Knepley PetscCall(DMDestroy(&cdmNew)); 11113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11126858538eSMatthew G. Knepley } 11136858538eSMatthew G. Knepley 11146858538eSMatthew G. Knepley /*@ 11150b3275a6SBarry Smith DMLocatePoints - Locate the points in `v` in the mesh and return a `PetscSF` of the containing cells 11166858538eSMatthew G. Knepley 111720f4b53cSBarry Smith Collective 11186858538eSMatthew G. Knepley 11196858538eSMatthew G. Knepley Input Parameters: 1120dce8aebaSBarry Smith + dm - The `DM` 1121dce8aebaSBarry Smith - ltype - The type of point location, e.g. `DM_POINTLOCATION_NONE` or `DM_POINTLOCATION_NEAREST` 11226858538eSMatthew G. Knepley 11236858538eSMatthew G. Knepley Input/Output Parameters: 1124dce8aebaSBarry Smith + v - The `Vec` of points, on output contains the nearest mesh points to the given points if `DM_POINTLOCATION_NEAREST` is used 112520f4b53cSBarry Smith - cellSF - Points to either `NULL`, or a `PetscSF` with guesses for which cells contain each point; 11260b3275a6SBarry Smith on output, the `PetscSF` containing the MPI ranks and local indices of the containing points 11276858538eSMatthew G. Knepley 11286858538eSMatthew G. Knepley Level: developer 11296858538eSMatthew G. Knepley 11306858538eSMatthew G. Knepley Notes: 11310b3275a6SBarry Smith To do a search of the local cells of the mesh, `v` should have `PETSC_COMM_SELF` as its communicator. 11320b3275a6SBarry Smith To do a search of all the cells in the distributed mesh, `v` should have the same MPI communicator as `dm`. 11336858538eSMatthew G. Knepley 1134d8206211SMatthew G. Knepley Points will only be located in owned cells, not overlap cells arising from `DMPlexDistribute()` or other overlapping distributions. 1135d8206211SMatthew G. Knepley 113620f4b53cSBarry Smith If *cellSF is `NULL` on input, a `PetscSF` will be created. 113720f4b53cSBarry Smith If *cellSF is not `NULL` on input, it should point to an existing `PetscSF`, whose graph will be used as initial guesses. 11386858538eSMatthew G. Knepley 11396858538eSMatthew G. Knepley An array that maps each point to its containing cell can be obtained with 1140dce8aebaSBarry Smith .vb 1141dce8aebaSBarry Smith const PetscSFNode *cells; 1142dce8aebaSBarry Smith PetscInt nFound; 1143dce8aebaSBarry Smith const PetscInt *found; 11446858538eSMatthew G. Knepley 1145dce8aebaSBarry Smith PetscSFGetGraph(cellSF,NULL,&nFound,&found,&cells); 1146dce8aebaSBarry Smith .ve 11476858538eSMatthew G. Knepley 11480b3275a6SBarry Smith Where cells[i].rank is the MPI rank of the process owning the cell containing point found[i] (or i if found == NULL), and cells[i].index is 11490b3275a6SBarry Smith the index of the cell in its MPI process' local numbering. This rank is in the communicator for `v`, so if `v` is on `PETSC_COMM_SELF` then the rank will always be 0. 11506858538eSMatthew G. Knepley 1151dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMPointLocationType` 11526858538eSMatthew G. Knepley @*/ 1153d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocatePoints(DM dm, Vec v, DMPointLocationType ltype, PetscSF *cellSF) 1154d71ae5a4SJacob Faibussowitsch { 11556858538eSMatthew G. Knepley PetscFunctionBegin; 11566858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 11576858538eSMatthew G. Knepley PetscValidHeaderSpecific(v, VEC_CLASSID, 2); 11584f572ea9SToby Isaac PetscAssertPointer(cellSF, 4); 11596858538eSMatthew G. Knepley if (*cellSF) { 11606858538eSMatthew G. Knepley PetscMPIInt result; 11616858538eSMatthew G. Knepley 11626858538eSMatthew G. Knepley PetscValidHeaderSpecific(*cellSF, PETSCSF_CLASSID, 4); 11636858538eSMatthew G. Knepley PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)v), PetscObjectComm((PetscObject)*cellSF), &result)); 11646858538eSMatthew 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"); 11656858538eSMatthew G. Knepley } else { 11666858538eSMatthew G. Knepley PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)v), cellSF)); 11676858538eSMatthew G. Knepley } 11686858538eSMatthew G. Knepley PetscCall(PetscLogEventBegin(DM_LocatePoints, dm, 0, 0, 0)); 1169dbbe0bcdSBarry Smith PetscUseTypeMethod(dm, locatepoints, v, ltype, *cellSF); 11706858538eSMatthew G. Knepley PetscCall(PetscLogEventEnd(DM_LocatePoints, dm, 0, 0, 0)); 11713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11726858538eSMatthew G. Knepley } 1173