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) { 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) { CHKERRQ(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)); } if (dim == 3) { CHKERRQ(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(DMDASetElementType(celldm,DMDA_ELEMENT_Q1)); CHKERRQ(DMSetFromOptions(celldm)); CHKERRQ(DMSetUp(celldm)); CHKERRQ(DMDASetUniformCoordinates(celldm,0.0,2.0,0.0,1.0,0.0,1.5)); /* Create the DMSwarm */ CHKERRQ(DMCreate(PETSC_COMM_WORLD,&swarm)); CHKERRQ(PetscObjectSetName((PetscObject)swarm,"Swarm")); CHKERRQ(DMSetType(swarm,DMSWARM)); CHKERRQ(DMSetDimension(swarm,dim)); /* Configure swarm to be of type PIC */ CHKERRQ(DMSwarmSetType(swarm,DMSWARM_PIC)); CHKERRQ(DMSwarmSetCellDM(swarm,celldm)); /* Register two scalar fields within the DMSwarm */ CHKERRQ(DMSwarmRegisterPetscDatatypeField(swarm,"viscosity",1,PETSC_DOUBLE)); CHKERRQ(DMSwarmRegisterPetscDatatypeField(swarm,"density",1,PETSC_DOUBLE)); CHKERRQ(DMSwarmFinalizeFieldRegister(swarm)); /* Set initial local sizes of the DMSwarm with a buffer length of zero */ CHKERRQ(DMSwarmSetLocalSizes(swarm,4,0)); /* Insert swarm coordinates cell-wise */ CHKERRQ(DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_REGULAR,3)); 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; CHKERRQ(DMSwarmSetPointsUniformCoordinates(swarm,min,max,ndir,ADD_VALUES)); /* This should be dispatched from a regular DMView() */ CHKERRQ(DMSwarmViewXDMF(swarm,"ex20.xmf")); CHKERRQ(DMView(celldm,PETSC_VIEWER_STDOUT_WORLD)); CHKERRQ(DMView(swarm,PETSC_VIEWER_STDOUT_WORLD)); { PetscInt npoints,*list; PetscMPIInt rank; CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); CHKERRQ(DMSwarmSortGetAccess(swarm)); CHKERRQ(DMSwarmSortGetNumberOfPointsPerCell(swarm,0,&npoints)); CHKERRQ(DMSwarmSortGetPointsPerCell(swarm,rank,&npoints,&list)); CHKERRQ(PetscFree(list)); CHKERRQ(DMSwarmSortRestoreAccess(swarm)); } CHKERRQ(DMSwarmMigrate(swarm,PETSC_FALSE)); CHKERRQ(DMDestroy(&celldm)); CHKERRQ(DMDestroy(&swarm)); PetscFunctionReturn(0); } PetscErrorCode pic_insert_DMPLEX_with_cell_list(PetscInt dim) { 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]; CHKERRQ(PetscMalloc1(n_tricells*3,&tricells)); CHKERRQ(PetscMalloc1(nx[0]*nx[1]*2,&trivert)); /* verts */ cnt = 0; dx = 2.0/((PetscReal)cells_per_dim[0]); dy = 1.0/((PetscReal)cells_per_dim[1]); for (jj=0; jj