static char help[] = "Tests DMSwarm with DMShell\n\n"; #include #include #include #include #include #include PetscErrorCode _DMLocatePoints_DMDARegular_IS(DM dm,Vec pos,IS *iscell) { PetscInt p,n,bs,npoints,si,sj,milocal,mjlocal,mx,my; DM dmregular; PetscInt *cellidx; const PetscScalar *coor; PetscReal dx,dy; PetscMPIInt rank; PetscFunctionBegin; CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); CHKERRQ(VecGetLocalSize(pos,&n)); CHKERRQ(VecGetBlockSize(pos,&bs)); npoints = n/bs; CHKERRQ(PetscMalloc1(npoints,&cellidx)); CHKERRQ(DMGetApplicationContext(dm,&dmregular)); CHKERRQ(DMDAGetCorners(dmregular,&si,&sj,NULL,&milocal,&mjlocal,NULL)); CHKERRQ(DMDAGetInfo(dmregular,NULL,&mx,&my,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL)); dx = 2.0/((PetscReal)mx); dy = 2.0/((PetscReal)my); CHKERRQ(VecGetArrayRead(pos,&coor)); for (p=0; p= sj) && (mj < sj + mjlocal)) { if ((mi >= si) && (mi < si + milocal)) { cellidx[p] = (mi-si) + (mj-sj) * milocal; } } if (coorx < -1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND; if (coorx > 1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND; if (coory < -1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND; if (coory > 1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND; } CHKERRQ(VecRestoreArrayRead(pos,&coor)); CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,npoints,cellidx,PETSC_OWN_POINTER,iscell)); PetscFunctionReturn(0); } PetscErrorCode DMLocatePoints_DMDARegular(DM dm,Vec pos,DMPointLocationType ltype, PetscSF cellSF) { IS iscell; PetscSFNode *cells; PetscInt p,bs,npoints,nfound; const PetscInt *boxCells; PetscFunctionBegin; CHKERRQ(_DMLocatePoints_DMDARegular_IS(dm,pos,&iscell)); CHKERRQ(VecGetLocalSize(pos,&npoints)); CHKERRQ(VecGetBlockSize(pos,&bs)); npoints = npoints / bs; CHKERRQ(PetscMalloc1(npoints, &cells)); CHKERRQ(ISGetIndices(iscell, &boxCells)); for (p=0; pops->locatepoints = DMLocatePoints_DMDARegular; dmcell->ops->getneighbors = DMGetNeighbors_DMDARegular; /* Create the swarm */ CHKERRQ(DMCreate(PETSC_COMM_WORLD,&dms)); CHKERRQ(DMSetType(dms,DMSWARM)); CHKERRQ(DMSetDimension(dms,2)); CHKERRQ(DMSwarmSetType(dms,DMSWARM_PIC)); CHKERRQ(DMSwarmSetCellDM(dms,dmcell)); CHKERRQ(DMSwarmRegisterPetscDatatypeField(dms,"itag",1,PETSC_INT)); CHKERRQ(DMSwarmFinalizeFieldRegister(dms)); { PetscInt si,sj,milocal,mjlocal; const PetscScalar *LA_coors; Vec coors; PetscInt cnt; CHKERRQ(DMDAGetCorners(dmregular,&si,&sj,NULL,&milocal,&mjlocal,NULL)); CHKERRQ(DMGetCoordinates(dmregular,&coors)); CHKERRQ(VecGetArrayRead(coors,&LA_coors)); CHKERRQ(DMSwarmSetLocalSizes(dms,milocal*mjlocal,4)); CHKERRQ(DMSwarmGetLocalSize(dms,&nlocal)); CHKERRQ(DMSwarmGetField(dms,DMSwarmPICField_coor,&bs,NULL,(void**)&array)); cnt = 0; CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&rand)); CHKERRQ(PetscRandomSetInterval(rand,-1.0,1.0)); for (p=0; p