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 13c4762a1bSJed Brown PetscErrorCode pic_insert_DMDA(PetscInt dim) 14c4762a1bSJed Brown { 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; 24c4762a1bSJed Brown if (dim == 2) { 25*9566063dSJacob Faibussowitsch 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)); 26c4762a1bSJed Brown } 27c4762a1bSJed Brown if (dim == 3) { 28*9566063dSJacob Faibussowitsch 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)); 29c4762a1bSJed Brown } 30c4762a1bSJed Brown 31*9566063dSJacob Faibussowitsch PetscCall(DMDASetElementType(celldm,DMDA_ELEMENT_Q1)); 32*9566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(celldm)); 33*9566063dSJacob Faibussowitsch PetscCall(DMSetUp(celldm)); 34c4762a1bSJed Brown 35*9566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(celldm,0.0,2.0,0.0,1.0,0.0,1.5)); 36c4762a1bSJed Brown 37c4762a1bSJed Brown /* Create the DMSwarm */ 38*9566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD,&swarm)); 39*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)swarm,"Swarm")); 40*9566063dSJacob Faibussowitsch PetscCall(DMSetType(swarm,DMSWARM)); 41*9566063dSJacob Faibussowitsch PetscCall(DMSetDimension(swarm,dim)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Configure swarm to be of type PIC */ 44*9566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(swarm,DMSWARM_PIC)); 45*9566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(swarm,celldm)); 46c4762a1bSJed Brown 47c4762a1bSJed Brown /* Register two scalar fields within the DMSwarm */ 48*9566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(swarm,"viscosity",1,PETSC_DOUBLE)); 49*9566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(swarm,"density",1,PETSC_DOUBLE)); 50*9566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(swarm)); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Set initial local sizes of the DMSwarm with a buffer length of zero */ 53*9566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(swarm,4,0)); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Insert swarm coordinates cell-wise */ 56*9566063dSJacob Faibussowitsch PetscCall(DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_REGULAR,3)); 57c4762a1bSJed Brown min[0] = 0.5; max[0] = 0.7; 58c4762a1bSJed Brown min[1] = 0.5; max[1] = 0.8; 59c4762a1bSJed Brown min[2] = 0.5; max[2] = 0.9; 60c4762a1bSJed Brown ndir[0] = ndir[1] = ndir[2] = 30; 61*9566063dSJacob Faibussowitsch PetscCall(DMSwarmSetPointsUniformCoordinates(swarm,min,max,ndir,ADD_VALUES)); 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* This should be dispatched from a regular DMView() */ 64*9566063dSJacob Faibussowitsch PetscCall(DMSwarmViewXDMF(swarm,"ex20.xmf")); 65*9566063dSJacob Faibussowitsch PetscCall(DMView(celldm,PETSC_VIEWER_STDOUT_WORLD)); 66*9566063dSJacob Faibussowitsch PetscCall(DMView(swarm,PETSC_VIEWER_STDOUT_WORLD)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown { 69c4762a1bSJed Brown PetscInt npoints,*list; 70c4762a1bSJed Brown PetscMPIInt rank; 71c4762a1bSJed Brown 72*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 73*9566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetAccess(swarm)); 74*9566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetNumberOfPointsPerCell(swarm,0,&npoints)); 75*9566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetPointsPerCell(swarm,rank,&npoints,&list)); 76*9566063dSJacob Faibussowitsch PetscCall(PetscFree(list)); 77*9566063dSJacob Faibussowitsch PetscCall(DMSwarmSortRestoreAccess(swarm)); 78c4762a1bSJed Brown } 79*9566063dSJacob Faibussowitsch PetscCall(DMSwarmMigrate(swarm,PETSC_FALSE)); 80*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&celldm)); 81*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&swarm)); 82c4762a1bSJed Brown PetscFunctionReturn(0); 83c4762a1bSJed Brown } 84c4762a1bSJed Brown 85c4762a1bSJed Brown PetscErrorCode pic_insert_DMPLEX_with_cell_list(PetscInt dim) 86c4762a1bSJed Brown { 87c4762a1bSJed Brown DM celldm = NULL,swarm,distributedMesh = NULL; 88c4762a1bSJed Brown const char *fieldnames[] = {"viscosity"}; 89c4762a1bSJed Brown 90c4762a1bSJed Brown PetscFunctionBegin; 91c4762a1bSJed Brown /* Create the background cell DM */ 92c4762a1bSJed Brown if (dim == 2) { 93c4762a1bSJed Brown PetscInt cells_per_dim[2],nx[2]; 94c4762a1bSJed Brown PetscInt n_tricells; 95c4762a1bSJed Brown PetscInt n_trivert; 96a4a685f2SJacob Faibussowitsch PetscInt *tricells; 97a4a685f2SJacob Faibussowitsch PetscReal *trivert,dx,dy; 98c4762a1bSJed Brown PetscInt ii,jj,cnt; 99c4762a1bSJed Brown 100c4762a1bSJed Brown cells_per_dim[0] = 4; 101c4762a1bSJed Brown cells_per_dim[1] = 4; 102c4762a1bSJed Brown n_tricells = cells_per_dim[0] * cells_per_dim[1] * 2; 103c4762a1bSJed Brown nx[0] = cells_per_dim[0] + 1; 104c4762a1bSJed Brown nx[1] = cells_per_dim[1] + 1; 105c4762a1bSJed Brown n_trivert = nx[0] * nx[1]; 106c4762a1bSJed Brown 107*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_tricells*3,&tricells)); 108*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nx[0]*nx[1]*2,&trivert)); 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* verts */ 111c4762a1bSJed Brown cnt = 0; 112c4762a1bSJed Brown dx = 2.0/((PetscReal)cells_per_dim[0]); 113c4762a1bSJed Brown dy = 1.0/((PetscReal)cells_per_dim[1]); 114c4762a1bSJed Brown for (jj=0; jj<nx[1]; jj++) { 115c4762a1bSJed Brown for (ii=0; ii<nx[0]; ii++) { 116c4762a1bSJed Brown trivert[2*cnt+0] = 0.0 + ii * dx; 117c4762a1bSJed Brown trivert[2*cnt+1] = 0.0 + jj * dy; 118c4762a1bSJed Brown cnt++; 119c4762a1bSJed Brown } 120c4762a1bSJed Brown } 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* connectivity */ 123c4762a1bSJed Brown cnt = 0; 124c4762a1bSJed Brown for (jj=0; jj<cells_per_dim[1]; jj++) { 125c4762a1bSJed Brown for (ii=0; ii<cells_per_dim[0]; ii++) { 126c4762a1bSJed Brown PetscInt idx,idx0,idx1,idx2,idx3; 127c4762a1bSJed Brown 128c4762a1bSJed Brown idx = (ii) + (jj) * nx[0]; 129c4762a1bSJed Brown idx0 = idx; 130c4762a1bSJed Brown idx1 = idx0 + 1; 131c4762a1bSJed Brown idx2 = idx1 + nx[0]; 132c4762a1bSJed Brown idx3 = idx0 + nx[0]; 133c4762a1bSJed Brown 134c4762a1bSJed Brown tricells[3*cnt+0] = idx0; 135c4762a1bSJed Brown tricells[3*cnt+1] = idx1; 136c4762a1bSJed Brown tricells[3*cnt+2] = idx2; 137c4762a1bSJed Brown cnt++; 138c4762a1bSJed Brown 139c4762a1bSJed Brown tricells[3*cnt+0] = idx0; 140c4762a1bSJed Brown tricells[3*cnt+1] = idx2; 141c4762a1bSJed Brown tricells[3*cnt+2] = idx3; 142c4762a1bSJed Brown cnt++; 143c4762a1bSJed Brown } 144c4762a1bSJed Brown } 145*9566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD,dim,n_tricells,n_trivert,3,PETSC_TRUE,tricells,dim,trivert,&celldm)); 146*9566063dSJacob Faibussowitsch PetscCall(PetscFree(trivert)); 147*9566063dSJacob Faibussowitsch PetscCall(PetscFree(tricells)); 148c4762a1bSJed Brown } 1492c71b3e2SJacob Faibussowitsch PetscCheckFalse(dim == 3,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only 2D PLEX example supported"); 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* Distribute mesh over processes */ 152*9566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(celldm,0,NULL,&distributedMesh)); 153c4762a1bSJed Brown if (distributedMesh) { 154*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&celldm)); 155c4762a1bSJed Brown celldm = distributedMesh; 156c4762a1bSJed Brown } 157*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)celldm,"Cells")); 158*9566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(celldm)); 159c4762a1bSJed Brown { 160c4762a1bSJed Brown PetscInt numComp[] = {1}; 161c4762a1bSJed Brown PetscInt numDof[] = {1,0,0}; /* vert, edge, cell */ 162c4762a1bSJed Brown PetscInt numBC = 0; 163c4762a1bSJed Brown PetscSection section; 164c4762a1bSJed Brown 165*9566063dSJacob Faibussowitsch PetscCall(DMPlexCreateSection(celldm,NULL,numComp,numDof,numBC,NULL,NULL,NULL,NULL,§ion)); 166*9566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(celldm,section)); 167*9566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ion)); 168c4762a1bSJed Brown } 169*9566063dSJacob Faibussowitsch PetscCall(DMSetUp(celldm)); 170c4762a1bSJed Brown { 171c4762a1bSJed Brown PetscViewer viewer; 172c4762a1bSJed Brown 173*9566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PETSC_COMM_WORLD,&viewer)); 174*9566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer,PETSCVIEWERVTK)); 175*9566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE)); 176*9566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(viewer,"ex20plex.vtk")); 177*9566063dSJacob Faibussowitsch PetscCall(DMView(celldm,viewer)); 178*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 179c4762a1bSJed Brown } 180c4762a1bSJed Brown 181c4762a1bSJed Brown /* Create the DMSwarm */ 182*9566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD,&swarm)); 183*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)swarm,"Swarm")); 184*9566063dSJacob Faibussowitsch PetscCall(DMSetType(swarm,DMSWARM)); 185*9566063dSJacob Faibussowitsch PetscCall(DMSetDimension(swarm,dim)); 186c4762a1bSJed Brown 187*9566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(swarm,DMSWARM_PIC)); 188*9566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(swarm,celldm)); 189c4762a1bSJed Brown 190c4762a1bSJed Brown /* Register two scalar fields within the DMSwarm */ 191*9566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(swarm,"viscosity",1,PETSC_DOUBLE)); 192*9566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(swarm,"density",1,PETSC_DOUBLE)); 193*9566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(swarm)); 194c4762a1bSJed Brown 195c4762a1bSJed Brown /* Set initial local sizes of the DMSwarm with a buffer length of zero */ 196*9566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(swarm,4,0)); 197c4762a1bSJed Brown 198c4762a1bSJed Brown /* Insert swarm coordinates cell-wise */ 199*9566063dSJacob Faibussowitsch PetscCall(DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_SUBDIVISION,2)); 200*9566063dSJacob Faibussowitsch PetscCall(DMSwarmViewFieldsXDMF(swarm,"ex20.xmf",1,fieldnames)); 201*9566063dSJacob Faibussowitsch PetscCall(DMView(celldm,PETSC_VIEWER_STDOUT_WORLD)); 202*9566063dSJacob Faibussowitsch PetscCall(DMView(swarm,PETSC_VIEWER_STDOUT_WORLD)); 203*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&celldm)); 204*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&swarm)); 205c4762a1bSJed Brown PetscFunctionReturn(0); 206c4762a1bSJed Brown } 207c4762a1bSJed Brown 208c4762a1bSJed Brown PetscErrorCode pic_insert_DMPLEX(PetscBool is_simplex,PetscInt dim) 209c4762a1bSJed Brown { 210c4762a1bSJed Brown DM celldm,swarm,distributedMesh = NULL; 211c4762a1bSJed Brown const char *fieldnames[] = {"viscosity","DMSwarm_rank"}; 212c4762a1bSJed Brown 213c4762a1bSJed Brown PetscFunctionBegin; 214c4762a1bSJed Brown 215c4762a1bSJed Brown /* Create the background cell DM */ 216c4762a1bSJed Brown { 217c4762a1bSJed Brown PetscInt faces[3] = {4, 2, 4}; 218*9566063dSJacob Faibussowitsch PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, is_simplex, faces, NULL, NULL, NULL, PETSC_TRUE, &celldm)); 219c4762a1bSJed Brown } 220c4762a1bSJed Brown 221c4762a1bSJed Brown /* Distribute mesh over processes */ 222*9566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(celldm,0,NULL,&distributedMesh)); 223c4762a1bSJed Brown if (distributedMesh) { 224*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&celldm)); 225c4762a1bSJed Brown celldm = distributedMesh; 226c4762a1bSJed Brown } 227*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)celldm,"Cells")); 228*9566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(celldm)); 229c4762a1bSJed Brown { 230c4762a1bSJed Brown PetscInt numComp[] = {1}; 231c4762a1bSJed Brown PetscInt numDof[] = {1,0,0}; /* vert, edge, cell */ 232c4762a1bSJed Brown PetscInt numBC = 0; 233c4762a1bSJed Brown PetscSection section; 234c4762a1bSJed Brown 235*9566063dSJacob Faibussowitsch PetscCall(DMPlexCreateSection(celldm,NULL,numComp,numDof,numBC,NULL,NULL,NULL,NULL,§ion)); 236*9566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(celldm,section)); 237*9566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ion)); 238c4762a1bSJed Brown } 239*9566063dSJacob Faibussowitsch PetscCall(DMSetUp(celldm)); 240c4762a1bSJed Brown { 241c4762a1bSJed Brown PetscViewer viewer; 242c4762a1bSJed Brown 243*9566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PETSC_COMM_WORLD,&viewer)); 244*9566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer,PETSCVIEWERVTK)); 245*9566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE)); 246*9566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(viewer,"ex20plex.vtk")); 247*9566063dSJacob Faibussowitsch PetscCall(DMView(celldm,viewer)); 248*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 249c4762a1bSJed Brown } 250c4762a1bSJed Brown 251*9566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD,&swarm)); 252*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)swarm,"Swarm")); 253*9566063dSJacob Faibussowitsch PetscCall(DMSetType(swarm,DMSWARM)); 254*9566063dSJacob Faibussowitsch PetscCall(DMSetDimension(swarm,dim)); 255c4762a1bSJed Brown 256*9566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(swarm,DMSWARM_PIC)); 257*9566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(swarm,celldm)); 258c4762a1bSJed Brown 259c4762a1bSJed Brown /* Register two scalar fields within the DMSwarm */ 260*9566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(swarm,"viscosity",1,PETSC_DOUBLE)); 261*9566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(swarm,"density",1,PETSC_DOUBLE)); 262*9566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(swarm)); 263c4762a1bSJed Brown 264c4762a1bSJed Brown /* Set initial local sizes of the DMSwarm with a buffer length of zero */ 265*9566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(swarm,4,0)); 266c4762a1bSJed Brown 267c4762a1bSJed Brown /* Insert swarm coordinates cell-wise */ 268*9566063dSJacob Faibussowitsch PetscCall(DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_GAUSS,3)); 269*9566063dSJacob Faibussowitsch PetscCall(DMSwarmViewFieldsXDMF(swarm,"ex20.xmf",2,fieldnames)); 270*9566063dSJacob Faibussowitsch PetscCall(DMView(celldm,PETSC_VIEWER_STDOUT_WORLD)); 271*9566063dSJacob Faibussowitsch PetscCall(DMView(swarm,PETSC_VIEWER_STDOUT_WORLD)); 272*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&celldm)); 273*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&swarm)); 274c4762a1bSJed Brown PetscFunctionReturn(0); 275c4762a1bSJed Brown } 276c4762a1bSJed Brown 277c4762a1bSJed Brown int main(int argc,char **args) 278c4762a1bSJed Brown { 279c4762a1bSJed Brown PetscInt mode = 0; 280c4762a1bSJed Brown PetscInt dim = 2; 281c4762a1bSJed Brown 282*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 283*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mode",&mode,NULL)); 284*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL)); 285c4762a1bSJed Brown switch (mode) { 286c4762a1bSJed Brown case 0: 287*9566063dSJacob Faibussowitsch PetscCall(pic_insert_DMDA(dim)); 288c4762a1bSJed Brown break; 289c4762a1bSJed Brown case 1: 290c4762a1bSJed Brown /* tri / tet */ 291*9566063dSJacob Faibussowitsch PetscCall(pic_insert_DMPLEX(PETSC_TRUE,dim)); 292c4762a1bSJed Brown break; 293c4762a1bSJed Brown case 2: 294c4762a1bSJed Brown /* quad / hex */ 295*9566063dSJacob Faibussowitsch PetscCall(pic_insert_DMPLEX(PETSC_FALSE,dim)); 296c4762a1bSJed Brown break; 297c4762a1bSJed Brown default: 298*9566063dSJacob Faibussowitsch PetscCall(pic_insert_DMDA(dim)); 299c4762a1bSJed Brown break; 300c4762a1bSJed Brown } 301*9566063dSJacob 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