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