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