xref: /petsc/src/dm/impls/swarm/tests/ex3.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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;
1730602db0SMatthew G. Knepley 
18d0c080abSJoseph Pusztay   ierr = PetscOptionsBegin(comm, "", "CellSwarm Options", "DMSWARM");CHKERRQ(ierr);
195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-particles_per_cell", "Number of particles per cell", "ex3.c", options->particlesPerCell, &options->particlesPerCell, NULL));
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   PetscFunctionBeginUser;
275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, dm));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm, DMPLEX));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view"));
31d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
32d0c080abSJoseph Pusztay }
33d0c080abSJoseph Pusztay 
34d0c080abSJoseph Pusztay static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
35d0c080abSJoseph Pusztay {
36d0c080abSJoseph Pusztay   PetscInt      *cellid;
37d0c080abSJoseph Pusztay   PetscInt       dim, cStart, cEnd, c, Np = user->particlesPerCell, p;
38d0c080abSJoseph Pusztay 
39d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
405f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(PetscObjectComm((PetscObject) dm), sw));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*sw, DMSWARM));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetDimension(*sw, dim));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetType(*sw, DMSWARM_PIC));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetCellDM(*sw, dm));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", 2, PETSC_REAL));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmFinalizeFieldRegister(*sw));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*sw));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid));
52d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
53d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
54d0c080abSJoseph Pusztay       const PetscInt n = c*Np + p;
55d0c080abSJoseph Pusztay 
56d0c080abSJoseph Pusztay       cellid[n] = c;
57d0c080abSJoseph Pusztay     }
58d0c080abSJoseph Pusztay   }
595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) *sw, "Particles"));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*sw, NULL, "-sw_view"));
62d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
63d0c080abSJoseph Pusztay }
64d0c080abSJoseph Pusztay 
65d0c080abSJoseph Pusztay int main(int argc,char **argv)
66d0c080abSJoseph Pusztay {
67d0c080abSJoseph Pusztay   DM             dm, sw, cellsw; /* Mesh and particle managers */
68d0c080abSJoseph Pusztay   MPI_Comm       comm;
69d0c080abSJoseph Pusztay   AppCtx         user;
70d0c080abSJoseph Pusztay 
71*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL, help));
72d0c080abSJoseph Pusztay   comm = PETSC_COMM_WORLD;
735f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(comm, &user));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(comm, &dm, &user));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateParticles(dm, &sw, &user));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetApplicationContext(sw, &user));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, &cellsw));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetCellSwarm(sw, 1, cellsw));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(cellsw, NULL, "-subswarm_view"));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreCellSwarm(sw, 1, cellsw));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&sw));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&cellsw));
84*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
85*b122ec5aSJacob Faibussowitsch   return 0;
86d0c080abSJoseph Pusztay }
87d0c080abSJoseph Pusztay 
88d0c080abSJoseph Pusztay /*TEST
89d0c080abSJoseph Pusztay    build:
90d0c080abSJoseph Pusztay      requires: triangle !single !complex
91d0c080abSJoseph Pusztay    test:
92d0c080abSJoseph Pusztay      suffix: 1
93d0c080abSJoseph Pusztay      args: -particles_per_cell 2 -dm_plex_box_faces 2,1 -dm_view -sw_view -subswarm_view
94d0c080abSJoseph Pusztay TEST*/
95