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) { 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; PetscCall(PetscOptionsGetBool(NULL,NULL,"-view",&view,NULL)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-nt",&nt,NULL)); /* Create the background cell DM */ if (meshtype == 0) { /* DA */ PetscInt nxy; PetscInt dof = 1; PetscInt stencil_width = 1; PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Mesh type: DMDA\n")); nxy = 33; PetscCall(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)); PetscCall(DMDASetElementType(celldm,DMDA_ELEMENT_Q1)); PetscCall(DMSetFromOptions(celldm)); PetscCall(DMSetUp(celldm)); PetscCall(DMDASetUniformCoordinates(celldm,0.0,1.0,0.0,1.0,0.0,1.5)); minradius = 1.0/((PetscReal)(nxy-1)); PetscCall(PetscPrintf(PETSC_COMM_WORLD,"DA(minradius) %1.4e\n",(double)minradius)); } 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; PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Mesh type: DMPLEX\n")); PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_TRUE, faces, NULL, NULL, PETSC_TRUE, &celldm)); /* Distribute mesh over processes */ PetscCall(DMPlexDistribute(celldm,0,NULL,&distributedMesh)); if (distributedMesh) { PetscCall(DMDestroy(&celldm)); celldm = distributedMesh; } PetscCall(DMSetFromOptions(celldm)); PetscCall(DMPlexCreateSection(celldm,NULL,numComp,numDof,numBC,NULL,NULL,NULL,NULL,§ion)); PetscCall(DMSetLocalSection(celldm,section)); PetscCall(DMSetUp(celldm)); /* Calling DMPlexComputeGeometryFVM() generates the value returned by DMPlexGetMinRadius() */ PetscCall(DMPlexComputeGeometryFVM(celldm,&cellgeom,&facegeom)); PetscCall(DMPlexGetMinRadius(celldm,&minradius)); PetscCall(PetscPrintf(PETSC_COMM_WORLD,"PLEX(minradius) %1.4e\n",(double)minradius)); PetscCall(VecDestroy(&cellgeom)); PetscCall(VecDestroy(&facegeom)); PetscCall(PetscSectionDestroy(§ion)); } /* Create the DMSwarm */ PetscCall(DMCreate(PETSC_COMM_WORLD,&swarm)); PetscCall(DMSetType(swarm,DMSWARM)); PetscCall(DMSetDimension(swarm,dim)); /* Configure swarm to be of type PIC */ PetscCall(DMSwarmSetType(swarm,DMSWARM_PIC)); PetscCall(DMSwarmSetCellDM(swarm,celldm)); /* Register two scalar fields within the DMSwarm */ PetscCall(DMSwarmRegisterPetscDatatypeField(swarm,"phi",1,PETSC_REAL)); PetscCall(DMSwarmRegisterPetscDatatypeField(swarm,"region",1,PETSC_REAL)); PetscCall(DMSwarmFinalizeFieldRegister(swarm)); /* Set initial local sizes of the DMSwarm with a buffer length of zero */ PetscCall(DMSwarmSetLocalSizes(swarm,4,0)); /* Insert swarm coordinates cell-wise */ /*PetscCall(DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_REGULAR,ppcell));*/ PetscCall(DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_SUBDIVISION,ppcell)); /* Define initial conditions for th swarm fields "phi" and "region" */ { PetscReal *s_coor,*s_phi,*s_region; PetscInt npoints,p; PetscCall(DMSwarmGetLocalSize(swarm,&npoints)); PetscCall(DMSwarmGetField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&s_coor)); PetscCall(DMSwarmGetField(swarm,"phi",NULL,NULL,(void**)&s_phi)); PetscCall(DMSwarmGetField(swarm,"region",NULL,NULL,(void**)&s_region)); for (p=0; p must be 0 or 1"); PetscCall(pic_advect(ppcell,meshtype)); PetscCall(PetscFinalize()); return 0; }