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