1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "DMSwarm-PIC demonstator of inserting points into a cell DM \n\ 3c4762a1bSJed Brown Options: \n\ 4c4762a1bSJed Brown -mode {0,1} : 0 ==> DMDA, 1 ==> DMPLEX cell DM \n\ 5c4762a1bSJed Brown -dim {2,3} : spatial dimension\n"; 6c4762a1bSJed Brown 7c4762a1bSJed Brown #include <petsc.h> 8c4762a1bSJed Brown #include <petscdm.h> 9c4762a1bSJed Brown #include <petscdmda.h> 10c4762a1bSJed Brown #include <petscdmplex.h> 11c4762a1bSJed Brown #include <petscdmswarm.h> 12c4762a1bSJed Brown 13d71ae5a4SJacob Faibussowitsch PetscErrorCode pic_insert_DMDA(PetscInt dim) 14d71ae5a4SJacob Faibussowitsch { 15c4762a1bSJed Brown DM celldm = NULL, swarm; 16c4762a1bSJed Brown PetscInt dof, stencil_width; 17c4762a1bSJed Brown PetscReal min[3], max[3]; 18c4762a1bSJed Brown PetscInt ndir[3]; 19c4762a1bSJed Brown 20c4762a1bSJed Brown PetscFunctionBegin; 21c4762a1bSJed Brown /* Create the background cell DM */ 22c4762a1bSJed Brown dof = 1; 23c4762a1bSJed Brown stencil_width = 1; 2448a46eb9SPierre Jolivet if (dim == 2) PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 25, 13, PETSC_DECIDE, PETSC_DECIDE, dof, stencil_width, NULL, NULL, &celldm)); 2548a46eb9SPierre Jolivet if (dim == 3) PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 25, 13, 19, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, dof, stencil_width, NULL, NULL, NULL, &celldm)); 26c4762a1bSJed Brown 279566063dSJacob Faibussowitsch PetscCall(DMDASetElementType(celldm, DMDA_ELEMENT_Q1)); 289566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(celldm)); 299566063dSJacob Faibussowitsch PetscCall(DMSetUp(celldm)); 30c4762a1bSJed Brown 319566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(celldm, 0.0, 2.0, 0.0, 1.0, 0.0, 1.5)); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* Create the DMSwarm */ 349566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &swarm)); 359566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)swarm, "Swarm")); 369566063dSJacob Faibussowitsch PetscCall(DMSetType(swarm, DMSWARM)); 379566063dSJacob Faibussowitsch PetscCall(DMSetDimension(swarm, dim)); 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* Configure swarm to be of type PIC */ 409566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(swarm, DMSWARM_PIC)); 419566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(swarm, celldm)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Register two scalar fields within the DMSwarm */ 449566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "viscosity", 1, PETSC_DOUBLE)); 459566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "density", 1, PETSC_DOUBLE)); 469566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(swarm)); 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* Set initial local sizes of the DMSwarm with a buffer length of zero */ 499566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(swarm, 4, 0)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* Insert swarm coordinates cell-wise */ 529566063dSJacob Faibussowitsch PetscCall(DMSwarmInsertPointsUsingCellDM(swarm, DMSWARMPIC_LAYOUT_REGULAR, 3)); 539371c9d4SSatish Balay min[0] = 0.5; 549371c9d4SSatish Balay max[0] = 0.7; 559371c9d4SSatish Balay min[1] = 0.5; 569371c9d4SSatish Balay max[1] = 0.8; 579371c9d4SSatish Balay min[2] = 0.5; 589371c9d4SSatish Balay max[2] = 0.9; 59c4762a1bSJed Brown ndir[0] = ndir[1] = ndir[2] = 30; 609566063dSJacob Faibussowitsch PetscCall(DMSwarmSetPointsUniformCoordinates(swarm, min, max, ndir, ADD_VALUES)); 61c4762a1bSJed Brown 62c4762a1bSJed Brown /* This should be dispatched from a regular DMView() */ 639566063dSJacob Faibussowitsch PetscCall(DMSwarmViewXDMF(swarm, "ex20.xmf")); 649566063dSJacob Faibussowitsch PetscCall(DMView(celldm, PETSC_VIEWER_STDOUT_WORLD)); 659566063dSJacob Faibussowitsch PetscCall(DMView(swarm, PETSC_VIEWER_STDOUT_WORLD)); 66c4762a1bSJed Brown 67c4762a1bSJed Brown { 68c4762a1bSJed Brown PetscInt npoints, *list; 69c4762a1bSJed Brown PetscMPIInt rank; 70c4762a1bSJed Brown 719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 729566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetAccess(swarm)); 739566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetNumberOfPointsPerCell(swarm, 0, &npoints)); 749566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetPointsPerCell(swarm, rank, &npoints, &list)); 759566063dSJacob Faibussowitsch PetscCall(PetscFree(list)); 769566063dSJacob Faibussowitsch PetscCall(DMSwarmSortRestoreAccess(swarm)); 77c4762a1bSJed Brown } 789566063dSJacob Faibussowitsch PetscCall(DMSwarmMigrate(swarm, PETSC_FALSE)); 799566063dSJacob Faibussowitsch PetscCall(DMDestroy(&celldm)); 809566063dSJacob Faibussowitsch PetscCall(DMDestroy(&swarm)); 81*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 82c4762a1bSJed Brown } 83c4762a1bSJed Brown 84d71ae5a4SJacob Faibussowitsch PetscErrorCode pic_insert_DMPLEX_with_cell_list(PetscInt dim) 85d71ae5a4SJacob Faibussowitsch { 86c4762a1bSJed Brown DM celldm = NULL, swarm, distributedMesh = NULL; 87c4762a1bSJed Brown const char *fieldnames[] = {"viscosity"}; 88c4762a1bSJed Brown 89c4762a1bSJed Brown PetscFunctionBegin; 90c4762a1bSJed Brown /* Create the background cell DM */ 91c4762a1bSJed Brown if (dim == 2) { 92c4762a1bSJed Brown PetscInt cells_per_dim[2], nx[2]; 93c4762a1bSJed Brown PetscInt n_tricells; 94c4762a1bSJed Brown PetscInt n_trivert; 95a4a685f2SJacob Faibussowitsch PetscInt *tricells; 96a4a685f2SJacob Faibussowitsch PetscReal *trivert, dx, dy; 97c4762a1bSJed Brown PetscInt ii, jj, cnt; 98c4762a1bSJed Brown 99c4762a1bSJed Brown cells_per_dim[0] = 4; 100c4762a1bSJed Brown cells_per_dim[1] = 4; 101c4762a1bSJed Brown n_tricells = cells_per_dim[0] * cells_per_dim[1] * 2; 102c4762a1bSJed Brown nx[0] = cells_per_dim[0] + 1; 103c4762a1bSJed Brown nx[1] = cells_per_dim[1] + 1; 104c4762a1bSJed Brown n_trivert = nx[0] * nx[1]; 105c4762a1bSJed Brown 1069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_tricells * 3, &tricells)); 1079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nx[0] * nx[1] * 2, &trivert)); 108c4762a1bSJed Brown 109c4762a1bSJed Brown /* verts */ 110c4762a1bSJed Brown cnt = 0; 111c4762a1bSJed Brown dx = 2.0 / ((PetscReal)cells_per_dim[0]); 112c4762a1bSJed Brown dy = 1.0 / ((PetscReal)cells_per_dim[1]); 113c4762a1bSJed Brown for (jj = 0; jj < nx[1]; jj++) { 114c4762a1bSJed Brown for (ii = 0; ii < nx[0]; ii++) { 115c4762a1bSJed Brown trivert[2 * cnt + 0] = 0.0 + ii * dx; 116c4762a1bSJed Brown trivert[2 * cnt + 1] = 0.0 + jj * dy; 117c4762a1bSJed Brown cnt++; 118c4762a1bSJed Brown } 119c4762a1bSJed Brown } 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* connectivity */ 122c4762a1bSJed Brown cnt = 0; 123c4762a1bSJed Brown for (jj = 0; jj < cells_per_dim[1]; jj++) { 124c4762a1bSJed Brown for (ii = 0; ii < cells_per_dim[0]; ii++) { 125c4762a1bSJed Brown PetscInt idx, idx0, idx1, idx2, idx3; 126c4762a1bSJed Brown 127c4762a1bSJed Brown idx = (ii) + (jj)*nx[0]; 128c4762a1bSJed Brown idx0 = idx; 129c4762a1bSJed Brown idx1 = idx0 + 1; 130c4762a1bSJed Brown idx2 = idx1 + nx[0]; 131c4762a1bSJed Brown idx3 = idx0 + nx[0]; 132c4762a1bSJed Brown 133c4762a1bSJed Brown tricells[3 * cnt + 0] = idx0; 134c4762a1bSJed Brown tricells[3 * cnt + 1] = idx1; 135c4762a1bSJed Brown tricells[3 * cnt + 2] = idx2; 136c4762a1bSJed Brown cnt++; 137c4762a1bSJed Brown 138c4762a1bSJed Brown tricells[3 * cnt + 0] = idx0; 139c4762a1bSJed Brown tricells[3 * cnt + 1] = idx2; 140c4762a1bSJed Brown tricells[3 * cnt + 2] = idx3; 141c4762a1bSJed Brown cnt++; 142c4762a1bSJed Brown } 143c4762a1bSJed Brown } 1449566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, dim, n_tricells, n_trivert, 3, PETSC_TRUE, tricells, dim, trivert, &celldm)); 1459566063dSJacob Faibussowitsch PetscCall(PetscFree(trivert)); 1469566063dSJacob Faibussowitsch PetscCall(PetscFree(tricells)); 147c4762a1bSJed Brown } 14808401ef6SPierre Jolivet PetscCheck(dim != 3, PETSC_COMM_WORLD, PETSC_ERR_SUP, "Only 2D PLEX example supported"); 149c4762a1bSJed Brown 150c4762a1bSJed Brown /* Distribute mesh over processes */ 1519566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(celldm, 0, NULL, &distributedMesh)); 152c4762a1bSJed Brown if (distributedMesh) { 1539566063dSJacob Faibussowitsch PetscCall(DMDestroy(&celldm)); 154c4762a1bSJed Brown celldm = distributedMesh; 155c4762a1bSJed Brown } 1569566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)celldm, "Cells")); 1579566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(celldm)); 158c4762a1bSJed Brown { 159c4762a1bSJed Brown PetscInt numComp[] = {1}; 160c4762a1bSJed Brown PetscInt numDof[] = {1, 0, 0}; /* vert, edge, cell */ 161c4762a1bSJed Brown PetscInt numBC = 0; 162c4762a1bSJed Brown PetscSection section; 163c4762a1bSJed Brown 1649566063dSJacob Faibussowitsch PetscCall(DMPlexCreateSection(celldm, NULL, numComp, numDof, numBC, NULL, NULL, NULL, NULL, §ion)); 1659566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(celldm, section)); 1669566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ion)); 167c4762a1bSJed Brown } 1689566063dSJacob Faibussowitsch PetscCall(DMSetUp(celldm)); 169c4762a1bSJed Brown { 170c4762a1bSJed Brown PetscViewer viewer; 171c4762a1bSJed Brown 1729566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer)); 1739566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer, PETSCVIEWERVTK)); 1749566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_WRITE)); 1759566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(viewer, "ex20plex.vtk")); 1769566063dSJacob Faibussowitsch PetscCall(DMView(celldm, viewer)); 1779566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 178c4762a1bSJed Brown } 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* Create the DMSwarm */ 1819566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &swarm)); 1829566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)swarm, "Swarm")); 1839566063dSJacob Faibussowitsch PetscCall(DMSetType(swarm, DMSWARM)); 1849566063dSJacob Faibussowitsch PetscCall(DMSetDimension(swarm, dim)); 185c4762a1bSJed Brown 1869566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(swarm, DMSWARM_PIC)); 1879566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(swarm, celldm)); 188c4762a1bSJed Brown 189c4762a1bSJed Brown /* Register two scalar fields within the DMSwarm */ 1909566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "viscosity", 1, PETSC_DOUBLE)); 1919566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "density", 1, PETSC_DOUBLE)); 1929566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(swarm)); 193c4762a1bSJed Brown 194c4762a1bSJed Brown /* Set initial local sizes of the DMSwarm with a buffer length of zero */ 1959566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(swarm, 4, 0)); 196c4762a1bSJed Brown 197c4762a1bSJed Brown /* Insert swarm coordinates cell-wise */ 1989566063dSJacob Faibussowitsch PetscCall(DMSwarmInsertPointsUsingCellDM(swarm, DMSWARMPIC_LAYOUT_SUBDIVISION, 2)); 1999566063dSJacob Faibussowitsch PetscCall(DMSwarmViewFieldsXDMF(swarm, "ex20.xmf", 1, fieldnames)); 2009566063dSJacob Faibussowitsch PetscCall(DMView(celldm, PETSC_VIEWER_STDOUT_WORLD)); 2019566063dSJacob Faibussowitsch PetscCall(DMView(swarm, PETSC_VIEWER_STDOUT_WORLD)); 2029566063dSJacob Faibussowitsch PetscCall(DMDestroy(&celldm)); 2039566063dSJacob Faibussowitsch PetscCall(DMDestroy(&swarm)); 204*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 205c4762a1bSJed Brown } 206c4762a1bSJed Brown 207d71ae5a4SJacob Faibussowitsch PetscErrorCode pic_insert_DMPLEX(PetscBool is_simplex, PetscInt dim) 208d71ae5a4SJacob Faibussowitsch { 209c4762a1bSJed Brown DM celldm, swarm, distributedMesh = NULL; 210c4762a1bSJed Brown const char *fieldnames[] = {"viscosity", "DMSwarm_rank"}; 211c4762a1bSJed Brown 212c4762a1bSJed Brown PetscFunctionBegin; 213c4762a1bSJed Brown 214c4762a1bSJed Brown /* Create the background cell DM */ 215c4762a1bSJed Brown { 216c4762a1bSJed Brown PetscInt faces[3] = {4, 2, 4}; 2179566063dSJacob Faibussowitsch PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, is_simplex, faces, NULL, NULL, NULL, PETSC_TRUE, &celldm)); 218c4762a1bSJed Brown } 219c4762a1bSJed Brown 220c4762a1bSJed Brown /* Distribute mesh over processes */ 2219566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(celldm, 0, NULL, &distributedMesh)); 222c4762a1bSJed Brown if (distributedMesh) { 2239566063dSJacob Faibussowitsch PetscCall(DMDestroy(&celldm)); 224c4762a1bSJed Brown celldm = distributedMesh; 225c4762a1bSJed Brown } 2269566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)celldm, "Cells")); 2279566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(celldm)); 228c4762a1bSJed Brown { 229c4762a1bSJed Brown PetscInt numComp[] = {1}; 230c4762a1bSJed Brown PetscInt numDof[] = {1, 0, 0}; /* vert, edge, cell */ 231c4762a1bSJed Brown PetscInt numBC = 0; 232c4762a1bSJed Brown PetscSection section; 233c4762a1bSJed Brown 2349566063dSJacob Faibussowitsch PetscCall(DMPlexCreateSection(celldm, NULL, numComp, numDof, numBC, NULL, NULL, NULL, NULL, §ion)); 2359566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(celldm, section)); 2369566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ion)); 237c4762a1bSJed Brown } 2389566063dSJacob Faibussowitsch PetscCall(DMSetUp(celldm)); 239c4762a1bSJed Brown { 240c4762a1bSJed Brown PetscViewer viewer; 241c4762a1bSJed Brown 2429566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer)); 2439566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer, PETSCVIEWERVTK)); 2449566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_WRITE)); 2459566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(viewer, "ex20plex.vtk")); 2469566063dSJacob Faibussowitsch PetscCall(DMView(celldm, viewer)); 2479566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 248c4762a1bSJed Brown } 249c4762a1bSJed Brown 2509566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &swarm)); 2519566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)swarm, "Swarm")); 2529566063dSJacob Faibussowitsch PetscCall(DMSetType(swarm, DMSWARM)); 2539566063dSJacob Faibussowitsch PetscCall(DMSetDimension(swarm, dim)); 254c4762a1bSJed Brown 2559566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(swarm, DMSWARM_PIC)); 2569566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(swarm, celldm)); 257c4762a1bSJed Brown 258c4762a1bSJed Brown /* Register two scalar fields within the DMSwarm */ 2599566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "viscosity", 1, PETSC_DOUBLE)); 2609566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "density", 1, PETSC_DOUBLE)); 2619566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(swarm)); 262c4762a1bSJed Brown 263c4762a1bSJed Brown /* Set initial local sizes of the DMSwarm with a buffer length of zero */ 2649566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(swarm, 4, 0)); 265c4762a1bSJed Brown 266c4762a1bSJed Brown /* Insert swarm coordinates cell-wise */ 2679566063dSJacob Faibussowitsch PetscCall(DMSwarmInsertPointsUsingCellDM(swarm, DMSWARMPIC_LAYOUT_GAUSS, 3)); 2689566063dSJacob Faibussowitsch PetscCall(DMSwarmViewFieldsXDMF(swarm, "ex20.xmf", 2, fieldnames)); 2699566063dSJacob Faibussowitsch PetscCall(DMView(celldm, PETSC_VIEWER_STDOUT_WORLD)); 2709566063dSJacob Faibussowitsch PetscCall(DMView(swarm, PETSC_VIEWER_STDOUT_WORLD)); 2719566063dSJacob Faibussowitsch PetscCall(DMDestroy(&celldm)); 2729566063dSJacob Faibussowitsch PetscCall(DMDestroy(&swarm)); 273*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 274c4762a1bSJed Brown } 275c4762a1bSJed Brown 276d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 277d71ae5a4SJacob Faibussowitsch { 278c4762a1bSJed Brown PetscInt mode = 0; 279c4762a1bSJed Brown PetscInt dim = 2; 280c4762a1bSJed Brown 281327415f7SBarry Smith PetscFunctionBeginUser; 2829566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 2839566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mode", &mode, NULL)); 2849566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL)); 285c4762a1bSJed Brown switch (mode) { 286d71ae5a4SJacob Faibussowitsch case 0: 287d71ae5a4SJacob Faibussowitsch PetscCall(pic_insert_DMDA(dim)); 288d71ae5a4SJacob Faibussowitsch break; 289c4762a1bSJed Brown case 1: 290c4762a1bSJed Brown /* tri / tet */ 2919566063dSJacob Faibussowitsch PetscCall(pic_insert_DMPLEX(PETSC_TRUE, dim)); 292c4762a1bSJed Brown break; 293c4762a1bSJed Brown case 2: 294c4762a1bSJed Brown /* quad / hex */ 2959566063dSJacob Faibussowitsch PetscCall(pic_insert_DMPLEX(PETSC_FALSE, dim)); 296c4762a1bSJed Brown break; 297d71ae5a4SJacob Faibussowitsch default: 298d71ae5a4SJacob Faibussowitsch PetscCall(pic_insert_DMDA(dim)); 299d71ae5a4SJacob Faibussowitsch break; 300c4762a1bSJed Brown } 3019566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 302b122ec5aSJacob Faibussowitsch return 0; 303c4762a1bSJed Brown } 304c4762a1bSJed Brown 305c4762a1bSJed Brown /*TEST 306c4762a1bSJed Brown 307c4762a1bSJed Brown test: 308c4762a1bSJed Brown args: 309c4762a1bSJed Brown requires: !complex double 3107ed4f029SJed Brown filter: grep -v atomic 311c4762a1bSJed Brown filter_output: grep -v atomic 312c4762a1bSJed Brown 313c4762a1bSJed Brown test: 314c4762a1bSJed Brown suffix: 2 315c4762a1bSJed Brown requires: triangle double !complex 316c4762a1bSJed Brown args: -mode 1 3177ed4f029SJed Brown filter: grep -v atomic 318c4762a1bSJed Brown filter_output: grep -v atomic 319c4762a1bSJed Brown 320c4762a1bSJed Brown TEST*/ 321