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
ProcessOptions(MPI_Comm comm,AppCtx * options)11d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
12d71ae5a4SJacob Faibussowitsch {
13d0c080abSJoseph Pusztay PetscFunctionBeginUser;
14d0c080abSJoseph Pusztay options->particlesPerCell = 1;
1530602db0SMatthew G. Knepley
16d0609cedSBarry Smith PetscOptionsBegin(comm, "", "CellSwarm Options", "DMSWARM");
179566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-particles_per_cell", "Number of particles per cell", "ex3.c", options->particlesPerCell, &options->particlesPerCell, NULL));
18d0609cedSBarry Smith PetscOptionsEnd();
193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
20d0c080abSJoseph Pusztay }
21d0c080abSJoseph Pusztay
CreateMesh(MPI_Comm comm,DM * dm,AppCtx * user)22d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
23d71ae5a4SJacob Faibussowitsch {
24d0c080abSJoseph Pusztay PetscFunctionBeginUser;
259566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm));
269566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX));
279566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm));
289566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
30d0c080abSJoseph Pusztay }
31d0c080abSJoseph Pusztay
CreateParticles(DM dm,DM * sw,AppCtx * user)32d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
33d71ae5a4SJacob Faibussowitsch {
34*19307e5cSMatthew G. Knepley DMSwarmCellDM celldm;
35*19307e5cSMatthew G. Knepley PetscInt *swarm_cellid;
36d0c080abSJoseph Pusztay PetscInt dim, cStart, cEnd, c, Np = user->particlesPerCell, p;
37*19307e5cSMatthew G. Knepley const char *cellid;
38d0c080abSJoseph Pusztay
39d0c080abSJoseph Pusztay PetscFunctionBeginUser;
409566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
419566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw));
429566063dSJacob Faibussowitsch PetscCall(DMSetType(*sw, DMSWARM));
439566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*sw, dim));
449566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
459566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(*sw, dm));
469566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", 2, PETSC_REAL));
479566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(*sw));
48*19307e5cSMatthew G. Knepley PetscCall(DMSwarmGetCellDMActive(*sw, &celldm));
49*19307e5cSMatthew G. Knepley PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid));
509566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
519566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0));
529566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*sw));
53*19307e5cSMatthew G. Knepley PetscCall(DMSwarmGetField(*sw, cellid, NULL, NULL, (void **)&swarm_cellid));
54d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) {
55d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) {
56d0c080abSJoseph Pusztay const PetscInt n = c * Np + p;
57d0c080abSJoseph Pusztay
58*19307e5cSMatthew G. Knepley swarm_cellid[n] = c;
59d0c080abSJoseph Pusztay }
60d0c080abSJoseph Pusztay }
61*19307e5cSMatthew G. Knepley PetscCall(DMSwarmRestoreField(*sw, cellid, NULL, NULL, (void **)&swarm_cellid));
629566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles"));
639566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
65d0c080abSJoseph Pusztay }
66d0c080abSJoseph Pusztay
main(int argc,char ** argv)67d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
68d71ae5a4SJacob Faibussowitsch {
69d0c080abSJoseph Pusztay DM dm, sw, cellsw; /* Mesh and particle managers */
70d0c080abSJoseph Pusztay MPI_Comm comm;
71d0c080abSJoseph Pusztay AppCtx user;
72d0c080abSJoseph Pusztay
73327415f7SBarry Smith PetscFunctionBeginUser;
749566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
75d0c080abSJoseph Pusztay comm = PETSC_COMM_WORLD;
769566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &user));
779566063dSJacob Faibussowitsch PetscCall(CreateMesh(comm, &dm, &user));
789566063dSJacob Faibussowitsch PetscCall(CreateParticles(dm, &sw, &user));
799566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(sw, &user));
809566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, &cellsw));
816a5217c0SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)cellsw, "SubParticles"));
829566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellSwarm(sw, 1, cellsw));
839566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(cellsw, NULL, "-subswarm_view"));
849566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreCellSwarm(sw, 1, cellsw));
859566063dSJacob Faibussowitsch PetscCall(DMDestroy(&sw));
869566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
879566063dSJacob Faibussowitsch PetscCall(DMDestroy(&cellsw));
889566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
89b122ec5aSJacob Faibussowitsch return 0;
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