1c4762a1bSJed Brown static char help[] = "Tests for point location\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscsf.h> 4c4762a1bSJed Brown #include <petscdmplex.h> 5c4762a1bSJed Brown 6*65f07beaSMatthew G. Knepley /* To inspect the location process, use 7*65f07beaSMatthew G. Knepley 8*65f07beaSMatthew G. Knepley -dm_plex_print_locate 5 -info :dm 9*65f07beaSMatthew G. Knepley */ 10*65f07beaSMatthew G. Knepley 11274abdb2SMatthew G. Knepley typedef struct { 12274abdb2SMatthew G. Knepley PetscBool centroids; 13274abdb2SMatthew G. Knepley PetscBool custom; 14274abdb2SMatthew G. Knepley } AppCtx; 15274abdb2SMatthew G. Knepley 16274abdb2SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 17274abdb2SMatthew G. Knepley { 18274abdb2SMatthew G. Knepley PetscFunctionBeginUser; 19274abdb2SMatthew G. Knepley options->centroids = PETSC_TRUE; 20274abdb2SMatthew G. Knepley options->custom = PETSC_FALSE; 21274abdb2SMatthew G. Knepley 22274abdb2SMatthew G. Knepley PetscOptionsBegin(comm, "", "Point Location Options", "DMPLEX"); 23274abdb2SMatthew G. Knepley PetscCall(PetscOptionsBool("-centroids", "Locate cell centroids", "ex17.c", options->centroids, &options->centroids, NULL)); 24274abdb2SMatthew G. Knepley PetscCall(PetscOptionsBool("-custom", "Locate user-defined points", "ex17.c", options->custom, &options->custom, NULL)); 25274abdb2SMatthew G. Knepley PetscOptionsEnd(); 263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 27274abdb2SMatthew G. Knepley } 28274abdb2SMatthew G. Knepley 29d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) 30d71ae5a4SJacob Faibussowitsch { 31c4762a1bSJed Brown PetscFunctionBeginUser; 329566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 339566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 349566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 359566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 37c4762a1bSJed Brown } 38c4762a1bSJed Brown 39274abdb2SMatthew G. Knepley static PetscErrorCode TestCentroidLocation(DM dm, AppCtx *user) 40d71ae5a4SJacob Faibussowitsch { 413285882aSMatthew G. Knepley Vec points; 423285882aSMatthew G. Knepley PetscSF cellSF = NULL; 433285882aSMatthew G. Knepley const PetscSFNode *cells; 443285882aSMatthew G. Knepley PetscScalar *a; 453285882aSMatthew G. Knepley PetscInt cdim, n; 46c4762a1bSJed Brown PetscInt cStart, cEnd, c; 47c4762a1bSJed Brown 48c4762a1bSJed Brown PetscFunctionBeginUser; 493ba16761SJacob Faibussowitsch if (!user->centroids) PetscFunctionReturn(PETSC_SUCCESS); 509566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 518fb5bd83SMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(dm)); 529566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 53c4762a1bSJed Brown /* Locate all centroids */ 549566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, (cEnd - cStart) * cdim, &points)); 559566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(points, cdim)); 569566063dSJacob Faibussowitsch PetscCall(VecGetArray(points, &a)); 57c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 58c4762a1bSJed Brown PetscReal centroid[3]; 593285882aSMatthew G. Knepley PetscInt off = (c - cStart) * cdim, d; 60c4762a1bSJed Brown 619566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL)); 623285882aSMatthew G. Knepley for (d = 0; d < cdim; ++d) a[off + d] = centroid[d]; 63c4762a1bSJed Brown } 649566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(points, &a)); 659566063dSJacob Faibussowitsch PetscCall(DMLocatePoints(dm, points, DM_POINTLOCATION_NONE, &cellSF)); 669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&points)); 679566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(cellSF, NULL, &n, NULL, &cells)); 683285882aSMatthew G. Knepley if (n != (cEnd - cStart)) { 693285882aSMatthew G. Knepley for (c = 0; c < n; ++c) { 7063a3b9bcSJacob Faibussowitsch if (cells[c].index != c + cStart) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Could not locate centroid of cell %" PetscInt_FMT ", error %" PetscInt_FMT "\n", c + cStart, cells[c].index)); 713285882aSMatthew G. Knepley } 7263a3b9bcSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Located %" PetscInt_FMT " points instead of %" PetscInt_FMT, n, cEnd - cStart); 733285882aSMatthew G. Knepley } 74ad540459SPierre Jolivet for (c = cStart; c < cEnd; ++c) PetscCheck(cells[c - cStart].index == c, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not locate centroid of cell %" PetscInt_FMT ", instead found %" PetscInt_FMT, c, cells[c - cStart].index); 759566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&cellSF)); 763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 77c4762a1bSJed Brown } 78c4762a1bSJed Brown 79274abdb2SMatthew G. Knepley static PetscErrorCode TestCustomLocation(DM dm, AppCtx *user) 80274abdb2SMatthew G. Knepley { 81274abdb2SMatthew G. Knepley PetscSF cellSF = NULL; 82274abdb2SMatthew G. Knepley const PetscSFNode *cells; 83274abdb2SMatthew G. Knepley const PetscInt *found; 84274abdb2SMatthew G. Knepley Vec points; 85274abdb2SMatthew G. Knepley PetscScalar coords[2] = {0.5, 0.5}; 86274abdb2SMatthew G. Knepley PetscInt cdim, Np = 1, Nfd; 87274abdb2SMatthew G. Knepley PetscMPIInt rank; 88274abdb2SMatthew G. Knepley MPI_Comm comm; 89274abdb2SMatthew G. Knepley 90274abdb2SMatthew G. Knepley PetscFunctionBeginUser; 913ba16761SJacob Faibussowitsch if (!user->custom) PetscFunctionReturn(PETSC_SUCCESS); 92274abdb2SMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &cdim)); 93274abdb2SMatthew G. Knepley 94274abdb2SMatthew G. Knepley // Locate serially on each process 95274abdb2SMatthew G. Knepley PetscCall(VecCreate(PETSC_COMM_SELF, &points)); 96274abdb2SMatthew G. Knepley PetscCall(VecSetBlockSize(points, cdim)); 97274abdb2SMatthew G. Knepley PetscCall(VecSetSizes(points, Np * cdim, PETSC_DETERMINE)); 98274abdb2SMatthew G. Knepley PetscCall(VecSetFromOptions(points)); 99274abdb2SMatthew G. Knepley for (PetscInt p = 0; p < Np; ++p) { 100274abdb2SMatthew G. Knepley const PetscInt idx[2] = {p * cdim, p * cdim + 1}; 101274abdb2SMatthew G. Knepley PetscCall(VecSetValues(points, cdim, idx, coords, INSERT_VALUES)); 102274abdb2SMatthew G. Knepley } 103274abdb2SMatthew G. Knepley PetscCall(VecAssemblyBegin(points)); 104274abdb2SMatthew G. Knepley PetscCall(VecAssemblyEnd(points)); 105274abdb2SMatthew G. Knepley 106274abdb2SMatthew G. Knepley PetscCall(DMLocatePoints(dm, points, DM_POINTLOCATION_NONE, &cellSF)); 107274abdb2SMatthew G. Knepley 108274abdb2SMatthew G. Knepley PetscCall(PetscSFGetGraph(cellSF, NULL, &Nfd, &found, &cells)); 109274abdb2SMatthew G. Knepley PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 110274abdb2SMatthew G. Knepley PetscCallMPI(MPI_Comm_rank(comm, &rank)); 111274abdb2SMatthew G. Knepley PetscCall(PetscSynchronizedPrintf(comm, "[%d] Found %" PetscInt_FMT " particles\n", rank, Nfd)); 112274abdb2SMatthew G. Knepley for (PetscInt p = 0; p < Nfd; ++p) { 113274abdb2SMatthew G. Knepley const PetscInt point = found ? found[p] : p; 114274abdb2SMatthew G. Knepley const PetscScalar *array; 115274abdb2SMatthew G. Knepley PetscScalar *ccoords = NULL; 116274abdb2SMatthew G. Knepley PetscInt numCoords; 117274abdb2SMatthew G. Knepley PetscBool isDG; 118274abdb2SMatthew G. Knepley 119274abdb2SMatthew G. Knepley // Since the v comm is SELF, rank is always 0 120274abdb2SMatthew G. Knepley PetscCall(PetscSynchronizedPrintf(comm, " point %" PetscInt_FMT " cell %" PetscInt_FMT "\n", point, cells[p].index)); 121274abdb2SMatthew G. Knepley PetscCall(DMPlexGetCellCoordinates(dm, cells[p].index, &isDG, &numCoords, &array, &ccoords)); 122274abdb2SMatthew G. Knepley for (PetscInt c = 0; c < numCoords / cdim; ++c) { 123274abdb2SMatthew G. Knepley PetscCall(PetscSynchronizedPrintf(comm, " ")); 124274abdb2SMatthew G. Knepley for (PetscInt d = 0; d < cdim; ++d) PetscCall(PetscSynchronizedPrintf(comm, " %g", (double)PetscRealPart(ccoords[c * cdim + d]))); 125274abdb2SMatthew G. Knepley PetscCall(PetscSynchronizedPrintf(comm, "\n")); 126274abdb2SMatthew G. Knepley } 127274abdb2SMatthew G. Knepley PetscCall(DMPlexRestoreCellCoordinates(dm, cells[p].index, &isDG, &numCoords, &array, &ccoords)); 128274abdb2SMatthew G. Knepley } 129274abdb2SMatthew G. Knepley PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT)); 130274abdb2SMatthew G. Knepley 131274abdb2SMatthew G. Knepley PetscCall(PetscSFDestroy(&cellSF)); 132274abdb2SMatthew G. Knepley PetscCall(VecDestroy(&points)); 1333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 134274abdb2SMatthew G. Knepley } 135274abdb2SMatthew G. Knepley 136d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 137d71ae5a4SJacob Faibussowitsch { 138c4762a1bSJed Brown DM dm; 139274abdb2SMatthew G. Knepley AppCtx user; 140c4762a1bSJed Brown 141327415f7SBarry Smith PetscFunctionBeginUser; 1429566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 143274abdb2SMatthew G. Knepley PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 1449566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm)); 145274abdb2SMatthew G. Knepley PetscCall(TestCentroidLocation(dm, &user)); 146274abdb2SMatthew G. Knepley PetscCall(TestCustomLocation(dm, &user)); 1479566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1489566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 149b122ec5aSJacob Faibussowitsch return 0; 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152c4762a1bSJed Brown /*TEST 153c4762a1bSJed Brown 154b26b5bf9SMatthew G. Knepley testset: 155b26b5bf9SMatthew G. Knepley args: -dm_plex_dim 1 -dm_plex_box_faces 10 156b26b5bf9SMatthew G. Knepley 157c4762a1bSJed Brown test: 1581af33867SMatthew G. Knepley suffix: seg 159b26b5bf9SMatthew G. Knepley 160b26b5bf9SMatthew G. Knepley test: 161b26b5bf9SMatthew G. Knepley suffix: seg_hash 162b26b5bf9SMatthew G. Knepley args: -dm_refine 2 -dm_plex_hash_location 1631af33867SMatthew G. Knepley 1643285882aSMatthew G. Knepley testset: 1653285882aSMatthew G. Knepley args: -dm_plex_box_faces 5,5 1663285882aSMatthew G. Knepley 1671af33867SMatthew G. Knepley test: 1681af33867SMatthew G. Knepley suffix: tri 169c4762a1bSJed Brown requires: triangle 1703285882aSMatthew G. Knepley 1713285882aSMatthew G. Knepley test: 1723285882aSMatthew G. Knepley suffix: tri_hash 1733285882aSMatthew G. Knepley requires: triangle 1743285882aSMatthew G. Knepley args: -dm_refine 2 -dm_plex_hash_location 1751af33867SMatthew G. Knepley 1761af33867SMatthew G. Knepley test: 1771af33867SMatthew G. Knepley suffix: quad 1783285882aSMatthew G. Knepley args: -dm_plex_simplex 0 1793285882aSMatthew G. Knepley 1803285882aSMatthew G. Knepley test: 181dd301514SZach Atkins suffix: quad_order_2 182dd301514SZach Atkins args: -dm_plex_simplex 0 -dm_coord_petscspace_degree 2 183dd301514SZach Atkins 184dd301514SZach Atkins test: 1853285882aSMatthew G. Knepley suffix: quad_hash 1863285882aSMatthew G. Knepley args: -dm_plex_simplex 0 -dm_refine 2 -dm_plex_hash_location 1873285882aSMatthew G. Knepley 1883285882aSMatthew G. Knepley testset: 1893285882aSMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 1901af33867SMatthew G. Knepley 1911af33867SMatthew G. Knepley test: 1921af33867SMatthew G. Knepley suffix: tet 1931af33867SMatthew G. Knepley requires: ctetgen 1943285882aSMatthew G. Knepley 1953285882aSMatthew G. Knepley test: 1963285882aSMatthew G. Knepley suffix: tet_hash 1973285882aSMatthew G. Knepley requires: ctetgen 1983285882aSMatthew G. Knepley args: -dm_refine 1 -dm_plex_hash_location 1991af33867SMatthew G. Knepley 2001af33867SMatthew G. Knepley test: 2011af33867SMatthew G. Knepley suffix: hex 2023285882aSMatthew G. Knepley args: -dm_plex_simplex 0 2033285882aSMatthew G. Knepley 2043285882aSMatthew G. Knepley test: 2053285882aSMatthew G. Knepley suffix: hex_hash 2063285882aSMatthew G. Knepley args: -dm_plex_simplex 0 -dm_refine 1 -dm_plex_hash_location 207c4762a1bSJed Brown 208dd301514SZach Atkins test: 209dd301514SZach Atkins suffix: hex_order_2 210dd301514SZach Atkins args: -dm_plex_simplex 0 -dm_refine 1 -dm_coord_petscspace_degree 2 211dd301514SZach Atkins nsize: 2 212dd301514SZach Atkins 213dd301514SZach Atkins test: 214dd301514SZach Atkins suffix: hex_order_3 215dd301514SZach Atkins args: -dm_plex_simplex 0 -dm_coord_petscspace_degree 3 216dd301514SZach Atkins 217274abdb2SMatthew G. Knepley testset: 218274abdb2SMatthew G. Knepley args: -centroids 0 -custom \ 219274abdb2SMatthew G. Knepley -dm_plex_simplex 0 -dm_plex_box_faces 21,21 -dm_distribute_overlap 4 -petscpartitioner_type simple 220274abdb2SMatthew G. Knepley nsize: 2 221274abdb2SMatthew G. Knepley 222274abdb2SMatthew G. Knepley test: 223274abdb2SMatthew G. Knepley suffix: quad_overlap 224274abdb2SMatthew G. Knepley args: -dm_plex_hash_location {{0 1}} 225274abdb2SMatthew G. Knepley 226*65f07beaSMatthew G. Knepley # Test location on a Monge Manifold 227*65f07beaSMatthew G. Knepley testset: 228*65f07beaSMatthew G. Knepley args: -dm_refine 3 -dm_coord_space 0 \ 229*65f07beaSMatthew G. Knepley -dm_plex_option_phases proj_ -cdm_proj_dm_plex_coordinate_dim 3 -proj_dm_coord_space \ 230*65f07beaSMatthew G. Knepley -proj_dm_coord_remap -proj_dm_coord_map sinusoid -proj_dm_coord_map_params 0.1,1.,1. 231*65f07beaSMatthew G. Knepley 232*65f07beaSMatthew G. Knepley test: 233*65f07beaSMatthew G. Knepley requires: triangle 234*65f07beaSMatthew G. Knepley suffix: tri_monge 235*65f07beaSMatthew G. Knepley 236*65f07beaSMatthew G. Knepley test: 237*65f07beaSMatthew G. Knepley suffix: quad_monge 238*65f07beaSMatthew G. Knepley args: -dm_plex_simplex 0 239*65f07beaSMatthew G. Knepley 240c4762a1bSJed Brown TEST*/ 241