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