1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests DMSwarm\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscdm.h> 5c4762a1bSJed Brown #include <petscdmda.h> 6c4762a1bSJed Brown #include <petscdmswarm.h> 7c4762a1bSJed Brown 8c4762a1bSJed Brown PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM,PetscErrorCode (*)(DM,void*,PetscInt*,PetscInt**),size_t,void*,PetscInt*); 9c4762a1bSJed Brown PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM,PetscInt*); 10c4762a1bSJed Brown 11c4762a1bSJed Brown PetscErrorCode ex1_1(void) 12c4762a1bSJed Brown { 13c4762a1bSJed Brown DM dms; 14c4762a1bSJed Brown Vec x; 15c4762a1bSJed Brown PetscMPIInt rank,size; 16c4762a1bSJed Brown PetscInt p; 17c4762a1bSJed Brown 18362febeeSStefano Zampini PetscFunctionBegin; 199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 21*08401ef6SPierre Jolivet PetscCheck(!(size > 1) || !(size != 4),PETSC_COMM_WORLD,PETSC_ERR_SUP,"Must be run wuth 4 MPI ranks"); 22c4762a1bSJed Brown 239566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD,&dms)); 249566063dSJacob Faibussowitsch PetscCall(DMSetType(dms,DMSWARM)); 25c4762a1bSJed Brown 269566063dSJacob Faibussowitsch PetscCall(DMSwarmInitializeFieldRegister(dms)); 279566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL)); 289566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"strain",1,PETSC_REAL)); 299566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(dms)); 309566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dms,5+rank,4)); 319566063dSJacob Faibussowitsch PetscCall(DMView(dms,PETSC_VIEWER_STDOUT_WORLD)); 32c4762a1bSJed Brown 33c4762a1bSJed Brown { 34c4762a1bSJed Brown PetscReal *array; 359566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms,"viscosity",NULL,NULL,(void**)&array)); 36c4762a1bSJed Brown for (p=0; p<5+rank; p++) { 37c4762a1bSJed Brown array[p] = 11.1 + p*0.1 + rank*100.0; 38c4762a1bSJed Brown } 399566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms,"viscosity",NULL,NULL,(void**)&array)); 40c4762a1bSJed Brown } 41c4762a1bSJed Brown 42c4762a1bSJed Brown { 43c4762a1bSJed Brown PetscReal *array; 449566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms,"strain",NULL,NULL,(void**)&array)); 45c4762a1bSJed Brown for (p=0; p<5+rank; p++) { 46c4762a1bSJed Brown array[p] = 2.0e-2 + p*0.001 + rank*1.0; 47c4762a1bSJed Brown } 489566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms,"strain",NULL,NULL,(void**)&array)); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown 519566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x)); 529566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x)); 53c4762a1bSJed Brown 549566063dSJacob Faibussowitsch PetscCall(DMSwarmVectorDefineField(dms,"strain")); 559566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dms,&x)); 569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 57c4762a1bSJed Brown 58c4762a1bSJed Brown { 59c4762a1bSJed Brown PetscInt *rankval; 60c4762a1bSJed Brown PetscInt npoints[2],npoints_orig[2]; 61c4762a1bSJed Brown 629566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dms,&npoints_orig[0])); 639566063dSJacob Faibussowitsch PetscCall(DMSwarmGetSize(dms,&npoints_orig[1])); 649566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval)); 65c4762a1bSJed Brown if ((rank == 0) && (size > 1)) { 66c4762a1bSJed Brown rankval[0] = 1; 67c4762a1bSJed Brown rankval[3] = 1; 68c4762a1bSJed Brown } 69c4762a1bSJed Brown if (rank == 3) { 70c4762a1bSJed Brown rankval[2] = 1; 71c4762a1bSJed Brown } 729566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval)); 739566063dSJacob Faibussowitsch PetscCall(DMSwarmMigrate(dms,PETSC_TRUE)); 749566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dms,&npoints[0])); 759566063dSJacob Faibussowitsch PetscCall(DMSwarmGetSize(dms,&npoints[1])); 769566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 779566063dSJacob Faibussowitsch PetscCall(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])); 789566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 799566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 80c4762a1bSJed Brown } 81c4762a1bSJed Brown { 829566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x)); 839566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 849566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x)); 85c4762a1bSJed Brown } 86c4762a1bSJed Brown { 879566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(dms,"strain",&x)); 889566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 899566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(dms,"strain",&x)); 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 929566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dms)); 93c4762a1bSJed Brown PetscFunctionReturn(0); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 96c4762a1bSJed Brown PetscErrorCode ex1_2(void) 97c4762a1bSJed Brown { 98c4762a1bSJed Brown DM dms; 99c4762a1bSJed Brown Vec x; 100c4762a1bSJed Brown PetscMPIInt rank,size; 101c4762a1bSJed Brown PetscInt p; 102c4762a1bSJed Brown 103362febeeSStefano Zampini PetscFunctionBegin; 1049566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 1059566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 1069566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD,&dms)); 1079566063dSJacob Faibussowitsch PetscCall(DMSetType(dms,DMSWARM)); 1089566063dSJacob Faibussowitsch PetscCall(DMSwarmInitializeFieldRegister(dms)); 109c4762a1bSJed Brown 1109566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL)); 1119566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"strain",1,PETSC_REAL)); 1129566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(dms)); 1139566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dms,5+rank,4)); 1149566063dSJacob Faibussowitsch PetscCall(DMView(dms,PETSC_VIEWER_STDOUT_WORLD)); 115c4762a1bSJed Brown { 116c4762a1bSJed Brown PetscReal *array; 1179566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms,"viscosity",NULL,NULL,(void**)&array)); 118c4762a1bSJed Brown for (p=0; p<5+rank; p++) { 119c4762a1bSJed Brown array[p] = 11.1 + p*0.1 + rank*100.0; 120c4762a1bSJed Brown } 1219566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms,"viscosity",NULL,NULL,(void**)&array)); 122c4762a1bSJed Brown } 123c4762a1bSJed Brown { 124c4762a1bSJed Brown PetscReal *array; 1259566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms,"strain",NULL,NULL,(void**)&array)); 126c4762a1bSJed Brown for (p=0; p<5+rank; p++) { 127c4762a1bSJed Brown array[p] = 2.0e-2 + p*0.001 + rank*1.0; 128c4762a1bSJed Brown } 1299566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms,"strain",NULL,NULL,(void**)&array)); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown { 132c4762a1bSJed Brown PetscInt *rankval; 133c4762a1bSJed Brown PetscInt npoints[2],npoints_orig[2]; 134c4762a1bSJed Brown 1359566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dms,&npoints_orig[0])); 1369566063dSJacob Faibussowitsch PetscCall(DMSwarmGetSize(dms,&npoints_orig[1])); 1379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 1389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] before(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1])); 1399566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 1409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 141c4762a1bSJed Brown 1429566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval)); 143c4762a1bSJed Brown 144c4762a1bSJed Brown if (rank == 1) { 145c4762a1bSJed Brown rankval[0] = -1; 146c4762a1bSJed Brown } 147c4762a1bSJed Brown if (rank == 2) { 148c4762a1bSJed Brown rankval[1] = -1; 149c4762a1bSJed Brown } 150c4762a1bSJed Brown if (rank == 3) { 151c4762a1bSJed Brown rankval[3] = -1; 152c4762a1bSJed Brown rankval[4] = -1; 153c4762a1bSJed Brown } 1549566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval)); 1559566063dSJacob Faibussowitsch PetscCall(DMSwarmCollectViewCreate(dms)); 1569566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dms,&npoints[0])); 1579566063dSJacob Faibussowitsch PetscCall(DMSwarmGetSize(dms,&npoints[1])); 1589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 1599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after(%D,%D)\n",rank,npoints[0],npoints[1])); 1609566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 1619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 162c4762a1bSJed Brown 1639566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x)); 1649566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 1659566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x)); 166c4762a1bSJed Brown 1679566063dSJacob Faibussowitsch PetscCall(DMSwarmCollectViewDestroy(dms)); 1689566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dms,&npoints[0])); 1699566063dSJacob Faibussowitsch PetscCall(DMSwarmGetSize(dms,&npoints[1])); 1709566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 1719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after_v(%D,%D)\n",rank,npoints[0],npoints[1])); 1729566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 1739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 174c4762a1bSJed Brown 1759566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x)); 1769566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 1779566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x)); 178c4762a1bSJed Brown } 1799566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dms)); 180c4762a1bSJed Brown PetscFunctionReturn(0); 181c4762a1bSJed Brown } 182c4762a1bSJed Brown 183c4762a1bSJed Brown /* 184c4762a1bSJed Brown splot "c-rank0.gp","c-rank1.gp","c-rank2.gp","c-rank3.gp" 185c4762a1bSJed Brown */ 186c4762a1bSJed Brown PetscErrorCode ex1_3(void) 187c4762a1bSJed Brown { 188c4762a1bSJed Brown DM dms; 189c4762a1bSJed Brown PetscMPIInt rank,size; 190c4762a1bSJed Brown PetscInt is,js,ni,nj,overlap; 191c4762a1bSJed Brown DM dmcell; 192c4762a1bSJed Brown 193362febeeSStefano Zampini PetscFunctionBegin; 1949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 1959566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 196c4762a1bSJed Brown overlap = 2; 1979566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,13,13,PETSC_DECIDE,PETSC_DECIDE,1,overlap,NULL,NULL,&dmcell)); 1989566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dmcell)); 1999566063dSJacob Faibussowitsch PetscCall(DMSetUp(dmcell)); 2009566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dmcell,-1.0,1.0,-1.0,1.0,-1.0,1.0)); 2019566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dmcell,&is,&js,NULL,&ni,&nj,NULL)); 2029566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD,&dms)); 2039566063dSJacob Faibussowitsch PetscCall(DMSetType(dms,DMSWARM)); 2049566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(dms,dmcell)); 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* load in data types */ 2079566063dSJacob Faibussowitsch PetscCall(DMSwarmInitializeFieldRegister(dms)); 2089566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL)); 2099566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"coorx",1,PETSC_REAL)); 2109566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"coory",1,PETSC_REAL)); 2119566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(dms)); 2129566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dms,ni*nj*4,4)); 2139566063dSJacob Faibussowitsch PetscCall(DMView(dms,PETSC_VIEWER_STDOUT_WORLD)); 214c4762a1bSJed Brown 215c4762a1bSJed Brown /* set values within the swarm */ 216c4762a1bSJed Brown { 217c4762a1bSJed Brown PetscReal *array_x,*array_y; 218c4762a1bSJed Brown PetscInt npoints,i,j,cnt; 219c4762a1bSJed Brown DMDACoor2d **LA_coor; 220c4762a1bSJed Brown Vec coor; 221c4762a1bSJed Brown DM dmcellcdm; 222c4762a1bSJed Brown 2239566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dmcell,&dmcellcdm)); 2249566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dmcell,&coor)); 2259566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dmcellcdm,coor,&LA_coor)); 2269566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dms,&npoints)); 2279566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x)); 2289566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y)); 229c4762a1bSJed Brown cnt = 0; 230c4762a1bSJed Brown for (j=js; j<js+nj; j++) { 231c4762a1bSJed Brown for (i=is; i<is+ni; i++) { 232c4762a1bSJed Brown PetscReal xp,yp; 233c4762a1bSJed Brown xp = PetscRealPart(LA_coor[j][i].x); 234c4762a1bSJed Brown yp = PetscRealPart(LA_coor[j][i].y); 235c4762a1bSJed 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; } 236c4762a1bSJed 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; } 237c4762a1bSJed 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; } 238c4762a1bSJed 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; } 239c4762a1bSJed Brown 240c4762a1bSJed 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; } 241c4762a1bSJed 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; } 242c4762a1bSJed 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; } 243c4762a1bSJed 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; } 244c4762a1bSJed Brown cnt++; 245c4762a1bSJed Brown } 246c4762a1bSJed Brown } 2479566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y)); 2489566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x)); 2499566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dmcellcdm,coor,&LA_coor)); 250c4762a1bSJed Brown } 251c4762a1bSJed Brown { 252c4762a1bSJed Brown PetscInt npoints[2],npoints_orig[2],ng; 253c4762a1bSJed Brown 2549566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dms,&npoints_orig[0])); 2559566063dSJacob Faibussowitsch PetscCall(DMSwarmGetSize(dms,&npoints_orig[1])); 2569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 2579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] before(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1])); 2589566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 2599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 2609566063dSJacob Faibussowitsch PetscCall(DMSwarmCollect_DMDABoundingBox(dms,&ng)); 261c4762a1bSJed Brown 2629566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dms,&npoints[0])); 2639566063dSJacob Faibussowitsch PetscCall(DMSwarmGetSize(dms,&npoints[1])); 2649566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 2659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after(%D,%D)\n",rank,npoints[0],npoints[1])); 2669566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 2679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 268c4762a1bSJed Brown } 269c4762a1bSJed Brown { 270c4762a1bSJed Brown PetscReal *array_x,*array_y; 271c4762a1bSJed Brown PetscInt npoints,p; 272c4762a1bSJed Brown FILE *fp = NULL; 273c4762a1bSJed Brown char name[PETSC_MAX_PATH_LEN]; 274c4762a1bSJed Brown 2759566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"c-rank%d.gp",rank)); 276c4762a1bSJed Brown fp = fopen(name,"w"); 27728b400f6SJacob Faibussowitsch PetscCheck(fp,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file %s",name); 2789566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dms,&npoints)); 2799566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x)); 2809566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y)); 281c4762a1bSJed Brown for (p=0; p<npoints; p++) { 282c4762a1bSJed Brown fprintf(fp,"%+1.4e %+1.4e %1.4e\n",array_x[p],array_y[p],(double)rank); 283c4762a1bSJed Brown } 2849566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y)); 2859566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x)); 286c4762a1bSJed Brown fclose(fp); 287c4762a1bSJed Brown } 2889566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmcell)); 2899566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dms)); 290c4762a1bSJed Brown PetscFunctionReturn(0); 291c4762a1bSJed Brown } 292c4762a1bSJed Brown 293c4762a1bSJed Brown typedef struct { 294c4762a1bSJed Brown PetscReal cx[2]; 295c4762a1bSJed Brown PetscReal radius; 296c4762a1bSJed Brown } CollectZoneCtx; 297c4762a1bSJed Brown 298c4762a1bSJed Brown PetscErrorCode collect_zone(DM dm,void *ctx,PetscInt *nfound,PetscInt **foundlist) 299c4762a1bSJed Brown { 300c4762a1bSJed Brown CollectZoneCtx *zone = (CollectZoneCtx*)ctx; 301c4762a1bSJed Brown PetscInt p,npoints; 302c4762a1bSJed Brown PetscReal *array_x,*array_y,r2; 303c4762a1bSJed Brown PetscInt p2collect,*plist; 304c4762a1bSJed Brown PetscMPIInt rank; 305c4762a1bSJed Brown 306362febeeSStefano Zampini PetscFunctionBegin; 3079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 3089566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dm,&npoints)); 3099566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x)); 3109566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y)); 311c4762a1bSJed Brown 312c4762a1bSJed Brown r2 = zone->radius * zone->radius; 313c4762a1bSJed Brown p2collect = 0; 314c4762a1bSJed Brown for (p=0; p<npoints; p++) { 315c4762a1bSJed Brown PetscReal sep2; 316c4762a1bSJed Brown 317c4762a1bSJed Brown sep2 = (array_x[p] - zone->cx[0])*(array_x[p] - zone->cx[0]); 318c4762a1bSJed Brown sep2 += (array_y[p] - zone->cx[1])*(array_y[p] - zone->cx[1]); 319c4762a1bSJed Brown if (sep2 < r2) { 320c4762a1bSJed Brown p2collect++; 321c4762a1bSJed Brown } 322c4762a1bSJed Brown } 323c4762a1bSJed Brown 3249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(p2collect+1,&plist)); 325c4762a1bSJed Brown p2collect = 0; 326c4762a1bSJed Brown for (p=0; p<npoints; p++) { 327c4762a1bSJed Brown PetscReal sep2; 328c4762a1bSJed Brown 329c4762a1bSJed Brown sep2 = (array_x[p] - zone->cx[0])*(array_x[p] - zone->cx[0]); 330c4762a1bSJed Brown sep2 += (array_y[p] - zone->cx[1])*(array_y[p] - zone->cx[1]); 331c4762a1bSJed Brown if (sep2 < r2) { 332c4762a1bSJed Brown plist[p2collect] = p; 333c4762a1bSJed Brown p2collect++; 334c4762a1bSJed Brown } 335c4762a1bSJed Brown } 3369566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y)); 3379566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x)); 338c4762a1bSJed Brown 339c4762a1bSJed Brown *nfound = p2collect; 340c4762a1bSJed Brown *foundlist = plist; 341c4762a1bSJed Brown PetscFunctionReturn(0); 342c4762a1bSJed Brown } 343c4762a1bSJed Brown 344c4762a1bSJed Brown PetscErrorCode ex1_4(void) 345c4762a1bSJed Brown { 346c4762a1bSJed Brown DM dms; 347c4762a1bSJed Brown PetscMPIInt rank,size; 348c4762a1bSJed Brown PetscInt is,js,ni,nj,overlap,nn; 349c4762a1bSJed Brown DM dmcell; 350c4762a1bSJed Brown CollectZoneCtx *zone; 351c4762a1bSJed Brown PetscReal dx; 352c4762a1bSJed Brown 353362febeeSStefano Zampini PetscFunctionBegin; 3549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 3559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 356c4762a1bSJed Brown nn = 101; 357c4762a1bSJed Brown dx = 2.0/ (PetscReal)(nn-1); 358c4762a1bSJed Brown overlap = 0; 3599566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,nn,nn,PETSC_DECIDE,PETSC_DECIDE,1,overlap,NULL,NULL,&dmcell)); 3609566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dmcell)); 3619566063dSJacob Faibussowitsch PetscCall(DMSetUp(dmcell)); 3629566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dmcell,-1.0,1.0,-1.0,1.0,-1.0,1.0)); 3639566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dmcell,&is,&js,NULL,&ni,&nj,NULL)); 3649566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD,&dms)); 3659566063dSJacob Faibussowitsch PetscCall(DMSetType(dms,DMSWARM)); 366c4762a1bSJed Brown 367c4762a1bSJed Brown /* load in data types */ 3689566063dSJacob Faibussowitsch PetscCall(DMSwarmInitializeFieldRegister(dms)); 3699566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL)); 3709566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"coorx",1,PETSC_REAL)); 3719566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"coory",1,PETSC_REAL)); 3729566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(dms)); 3739566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dms,ni*nj*4,4)); 3749566063dSJacob Faibussowitsch PetscCall(DMView(dms,PETSC_VIEWER_STDOUT_WORLD)); 375c4762a1bSJed Brown 376c4762a1bSJed Brown /* set values within the swarm */ 377c4762a1bSJed Brown { 378c4762a1bSJed Brown PetscReal *array_x,*array_y; 379c4762a1bSJed Brown PetscInt npoints,i,j,cnt; 380c4762a1bSJed Brown DMDACoor2d **LA_coor; 381c4762a1bSJed Brown Vec coor; 382c4762a1bSJed Brown DM dmcellcdm; 383c4762a1bSJed Brown 3849566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dmcell,&dmcellcdm)); 3859566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dmcell,&coor)); 3869566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dmcellcdm,coor,&LA_coor)); 3879566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dms,&npoints)); 3889566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x)); 3899566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y)); 390c4762a1bSJed Brown cnt = 0; 391c4762a1bSJed Brown for (j=js; j<js+nj; j++) { 392c4762a1bSJed Brown for (i=is; i<is+ni; i++) { 393c4762a1bSJed Brown PetscReal xp,yp; 394c4762a1bSJed Brown 395c4762a1bSJed Brown xp = PetscRealPart(LA_coor[j][i].x); 396c4762a1bSJed Brown yp = PetscRealPart(LA_coor[j][i].y); 397c4762a1bSJed 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; }*/ 398c4762a1bSJed 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; }*/ 399c4762a1bSJed 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; }*/ 400c4762a1bSJed 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; }*/ 401c4762a1bSJed 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; }*/ 402c4762a1bSJed 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; }*/ 403c4762a1bSJed 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; }*/ 404c4762a1bSJed 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; }*/ 405c4762a1bSJed Brown cnt++; 406c4762a1bSJed Brown } 407c4762a1bSJed Brown } 4089566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y)); 4099566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x)); 4109566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dmcellcdm,coor,&LA_coor)); 411c4762a1bSJed Brown } 4129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1,&zone)); 413c4762a1bSJed Brown if (size == 4) { 414c4762a1bSJed Brown if (rank == 0) { 415c4762a1bSJed Brown zone->cx[0] = 0.5; 416c4762a1bSJed Brown zone->cx[1] = 0.5; 417c4762a1bSJed Brown zone->radius = 0.3; 418c4762a1bSJed Brown } 419c4762a1bSJed Brown if (rank == 1) { 420c4762a1bSJed Brown zone->cx[0] = -0.5; 421c4762a1bSJed Brown zone->cx[1] = 0.5; 422c4762a1bSJed Brown zone->radius = 0.25; 423c4762a1bSJed Brown } 424c4762a1bSJed Brown if (rank == 2) { 425c4762a1bSJed Brown zone->cx[0] = 0.5; 426c4762a1bSJed Brown zone->cx[1] = -0.5; 427c4762a1bSJed Brown zone->radius = 0.2; 428c4762a1bSJed Brown } 429c4762a1bSJed Brown if (rank == 3) { 430c4762a1bSJed Brown zone->cx[0] = -0.5; 431c4762a1bSJed Brown zone->cx[1] = -0.5; 432c4762a1bSJed Brown zone->radius = 0.1; 433c4762a1bSJed Brown } 434c4762a1bSJed Brown } else { 435c4762a1bSJed Brown if (rank == 0) { 436c4762a1bSJed Brown zone->cx[0] = 0.5; 437c4762a1bSJed Brown zone->cx[1] = 0.5; 438c4762a1bSJed Brown zone->radius = 0.8; 439c4762a1bSJed Brown } else { 440c4762a1bSJed Brown zone->cx[0] = 10.0; 441c4762a1bSJed Brown zone->cx[1] = 10.0; 442c4762a1bSJed Brown zone->radius = 0.0; 443c4762a1bSJed Brown } 444c4762a1bSJed Brown } 445c4762a1bSJed Brown { 446c4762a1bSJed Brown PetscInt npoints[2],npoints_orig[2],ng; 447c4762a1bSJed Brown 4489566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dms,&npoints_orig[0])); 4499566063dSJacob Faibussowitsch PetscCall(DMSwarmGetSize(dms,&npoints_orig[1])); 4509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 4519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] before(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1])); 4529566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 4539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 4549566063dSJacob Faibussowitsch PetscCall(DMSwarmCollect_General(dms,collect_zone,sizeof(CollectZoneCtx),zone,&ng)); 4559566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dms,&npoints[0])); 4569566063dSJacob Faibussowitsch PetscCall(DMSwarmGetSize(dms,&npoints[1])); 4579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 4589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after(%D,%D)\n",rank,npoints[0],npoints[1])); 4599566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 4609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 461c4762a1bSJed Brown } 462c4762a1bSJed Brown { 463c4762a1bSJed Brown PetscReal *array_x,*array_y; 464c4762a1bSJed Brown PetscInt npoints,p; 465c4762a1bSJed Brown FILE *fp = NULL; 466c4762a1bSJed Brown char name[PETSC_MAX_PATH_LEN]; 467c4762a1bSJed Brown 4689566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"c-rank%d.gp",rank)); 469c4762a1bSJed Brown fp = fopen(name,"w"); 47028b400f6SJacob Faibussowitsch PetscCheck(fp,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file %s",name); 4719566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dms,&npoints)); 4729566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x)); 4739566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y)); 474c4762a1bSJed Brown for (p=0; p<npoints; p++) { 475c4762a1bSJed Brown fprintf(fp,"%+1.4e %+1.4e %1.4e\n",array_x[p],array_y[p],(double)rank); 476c4762a1bSJed Brown } 4779566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y)); 4789566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x)); 479c4762a1bSJed Brown fclose(fp); 480c4762a1bSJed Brown } 4819566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmcell)); 4829566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dms)); 4839566063dSJacob Faibussowitsch PetscCall(PetscFree(zone)); 484c4762a1bSJed Brown PetscFunctionReturn(0); 485c4762a1bSJed Brown } 486c4762a1bSJed Brown 487c4762a1bSJed Brown int main(int argc,char **argv) 488c4762a1bSJed Brown { 489c4762a1bSJed Brown PetscInt test_mode = 4; 490c4762a1bSJed Brown 4919566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 4929566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-test_mode",&test_mode,NULL)); 493c4762a1bSJed Brown if (test_mode == 1) { 4949566063dSJacob Faibussowitsch PetscCall(ex1_1()); 495c4762a1bSJed Brown } else if (test_mode == 2) { 4969566063dSJacob Faibussowitsch PetscCall(ex1_2()); 497c4762a1bSJed Brown } else if (test_mode == 3) { 4989566063dSJacob Faibussowitsch PetscCall(ex1_3()); 499c4762a1bSJed Brown } else if (test_mode == 4) { 5009566063dSJacob Faibussowitsch PetscCall(ex1_4()); 501c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown test_mode value, should be 1,2,3,4"); 5029566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 503b122ec5aSJacob Faibussowitsch return 0; 504c4762a1bSJed Brown } 505c4762a1bSJed Brown 506c4762a1bSJed Brown /*TEST 507c4762a1bSJed Brown 508c4762a1bSJed Brown build: 509c4762a1bSJed Brown requires: !complex double 510c4762a1bSJed Brown 511c4762a1bSJed Brown test: 512c4762a1bSJed Brown args: -test_mode 1 513c4762a1bSJed Brown filter: grep -v atomic 514c4762a1bSJed Brown filter_output: grep -v atomic 515c4762a1bSJed Brown 516c4762a1bSJed Brown test: 517c4762a1bSJed Brown suffix: 2 518c4762a1bSJed Brown args: -test_mode 2 519c4762a1bSJed Brown filter: grep -v atomic 520c4762a1bSJed Brown filter_output: grep -v atomic 521c4762a1bSJed Brown 522c4762a1bSJed Brown test: 523c4762a1bSJed Brown suffix: 3 524c4762a1bSJed Brown args: -test_mode 3 525c4762a1bSJed Brown filter: grep -v atomic 526c4762a1bSJed Brown filter_output: grep -v atomic 527c4762a1bSJed Brown TODO: broken 528c4762a1bSJed Brown 529c4762a1bSJed Brown test: 530c4762a1bSJed Brown suffix: 4 531c4762a1bSJed Brown args: -test_mode 4 532c4762a1bSJed Brown filter: grep -v atomic 533c4762a1bSJed Brown filter_output: grep -v atomic 534c4762a1bSJed Brown 535c4762a1bSJed Brown test: 536c4762a1bSJed Brown suffix: 5 537c4762a1bSJed Brown nsize: 4 538c4762a1bSJed Brown args: -test_mode 1 539c4762a1bSJed Brown filter: grep -v atomic 540c4762a1bSJed Brown filter_output: grep -v atomic 541c4762a1bSJed Brown 542c4762a1bSJed Brown test: 543c4762a1bSJed Brown suffix: 6 544c4762a1bSJed Brown nsize: 4 545c4762a1bSJed Brown args: -test_mode 2 546c4762a1bSJed Brown filter: grep -v atomic 547c4762a1bSJed Brown filter_output: grep -v atomic 548c4762a1bSJed Brown 549c4762a1bSJed Brown test: 550c4762a1bSJed Brown suffix: 7 551c4762a1bSJed Brown nsize: 4 552c4762a1bSJed Brown args: -test_mode 3 553c4762a1bSJed Brown filter: grep -v atomic 554c4762a1bSJed Brown filter_output: grep -v atomic 555c4762a1bSJed Brown TODO: broken 556c4762a1bSJed Brown 557c4762a1bSJed Brown test: 558c4762a1bSJed Brown suffix: 8 559c4762a1bSJed Brown nsize: 4 560c4762a1bSJed Brown args: -test_mode 4 561c4762a1bSJed Brown filter: grep -v atomic 562c4762a1bSJed Brown filter_output: grep -v atomic 563c4762a1bSJed Brown 564c4762a1bSJed Brown TEST*/ 565