xref: /petsc/src/dm/impls/swarm/tests/ex3.c (revision d0c080ab42685f06eb6648c3b073f69b5fa02d59)
1*d0c080abSJoseph 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";
2*d0c080abSJoseph Pusztay 
3*d0c080abSJoseph Pusztay #include <petscdmplex.h>
4*d0c080abSJoseph Pusztay #include <petscdmswarm.h>
5*d0c080abSJoseph Pusztay #include <petscts.h>
6*d0c080abSJoseph Pusztay 
7*d0c080abSJoseph Pusztay typedef struct {
8*d0c080abSJoseph Pusztay   PetscInt  dim;                          /* The topological mesh dimension */
9*d0c080abSJoseph Pusztay   PetscBool simplex;                      /* Flag for simplices or tensor cells */
10*d0c080abSJoseph Pusztay   char      filename[PETSC_MAX_PATH_LEN]; /* Name of the mesh filename if any */
11*d0c080abSJoseph Pusztay   PetscInt  particlesPerCell;             /* The number of partices per cell */
12*d0c080abSJoseph Pusztay } AppCtx;
13*d0c080abSJoseph Pusztay 
14*d0c080abSJoseph Pusztay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
15*d0c080abSJoseph Pusztay {
16*d0c080abSJoseph Pusztay   PetscErrorCode ierr;
17*d0c080abSJoseph Pusztay 
18*d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
19*d0c080abSJoseph Pusztay   options->dim              = 2;
20*d0c080abSJoseph Pusztay   options->simplex          = PETSC_TRUE;
21*d0c080abSJoseph Pusztay   options->particlesPerCell = 1;
22*d0c080abSJoseph Pusztay   ierr = PetscStrcpy(options->filename, "");CHKERRQ(ierr);
23*d0c080abSJoseph Pusztay   ierr = PetscOptionsBegin(comm, "", "CellSwarm Options", "DMSWARM");CHKERRQ(ierr);
24*d0c080abSJoseph Pusztay   ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex3.c", options->dim, &options->dim, NULL);CHKERRQ(ierr);
25*d0c080abSJoseph Pusztay   ierr = PetscOptionsBool("-simplex", "The flag for simplices or tensor cells", "ex3.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
26*d0c080abSJoseph Pusztay   ierr = PetscOptionsString("-mesh", "Name of the mesh filename if any", "ex3.c", options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr);
27*d0c080abSJoseph Pusztay   ierr = PetscOptionsInt("-particles_per_cell", "Number of particles per cell", "ex3.c", options->particlesPerCell, &options->particlesPerCell, NULL);CHKERRQ(ierr);
28*d0c080abSJoseph Pusztay   ierr = PetscOptionsEnd();CHKERRQ(ierr);
29*d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
30*d0c080abSJoseph Pusztay }
31*d0c080abSJoseph Pusztay 
32*d0c080abSJoseph Pusztay static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
33*d0c080abSJoseph Pusztay {
34*d0c080abSJoseph Pusztay   PetscBool      flg;
35*d0c080abSJoseph Pusztay   PetscErrorCode ierr;
36*d0c080abSJoseph Pusztay 
37*d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
38*d0c080abSJoseph Pusztay   ierr = PetscStrcmp(user->filename, "", &flg);CHKERRQ(ierr);
39*d0c080abSJoseph Pusztay   if (flg) {
40*d0c080abSJoseph Pusztay     ierr = DMPlexCreateBoxMesh(comm, user->dim, user->simplex, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
41*d0c080abSJoseph Pusztay   } else {
42*d0c080abSJoseph Pusztay     ierr = DMPlexCreateFromFile(comm, user->filename, PETSC_TRUE, dm);CHKERRQ(ierr);
43*d0c080abSJoseph Pusztay     ierr = DMGetDimension(*dm, &user->dim);CHKERRQ(ierr);
44*d0c080abSJoseph Pusztay   }
45*d0c080abSJoseph Pusztay   {
46*d0c080abSJoseph Pusztay     DM distributedMesh = NULL;
47*d0c080abSJoseph Pusztay 
48*d0c080abSJoseph Pusztay     ierr = DMPlexDistribute(*dm, 0, NULL, &distributedMesh);CHKERRQ(ierr);
49*d0c080abSJoseph Pusztay     if (distributedMesh) {
50*d0c080abSJoseph Pusztay       ierr = DMDestroy(dm);CHKERRQ(ierr);
51*d0c080abSJoseph Pusztay       *dm  = distributedMesh;
52*d0c080abSJoseph Pusztay     }
53*d0c080abSJoseph Pusztay   }
54*d0c080abSJoseph Pusztay   ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr); /* needed for periodic */
55*d0c080abSJoseph Pusztay   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
56*d0c080abSJoseph Pusztay   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
57*d0c080abSJoseph Pusztay   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
58*d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
59*d0c080abSJoseph Pusztay }
60*d0c080abSJoseph Pusztay 
61*d0c080abSJoseph Pusztay static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
62*d0c080abSJoseph Pusztay {
63*d0c080abSJoseph Pusztay   PetscInt      *cellid;
64*d0c080abSJoseph Pusztay   PetscInt       dim, cStart, cEnd, c, Np = user->particlesPerCell, p;
65*d0c080abSJoseph Pusztay   PetscErrorCode ierr;
66*d0c080abSJoseph Pusztay 
67*d0c080abSJoseph Pusztay   PetscFunctionBeginUser;
68*d0c080abSJoseph Pusztay   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
69*d0c080abSJoseph Pusztay   ierr = DMCreate(PetscObjectComm((PetscObject) dm), sw);CHKERRQ(ierr);
70*d0c080abSJoseph Pusztay   ierr = DMSetType(*sw, DMSWARM);CHKERRQ(ierr);
71*d0c080abSJoseph Pusztay   ierr = DMSetDimension(*sw, dim);CHKERRQ(ierr);
72*d0c080abSJoseph Pusztay   ierr = DMSwarmSetType(*sw, DMSWARM_PIC);CHKERRQ(ierr);
73*d0c080abSJoseph Pusztay   ierr = DMSwarmSetCellDM(*sw, dm);CHKERRQ(ierr);
74*d0c080abSJoseph Pusztay   ierr = DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", 2, PETSC_REAL);CHKERRQ(ierr);
75*d0c080abSJoseph Pusztay   ierr = DMSwarmFinalizeFieldRegister(*sw);CHKERRQ(ierr);
76*d0c080abSJoseph Pusztay   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
77*d0c080abSJoseph Pusztay   ierr = DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0);CHKERRQ(ierr);
78*d0c080abSJoseph Pusztay   ierr = DMSetFromOptions(*sw);CHKERRQ(ierr);
79*d0c080abSJoseph Pusztay   ierr = DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr);
80*d0c080abSJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
81*d0c080abSJoseph Pusztay     for (p = 0; p < Np; ++p) {
82*d0c080abSJoseph Pusztay       const PetscInt n = c*Np + p;
83*d0c080abSJoseph Pusztay 
84*d0c080abSJoseph Pusztay       cellid[n] = c;
85*d0c080abSJoseph Pusztay     }
86*d0c080abSJoseph Pusztay   }
87*d0c080abSJoseph Pusztay   ierr = DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr);
88*d0c080abSJoseph Pusztay   ierr = PetscObjectSetName((PetscObject) *sw, "Particles");CHKERRQ(ierr);
89*d0c080abSJoseph Pusztay   ierr = DMViewFromOptions(*sw, NULL, "-sw_view");CHKERRQ(ierr);
90*d0c080abSJoseph Pusztay   PetscFunctionReturn(0);
91*d0c080abSJoseph Pusztay }
92*d0c080abSJoseph Pusztay 
93*d0c080abSJoseph Pusztay int main(int argc,char **argv)
94*d0c080abSJoseph Pusztay {
95*d0c080abSJoseph Pusztay   DM             dm, sw, cellsw; /* Mesh and particle managers */
96*d0c080abSJoseph Pusztay   MPI_Comm       comm;
97*d0c080abSJoseph Pusztay   AppCtx         user;
98*d0c080abSJoseph Pusztay   PetscErrorCode ierr;
99*d0c080abSJoseph Pusztay 
100*d0c080abSJoseph Pusztay   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
101*d0c080abSJoseph Pusztay   comm = PETSC_COMM_WORLD;
102*d0c080abSJoseph Pusztay   ierr = ProcessOptions(comm, &user);CHKERRQ(ierr);
103*d0c080abSJoseph Pusztay   ierr = CreateMesh(comm, &dm, &user);CHKERRQ(ierr);
104*d0c080abSJoseph Pusztay   ierr = CreateParticles(dm, &sw, &user);CHKERRQ(ierr);
105*d0c080abSJoseph Pusztay   ierr = DMSetApplicationContext(sw, &user);CHKERRQ(ierr);
106*d0c080abSJoseph Pusztay   ierr = DMCreate(comm, &cellsw);CHKERRQ(ierr);
107*d0c080abSJoseph Pusztay   ierr = DMSwarmGetCellSwarm(sw, 1, cellsw);CHKERRQ(ierr);
108*d0c080abSJoseph Pusztay   ierr = DMViewFromOptions(cellsw, NULL, "-subswarm_view");CHKERRQ(ierr);
109*d0c080abSJoseph Pusztay   ierr = DMSwarmRestoreCellSwarm(sw, 1, cellsw);CHKERRQ(ierr);
110*d0c080abSJoseph Pusztay   ierr = DMDestroy(&sw);CHKERRQ(ierr);
111*d0c080abSJoseph Pusztay   ierr = DMDestroy(&dm);CHKERRQ(ierr);
112*d0c080abSJoseph Pusztay   ierr = DMDestroy(&cellsw);
113*d0c080abSJoseph Pusztay   ierr = PetscFinalize();
114*d0c080abSJoseph Pusztay   return ierr;
115*d0c080abSJoseph Pusztay }
116*d0c080abSJoseph Pusztay 
117*d0c080abSJoseph Pusztay /*TEST
118*d0c080abSJoseph Pusztay    build:
119*d0c080abSJoseph Pusztay      requires: triangle !single !complex
120*d0c080abSJoseph Pusztay    test:
121*d0c080abSJoseph Pusztay      suffix: 1
122*d0c080abSJoseph Pusztay      args: -particles_per_cell 2 -dm_plex_box_faces 2,1 -dm_view -sw_view -subswarm_view
123*d0c080abSJoseph Pusztay TEST*/
124