xref: /petsc/src/dm/impls/swarm/tests/ex3.c (revision 30602db00db74b7e41a0c75e517aefe6711423f0)
1d0c080abSJoseph Pusztay static char help[] = "Example usage of extracting single cells with their associated fields from a swarm and putting it in a new swarm object\n";
2d0c080abSJoseph Pusztay 
3d0c080abSJoseph Pusztay #include <petscdmplex.h>
4d0c080abSJoseph Pusztay #include <petscdmswarm.h>
5d0c080abSJoseph Pusztay #include <petscts.h>
6d0c080abSJoseph Pusztay 
7d0c080abSJoseph Pusztay typedef struct {
8d0c080abSJoseph Pusztay   PetscInt particlesPerCell; /* The number of partices per cell */
9d0c080abSJoseph Pusztay } AppCtx;
10d0c080abSJoseph Pusztay 
11d0c080abSJoseph Pusztay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
12d0c080abSJoseph Pusztay {
13d0c080abSJoseph Pusztay   PetscErrorCode ierr;
14d0c080abSJoseph Pusztay 
15d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
16d0c080abSJoseph Pusztay   options->particlesPerCell = 1;
17*30602db0SMatthew G. Knepley 
18d0c080abSJoseph Pusztay   ierr = PetscOptionsBegin(comm, "", "CellSwarm Options", "DMSWARM");CHKERRQ(ierr);
19d0c080abSJoseph Pusztay   ierr = PetscOptionsInt("-particles_per_cell", "Number of particles per cell", "ex3.c", options->particlesPerCell, &options->particlesPerCell, NULL);CHKERRQ(ierr);
20d0c080abSJoseph Pusztay   ierr = PetscOptionsEnd();CHKERRQ(ierr);
21d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
22d0c080abSJoseph Pusztay }
23d0c080abSJoseph Pusztay 
24d0c080abSJoseph Pusztay static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
25d0c080abSJoseph Pusztay {
26d0c080abSJoseph Pusztay   PetscErrorCode ierr;
27d0c080abSJoseph Pusztay 
28d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
29*30602db0SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
30*30602db0SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
31d0c080abSJoseph Pusztay   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
32d0c080abSJoseph Pusztay   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
33d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
34d0c080abSJoseph Pusztay }
35d0c080abSJoseph Pusztay 
36d0c080abSJoseph Pusztay static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
37d0c080abSJoseph Pusztay {
38d0c080abSJoseph Pusztay   PetscInt      *cellid;
39d0c080abSJoseph Pusztay   PetscInt       dim, cStart, cEnd, c, Np = user->particlesPerCell, p;
40d0c080abSJoseph Pusztay   PetscErrorCode ierr;
41d0c080abSJoseph Pusztay 
42d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
43d0c080abSJoseph Pusztay   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
44d0c080abSJoseph Pusztay   ierr = DMCreate(PetscObjectComm((PetscObject) dm), sw);CHKERRQ(ierr);
45d0c080abSJoseph Pusztay   ierr = DMSetType(*sw, DMSWARM);CHKERRQ(ierr);
46d0c080abSJoseph Pusztay   ierr = DMSetDimension(*sw, dim);CHKERRQ(ierr);
47d0c080abSJoseph Pusztay   ierr = DMSwarmSetType(*sw, DMSWARM_PIC);CHKERRQ(ierr);
48d0c080abSJoseph Pusztay   ierr = DMSwarmSetCellDM(*sw, dm);CHKERRQ(ierr);
49d0c080abSJoseph Pusztay   ierr = DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", 2, PETSC_REAL);CHKERRQ(ierr);
50d0c080abSJoseph Pusztay   ierr = DMSwarmFinalizeFieldRegister(*sw);CHKERRQ(ierr);
51d0c080abSJoseph Pusztay   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
52d0c080abSJoseph Pusztay   ierr = DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0);CHKERRQ(ierr);
53d0c080abSJoseph Pusztay   ierr = DMSetFromOptions(*sw);CHKERRQ(ierr);
54d0c080abSJoseph Pusztay   ierr = DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr);
55d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
56d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
57d0c080abSJoseph Pusztay       const PetscInt n = c*Np + p;
58d0c080abSJoseph Pusztay 
59d0c080abSJoseph Pusztay       cellid[n] = c;
60d0c080abSJoseph Pusztay     }
61d0c080abSJoseph Pusztay   }
62d0c080abSJoseph Pusztay   ierr = DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr);
63d0c080abSJoseph Pusztay   ierr = PetscObjectSetName((PetscObject) *sw, "Particles");CHKERRQ(ierr);
64d0c080abSJoseph Pusztay   ierr = DMViewFromOptions(*sw, NULL, "-sw_view");CHKERRQ(ierr);
65d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
66d0c080abSJoseph Pusztay }
67d0c080abSJoseph Pusztay 
68d0c080abSJoseph Pusztay int main(int argc,char **argv)
69d0c080abSJoseph Pusztay {
70d0c080abSJoseph Pusztay   DM             dm, sw, cellsw; /* Mesh and particle managers */
71d0c080abSJoseph Pusztay   MPI_Comm       comm;
72d0c080abSJoseph Pusztay   AppCtx         user;
73d0c080abSJoseph Pusztay   PetscErrorCode ierr;
74d0c080abSJoseph Pusztay 
75d0c080abSJoseph Pusztay   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
76d0c080abSJoseph Pusztay   comm = PETSC_COMM_WORLD;
77d0c080abSJoseph Pusztay   ierr = ProcessOptions(comm, &user);CHKERRQ(ierr);
78d0c080abSJoseph Pusztay   ierr = CreateMesh(comm, &dm, &user);CHKERRQ(ierr);
79d0c080abSJoseph Pusztay   ierr = CreateParticles(dm, &sw, &user);CHKERRQ(ierr);
80d0c080abSJoseph Pusztay   ierr = DMSetApplicationContext(sw, &user);CHKERRQ(ierr);
81d0c080abSJoseph Pusztay   ierr = DMCreate(comm, &cellsw);CHKERRQ(ierr);
82d0c080abSJoseph Pusztay   ierr = DMSwarmGetCellSwarm(sw, 1, cellsw);CHKERRQ(ierr);
83d0c080abSJoseph Pusztay   ierr = DMViewFromOptions(cellsw, NULL, "-subswarm_view");CHKERRQ(ierr);
84d0c080abSJoseph Pusztay   ierr = DMSwarmRestoreCellSwarm(sw, 1, cellsw);CHKERRQ(ierr);
85d0c080abSJoseph Pusztay   ierr = DMDestroy(&sw);CHKERRQ(ierr);
86d0c080abSJoseph Pusztay   ierr = DMDestroy(&dm);CHKERRQ(ierr);
87d0c080abSJoseph Pusztay   ierr = DMDestroy(&cellsw);
88d0c080abSJoseph Pusztay   ierr = PetscFinalize();
89d0c080abSJoseph Pusztay   return ierr;
90d0c080abSJoseph Pusztay }
91d0c080abSJoseph Pusztay 
92d0c080abSJoseph Pusztay /*TEST
93d0c080abSJoseph Pusztay    build:
94d0c080abSJoseph Pusztay      requires: triangle !single !complex
95d0c080abSJoseph Pusztay    test:
96d0c080abSJoseph Pusztay      suffix: 1
97d0c080abSJoseph Pusztay      args: -particles_per_cell 2 -dm_plex_box_faces 2,1 -dm_view -sw_view -subswarm_view
98d0c080abSJoseph Pusztay TEST*/
99