1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "DMSwarm-PIC demonstator of advecting points within cell DM defined by a DA or PLEX object \n\ 3c4762a1bSJed Brown Options: \n\ 4c4762a1bSJed Brown -ppcell : Number of times to sub-divide the reference cell when layout the initial particle coordinates \n\ 5c4762a1bSJed Brown -meshtype : 0 ==> DA , 1 ==> PLEX \n\ 6c4762a1bSJed Brown -nt : Number of timestep to perform \n\ 7c4762a1bSJed Brown -view : Write out initial condition and time dependent data \n"; 8c4762a1bSJed Brown 9c4762a1bSJed Brown #include <petsc.h> 10c4762a1bSJed Brown #include <petscdm.h> 11c4762a1bSJed Brown #include <petscdmda.h> 12c4762a1bSJed Brown #include <petscdmswarm.h> 13c4762a1bSJed Brown 14c4762a1bSJed Brown PetscErrorCode pic_advect(PetscInt ppcell,PetscInt meshtype) 15c4762a1bSJed Brown { 16c4762a1bSJed Brown PetscErrorCode ierr; 17c4762a1bSJed Brown const PetscInt dim = 2; 18c4762a1bSJed Brown DM celldm,swarm; 19c4762a1bSJed Brown PetscInt tk,nt = 200; 20c4762a1bSJed Brown PetscBool view = PETSC_FALSE; 21c4762a1bSJed Brown Vec *pfields; 22c4762a1bSJed Brown PetscReal minradius; 23c4762a1bSJed Brown PetscReal dt; 24c4762a1bSJed Brown PetscReal vel[] = { 1.0, 0.16 }; 25c4762a1bSJed Brown const char *fieldnames[] = { "phi" }; 26c4762a1bSJed Brown PetscViewer viewer; 27c4762a1bSJed Brown 28c4762a1bSJed Brown PetscFunctionBegin; 295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-view",&view,NULL)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nt",&nt,NULL)); 31c4762a1bSJed Brown 32c4762a1bSJed Brown /* Create the background cell DM */ 33c4762a1bSJed Brown if (meshtype == 0) { /* DA */ 34c4762a1bSJed Brown PetscInt nxy; 35c4762a1bSJed Brown PetscInt dof = 1; 36c4762a1bSJed Brown PetscInt stencil_width = 1; 37c4762a1bSJed Brown 385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Mesh type: DMDA\n")); 39c4762a1bSJed Brown nxy = 33; 405f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 41c4762a1bSJed Brown 425f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetElementType(celldm,DMDA_ELEMENT_Q1)); 43c4762a1bSJed Brown 445f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(celldm)); 45c4762a1bSJed Brown 465f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(celldm)); 47c4762a1bSJed Brown 485f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetUniformCoordinates(celldm,0.0,1.0,0.0,1.0,0.0,1.5)); 49c4762a1bSJed Brown 50c4762a1bSJed Brown minradius = 1.0/((PetscReal)(nxy-1)); 51c4762a1bSJed Brown 525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"DA(minradius) %1.4e\n",(double)minradius)); 53c4762a1bSJed Brown } 54c4762a1bSJed Brown 55c4762a1bSJed Brown if (meshtype == 1){ /* PLEX */ 56c4762a1bSJed Brown DM distributedMesh = NULL; 57c4762a1bSJed Brown PetscInt numComp[] = {1}; 58c4762a1bSJed Brown PetscInt numDof[] = {1,0,0}; /* vert, edge, cell */ 59c4762a1bSJed Brown PetscInt faces[] = {1,1,1}; 60c4762a1bSJed Brown PetscInt numBC = 0; 61c4762a1bSJed Brown PetscSection section; 62c4762a1bSJed Brown Vec cellgeom = NULL; 63c4762a1bSJed Brown Vec facegeom = NULL; 64c4762a1bSJed Brown 655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Mesh type: DMPLEX\n")); 665f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_TRUE, faces, NULL, NULL, PETSC_TRUE, &celldm)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* Distribute mesh over processes */ 695f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistribute(celldm,0,NULL,&distributedMesh)); 70c4762a1bSJed Brown if (distributedMesh) { 715f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&celldm)); 72c4762a1bSJed Brown celldm = distributedMesh; 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 755f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(celldm)); 76c4762a1bSJed Brown 775f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateSection(celldm,NULL,numComp,numDof,numBC,NULL,NULL,NULL,NULL,§ion)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetLocalSection(celldm,section)); 79c4762a1bSJed Brown 805f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(celldm)); 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* Calling DMPlexComputeGeometryFVM() generates the value returned by DMPlexGetMinRadius() */ 835f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexComputeGeometryFVM(celldm,&cellgeom,&facegeom)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetMinRadius(celldm,&minradius)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"PLEX(minradius) %1.4e\n",(double)minradius)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&cellgeom)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&facegeom)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionDestroy(§ion)); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* Create the DMSwarm */ 925f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(PETSC_COMM_WORLD,&swarm)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(swarm,DMSWARM)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(swarm,dim)); 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* Configure swarm to be of type PIC */ 975f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetType(swarm,DMSWARM_PIC)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetCellDM(swarm,celldm)); 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* Register two scalar fields within the DMSwarm */ 1015f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRegisterPetscDatatypeField(swarm,"phi",1,PETSC_REAL)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRegisterPetscDatatypeField(swarm,"region",1,PETSC_REAL)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmFinalizeFieldRegister(swarm)); 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* Set initial local sizes of the DMSwarm with a buffer length of zero */ 1065f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetLocalSizes(swarm,4,0)); 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* Insert swarm coordinates cell-wise */ 1095f80ce2aSJacob Faibussowitsch /*CHKERRQ(DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_REGULAR,ppcell));*/ 1105f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_SUBDIVISION,ppcell)); 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* Define initial conditions for th swarm fields "phi" and "region" */ 113c4762a1bSJed Brown { 114c4762a1bSJed Brown PetscReal *s_coor,*s_phi,*s_region; 115c4762a1bSJed Brown PetscInt npoints,p; 116c4762a1bSJed Brown 1175f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetLocalSize(swarm,&npoints)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&s_coor)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(swarm,"phi",NULL,NULL,(void**)&s_phi)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(swarm,"region",NULL,NULL,(void**)&s_region)); 121c4762a1bSJed Brown for (p=0; p<npoints; p++) { 122c4762a1bSJed Brown PetscReal pos[2]; 123c4762a1bSJed Brown pos[0] = s_coor[2*p+0]; 124c4762a1bSJed Brown pos[1] = s_coor[2*p+1]; 125c4762a1bSJed Brown 126c4762a1bSJed Brown s_region[p] = 1.0; 127c4762a1bSJed Brown s_phi[p] = 1.0 + PetscExpReal(-200.0*((pos[0]-0.5)*(pos[0]-0.5) + (pos[1]-0.5)*(pos[1]-0.5))); 128c4762a1bSJed Brown } 1295f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(swarm,"region",NULL,NULL,(void**)&s_region)); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(swarm,"phi",NULL,NULL,(void**)&s_phi)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&s_coor)); 132c4762a1bSJed Brown } 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* Project initial value of phi onto the mesh */ 1355f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmProjectFields(swarm,1,fieldnames,&pfields,PETSC_FALSE)); 136c4762a1bSJed Brown 137c4762a1bSJed Brown if (view) { 138c4762a1bSJed Brown /* View swarm all swarm fields using data type PETSC_REAL */ 1395f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmViewXDMF(swarm,"ic_dms.xmf")); 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* View projected swarm field "phi" */ 1425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerCreate(PETSC_COMM_WORLD,&viewer)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerSetType(viewer,PETSCVIEWERVTK)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE)); 145c4762a1bSJed Brown if (meshtype == 0) { /* DA */ 1465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFileSetName(viewer,"ic_dmda.vts")); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(pfields[0],viewer)); 148c4762a1bSJed Brown } 149c4762a1bSJed Brown if (meshtype == 1) { /* PLEX */ 1505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFileSetName(viewer,"ic_dmplex.vtk")); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(DMView(celldm,viewer)); 1525f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(pfields[0],viewer)); 153c4762a1bSJed Brown } 1545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 155c4762a1bSJed Brown } 156c4762a1bSJed Brown 1575f80ce2aSJacob Faibussowitsch CHKERRQ(DMView(celldm,PETSC_VIEWER_STDOUT_WORLD)); 1585f80ce2aSJacob Faibussowitsch CHKERRQ(DMView(swarm,PETSC_VIEWER_STDOUT_WORLD)); 159c4762a1bSJed Brown 160c4762a1bSJed Brown dt = 0.5 * minradius / PetscSqrtReal(vel[0]*vel[0] + vel[1]*vel[1]); 161c4762a1bSJed Brown for (tk=1; tk<=nt; tk++) { 162c4762a1bSJed Brown PetscReal *s_coor; 163c4762a1bSJed Brown PetscInt npoints,p; 164c4762a1bSJed Brown 1655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[step %D]\n",tk)); 166c4762a1bSJed Brown /* advect with analytic prescribed (constant) velocity field */ 1675f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetLocalSize(swarm,&npoints)); 1685f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&s_coor)); 169c4762a1bSJed Brown for (p=0; p<npoints; p++) { 170c4762a1bSJed Brown s_coor[2*p+0] += dt * vel[0]; 171c4762a1bSJed Brown s_coor[2*p+1] += dt * vel[1]; 172c4762a1bSJed Brown } 1735f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&s_coor)); 174c4762a1bSJed Brown 1755f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmMigrate(swarm,PETSC_TRUE)); 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* Ad-hoc cell filling algorithm */ 178c4762a1bSJed Brown /* 179c4762a1bSJed Brown The injection frequency is chosen for default DA case. 180c4762a1bSJed Brown They will likely not be appropriate for the general case using an unstructure PLEX mesh. 181c4762a1bSJed Brown */ 182c4762a1bSJed Brown if (tk%10 == 0) { 183c4762a1bSJed Brown PetscReal dx = 1.0/32.0; 184c4762a1bSJed Brown PetscInt npoints_dir_x[] = { 32, 1 }; 185c4762a1bSJed Brown PetscReal min[2],max[2]; 186c4762a1bSJed Brown 187c4762a1bSJed Brown min[0] = 0.5 * dx; max[0] = 0.5 * dx + 31.0 * dx; 188c4762a1bSJed Brown min[1] = 0.5 * dx; max[1] = 0.5 * dx; 1895f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetPointsUniformCoordinates(swarm,min,max,npoints_dir_x,ADD_VALUES)); 190c4762a1bSJed Brown } 191c4762a1bSJed Brown if (tk%2 == 0) { 192c4762a1bSJed Brown PetscReal dx = 1.0/32.0; 193c4762a1bSJed Brown PetscInt npoints_dir_y[] = { 2, 31 }; 194c4762a1bSJed Brown PetscReal min[2],max[2]; 195c4762a1bSJed Brown 196c4762a1bSJed Brown min[0] = 0.05 * dx; max[0] = 0.5 * dx; 197c4762a1bSJed Brown min[1] = 0.5 * dx; max[1] = 0.5 * dx + 31.0 * dx; 1985f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetPointsUniformCoordinates(swarm,min,max,npoints_dir_y,ADD_VALUES)); 199c4762a1bSJed Brown } 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* Project swarm field "phi" onto the cell DM */ 2025f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmProjectFields(swarm,1,fieldnames,&pfields,PETSC_TRUE)); 203c4762a1bSJed Brown 204c4762a1bSJed Brown if (view) { 205c4762a1bSJed Brown PetscViewer viewer; 206c4762a1bSJed Brown char fname[PETSC_MAX_PATH_LEN]; 207c4762a1bSJed Brown 208c4762a1bSJed Brown /* View swarm fields */ 2095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(fname,PETSC_MAX_PATH_LEN-1,"step%.4D_dms.xmf",tk)); 2105f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmViewXDMF(swarm,fname)); 211c4762a1bSJed Brown 212c4762a1bSJed Brown /* View projected field */ 2135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerCreate(PETSC_COMM_WORLD,&viewer)); 2145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerSetType(viewer,PETSCVIEWERVTK)); 2155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE)); 216c4762a1bSJed Brown 217c4762a1bSJed Brown if (meshtype == 0) { /* DA */ 2185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(fname,PETSC_MAX_PATH_LEN-1,"step%.4D_dmda.vts",tk)); 2195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFileSetName(viewer,fname)); 2205f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(pfields[0],viewer)); 221c4762a1bSJed Brown } 222c4762a1bSJed Brown if (meshtype == 1) { /* PLEX */ 2235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(fname,PETSC_MAX_PATH_LEN-1,"step%.4D_dmplex.vtk",tk)); 2245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFileSetName(viewer,fname)); 2255f80ce2aSJacob Faibussowitsch CHKERRQ(DMView(celldm,viewer)); 2265f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(pfields[0],viewer)); 227c4762a1bSJed Brown } 2285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 229c4762a1bSJed Brown } 230c4762a1bSJed Brown 231c4762a1bSJed Brown } 2325f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&pfields[0])); 2335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pfields)); 2345f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&celldm)); 2355f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&swarm)); 236c4762a1bSJed Brown 237c4762a1bSJed Brown PetscFunctionReturn(0); 238c4762a1bSJed Brown } 239c4762a1bSJed Brown 240c4762a1bSJed Brown int main(int argc,char **args) 241c4762a1bSJed Brown { 242c4762a1bSJed Brown PetscErrorCode ierr; 243c4762a1bSJed Brown PetscInt ppcell = 1; 244c4762a1bSJed Brown PetscInt meshtype = 0; 245c4762a1bSJed Brown 246*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 2475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ppcell",&ppcell,NULL)); 2485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-meshtype",&meshtype,NULL)); 2492c71b3e2SJacob Faibussowitsch PetscCheckFalse(meshtype > 1,PETSC_COMM_WORLD,PETSC_ERR_USER,"-meshtype <value> must be 0 or 1"); 250c4762a1bSJed Brown 2515f80ce2aSJacob Faibussowitsch CHKERRQ(pic_advect(ppcell,meshtype)); 252c4762a1bSJed Brown 253*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 254*b122ec5aSJacob Faibussowitsch return 0; 255c4762a1bSJed Brown } 256