static char help[] = "DMSwarm-PIC demonstator of inserting points into a cell DM \n\ Options: \n\ -mode {0,1} : 0 ==> DMDA, 1 ==> DMPLEX cell DM \n\ -dim {2,3} : spatial dimension\n"; #include #include #include #include #include PetscErrorCode pic_insert_DMDA(PetscInt dim) { PetscErrorCode ierr; DM celldm = NULL,swarm; PetscInt dof,stencil_width; PetscReal min[3],max[3]; PetscInt ndir[3]; PetscFunctionBegin; /* Create the background cell DM */ dof = 1; stencil_width = 1; if (dim == 2) { ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,25,13,PETSC_DECIDE,PETSC_DECIDE,dof,stencil_width,NULL,NULL,&celldm);CHKERRQ(ierr); } if (dim == 3) { ierr = DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,25,13,19,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,dof,stencil_width,NULL,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,2.0,0.0,1.0,0.0,1.5);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,"viscosity",1,PETSC_DOUBLE);CHKERRQ(ierr); ierr = DMSwarmRegisterPetscDatatypeField(swarm,"density",1,PETSC_DOUBLE);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,3);CHKERRQ(ierr); min[0] = 0.5; max[0] = 0.7; min[1] = 0.5; max[1] = 0.8; min[2] = 0.5; max[2] = 0.9; ndir[0] = ndir[1] = ndir[2] = 30; ierr = DMSwarmSetPointsUniformCoordinates(swarm,min,max,ndir,ADD_VALUES);CHKERRQ(ierr); /* This should be dispatched from a regular DMView() */ ierr = DMSwarmViewXDMF(swarm,"ex20.xmf");CHKERRQ(ierr); ierr = DMView(celldm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = DMView(swarm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); { PetscInt npoints,*list; PetscMPIInt rank; ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); ierr = DMSwarmSortGetAccess(swarm);CHKERRQ(ierr); ierr = DMSwarmSortGetNumberOfPointsPerCell(swarm,0,&npoints);CHKERRQ(ierr); ierr = DMSwarmSortGetPointsPerCell(swarm,rank,&npoints,&list);CHKERRQ(ierr); ierr = PetscFree(list);CHKERRQ(ierr); ierr = DMSwarmSortRestoreAccess(swarm);CHKERRQ(ierr); } ierr = DMSwarmMigrate(swarm,PETSC_FALSE);CHKERRQ(ierr); ierr = DMDestroy(&celldm);CHKERRQ(ierr); ierr = DMDestroy(&swarm);CHKERRQ(ierr); PetscFunctionReturn(0); } PetscErrorCode pic_insert_DMPLEX_with_cell_list(PetscInt dim) { PetscErrorCode ierr; DM celldm = NULL,swarm,distributedMesh = NULL; const char *fieldnames[] = {"viscosity"}; PetscFunctionBegin; /* Create the background cell DM */ if (dim == 2) { PetscInt cells_per_dim[2],nx[2]; PetscInt n_tricells; PetscInt n_trivert; PetscInt *tricells; PetscReal *trivert,dx,dy; PetscInt ii,jj,cnt; cells_per_dim[0] = 4; cells_per_dim[1] = 4; n_tricells = cells_per_dim[0] * cells_per_dim[1] * 2; nx[0] = cells_per_dim[0] + 1; nx[1] = cells_per_dim[1] + 1; n_trivert = nx[0] * nx[1]; ierr = PetscMalloc1(n_tricells*3,&tricells);CHKERRQ(ierr); ierr = PetscMalloc1(nx[0]*nx[1]*2,&trivert);CHKERRQ(ierr); /* verts */ cnt = 0; dx = 2.0/((PetscReal)cells_per_dim[0]); dy = 1.0/((PetscReal)cells_per_dim[1]); for (jj=0; jj