1c4762a1bSJed Brown static char help[] = "Test DMCreateLocalVector_Plex, DMPlexGetCellFields and DMPlexRestoreCellFields work properly for 0 fields/cells/DS dimension\n\n"; 2c4762a1bSJed Brown static char FILENAME[] = "ex25.c"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscdmplex.h> 5c4762a1bSJed Brown #include <petscds.h> 6c4762a1bSJed Brown #include <petscsnes.h> 7c4762a1bSJed Brown 8c4762a1bSJed Brown 9c4762a1bSJed Brown typedef struct { 10c4762a1bSJed Brown PetscInt test; 11c4762a1bSJed Brown PetscInt dim; /* The topological mesh dimension */ 12c4762a1bSJed Brown PetscInt faces[3]; /* Number of faces per dimension */ 13c4762a1bSJed Brown PetscBool simplex; /* Use simplices or hexes */ 14c4762a1bSJed Brown PetscBool interpolate; /* Interpolate mesh */ 15c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN]; /* Import mesh from file */ 16c4762a1bSJed Brown } AppCtx; 17c4762a1bSJed Brown 18c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 19c4762a1bSJed Brown { 20c4762a1bSJed Brown PetscInt dim; 21c4762a1bSJed Brown PetscErrorCode ierr; 22c4762a1bSJed Brown 23c4762a1bSJed Brown PetscFunctionBegin; 24c4762a1bSJed Brown options->test = 0; 25c4762a1bSJed Brown options->dim = 3; 26c4762a1bSJed Brown options->simplex = PETSC_TRUE; 27c4762a1bSJed Brown options->interpolate = PETSC_FALSE; 28c4762a1bSJed Brown options->filename[0] = '\0'; 29c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Zero-sized DMPlexGetCellFields Test Options", "DMPLEX");CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = PetscOptionsBoundedInt("-test", "Test to run", FILENAME, options->test, &options->test, NULL,0);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", FILENAME, options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr); 32c4762a1bSJed Brown if (options->dim > 3) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "dimension set to %d, must be <= 3", options->dim); 33c4762a1bSJed Brown ierr = PetscOptionsBool("-simplex", "Use simplices if true, otherwise hexes", FILENAME, options->simplex, &options->simplex, NULL);CHKERRQ(ierr); 34c4762a1bSJed Brown ierr = PetscOptionsBool("-interpolate", "Interpolate the mesh", FILENAME, options->interpolate, &options->interpolate, NULL);CHKERRQ(ierr); 35589a23caSBarry Smith ierr = PetscOptionsString("-filename", "The mesh file", FILENAME, options->filename, options->filename, sizeof(options->filename), NULL);CHKERRQ(ierr); 36c4762a1bSJed Brown options->faces[0] = 1; options->faces[1] = 1; options->faces[2] = 1; 37c4762a1bSJed Brown dim = options->dim; 38c4762a1bSJed Brown ierr = PetscOptionsIntArray("-faces", "Number of faces per dimension", FILENAME, options->faces, &dim, NULL);CHKERRQ(ierr); 39c4762a1bSJed Brown if (dim) options->dim = dim; 40c4762a1bSJed Brown ierr = PetscOptionsEnd(); 41c4762a1bSJed Brown PetscFunctionReturn(0); 42c4762a1bSJed Brown } 43c4762a1bSJed Brown 44c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *options, DM *dm) 45c4762a1bSJed Brown { 46c4762a1bSJed Brown PetscInt dim = options->dim; 47c4762a1bSJed Brown PetscInt *faces = options->faces; 48c4762a1bSJed Brown PetscBool simplex = options->simplex; 49c4762a1bSJed Brown PetscBool interpolate = options->interpolate; 50c4762a1bSJed Brown const char *filename = options->filename; 51c4762a1bSJed Brown size_t len; 52c4762a1bSJed Brown PetscMPIInt rank; 53c4762a1bSJed Brown PetscErrorCode ierr; 54c4762a1bSJed Brown 55c4762a1bSJed Brown PetscFunctionBegin; 56*ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 57c4762a1bSJed Brown ierr = PetscStrlen(filename, &len);CHKERRQ(ierr); 58c4762a1bSJed Brown if (len) { 59c4762a1bSJed Brown ierr = DMPlexCreateFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr); 60c4762a1bSJed Brown ierr = DMGetDimension(*dm, &options->dim);CHKERRQ(ierr); 61c4762a1bSJed Brown } else { 62c4762a1bSJed Brown ierr = DMPlexCreateBoxMesh(comm, dim, simplex, faces, NULL, NULL, NULL, interpolate, dm);CHKERRQ(ierr); 63c4762a1bSJed Brown } 64c4762a1bSJed Brown PetscFunctionReturn(0); 65c4762a1bSJed Brown } 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* no discretization is given so DMGetNumFields yields 0 */ 68c4762a1bSJed Brown static PetscErrorCode test0(DM dm, AppCtx *options) 69c4762a1bSJed Brown { 70c4762a1bSJed Brown Vec locX; 71c4762a1bSJed Brown PetscErrorCode ierr; 72c4762a1bSJed Brown 73c4762a1bSJed Brown PetscFunctionBegin; 74c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr); 75c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr); 76c4762a1bSJed Brown PetscFunctionReturn(0); 77c4762a1bSJed Brown } 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* no discretization is given so DMGetNumFields and PetscDSGetTotalDimension yield 0 */ 80c4762a1bSJed Brown static PetscErrorCode test1(DM dm, AppCtx *options) 81c4762a1bSJed Brown { 82c4762a1bSJed Brown IS cells; 83c4762a1bSJed Brown Vec locX, locX_t, locA; 84c4762a1bSJed Brown PetscScalar *u, *u_t, *a; 85c4762a1bSJed Brown PetscErrorCode ierr; 86c4762a1bSJed Brown 87c4762a1bSJed Brown PetscFunctionBegin; 88c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &cells);CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &locX_t);CHKERRQ(ierr); 91c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &locA);CHKERRQ(ierr); 92c4762a1bSJed Brown ierr = DMPlexGetCellFields( dm, cells, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr); 93c4762a1bSJed Brown ierr = DMPlexRestoreCellFields(dm, cells, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr); 94c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr); 95c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &locX_t);CHKERRQ(ierr); 96c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &locA);CHKERRQ(ierr); 97c4762a1bSJed Brown ierr = ISDestroy(&cells);CHKERRQ(ierr); 98c4762a1bSJed Brown PetscFunctionReturn(0); 99c4762a1bSJed Brown } 100c4762a1bSJed Brown 101c4762a1bSJed Brown /* no discretization is given so DMGetNumFields and PetscDSGetTotalDimension yield 0 */ 102c4762a1bSJed Brown static PetscErrorCode test2(DM dm, AppCtx *options) 103c4762a1bSJed Brown { 104c4762a1bSJed Brown IS cells; 105c4762a1bSJed Brown Vec locX, locX_t, locA; 106c4762a1bSJed Brown PetscScalar *u, *u_t, *a; 107c4762a1bSJed Brown PetscMPIInt rank; 108c4762a1bSJed Brown PetscErrorCode ierr; 109c4762a1bSJed Brown 110c4762a1bSJed Brown PetscFunctionBegin; 111*ffc4695bSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRMPI(ierr); 112c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF, rank ? 0 : 1, 0, 1, &cells);CHKERRQ(ierr); 113c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr); 114c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &locX_t);CHKERRQ(ierr); 115c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &locA);CHKERRQ(ierr); 116c4762a1bSJed Brown ierr = DMPlexGetCellFields( dm, cells, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr); 117c4762a1bSJed Brown ierr = DMPlexRestoreCellFields(dm, cells, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr); 118c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &locX_t);CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &locA);CHKERRQ(ierr); 121c4762a1bSJed Brown ierr = ISDestroy(&cells);CHKERRQ(ierr); 122c4762a1bSJed Brown PetscFunctionReturn(0); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown 125c4762a1bSJed Brown static PetscErrorCode test3(DM dm, AppCtx *options) 126c4762a1bSJed Brown { 127c4762a1bSJed Brown PetscDS ds; 128c4762a1bSJed Brown PetscFE fe; 129c4762a1bSJed Brown PetscErrorCode ierr; 130c4762a1bSJed Brown 131c4762a1bSJed Brown PetscFunctionBegin; 132c4762a1bSJed Brown ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 133c4762a1bSJed Brown ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), options->dim, 1, options->simplex, NULL, -1, &fe);CHKERRQ(ierr); 134c4762a1bSJed Brown ierr = PetscDSSetDiscretization(ds, 0, (PetscObject)fe);CHKERRQ(ierr); 135c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 136c4762a1bSJed Brown ierr = test1(dm, options);CHKERRQ(ierr); 137c4762a1bSJed Brown PetscFunctionReturn(0); 138c4762a1bSJed Brown } 139c4762a1bSJed Brown 140c4762a1bSJed Brown static PetscErrorCode test4(DM dm, AppCtx *options) 141c4762a1bSJed Brown { 142c4762a1bSJed Brown PetscFE fe; 143c4762a1bSJed Brown PetscErrorCode ierr; 144c4762a1bSJed Brown 145c4762a1bSJed Brown PetscFunctionBegin; 146c4762a1bSJed Brown ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), options->dim, 1, options->simplex, NULL, -1, &fe);CHKERRQ(ierr); 147c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject)fe);CHKERRQ(ierr); 148c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 149c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 150c4762a1bSJed Brown ierr = test2(dm, options);CHKERRQ(ierr); 151c4762a1bSJed Brown PetscFunctionReturn(0); 152c4762a1bSJed Brown } 153c4762a1bSJed Brown 154c4762a1bSJed Brown static PetscErrorCode test5(DM dm, AppCtx *options) 155c4762a1bSJed Brown { 156c4762a1bSJed Brown IS cells; 157c4762a1bSJed Brown Vec locX, locX_t, locA; 158c4762a1bSJed Brown PetscScalar *u, *u_t, *a; 159c4762a1bSJed Brown PetscErrorCode ierr; 160c4762a1bSJed Brown 161c4762a1bSJed Brown PetscFunctionBegin; 162c4762a1bSJed Brown locX_t = NULL; 163c4762a1bSJed Brown locA = NULL; 164c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &cells);CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = DMPlexGetCellFields( dm, cells, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = DMPlexRestoreCellFields(dm, cells, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr); 168c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr); 169c4762a1bSJed Brown ierr = ISDestroy(&cells);CHKERRQ(ierr); 170c4762a1bSJed Brown PetscFunctionReturn(0); 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 173c4762a1bSJed Brown static PetscErrorCode test6(DM dm, AppCtx *options) 174c4762a1bSJed Brown { 175c4762a1bSJed Brown IS cells; 176c4762a1bSJed Brown Vec locX, locX_t, locA; 177c4762a1bSJed Brown PetscScalar *u, *u_t, *a; 178c4762a1bSJed Brown PetscMPIInt rank; 179c4762a1bSJed Brown PetscErrorCode ierr; 180c4762a1bSJed Brown 181c4762a1bSJed Brown PetscFunctionBegin; 182*ffc4695bSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRMPI(ierr); 183c4762a1bSJed Brown locX_t = NULL; 184c4762a1bSJed Brown locA = NULL; 185c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF, rank ? 0 : 1, 0, 1, &cells);CHKERRQ(ierr); 186c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr); 187c4762a1bSJed Brown ierr = DMPlexGetCellFields( dm, cells, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr); 188c4762a1bSJed Brown ierr = DMPlexRestoreCellFields(dm, cells, locX, locX_t, locA, &u, &u_t, &a);CHKERRQ(ierr); 189c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr); 190c4762a1bSJed Brown ierr = ISDestroy(&cells);CHKERRQ(ierr); 191c4762a1bSJed Brown PetscFunctionReturn(0); 192c4762a1bSJed Brown } 193c4762a1bSJed Brown 194c4762a1bSJed Brown static PetscErrorCode test7(DM dm, AppCtx *options) 195c4762a1bSJed Brown { 196c4762a1bSJed Brown PetscFE fe; 197c4762a1bSJed Brown PetscErrorCode ierr; 198c4762a1bSJed Brown 199c4762a1bSJed Brown PetscFunctionBegin; 200c4762a1bSJed Brown ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), options->dim, 1, options->simplex, NULL, -1, &fe);CHKERRQ(ierr); 201c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject)fe);CHKERRQ(ierr); 202c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 203c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 204c4762a1bSJed Brown ierr = test5(dm, options);CHKERRQ(ierr); 205c4762a1bSJed Brown PetscFunctionReturn(0); 206c4762a1bSJed Brown } 207c4762a1bSJed Brown 208c4762a1bSJed Brown static PetscErrorCode test8(DM dm, AppCtx *options) 209c4762a1bSJed Brown { 210c4762a1bSJed Brown PetscFE fe; 211c4762a1bSJed Brown PetscErrorCode ierr; 212c4762a1bSJed Brown 213c4762a1bSJed Brown PetscFunctionBegin; 214c4762a1bSJed Brown ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), options->dim, 1, options->simplex, NULL, -1, &fe);CHKERRQ(ierr); 215c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject)fe);CHKERRQ(ierr); 216c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 217c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 218c4762a1bSJed Brown ierr = test6(dm, options);CHKERRQ(ierr); 219c4762a1bSJed Brown PetscFunctionReturn(0); 220c4762a1bSJed Brown } 221c4762a1bSJed Brown 222c4762a1bSJed Brown int main(int argc, char **argv) 223c4762a1bSJed Brown { 224c4762a1bSJed Brown MPI_Comm comm; 225c4762a1bSJed Brown DM dm; 226c4762a1bSJed Brown AppCtx options; 227c4762a1bSJed Brown PetscErrorCode ierr; 228c4762a1bSJed Brown 229c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 230c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 231c4762a1bSJed Brown ierr = ProcessOptions(comm, &options);CHKERRQ(ierr); 232c4762a1bSJed Brown ierr = CreateMesh(comm, &options, &dm);CHKERRQ(ierr); 233c4762a1bSJed Brown 234c4762a1bSJed Brown switch (options.test) { 235c4762a1bSJed Brown case 0: ierr = test0(dm, &options);CHKERRQ(ierr); break; 236c4762a1bSJed Brown case 1: ierr = test1(dm, &options);CHKERRQ(ierr); break; 237c4762a1bSJed Brown case 2: ierr = test2(dm, &options);CHKERRQ(ierr); break; 238c4762a1bSJed Brown case 3: ierr = test3(dm, &options);CHKERRQ(ierr); break; 239c4762a1bSJed Brown case 4: ierr = test4(dm, &options);CHKERRQ(ierr); break; 240c4762a1bSJed Brown case 5: ierr = test5(dm, &options);CHKERRQ(ierr); break; 241c4762a1bSJed Brown case 6: ierr = test6(dm, &options);CHKERRQ(ierr); break; 242c4762a1bSJed Brown case 7: ierr = test7(dm, &options);CHKERRQ(ierr); break; 243c4762a1bSJed Brown case 8: ierr = test8(dm, &options);CHKERRQ(ierr); break; 244c4762a1bSJed Brown default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No such test: %D", options.test); 245c4762a1bSJed Brown } 246c4762a1bSJed Brown 247c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 248c4762a1bSJed Brown ierr = PetscFinalize(); 249c4762a1bSJed Brown return ierr; 250c4762a1bSJed Brown } 251c4762a1bSJed Brown 252c4762a1bSJed Brown /*TEST 253c4762a1bSJed Brown 254c4762a1bSJed Brown test: 255c4762a1bSJed Brown suffix: 0 256c4762a1bSJed Brown requires: ctetgen 257c4762a1bSJed Brown args: -test 0 258c4762a1bSJed Brown test: 259c4762a1bSJed Brown suffix: 1 260c4762a1bSJed Brown requires: ctetgen 261c4762a1bSJed Brown args: -test 1 262c4762a1bSJed Brown test: 263c4762a1bSJed Brown suffix: 2 264c4762a1bSJed Brown requires: ctetgen 265c4762a1bSJed Brown args: -test 2 266c4762a1bSJed Brown test: 267c4762a1bSJed Brown suffix: 3 268c4762a1bSJed Brown requires: ctetgen 269c4762a1bSJed Brown args: -test 3 270c4762a1bSJed Brown test: 271c4762a1bSJed Brown suffix: 4 272c4762a1bSJed Brown requires: ctetgen 273c4762a1bSJed Brown args: -test 4 274c4762a1bSJed Brown test: 275c4762a1bSJed Brown suffix: 5 276c4762a1bSJed Brown requires: ctetgen 277c4762a1bSJed Brown args: -test 5 278c4762a1bSJed Brown test: 279c4762a1bSJed Brown suffix: 6 280c4762a1bSJed Brown requires: ctetgen 281c4762a1bSJed Brown args: -test 6 282c4762a1bSJed Brown test: 283c4762a1bSJed Brown suffix: 7 284c4762a1bSJed Brown requires: ctetgen 285c4762a1bSJed Brown args: -test 7 286c4762a1bSJed Brown test: 287c4762a1bSJed Brown suffix: 8 288c4762a1bSJed Brown requires: ctetgen 289c4762a1bSJed Brown args: -test 8 290c4762a1bSJed Brown test: 291c4762a1bSJed Brown suffix: 9 292c4762a1bSJed Brown requires: ctetgen 293c4762a1bSJed Brown nsize: 2 294c4762a1bSJed Brown args: -test 1 295c4762a1bSJed Brown test: 296c4762a1bSJed Brown suffix: 10 297c4762a1bSJed Brown requires: ctetgen 298c4762a1bSJed Brown nsize: 2 299c4762a1bSJed Brown args: -test 2 300c4762a1bSJed Brown test: 301c4762a1bSJed Brown suffix: 11 302c4762a1bSJed Brown requires: ctetgen 303c4762a1bSJed Brown nsize: 2 304c4762a1bSJed Brown args: -test 3 305c4762a1bSJed Brown test: 306c4762a1bSJed Brown suffix: 12 307c4762a1bSJed Brown requires: ctetgen 308c4762a1bSJed Brown nsize: 2 309c4762a1bSJed Brown args: -test 4 310c4762a1bSJed Brown 311c4762a1bSJed Brown TEST*/ 312c4762a1bSJed Brown 313