1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests DMSwarm\n\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown #include <petscdm.h> 5*c4762a1bSJed Brown #include <petscdmda.h> 6*c4762a1bSJed Brown #include <petscdmswarm.h> 7*c4762a1bSJed Brown 8*c4762a1bSJed Brown PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM,PetscErrorCode (*)(DM,void*,PetscInt*,PetscInt**),size_t,void*,PetscInt*); 9*c4762a1bSJed Brown PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM,PetscInt*); 10*c4762a1bSJed Brown 11*c4762a1bSJed Brown PetscErrorCode ex1_1(void) 12*c4762a1bSJed Brown { 13*c4762a1bSJed Brown DM dms; 14*c4762a1bSJed Brown PetscErrorCode ierr; 15*c4762a1bSJed Brown Vec x; 16*c4762a1bSJed Brown PetscMPIInt rank,size; 17*c4762a1bSJed Brown PetscInt p; 18*c4762a1bSJed Brown 19*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 20*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 21*c4762a1bSJed Brown if ((size > 1) && (size != 4)) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Must be run wuth 4 MPI ranks"); 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown ierr = DMCreate(PETSC_COMM_WORLD,&dms);CHKERRQ(ierr); 24*c4762a1bSJed Brown ierr = DMSetType(dms,DMSWARM);CHKERRQ(ierr); 25*c4762a1bSJed Brown 26*c4762a1bSJed Brown ierr = DMSwarmInitializeFieldRegister(dms);CHKERRQ(ierr); 27*c4762a1bSJed Brown ierr = DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);CHKERRQ(ierr); 28*c4762a1bSJed Brown ierr = DMSwarmRegisterPetscDatatypeField(dms,"strain",1,PETSC_REAL);CHKERRQ(ierr); 29*c4762a1bSJed Brown ierr = DMSwarmFinalizeFieldRegister(dms);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = DMSwarmSetLocalSizes(dms,5+rank,4);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = DMView(dms,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown { 34*c4762a1bSJed Brown PetscReal *array; 35*c4762a1bSJed Brown ierr = DMSwarmGetField(dms,"viscosity",NULL,NULL,(void**)&array);CHKERRQ(ierr); 36*c4762a1bSJed Brown for (p=0; p<5+rank; p++) { 37*c4762a1bSJed Brown array[p] = 11.1 + p*0.1 + rank*100.0; 38*c4762a1bSJed Brown } 39*c4762a1bSJed Brown ierr = DMSwarmRestoreField(dms,"viscosity",NULL,NULL,(void**)&array);CHKERRQ(ierr); 40*c4762a1bSJed Brown } 41*c4762a1bSJed Brown 42*c4762a1bSJed Brown { 43*c4762a1bSJed Brown PetscReal *array; 44*c4762a1bSJed Brown ierr = DMSwarmGetField(dms,"strain",NULL,NULL,(void**)&array);CHKERRQ(ierr); 45*c4762a1bSJed Brown for (p=0; p<5+rank; p++) { 46*c4762a1bSJed Brown array[p] = 2.0e-2 + p*0.001 + rank*1.0; 47*c4762a1bSJed Brown } 48*c4762a1bSJed Brown ierr = DMSwarmRestoreField(dms,"strain",NULL,NULL,(void**)&array);CHKERRQ(ierr); 49*c4762a1bSJed Brown } 50*c4762a1bSJed Brown 51*c4762a1bSJed Brown ierr = DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr); 52*c4762a1bSJed Brown ierr = DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr); 53*c4762a1bSJed Brown 54*c4762a1bSJed Brown ierr = DMSwarmVectorDefineField(dms,"strain");CHKERRQ(ierr); 55*c4762a1bSJed Brown ierr = DMCreateGlobalVector(dms,&x);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 57*c4762a1bSJed Brown 58*c4762a1bSJed Brown { 59*c4762a1bSJed Brown PetscInt *rankval; 60*c4762a1bSJed Brown PetscInt npoints[2],npoints_orig[2]; 61*c4762a1bSJed Brown 62*c4762a1bSJed Brown ierr = DMSwarmGetLocalSize(dms,&npoints_orig[0]);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = DMSwarmGetSize(dms,&npoints_orig[1]);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = DMSwarmGetField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 65*c4762a1bSJed Brown if ((rank == 0) && (size > 1)) { 66*c4762a1bSJed Brown rankval[0] = 1; 67*c4762a1bSJed Brown rankval[3] = 1; 68*c4762a1bSJed Brown } 69*c4762a1bSJed Brown if (rank == 3) { 70*c4762a1bSJed Brown rankval[2] = 1; 71*c4762a1bSJed Brown } 72*c4762a1bSJed Brown ierr = DMSwarmRestoreField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 73*c4762a1bSJed Brown ierr = DMSwarmMigrate(dms,PETSC_TRUE);CHKERRQ(ierr); 74*c4762a1bSJed Brown ierr = DMSwarmGetLocalSize(dms,&npoints[0]);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = DMSwarmGetSize(dms,&npoints[1]);CHKERRQ(ierr); 76*c4762a1bSJed Brown ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 77*c4762a1bSJed Brown ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] before(%D,%D) after(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1],npoints[0],npoints[1]);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 80*c4762a1bSJed Brown } 81*c4762a1bSJed Brown { 82*c4762a1bSJed Brown ierr = DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr); 83*c4762a1bSJed Brown ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 84*c4762a1bSJed Brown ierr = DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr); 85*c4762a1bSJed Brown } 86*c4762a1bSJed Brown { 87*c4762a1bSJed Brown ierr = DMSwarmCreateGlobalVectorFromField(dms,"strain",&x);CHKERRQ(ierr); 88*c4762a1bSJed Brown ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = DMSwarmDestroyGlobalVectorFromField(dms,"strain",&x);CHKERRQ(ierr); 90*c4762a1bSJed Brown } 91*c4762a1bSJed Brown 92*c4762a1bSJed Brown ierr = DMDestroy(&dms);CHKERRQ(ierr); 93*c4762a1bSJed Brown PetscFunctionReturn(0); 94*c4762a1bSJed Brown } 95*c4762a1bSJed Brown 96*c4762a1bSJed Brown PetscErrorCode ex1_2(void) 97*c4762a1bSJed Brown { 98*c4762a1bSJed Brown DM dms; 99*c4762a1bSJed Brown PetscErrorCode ierr; 100*c4762a1bSJed Brown Vec x; 101*c4762a1bSJed Brown PetscMPIInt rank,size; 102*c4762a1bSJed Brown PetscInt p; 103*c4762a1bSJed Brown 104*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 105*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 106*c4762a1bSJed Brown ierr = DMCreate(PETSC_COMM_WORLD,&dms);CHKERRQ(ierr); 107*c4762a1bSJed Brown ierr = DMSetType(dms,DMSWARM);CHKERRQ(ierr); 108*c4762a1bSJed Brown ierr = DMSwarmInitializeFieldRegister(dms);CHKERRQ(ierr); 109*c4762a1bSJed Brown 110*c4762a1bSJed Brown ierr = DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);CHKERRQ(ierr); 111*c4762a1bSJed Brown ierr = DMSwarmRegisterPetscDatatypeField(dms,"strain",1,PETSC_REAL);CHKERRQ(ierr); 112*c4762a1bSJed Brown ierr = DMSwarmFinalizeFieldRegister(dms);CHKERRQ(ierr); 113*c4762a1bSJed Brown ierr = DMSwarmSetLocalSizes(dms,5+rank,4);CHKERRQ(ierr); 114*c4762a1bSJed Brown ierr = DMView(dms,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 115*c4762a1bSJed Brown { 116*c4762a1bSJed Brown PetscReal *array; 117*c4762a1bSJed Brown ierr = DMSwarmGetField(dms,"viscosity",NULL,NULL,(void**)&array);CHKERRQ(ierr); 118*c4762a1bSJed Brown for (p=0; p<5+rank; p++) { 119*c4762a1bSJed Brown array[p] = 11.1 + p*0.1 + rank*100.0; 120*c4762a1bSJed Brown } 121*c4762a1bSJed Brown ierr = DMSwarmRestoreField(dms,"viscosity",NULL,NULL,(void**)&array);CHKERRQ(ierr); 122*c4762a1bSJed Brown } 123*c4762a1bSJed Brown { 124*c4762a1bSJed Brown PetscReal *array; 125*c4762a1bSJed Brown ierr = DMSwarmGetField(dms,"strain",NULL,NULL,(void**)&array);CHKERRQ(ierr); 126*c4762a1bSJed Brown for (p=0; p<5+rank; p++) { 127*c4762a1bSJed Brown array[p] = 2.0e-2 + p*0.001 + rank*1.0; 128*c4762a1bSJed Brown } 129*c4762a1bSJed Brown ierr = DMSwarmRestoreField(dms,"strain",NULL,NULL,(void**)&array);CHKERRQ(ierr); 130*c4762a1bSJed Brown } 131*c4762a1bSJed Brown { 132*c4762a1bSJed Brown PetscInt *rankval; 133*c4762a1bSJed Brown PetscInt npoints[2],npoints_orig[2]; 134*c4762a1bSJed Brown 135*c4762a1bSJed Brown ierr = DMSwarmGetLocalSize(dms,&npoints_orig[0]);CHKERRQ(ierr); 136*c4762a1bSJed Brown ierr = DMSwarmGetSize(dms,&npoints_orig[1]);CHKERRQ(ierr); 137*c4762a1bSJed Brown ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 138*c4762a1bSJed Brown ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] before(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1]);CHKERRQ(ierr); 139*c4762a1bSJed Brown ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 140*c4762a1bSJed Brown ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 141*c4762a1bSJed Brown 142*c4762a1bSJed Brown ierr = DMSwarmGetField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 143*c4762a1bSJed Brown 144*c4762a1bSJed Brown if (rank == 1) { 145*c4762a1bSJed Brown rankval[0] = -1; 146*c4762a1bSJed Brown } 147*c4762a1bSJed Brown if (rank == 2) { 148*c4762a1bSJed Brown rankval[1] = -1; 149*c4762a1bSJed Brown } 150*c4762a1bSJed Brown if (rank == 3) { 151*c4762a1bSJed Brown rankval[3] = -1; 152*c4762a1bSJed Brown rankval[4] = -1; 153*c4762a1bSJed Brown } 154*c4762a1bSJed Brown ierr = DMSwarmRestoreField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 155*c4762a1bSJed Brown ierr = DMSwarmCollectViewCreate(dms);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = DMSwarmGetLocalSize(dms,&npoints[0]);CHKERRQ(ierr); 157*c4762a1bSJed Brown ierr = DMSwarmGetSize(dms,&npoints[1]);CHKERRQ(ierr); 158*c4762a1bSJed Brown ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 159*c4762a1bSJed Brown ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after(%D,%D)\n",rank,npoints[0],npoints[1]);CHKERRQ(ierr); 160*c4762a1bSJed Brown ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 161*c4762a1bSJed Brown ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 162*c4762a1bSJed Brown 163*c4762a1bSJed Brown ierr = DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr); 164*c4762a1bSJed Brown ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 165*c4762a1bSJed Brown ierr = DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr); 166*c4762a1bSJed Brown 167*c4762a1bSJed Brown ierr = DMSwarmCollectViewDestroy(dms);CHKERRQ(ierr); 168*c4762a1bSJed Brown ierr = DMSwarmGetLocalSize(dms,&npoints[0]);CHKERRQ(ierr); 169*c4762a1bSJed Brown ierr = DMSwarmGetSize(dms,&npoints[1]);CHKERRQ(ierr); 170*c4762a1bSJed Brown ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 171*c4762a1bSJed Brown ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after_v(%D,%D)\n",rank,npoints[0],npoints[1]);CHKERRQ(ierr); 172*c4762a1bSJed Brown ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 173*c4762a1bSJed Brown ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 174*c4762a1bSJed Brown 175*c4762a1bSJed Brown ierr = DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr); 176*c4762a1bSJed Brown ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 177*c4762a1bSJed Brown ierr = DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr); 178*c4762a1bSJed Brown } 179*c4762a1bSJed Brown ierr = DMDestroy(&dms);CHKERRQ(ierr); 180*c4762a1bSJed Brown PetscFunctionReturn(0); 181*c4762a1bSJed Brown } 182*c4762a1bSJed Brown 183*c4762a1bSJed Brown /* 184*c4762a1bSJed Brown splot "c-rank0.gp","c-rank1.gp","c-rank2.gp","c-rank3.gp" 185*c4762a1bSJed Brown */ 186*c4762a1bSJed Brown PetscErrorCode ex1_3(void) 187*c4762a1bSJed Brown { 188*c4762a1bSJed Brown DM dms; 189*c4762a1bSJed Brown PetscErrorCode ierr; 190*c4762a1bSJed Brown PetscMPIInt rank,size; 191*c4762a1bSJed Brown PetscInt is,js,ni,nj,overlap; 192*c4762a1bSJed Brown DM dmcell; 193*c4762a1bSJed Brown 194*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 195*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 196*c4762a1bSJed Brown overlap = 2; 197*c4762a1bSJed Brown ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,13,13,PETSC_DECIDE,PETSC_DECIDE,1,overlap,NULL,NULL,&dmcell);CHKERRQ(ierr); 198*c4762a1bSJed Brown ierr = DMSetFromOptions(dmcell);CHKERRQ(ierr); 199*c4762a1bSJed Brown ierr = DMSetUp(dmcell);CHKERRQ(ierr); 200*c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(dmcell,-1.0,1.0,-1.0,1.0,-1.0,1.0);CHKERRQ(ierr); 201*c4762a1bSJed Brown ierr = DMDAGetCorners(dmcell,&is,&js,NULL,&ni,&nj,NULL);CHKERRQ(ierr); 202*c4762a1bSJed Brown ierr = DMCreate(PETSC_COMM_WORLD,&dms);CHKERRQ(ierr); 203*c4762a1bSJed Brown ierr = DMSetType(dms,DMSWARM);CHKERRQ(ierr); 204*c4762a1bSJed Brown ierr = DMSwarmSetCellDM(dms,dmcell);CHKERRQ(ierr); 205*c4762a1bSJed Brown 206*c4762a1bSJed Brown /* load in data types */ 207*c4762a1bSJed Brown ierr = DMSwarmInitializeFieldRegister(dms);CHKERRQ(ierr); 208*c4762a1bSJed Brown ierr = DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);CHKERRQ(ierr); 209*c4762a1bSJed Brown ierr = DMSwarmRegisterPetscDatatypeField(dms,"coorx",1,PETSC_REAL);CHKERRQ(ierr); 210*c4762a1bSJed Brown ierr = DMSwarmRegisterPetscDatatypeField(dms,"coory",1,PETSC_REAL);CHKERRQ(ierr); 211*c4762a1bSJed Brown ierr = DMSwarmFinalizeFieldRegister(dms);CHKERRQ(ierr); 212*c4762a1bSJed Brown ierr = DMSwarmSetLocalSizes(dms,ni*nj*4,4);CHKERRQ(ierr); 213*c4762a1bSJed Brown ierr = DMView(dms,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 214*c4762a1bSJed Brown 215*c4762a1bSJed Brown /* set values within the swarm */ 216*c4762a1bSJed Brown { 217*c4762a1bSJed Brown PetscReal *array_x,*array_y; 218*c4762a1bSJed Brown PetscInt npoints,i,j,cnt; 219*c4762a1bSJed Brown DMDACoor2d **LA_coor; 220*c4762a1bSJed Brown Vec coor; 221*c4762a1bSJed Brown DM dmcellcdm; 222*c4762a1bSJed Brown 223*c4762a1bSJed Brown ierr = DMGetCoordinateDM(dmcell,&dmcellcdm);CHKERRQ(ierr); 224*c4762a1bSJed Brown ierr = DMGetCoordinates(dmcell,&coor);CHKERRQ(ierr); 225*c4762a1bSJed Brown ierr = DMDAVecGetArray(dmcellcdm,coor,&LA_coor);CHKERRQ(ierr); 226*c4762a1bSJed Brown ierr = DMSwarmGetLocalSize(dms,&npoints);CHKERRQ(ierr); 227*c4762a1bSJed Brown ierr = DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 228*c4762a1bSJed Brown ierr = DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 229*c4762a1bSJed Brown cnt = 0; 230*c4762a1bSJed Brown for (j=js; j<js+nj; j++) { 231*c4762a1bSJed Brown for (i=is; i<is+ni; i++) { 232*c4762a1bSJed Brown PetscReal xp,yp; 233*c4762a1bSJed Brown xp = PetscRealPart(LA_coor[j][i].x); 234*c4762a1bSJed Brown yp = PetscRealPart(LA_coor[j][i].y); 235*c4762a1bSJed Brown array_x[4*cnt+0] = xp - 0.05; if (array_x[4*cnt+0] < -1.0) { array_x[4*cnt+0] = -1.0+1.0e-12; } 236*c4762a1bSJed Brown array_x[4*cnt+1] = xp + 0.05; if (array_x[4*cnt+1] > 1.0) { array_x[4*cnt+1] = 1.0-1.0e-12; } 237*c4762a1bSJed Brown array_x[4*cnt+2] = xp - 0.05; if (array_x[4*cnt+2] < -1.0) { array_x[4*cnt+2] = -1.0+1.0e-12; } 238*c4762a1bSJed Brown array_x[4*cnt+3] = xp + 0.05; if (array_x[4*cnt+3] > 1.0) { array_x[4*cnt+3] = 1.0-1.0e-12; } 239*c4762a1bSJed Brown 240*c4762a1bSJed Brown array_y[4*cnt+0] = yp - 0.05; if (array_y[4*cnt+0] < -1.0) { array_y[4*cnt+0] = -1.0+1.0e-12; } 241*c4762a1bSJed Brown array_y[4*cnt+1] = yp - 0.05; if (array_y[4*cnt+1] < -1.0) { array_y[4*cnt+1] = -1.0+1.0e-12; } 242*c4762a1bSJed Brown array_y[4*cnt+2] = yp + 0.05; if (array_y[4*cnt+2] > 1.0) { array_y[4*cnt+2] = 1.0-1.0e-12; } 243*c4762a1bSJed Brown array_y[4*cnt+3] = yp + 0.05; if (array_y[4*cnt+3] > 1.0) { array_y[4*cnt+3] = 1.0-1.0e-12; } 244*c4762a1bSJed Brown cnt++; 245*c4762a1bSJed Brown } 246*c4762a1bSJed Brown } 247*c4762a1bSJed Brown ierr = DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 248*c4762a1bSJed Brown ierr = DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 249*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(dmcellcdm,coor,&LA_coor);CHKERRQ(ierr); 250*c4762a1bSJed Brown } 251*c4762a1bSJed Brown { 252*c4762a1bSJed Brown PetscInt npoints[2],npoints_orig[2],ng; 253*c4762a1bSJed Brown 254*c4762a1bSJed Brown ierr = DMSwarmGetLocalSize(dms,&npoints_orig[0]);CHKERRQ(ierr); 255*c4762a1bSJed Brown ierr = DMSwarmGetSize(dms,&npoints_orig[1]);CHKERRQ(ierr); 256*c4762a1bSJed Brown ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 257*c4762a1bSJed Brown ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] before(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1]);CHKERRQ(ierr); 258*c4762a1bSJed Brown ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 259*c4762a1bSJed Brown ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 260*c4762a1bSJed Brown ierr = DMSwarmCollect_DMDABoundingBox(dms,&ng);CHKERRQ(ierr); 261*c4762a1bSJed Brown 262*c4762a1bSJed Brown ierr = DMSwarmGetLocalSize(dms,&npoints[0]);CHKERRQ(ierr); 263*c4762a1bSJed Brown ierr = DMSwarmGetSize(dms,&npoints[1]);CHKERRQ(ierr); 264*c4762a1bSJed Brown ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 265*c4762a1bSJed Brown ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after(%D,%D)\n",rank,npoints[0],npoints[1]);CHKERRQ(ierr); 266*c4762a1bSJed Brown ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 267*c4762a1bSJed Brown ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 268*c4762a1bSJed Brown } 269*c4762a1bSJed Brown { 270*c4762a1bSJed Brown PetscReal *array_x,*array_y; 271*c4762a1bSJed Brown PetscInt npoints,p; 272*c4762a1bSJed Brown FILE *fp = NULL; 273*c4762a1bSJed Brown char name[PETSC_MAX_PATH_LEN]; 274*c4762a1bSJed Brown 275*c4762a1bSJed Brown ierr = PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"c-rank%d.gp",rank);CHKERRQ(ierr); 276*c4762a1bSJed Brown fp = fopen(name,"w"); 277*c4762a1bSJed Brown if (!fp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file %s",name); 278*c4762a1bSJed Brown ierr = DMSwarmGetLocalSize(dms,&npoints);CHKERRQ(ierr); 279*c4762a1bSJed Brown ierr = DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 280*c4762a1bSJed Brown ierr = DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 281*c4762a1bSJed Brown for (p=0; p<npoints; p++) { 282*c4762a1bSJed Brown fprintf(fp,"%+1.4e %+1.4e %1.4e\n",array_x[p],array_y[p],(double)rank); 283*c4762a1bSJed Brown } 284*c4762a1bSJed Brown ierr = DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 285*c4762a1bSJed Brown ierr = DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 286*c4762a1bSJed Brown fclose(fp); 287*c4762a1bSJed Brown } 288*c4762a1bSJed Brown ierr = DMDestroy(&dmcell);CHKERRQ(ierr); 289*c4762a1bSJed Brown ierr = DMDestroy(&dms);CHKERRQ(ierr); 290*c4762a1bSJed Brown PetscFunctionReturn(0); 291*c4762a1bSJed Brown } 292*c4762a1bSJed Brown 293*c4762a1bSJed Brown typedef struct { 294*c4762a1bSJed Brown PetscReal cx[2]; 295*c4762a1bSJed Brown PetscReal radius; 296*c4762a1bSJed Brown } CollectZoneCtx; 297*c4762a1bSJed Brown 298*c4762a1bSJed Brown PetscErrorCode collect_zone(DM dm,void *ctx,PetscInt *nfound,PetscInt **foundlist) 299*c4762a1bSJed Brown { 300*c4762a1bSJed Brown CollectZoneCtx *zone = (CollectZoneCtx*)ctx; 301*c4762a1bSJed Brown PetscInt p,npoints; 302*c4762a1bSJed Brown PetscReal *array_x,*array_y,r2; 303*c4762a1bSJed Brown PetscInt p2collect,*plist; 304*c4762a1bSJed Brown PetscMPIInt rank; 305*c4762a1bSJed Brown PetscErrorCode ierr; 306*c4762a1bSJed Brown 307*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 308*c4762a1bSJed Brown ierr = DMSwarmGetLocalSize(dm,&npoints);CHKERRQ(ierr); 309*c4762a1bSJed Brown ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 310*c4762a1bSJed Brown ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 311*c4762a1bSJed Brown 312*c4762a1bSJed Brown r2 = zone->radius * zone->radius; 313*c4762a1bSJed Brown p2collect = 0; 314*c4762a1bSJed Brown for (p=0; p<npoints; p++) { 315*c4762a1bSJed Brown PetscReal sep2; 316*c4762a1bSJed Brown 317*c4762a1bSJed Brown sep2 = (array_x[p] - zone->cx[0])*(array_x[p] - zone->cx[0]); 318*c4762a1bSJed Brown sep2 += (array_y[p] - zone->cx[1])*(array_y[p] - zone->cx[1]); 319*c4762a1bSJed Brown if (sep2 < r2) { 320*c4762a1bSJed Brown p2collect++; 321*c4762a1bSJed Brown } 322*c4762a1bSJed Brown } 323*c4762a1bSJed Brown 324*c4762a1bSJed Brown ierr = PetscMalloc1(p2collect+1,&plist);CHKERRQ(ierr); 325*c4762a1bSJed Brown p2collect = 0; 326*c4762a1bSJed Brown for (p=0; p<npoints; p++) { 327*c4762a1bSJed Brown PetscReal sep2; 328*c4762a1bSJed Brown 329*c4762a1bSJed Brown sep2 = (array_x[p] - zone->cx[0])*(array_x[p] - zone->cx[0]); 330*c4762a1bSJed Brown sep2 += (array_y[p] - zone->cx[1])*(array_y[p] - zone->cx[1]); 331*c4762a1bSJed Brown if (sep2 < r2) { 332*c4762a1bSJed Brown plist[p2collect] = p; 333*c4762a1bSJed Brown p2collect++; 334*c4762a1bSJed Brown } 335*c4762a1bSJed Brown } 336*c4762a1bSJed Brown ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 337*c4762a1bSJed Brown ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 338*c4762a1bSJed Brown 339*c4762a1bSJed Brown *nfound = p2collect; 340*c4762a1bSJed Brown *foundlist = plist; 341*c4762a1bSJed Brown PetscFunctionReturn(0); 342*c4762a1bSJed Brown } 343*c4762a1bSJed Brown 344*c4762a1bSJed Brown PetscErrorCode ex1_4(void) 345*c4762a1bSJed Brown { 346*c4762a1bSJed Brown DM dms; 347*c4762a1bSJed Brown PetscErrorCode ierr; 348*c4762a1bSJed Brown PetscMPIInt rank,size; 349*c4762a1bSJed Brown PetscInt is,js,ni,nj,overlap,nn; 350*c4762a1bSJed Brown DM dmcell; 351*c4762a1bSJed Brown CollectZoneCtx *zone; 352*c4762a1bSJed Brown PetscReal dx; 353*c4762a1bSJed Brown 354*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 355*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 356*c4762a1bSJed Brown nn = 101; 357*c4762a1bSJed Brown dx = 2.0/ (PetscReal)(nn-1); 358*c4762a1bSJed Brown overlap = 0; 359*c4762a1bSJed Brown ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,nn,nn,PETSC_DECIDE,PETSC_DECIDE,1,overlap,NULL,NULL,&dmcell);CHKERRQ(ierr); 360*c4762a1bSJed Brown ierr = DMSetFromOptions(dmcell);CHKERRQ(ierr); 361*c4762a1bSJed Brown ierr = DMSetUp(dmcell);CHKERRQ(ierr); 362*c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(dmcell,-1.0,1.0,-1.0,1.0,-1.0,1.0);CHKERRQ(ierr); 363*c4762a1bSJed Brown ierr = DMDAGetCorners(dmcell,&is,&js,NULL,&ni,&nj,NULL);CHKERRQ(ierr); 364*c4762a1bSJed Brown ierr = DMCreate(PETSC_COMM_WORLD,&dms);CHKERRQ(ierr); 365*c4762a1bSJed Brown ierr = DMSetType(dms,DMSWARM);CHKERRQ(ierr); 366*c4762a1bSJed Brown 367*c4762a1bSJed Brown /* load in data types */ 368*c4762a1bSJed Brown ierr = DMSwarmInitializeFieldRegister(dms);CHKERRQ(ierr); 369*c4762a1bSJed Brown ierr = DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);CHKERRQ(ierr); 370*c4762a1bSJed Brown ierr = DMSwarmRegisterPetscDatatypeField(dms,"coorx",1,PETSC_REAL);CHKERRQ(ierr); 371*c4762a1bSJed Brown ierr = DMSwarmRegisterPetscDatatypeField(dms,"coory",1,PETSC_REAL);CHKERRQ(ierr); 372*c4762a1bSJed Brown ierr = DMSwarmFinalizeFieldRegister(dms);CHKERRQ(ierr); 373*c4762a1bSJed Brown ierr = DMSwarmSetLocalSizes(dms,ni*nj*4,4);CHKERRQ(ierr); 374*c4762a1bSJed Brown ierr = DMView(dms,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 375*c4762a1bSJed Brown 376*c4762a1bSJed Brown /* set values within the swarm */ 377*c4762a1bSJed Brown { 378*c4762a1bSJed Brown PetscReal *array_x,*array_y; 379*c4762a1bSJed Brown PetscInt npoints,i,j,cnt; 380*c4762a1bSJed Brown DMDACoor2d **LA_coor; 381*c4762a1bSJed Brown Vec coor; 382*c4762a1bSJed Brown DM dmcellcdm; 383*c4762a1bSJed Brown 384*c4762a1bSJed Brown ierr = DMGetCoordinateDM(dmcell,&dmcellcdm);CHKERRQ(ierr); 385*c4762a1bSJed Brown ierr = DMGetCoordinates(dmcell,&coor);CHKERRQ(ierr); 386*c4762a1bSJed Brown ierr = DMDAVecGetArray(dmcellcdm,coor,&LA_coor);CHKERRQ(ierr); 387*c4762a1bSJed Brown ierr = DMSwarmGetLocalSize(dms,&npoints);CHKERRQ(ierr); 388*c4762a1bSJed Brown ierr = DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 389*c4762a1bSJed Brown ierr = DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 390*c4762a1bSJed Brown cnt = 0; 391*c4762a1bSJed Brown for (j=js; j<js+nj; j++) { 392*c4762a1bSJed Brown for (i=is; i<is+ni; i++) { 393*c4762a1bSJed Brown PetscReal xp,yp; 394*c4762a1bSJed Brown 395*c4762a1bSJed Brown xp = PetscRealPart(LA_coor[j][i].x); 396*c4762a1bSJed Brown yp = PetscRealPart(LA_coor[j][i].y); 397*c4762a1bSJed Brown array_x[4*cnt+0] = xp - dx*0.1; /*if (array_x[4*cnt+0] < -1.0) { array_x[4*cnt+0] = -1.0+1.0e-12; }*/ 398*c4762a1bSJed Brown array_x[4*cnt+1] = xp + dx*0.1; /*if (array_x[4*cnt+1] > 1.0) { array_x[4*cnt+1] = 1.0-1.0e-12; }*/ 399*c4762a1bSJed Brown array_x[4*cnt+2] = xp - dx*0.1; /*if (array_x[4*cnt+2] < -1.0) { array_x[4*cnt+2] = -1.0+1.0e-12; }*/ 400*c4762a1bSJed Brown array_x[4*cnt+3] = xp + dx*0.1; /*if (array_x[4*cnt+3] > 1.0) { array_x[4*cnt+3] = 1.0-1.0e-12; }*/ 401*c4762a1bSJed Brown array_y[4*cnt+0] = yp - dx*0.1; /*if (array_y[4*cnt+0] < -1.0) { array_y[4*cnt+0] = -1.0+1.0e-12; }*/ 402*c4762a1bSJed Brown array_y[4*cnt+1] = yp - dx*0.1; /*if (array_y[4*cnt+1] < -1.0) { array_y[4*cnt+1] = -1.0+1.0e-12; }*/ 403*c4762a1bSJed Brown array_y[4*cnt+2] = yp + dx*0.1; /*if (array_y[4*cnt+2] > 1.0) { array_y[4*cnt+2] = 1.0-1.0e-12; }*/ 404*c4762a1bSJed Brown array_y[4*cnt+3] = yp + dx*0.1; /*if (array_y[4*cnt+3] > 1.0) { array_y[4*cnt+3] = 1.0-1.0e-12; }*/ 405*c4762a1bSJed Brown cnt++; 406*c4762a1bSJed Brown } 407*c4762a1bSJed Brown } 408*c4762a1bSJed Brown ierr = DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 409*c4762a1bSJed Brown ierr = DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 410*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(dmcellcdm,coor,&LA_coor);CHKERRQ(ierr); 411*c4762a1bSJed Brown } 412*c4762a1bSJed Brown ierr = PetscMalloc1(1,&zone);CHKERRQ(ierr); 413*c4762a1bSJed Brown if (size == 4) { 414*c4762a1bSJed Brown if (rank == 0) { 415*c4762a1bSJed Brown zone->cx[0] = 0.5; 416*c4762a1bSJed Brown zone->cx[1] = 0.5; 417*c4762a1bSJed Brown zone->radius = 0.3; 418*c4762a1bSJed Brown } 419*c4762a1bSJed Brown if (rank == 1) { 420*c4762a1bSJed Brown zone->cx[0] = -0.5; 421*c4762a1bSJed Brown zone->cx[1] = 0.5; 422*c4762a1bSJed Brown zone->radius = 0.25; 423*c4762a1bSJed Brown } 424*c4762a1bSJed Brown if (rank == 2) { 425*c4762a1bSJed Brown zone->cx[0] = 0.5; 426*c4762a1bSJed Brown zone->cx[1] = -0.5; 427*c4762a1bSJed Brown zone->radius = 0.2; 428*c4762a1bSJed Brown } 429*c4762a1bSJed Brown if (rank == 3) { 430*c4762a1bSJed Brown zone->cx[0] = -0.5; 431*c4762a1bSJed Brown zone->cx[1] = -0.5; 432*c4762a1bSJed Brown zone->radius = 0.1; 433*c4762a1bSJed Brown } 434*c4762a1bSJed Brown } else { 435*c4762a1bSJed Brown if (rank == 0) { 436*c4762a1bSJed Brown zone->cx[0] = 0.5; 437*c4762a1bSJed Brown zone->cx[1] = 0.5; 438*c4762a1bSJed Brown zone->radius = 0.8; 439*c4762a1bSJed Brown } else { 440*c4762a1bSJed Brown zone->cx[0] = 10.0; 441*c4762a1bSJed Brown zone->cx[1] = 10.0; 442*c4762a1bSJed Brown zone->radius = 0.0; 443*c4762a1bSJed Brown } 444*c4762a1bSJed Brown } 445*c4762a1bSJed Brown { 446*c4762a1bSJed Brown PetscInt npoints[2],npoints_orig[2],ng; 447*c4762a1bSJed Brown 448*c4762a1bSJed Brown ierr = DMSwarmGetLocalSize(dms,&npoints_orig[0]);CHKERRQ(ierr); 449*c4762a1bSJed Brown ierr = DMSwarmGetSize(dms,&npoints_orig[1]);CHKERRQ(ierr); 450*c4762a1bSJed Brown ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 451*c4762a1bSJed Brown ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] before(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1]);CHKERRQ(ierr); 452*c4762a1bSJed Brown ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 453*c4762a1bSJed Brown ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 454*c4762a1bSJed Brown ierr = DMSwarmCollect_General(dms,collect_zone,sizeof(CollectZoneCtx),zone,&ng);CHKERRQ(ierr); 455*c4762a1bSJed Brown ierr = DMSwarmGetLocalSize(dms,&npoints[0]);CHKERRQ(ierr); 456*c4762a1bSJed Brown ierr = DMSwarmGetSize(dms,&npoints[1]);CHKERRQ(ierr); 457*c4762a1bSJed Brown ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 458*c4762a1bSJed Brown ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after(%D,%D)\n",rank,npoints[0],npoints[1]);CHKERRQ(ierr); 459*c4762a1bSJed Brown ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 460*c4762a1bSJed Brown ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 461*c4762a1bSJed Brown } 462*c4762a1bSJed Brown { 463*c4762a1bSJed Brown PetscReal *array_x,*array_y; 464*c4762a1bSJed Brown PetscInt npoints,p; 465*c4762a1bSJed Brown FILE *fp = NULL; 466*c4762a1bSJed Brown char name[PETSC_MAX_PATH_LEN]; 467*c4762a1bSJed Brown 468*c4762a1bSJed Brown ierr = PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"c-rank%d.gp",rank);CHKERRQ(ierr); 469*c4762a1bSJed Brown fp = fopen(name,"w"); 470*c4762a1bSJed Brown if (!fp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file %s",name); 471*c4762a1bSJed Brown ierr = DMSwarmGetLocalSize(dms,&npoints);CHKERRQ(ierr); 472*c4762a1bSJed Brown ierr = DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 473*c4762a1bSJed Brown ierr = DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 474*c4762a1bSJed Brown for (p=0; p<npoints; p++) { 475*c4762a1bSJed Brown fprintf(fp,"%+1.4e %+1.4e %1.4e\n",array_x[p],array_y[p],(double)rank); 476*c4762a1bSJed Brown } 477*c4762a1bSJed Brown ierr = DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr); 478*c4762a1bSJed Brown ierr = DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr); 479*c4762a1bSJed Brown fclose(fp); 480*c4762a1bSJed Brown } 481*c4762a1bSJed Brown ierr = DMDestroy(&dmcell);CHKERRQ(ierr); 482*c4762a1bSJed Brown ierr = DMDestroy(&dms);CHKERRQ(ierr); 483*c4762a1bSJed Brown ierr = PetscFree(zone);CHKERRQ(ierr); 484*c4762a1bSJed Brown PetscFunctionReturn(0); 485*c4762a1bSJed Brown } 486*c4762a1bSJed Brown 487*c4762a1bSJed Brown int main(int argc,char **argv) 488*c4762a1bSJed Brown { 489*c4762a1bSJed Brown PetscErrorCode ierr; 490*c4762a1bSJed Brown PetscInt test_mode = 4; 491*c4762a1bSJed Brown 492*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 493*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-test_mode",&test_mode,NULL);CHKERRQ(ierr); 494*c4762a1bSJed Brown if (test_mode == 1) { 495*c4762a1bSJed Brown ierr = ex1_1();CHKERRQ(ierr); 496*c4762a1bSJed Brown } else if (test_mode == 2) { 497*c4762a1bSJed Brown ierr = ex1_2();CHKERRQ(ierr); 498*c4762a1bSJed Brown } else if (test_mode == 3) { 499*c4762a1bSJed Brown ierr = ex1_3();CHKERRQ(ierr); 500*c4762a1bSJed Brown } else if (test_mode == 4) { 501*c4762a1bSJed Brown ierr = ex1_4();CHKERRQ(ierr); 502*c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown test_mode value, should be 1,2,3,4"); 503*c4762a1bSJed Brown ierr = PetscFinalize(); 504*c4762a1bSJed Brown return ierr; 505*c4762a1bSJed Brown } 506*c4762a1bSJed Brown 507*c4762a1bSJed Brown /*TEST 508*c4762a1bSJed Brown 509*c4762a1bSJed Brown build: 510*c4762a1bSJed Brown requires: !complex double 511*c4762a1bSJed Brown 512*c4762a1bSJed Brown test: 513*c4762a1bSJed Brown args: -test_mode 1 514*c4762a1bSJed Brown filter: grep -v atomic 515*c4762a1bSJed Brown filter_output: grep -v atomic 516*c4762a1bSJed Brown 517*c4762a1bSJed Brown test: 518*c4762a1bSJed Brown suffix: 2 519*c4762a1bSJed Brown args: -test_mode 2 520*c4762a1bSJed Brown filter: grep -v atomic 521*c4762a1bSJed Brown filter_output: grep -v atomic 522*c4762a1bSJed Brown 523*c4762a1bSJed Brown test: 524*c4762a1bSJed Brown suffix: 3 525*c4762a1bSJed Brown args: -test_mode 3 526*c4762a1bSJed Brown filter: grep -v atomic 527*c4762a1bSJed Brown filter_output: grep -v atomic 528*c4762a1bSJed Brown TODO: broken 529*c4762a1bSJed Brown 530*c4762a1bSJed Brown test: 531*c4762a1bSJed Brown suffix: 4 532*c4762a1bSJed Brown args: -test_mode 4 533*c4762a1bSJed Brown filter: grep -v atomic 534*c4762a1bSJed Brown filter_output: grep -v atomic 535*c4762a1bSJed Brown 536*c4762a1bSJed Brown test: 537*c4762a1bSJed Brown suffix: 5 538*c4762a1bSJed Brown nsize: 4 539*c4762a1bSJed Brown args: -test_mode 1 540*c4762a1bSJed Brown filter: grep -v atomic 541*c4762a1bSJed Brown filter_output: grep -v atomic 542*c4762a1bSJed Brown 543*c4762a1bSJed Brown test: 544*c4762a1bSJed Brown suffix: 6 545*c4762a1bSJed Brown nsize: 4 546*c4762a1bSJed Brown args: -test_mode 2 547*c4762a1bSJed Brown filter: grep -v atomic 548*c4762a1bSJed Brown filter_output: grep -v atomic 549*c4762a1bSJed Brown 550*c4762a1bSJed Brown test: 551*c4762a1bSJed Brown suffix: 7 552*c4762a1bSJed Brown nsize: 4 553*c4762a1bSJed Brown args: -test_mode 3 554*c4762a1bSJed Brown filter: grep -v atomic 555*c4762a1bSJed Brown filter_output: grep -v atomic 556*c4762a1bSJed Brown TODO: broken 557*c4762a1bSJed Brown 558*c4762a1bSJed Brown test: 559*c4762a1bSJed Brown suffix: 8 560*c4762a1bSJed Brown nsize: 4 561*c4762a1bSJed Brown args: -test_mode 4 562*c4762a1bSJed Brown filter: grep -v atomic 563*c4762a1bSJed Brown filter_output: grep -v atomic 564*c4762a1bSJed Brown 565*c4762a1bSJed Brown TEST*/ 566