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 PetscFunctionBeginUser; 14d0c080abSJoseph Pusztay options->particlesPerCell = 1; 1530602db0SMatthew G. Knepley 16*d0609cedSBarry 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)); 18*d0609cedSBarry Smith PetscOptionsEnd(); 19d0c080abSJoseph Pusztay PetscFunctionReturn(0); 20d0c080abSJoseph Pusztay } 21d0c080abSJoseph Pusztay 22d0c080abSJoseph Pusztay static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user) 23d0c080abSJoseph Pusztay { 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")); 29d0c080abSJoseph Pusztay PetscFunctionReturn(0); 30d0c080abSJoseph Pusztay } 31d0c080abSJoseph Pusztay 32d0c080abSJoseph Pusztay static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user) 33d0c080abSJoseph Pusztay { 34d0c080abSJoseph Pusztay PetscInt *cellid; 35d0c080abSJoseph Pusztay PetscInt dim, cStart, cEnd, c, Np = user->particlesPerCell, p; 36d0c080abSJoseph Pusztay 37d0c080abSJoseph Pusztay PetscFunctionBeginUser; 389566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 399566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject) dm), sw)); 409566063dSJacob Faibussowitsch PetscCall(DMSetType(*sw, DMSWARM)); 419566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*sw, dim)); 429566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC)); 439566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(*sw, dm)); 449566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", 2, PETSC_REAL)); 459566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(*sw)); 469566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 479566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0)); 489566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*sw)); 499566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid)); 50d0c080abSJoseph Pusztay for (c = cStart; c < cEnd; ++c) { 51d0c080abSJoseph Pusztay for (p = 0; p < Np; ++p) { 52d0c080abSJoseph Pusztay const PetscInt n = c*Np + p; 53d0c080abSJoseph Pusztay 54d0c080abSJoseph Pusztay cellid[n] = c; 55d0c080abSJoseph Pusztay } 56d0c080abSJoseph Pusztay } 579566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid)); 589566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) *sw, "Particles")); 599566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view")); 60d0c080abSJoseph Pusztay PetscFunctionReturn(0); 61d0c080abSJoseph Pusztay } 62d0c080abSJoseph Pusztay 63d0c080abSJoseph Pusztay int main(int argc,char **argv) 64d0c080abSJoseph Pusztay { 65d0c080abSJoseph Pusztay DM dm, sw, cellsw; /* Mesh and particle managers */ 66d0c080abSJoseph Pusztay MPI_Comm comm; 67d0c080abSJoseph Pusztay AppCtx user; 68d0c080abSJoseph Pusztay 699566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 70d0c080abSJoseph Pusztay comm = PETSC_COMM_WORLD; 719566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &user)); 729566063dSJacob Faibussowitsch PetscCall(CreateMesh(comm, &dm, &user)); 739566063dSJacob Faibussowitsch PetscCall(CreateParticles(dm, &sw, &user)); 749566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(sw, &user)); 759566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, &cellsw)); 766a5217c0SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject) cellsw, "SubParticles")); 779566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellSwarm(sw, 1, cellsw)); 789566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(cellsw, NULL, "-subswarm_view")); 799566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreCellSwarm(sw, 1, cellsw)); 809566063dSJacob Faibussowitsch PetscCall(DMDestroy(&sw)); 819566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 829566063dSJacob Faibussowitsch PetscCall(DMDestroy(&cellsw)); 839566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 84b122ec5aSJacob Faibussowitsch return 0; 85d0c080abSJoseph Pusztay } 86d0c080abSJoseph Pusztay 87d0c080abSJoseph Pusztay /*TEST 88d0c080abSJoseph Pusztay build: 89d0c080abSJoseph Pusztay requires: triangle !single !complex 90d0c080abSJoseph Pusztay test: 91d0c080abSJoseph Pusztay suffix: 1 92d0c080abSJoseph Pusztay args: -particles_per_cell 2 -dm_plex_box_faces 2,1 -dm_view -sw_view -subswarm_view 93d0c080abSJoseph Pusztay TEST*/ 94