1 static char help[] = "Test point location in DA using DMSwarm\n"; 2 3 #include <petscdmda.h> 4 #include <petscdmswarm.h> 5 6 PetscErrorCode DMSwarmPrint(DM sw) 7 { 8 PetscReal *array; 9 PetscInt *pidArray, *cidArray; 10 PetscInt Np, bs; 11 PetscMPIInt rank; 12 13 PetscFunctionBeginUser; 14 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sw), &rank)); 15 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 16 PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, &bs, NULL, (void **)&array)); 17 PetscCall(DMSwarmGetField(sw, DMSwarmField_pid, &bs, NULL, (void **)&pidArray)); 18 PetscCall(DMSwarmGetField(sw, DMSwarmPICField_cellid, &bs, NULL, (void **)&cidArray)); 19 for (PetscInt p = 0; p < Np; ++p) { 20 const PetscReal th = PetscAtan2Real(array[2 * p + 1], array[2 * p]) / PETSC_PI; 21 const PetscReal r = PetscSqrtReal(array[2 * p + 1] * array[2 * p + 1] + array[2 * p] * array[2 * p]); 22 PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] p %" PetscInt_FMT " (%+1.4f,%+1.4f) r=%1.2f th=%1.3f*pi cellid=%" PetscInt_FMT "\n", rank, pidArray[p], (double)array[2 * p], (double)array[2 * p + 1], (double)r, (double)th, cidArray[p])); 23 } 24 PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, &bs, NULL, (void **)&array)); 25 PetscCall(DMSwarmRestoreField(sw, DMSwarmField_pid, &bs, NULL, (void **)&pidArray)); 26 PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_cellid, &bs, NULL, (void **)&pidArray)); 27 PetscFunctionReturn(PETSC_SUCCESS); 28 } 29 30 int main(int argc, char **argv) 31 { 32 DM dm, sw; 33 PetscDataType dtype; 34 PetscReal *coords, r, dr; 35 PetscInt Nx = 4, Ny = 4, Np = 8, bs; 36 PetscMPIInt rank, size; 37 38 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 39 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 40 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 41 42 PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_GHOSTED, DMDA_STENCIL_BOX, Nx, Ny, PETSC_DECIDE, PETSC_DECIDE, 1, 2, NULL, NULL, &dm)); 43 PetscCall(DMSetFromOptions(dm)); 44 PetscCall(DMSetUp(dm)); 45 PetscCall(DMDASetUniformCoordinates(dm, -1., 1., -1., 1., -1., 1.)); 46 PetscCall(DMViewFromOptions(dm, NULL, "-da_view")); 47 48 PetscCall(DMCreate(PETSC_COMM_WORLD, &sw)); 49 PetscCall(PetscObjectSetName((PetscObject)sw, "Particle Grid")); 50 PetscCall(DMSetType(sw, DMSWARM)); 51 PetscCall(DMSetDimension(sw, 2)); 52 PetscCall(DMSwarmSetType(sw, DMSWARM_PIC)); 53 PetscCall(DMSetFromOptions(sw)); 54 PetscCall(DMSwarmSetCellDM(sw, dm)); 55 PetscCall(DMSwarmInitializeFieldRegister(sw)); 56 PetscCall(DMSwarmRegisterPetscDatatypeField(sw, "u", 1, PETSC_SCALAR)); 57 PetscCall(DMSwarmFinalizeFieldRegister(sw)); 58 PetscCall(DMSwarmSetLocalSizes(sw, Np, 2)); 59 PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords)); 60 dr = 1.0 / (size + 1); 61 r = (rank + 1) * dr; 62 for (PetscInt p = 0; p < Np; ++p) { 63 const PetscReal th = (p + 0.5) * 2. * PETSC_PI / Np; 64 65 coords[p * 2 + 0] = r * PetscCosReal(th); 66 coords[p * 2 + 1] = r * PetscSinReal(th); 67 } 68 PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords)); 69 PetscCall(DMViewFromOptions(sw, NULL, "-swarm_view")); 70 PetscCall(DMSwarmPrint(sw)); 71 72 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n... calling DMSwarmMigrate ...\n\n")); 73 PetscCall(DMSwarmMigrate(sw, PETSC_TRUE)); 74 PetscCall(DMViewFromOptions(sw, NULL, "-swarm_view")); 75 PetscCall(DMSwarmPrint(sw)); 76 77 PetscCall(DMDestroy(&sw)); 78 PetscCall(DMDestroy(&dm)); 79 PetscCall(PetscFinalize()); 80 return 0; 81 } 82 83 /*TEST 84 85 # Swarm does not handle complex or quad 86 build: 87 requires: !complex double 88 89 test: 90 suffix: 0 91 92 TEST*/ 93