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); 19*5f80ce2aSJacob 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; 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 30*5f80ce2aSJacob 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; 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(PetscObjectComm((PetscObject) dm), sw)); 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*sw, DMSWARM)); 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(*sw, dim)); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetType(*sw, DMSWARM_PIC)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetCellDM(*sw, dm)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", 2, PETSC_REAL)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmFinalizeFieldRegister(*sw)); 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0)); 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*sw)); 51*5f80ce2aSJacob 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 } 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) *sw, "Particles")); 61*5f80ce2aSJacob 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 PetscErrorCode ierr; 71d0c080abSJoseph Pusztay 72d0c080abSJoseph Pusztay ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 73d0c080abSJoseph Pusztay comm = PETSC_COMM_WORLD; 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(comm, &user)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(comm, &dm, &user)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateParticles(dm, &sw, &user)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetApplicationContext(sw, &user)); 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, &cellsw)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetCellSwarm(sw, 1, cellsw)); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(cellsw, NULL, "-subswarm_view")); 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreCellSwarm(sw, 1, cellsw)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&sw)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&cellsw)); 85d0c080abSJoseph Pusztay ierr = PetscFinalize(); 86d0c080abSJoseph Pusztay return ierr; 87d0c080abSJoseph Pusztay } 88d0c080abSJoseph Pusztay 89d0c080abSJoseph Pusztay /*TEST 90d0c080abSJoseph Pusztay build: 91d0c080abSJoseph Pusztay requires: triangle !single !complex 92d0c080abSJoseph Pusztay test: 93d0c080abSJoseph Pusztay suffix: 1 94d0c080abSJoseph Pusztay args: -particles_per_cell 2 -dm_plex_box_faces 2,1 -dm_view -sw_view -subswarm_view 95d0c080abSJoseph Pusztay TEST*/ 96