static char help[] = "DMSwarm-PIC demonstator of advecting points within cell DM defined by a DA or PLEX object \n\ Options: \n\ -ppcell : Number of times to sub-divide the reference cell when layout the initial particle coordinates \n\ -meshtype : 0 ==> DA , 1 ==> PLEX \n\ -nt : Number of timestep to perform \n\ -view : Write out initial condition and time dependent data \n"; #include #include #include #include PetscErrorCode pic_advect(PetscInt ppcell,PetscInt meshtype) { PetscErrorCode ierr; const PetscInt dim = 2; DM celldm,swarm; PetscInt tk,nt = 200; PetscBool view = PETSC_FALSE; Vec *pfields; PetscReal minradius; PetscReal dt; PetscReal vel[] = { 1.0, 0.16 }; const char *fieldnames[] = { "phi" }; PetscViewer viewer; PetscFunctionBegin; ierr = PetscOptionsGetBool(NULL,NULL,"-view",&view,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-nt",&nt,NULL);CHKERRQ(ierr); /* Create the background cell DM */ if (meshtype == 0) { /* DA */ PetscInt nxy; PetscInt dof = 1; PetscInt stencil_width = 1; ierr = PetscPrintf(PETSC_COMM_WORLD,"Mesh type: DMDA\n");CHKERRQ(ierr); nxy = 33; ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,nxy,nxy,PETSC_DECIDE,PETSC_DECIDE,dof,stencil_width,NULL,NULL,&celldm);CHKERRQ(ierr); ierr = DMDASetElementType(celldm,DMDA_ELEMENT_Q1);CHKERRQ(ierr); ierr = DMSetFromOptions(celldm);CHKERRQ(ierr); ierr = DMSetUp(celldm);CHKERRQ(ierr); ierr = DMDASetUniformCoordinates(celldm,0.0,1.0,0.0,1.0,0.0,1.5);CHKERRQ(ierr); minradius = 1.0/((PetscReal)(nxy-1)); ierr = PetscPrintf(PETSC_COMM_WORLD,"DA(minradius) %1.4e\n",(double)minradius);CHKERRQ(ierr); } if (meshtype == 1){ /* PLEX */ DM distributedMesh = NULL; PetscInt numComp[] = {1}; PetscInt numDof[] = {1,0,0}; /* vert, edge, cell */ PetscInt faces[] = {1,1,1}; PetscInt numBC = 0; PetscSection section; Vec cellgeom = NULL; Vec facegeom = NULL; ierr = PetscPrintf(PETSC_COMM_WORLD,"Mesh type: DMPLEX\n");CHKERRQ(ierr); ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_TRUE, faces, NULL, NULL, PETSC_TRUE, &celldm);CHKERRQ(ierr); /* Distribute mesh over processes */ ierr = DMPlexDistribute(celldm,0,NULL,&distributedMesh);CHKERRQ(ierr); if (distributedMesh) { ierr = DMDestroy(&celldm);CHKERRQ(ierr); celldm = distributedMesh; } ierr = DMSetFromOptions(celldm);CHKERRQ(ierr); ierr = DMPlexCreateSection(celldm,NULL,numComp,numDof,numBC,NULL,NULL,NULL,NULL,§ion);CHKERRQ(ierr); ierr = DMSetLocalSection(celldm,section);CHKERRQ(ierr); ierr = DMSetUp(celldm);CHKERRQ(ierr); /* Calling DMPlexComputeGeometryFVM() generates the value returned by DMPlexGetMinRadius() */ ierr = DMPlexComputeGeometryFVM(celldm,&cellgeom,&facegeom);CHKERRQ(ierr); ierr = DMPlexGetMinRadius(celldm,&minradius);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_WORLD,"PLEX(minradius) %1.4e\n",(double)minradius);CHKERRQ(ierr); ierr = VecDestroy(&cellgeom);CHKERRQ(ierr); ierr = VecDestroy(&facegeom);CHKERRQ(ierr); ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); } /* Create the DMSwarm */ ierr = DMCreate(PETSC_COMM_WORLD,&swarm);CHKERRQ(ierr); ierr = DMSetType(swarm,DMSWARM);CHKERRQ(ierr); ierr = DMSetDimension(swarm,dim);CHKERRQ(ierr); /* Configure swarm to be of type PIC */ ierr = DMSwarmSetType(swarm,DMSWARM_PIC);CHKERRQ(ierr); ierr = DMSwarmSetCellDM(swarm,celldm);CHKERRQ(ierr); /* Register two scalar fields within the DMSwarm */ ierr = DMSwarmRegisterPetscDatatypeField(swarm,"phi",1,PETSC_REAL);CHKERRQ(ierr); ierr = DMSwarmRegisterPetscDatatypeField(swarm,"region",1,PETSC_REAL);CHKERRQ(ierr); ierr = DMSwarmFinalizeFieldRegister(swarm);CHKERRQ(ierr); /* Set initial local sizes of the DMSwarm with a buffer length of zero */ ierr = DMSwarmSetLocalSizes(swarm,4,0);CHKERRQ(ierr); /* Insert swarm coordinates cell-wise */ /*ierr = DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_REGULAR,ppcell);CHKERRQ(ierr);*/ ierr = DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_SUBDIVISION,ppcell);CHKERRQ(ierr); /* Define initial conditions for th swarm fields "phi" and "region" */ { PetscReal *s_coor,*s_phi,*s_region; PetscInt npoints,p; ierr = DMSwarmGetLocalSize(swarm,&npoints);CHKERRQ(ierr); ierr = DMSwarmGetField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&s_coor);CHKERRQ(ierr); ierr = DMSwarmGetField(swarm,"phi",NULL,NULL,(void**)&s_phi);CHKERRQ(ierr); ierr = DMSwarmGetField(swarm,"region",NULL,NULL,(void**)&s_region);CHKERRQ(ierr); for (p=0; p 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"-meshtype must be 0 or 1"); ierr = pic_advect(ppcell,meshtype);CHKERRQ(ierr); ierr = PetscFinalize(); return ierr; }