10e2ec84fSDave May 20e2ec84fSDave May #define PETSCDM_DLL 30e2ec84fSDave May #include <petsc/private/dmswarmimpl.h> /*I "petscdmswarm.h" I*/ 40e2ec84fSDave May #include <petscsf.h> 5b799feefSDave May #include <petscdmda.h> 6b799feefSDave May #include <petscdmplex.h> 70e2ec84fSDave May 80e2ec84fSDave May /* 90e2ec84fSDave May Error chceking macto to ensure the swarm type is correct and that a cell DM has been set 100e2ec84fSDave May */ 110e2ec84fSDave May #define DMSWARMPICVALID(dm) \ 120e2ec84fSDave May { \ 130e2ec84fSDave May DM_Swarm *_swarm = (DM_Swarm*)(dm)->data; \ 140e2ec84fSDave May if (_swarm->swarm_type != DMSWARM_PIC) SETERRQ(PetscObjectComm((PetscObject)(dm)),PETSC_ERR_SUP,"Only valid for DMSwarm-PIC. You must call DMSwarmSetType(dm,DMSWARM_PIC)"); \ 150e2ec84fSDave May else \ 160e2ec84fSDave May if (!_swarm->dmcell) SETERRQ(PetscObjectComm((PetscObject)(dm)),PETSC_ERR_SUP,"Only valid for DMSwarmPIC if the cell DM is set. You must call DMSwarmSetCellDM(dm,celldm)"); \ 170e2ec84fSDave May } 180e2ec84fSDave May 190e2ec84fSDave May /* Coordinate insertition/addition API */ 200e2ec84fSDave May /*@C 210e2ec84fSDave May DMSwarmSetPointsUniformCoordinates - Set point coordinates in a DMSwarm on a regular (ijk) grid 220e2ec84fSDave May 230e2ec84fSDave May Collective on DM 240e2ec84fSDave May 250e2ec84fSDave May Input parameters: 260e2ec84fSDave May + dm - the DMSwarm 270e2ec84fSDave May . min - minimum coordinate values in the x, y, z directions (array of length dim) 280e2ec84fSDave May . max - maximum coordinate values in the x, y, z directions (array of length dim) 290e2ec84fSDave May . npoints - number of points in each spatial direction (array of length dim) 300e2ec84fSDave May - mode - indicates whether to append points to the swarm (ADD_VALUES), or over-ride existing points (INSERT_VALUES) 310e2ec84fSDave May 320e2ec84fSDave May Level: beginner 330e2ec84fSDave May 340e2ec84fSDave May Notes: 350e2ec84fSDave May When using mode = INSERT_VALUES, this method will reset the number of particles in the DMSwarm 360e2ec84fSDave May to be npoints[0]*npoints[1] (2D) or npoints[0]*npoints[1]*npoints[2] (3D). When using mode = ADD_VALUES, 370e2ec84fSDave May new points will be appended to any already existing in the DMSwarm 380e2ec84fSDave May 390e2ec84fSDave May .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType 400e2ec84fSDave May @*/ 410e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmSetPointsUniformCoordinates(DM dm,PetscReal min[],PetscReal max[],PetscInt npoints[],InsertMode mode) 420e2ec84fSDave May { 430e2ec84fSDave May PetscErrorCode ierr; 440e2ec84fSDave May PetscReal gmin[] = {PETSC_MAX_REAL ,PETSC_MAX_REAL, PETSC_MAX_REAL}; 450e2ec84fSDave May PetscReal gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 460e2ec84fSDave May PetscInt i,j,k,N,bs,b,n_estimate,n_curr,n_new_est,p,n_found; 470e2ec84fSDave May Vec coorlocal; 480e2ec84fSDave May const PetscScalar *_coor; 490e2ec84fSDave May DM celldm; 500e2ec84fSDave May PetscReal dx[3]; 510e2ec84fSDave May Vec pos; 520e2ec84fSDave May PetscScalar *_pos; 530e2ec84fSDave May PetscReal *swarm_coor; 540e2ec84fSDave May PetscInt *swarm_cellid; 550e2ec84fSDave May PetscSF sfcell = NULL; 560e2ec84fSDave May const PetscSFNode *LA_sfcell; 570e2ec84fSDave May 580e2ec84fSDave May PetscFunctionBegin; 590e2ec84fSDave May DMSWARMPICVALID(dm); 600e2ec84fSDave May ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr); 610e2ec84fSDave May ierr = DMGetCoordinatesLocal(celldm,&coorlocal);CHKERRQ(ierr); 620e2ec84fSDave May ierr = VecGetSize(coorlocal,&N);CHKERRQ(ierr); 630e2ec84fSDave May ierr = VecGetBlockSize(coorlocal,&bs);CHKERRQ(ierr); 640e2ec84fSDave May N = N / bs; 650e2ec84fSDave May ierr = VecGetArrayRead(coorlocal,&_coor);CHKERRQ(ierr); 660e2ec84fSDave May for (i=0; i<N; i++) { 670e2ec84fSDave May for (b=0; b<bs; b++) { 680e2ec84fSDave May gmin[b] = PetscMin(gmin[b],_coor[bs*i+b]); 690e2ec84fSDave May gmax[b] = PetscMax(gmax[b],_coor[bs*i+b]); 700e2ec84fSDave May } 710e2ec84fSDave May } 720e2ec84fSDave May ierr = VecRestoreArrayRead(coorlocal,&_coor);CHKERRQ(ierr); 730e2ec84fSDave May 740e2ec84fSDave May for (b=0; b<bs; b++) { 75*71844bb8SDave May if (npoints[b] > 1) { 760e2ec84fSDave May dx[b] = (max[b] - min[b])/((PetscReal)(npoints[b]-1)); 77*71844bb8SDave May } else { 78*71844bb8SDave May dx[b] = 0.0; 79*71844bb8SDave May } 800e2ec84fSDave May } 810e2ec84fSDave May 820e2ec84fSDave May /* determine number of points living in the bounding box */ 830e2ec84fSDave May n_estimate = 0; 840e2ec84fSDave May if (bs == 2) { npoints[2] = 1; } 850e2ec84fSDave May for (k=0; k<npoints[2]; k++) { 860e2ec84fSDave May for (j=0; j<npoints[1]; j++) { 870e2ec84fSDave May for (i=0; i<npoints[0]; i++) { 880e2ec84fSDave May PetscReal xp[] = {0.0,0.0,0.0}; 890e2ec84fSDave May PetscInt ijk[3]; 900e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 910e2ec84fSDave May 920e2ec84fSDave May ijk[0] = i; 930e2ec84fSDave May ijk[1] = j; 940e2ec84fSDave May ijk[2] = k; 950e2ec84fSDave May for (b=0; b<bs; b++) { 960e2ec84fSDave May xp[b] = min[b] + ijk[b] * dx[b]; 970e2ec84fSDave May } 980e2ec84fSDave May for (b=0; b<bs; b++) { 990e2ec84fSDave May if (xp[b] < gmin[b]) { point_inside = PETSC_FALSE; } 1000e2ec84fSDave May if (xp[b] > gmax[b]) { point_inside = PETSC_FALSE; } 1010e2ec84fSDave May } 1020e2ec84fSDave May if (point_inside) { n_estimate++; } 1030e2ec84fSDave May } 1040e2ec84fSDave May } 1050e2ec84fSDave May } 1060e2ec84fSDave May 1070e2ec84fSDave May /* create candidate list */ 1080e2ec84fSDave May ierr = VecCreate(PetscObjectComm((PetscObject)dm),&pos);CHKERRQ(ierr); 1090e2ec84fSDave May ierr = VecSetSizes(pos,bs*n_estimate,PETSC_DECIDE);CHKERRQ(ierr); 1100e2ec84fSDave May ierr = VecSetBlockSize(pos,bs);CHKERRQ(ierr); 1110e2ec84fSDave May ierr = VecSetFromOptions(pos);CHKERRQ(ierr); 1120e2ec84fSDave May ierr = VecGetArray(pos,&_pos);CHKERRQ(ierr); 1130e2ec84fSDave May 1140e2ec84fSDave May n_estimate = 0; 1150e2ec84fSDave May for (k=0; k<npoints[2]; k++) { 1160e2ec84fSDave May for (j=0; j<npoints[1]; j++) { 1170e2ec84fSDave May for (i=0; i<npoints[0]; i++) { 1180e2ec84fSDave May PetscReal xp[] = {0.0,0.0,0.0}; 1190e2ec84fSDave May PetscInt ijk[3]; 1200e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 1210e2ec84fSDave May 1220e2ec84fSDave May ijk[0] = i; 1230e2ec84fSDave May ijk[1] = j; 1240e2ec84fSDave May ijk[2] = k; 1250e2ec84fSDave May for (b=0; b<bs; b++) { 1260e2ec84fSDave May xp[b] = min[b] + ijk[b] * dx[b]; 1270e2ec84fSDave May } 1280e2ec84fSDave May for (b=0; b<bs; b++) { 1290e2ec84fSDave May if (xp[b] < gmin[b]) { point_inside = PETSC_FALSE; } 1300e2ec84fSDave May if (xp[b] > gmax[b]) { point_inside = PETSC_FALSE; } 1310e2ec84fSDave May } 1320e2ec84fSDave May if (point_inside) { 1330e2ec84fSDave May for (b=0; b<bs; b++) { 1340e2ec84fSDave May _pos[bs*n_estimate+b] = xp[b]; 1350e2ec84fSDave May } 1360e2ec84fSDave May n_estimate++; 1370e2ec84fSDave May } 1380e2ec84fSDave May } 1390e2ec84fSDave May } 1400e2ec84fSDave May } 1410e2ec84fSDave May ierr = VecRestoreArray(pos,&_pos);CHKERRQ(ierr); 1420e2ec84fSDave May 1430e2ec84fSDave May /* locate points */ 1440e2ec84fSDave May ierr = DMLocatePoints(celldm,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr); 1450e2ec84fSDave May 1460e2ec84fSDave May ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr); 1470e2ec84fSDave May n_found = 0; 1480e2ec84fSDave May for (p=0; p<n_estimate; p++) { 1490e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 1500e2ec84fSDave May n_found++; 1510e2ec84fSDave May } 1520e2ec84fSDave May } 1530e2ec84fSDave May 1540e2ec84fSDave May /* adjust size */ 1550e2ec84fSDave May if (mode == ADD_VALUES) { 1560e2ec84fSDave May ierr = DMSwarmGetLocalSize(dm,&n_curr);CHKERRQ(ierr); 1570e2ec84fSDave May n_new_est = n_curr + n_found; 1580e2ec84fSDave May ierr = DMSwarmSetLocalSizes(dm,n_new_est,-1);CHKERRQ(ierr); 1590e2ec84fSDave May } 1600e2ec84fSDave May if (mode == INSERT_VALUES) { 1610e2ec84fSDave May n_curr = 0; 1620e2ec84fSDave May n_new_est = n_found; 1630e2ec84fSDave May ierr = DMSwarmSetLocalSizes(dm,n_new_est,-1);CHKERRQ(ierr); 1640e2ec84fSDave May } 1650e2ec84fSDave May 1660e2ec84fSDave May /* initialize new coords, cell owners, pid */ 1670e2ec84fSDave May ierr = VecGetArrayRead(pos,&_coor);CHKERRQ(ierr); 1680e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 1690e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 1700e2ec84fSDave May n_found = 0; 1710e2ec84fSDave May for (p=0; p<n_estimate; p++) { 1720e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 1730e2ec84fSDave May for (b=0; b<bs; b++) { 1740e2ec84fSDave May swarm_coor[bs*(n_curr + n_found) + b] = _coor[bs*p+b]; 1750e2ec84fSDave May } 1760e2ec84fSDave May swarm_cellid[n_curr + n_found] = LA_sfcell[p].index; 1770e2ec84fSDave May n_found++; 1780e2ec84fSDave May } 1790e2ec84fSDave May } 1800e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 1810e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 1820e2ec84fSDave May ierr = VecRestoreArrayRead(pos,&_coor);CHKERRQ(ierr); 1830e2ec84fSDave May 1840e2ec84fSDave May ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr); 1850e2ec84fSDave May ierr = VecDestroy(&pos);CHKERRQ(ierr); 1860e2ec84fSDave May 1870e2ec84fSDave May PetscFunctionReturn(0); 1880e2ec84fSDave May } 1890e2ec84fSDave May 1900e2ec84fSDave May /*@C 1910e2ec84fSDave May DMSwarmSetPointCoordinates - Set point coordinates in a DMSwarm from a user defined list 1920e2ec84fSDave May 1930e2ec84fSDave May Collective on DM 1940e2ec84fSDave May 1950e2ec84fSDave May Input parameters: 1960e2ec84fSDave May + dm - the DMSwarm 1970e2ec84fSDave May . npoints - the number of points to insert 1980e2ec84fSDave May . coor - the coordinate values 1990e2ec84fSDave May . redundant - if set to PETSC_TRUE, it is assumed that npoints and coor[] are only valid on rank 0 and should be broadcast to other ranks 2000e2ec84fSDave May - mode - indicates whether to append points to the swarm (ADD_VALUES), or over-ride existing points (INSERT_VALUES) 2010e2ec84fSDave May 2020e2ec84fSDave May Level: beginner 2030e2ec84fSDave May 2040e2ec84fSDave May Notes: 2050e2ec84fSDave May If the user has specified redundant = PETSC_FALSE, the cell DM will attempt to locate the coordinates provided by coor[] within 2060e2ec84fSDave May its sub-domain. If they any values within coor[] are not located in the sub-domain, they will be ignored and will not get 2070e2ec84fSDave May added to the DMSwarm. 2080e2ec84fSDave May 2090e2ec84fSDave May .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType, DMSwarmSetPointsUniformCoordinates() 2100e2ec84fSDave May @*/ 2110e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinates(DM dm,PetscInt npoints,PetscReal coor[],PetscBool redundant,InsertMode mode) 2120e2ec84fSDave May { 2130e2ec84fSDave May PetscErrorCode ierr; 2140e2ec84fSDave May PetscReal gmin[] = {PETSC_MAX_REAL ,PETSC_MAX_REAL, PETSC_MAX_REAL}; 2150e2ec84fSDave May PetscReal gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 2160e2ec84fSDave May PetscInt i,N,bs,b,n_estimate,n_curr,n_new_est,p,n_found; 2170e2ec84fSDave May Vec coorlocal; 2180e2ec84fSDave May const PetscScalar *_coor; 2190e2ec84fSDave May DM celldm; 2200e2ec84fSDave May Vec pos; 2210e2ec84fSDave May PetscScalar *_pos; 2220e2ec84fSDave May PetscReal *swarm_coor; 2230e2ec84fSDave May PetscInt *swarm_cellid; 2240e2ec84fSDave May PetscSF sfcell = NULL; 2250e2ec84fSDave May const PetscSFNode *LA_sfcell; 2260e2ec84fSDave May PetscReal *my_coor; 2270e2ec84fSDave May PetscInt my_npoints; 2280e2ec84fSDave May PetscMPIInt rank; 2290e2ec84fSDave May MPI_Comm comm; 2300e2ec84fSDave May 2310e2ec84fSDave May PetscFunctionBegin; 2320e2ec84fSDave May DMSWARMPICVALID(dm); 2330e2ec84fSDave May ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 2340e2ec84fSDave May ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 2350e2ec84fSDave May 2360e2ec84fSDave May ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr); 2370e2ec84fSDave May ierr = DMGetCoordinatesLocal(celldm,&coorlocal);CHKERRQ(ierr); 2380e2ec84fSDave May ierr = VecGetSize(coorlocal,&N);CHKERRQ(ierr); 2390e2ec84fSDave May ierr = VecGetBlockSize(coorlocal,&bs);CHKERRQ(ierr); 2400e2ec84fSDave May N = N / bs; 2410e2ec84fSDave May ierr = VecGetArrayRead(coorlocal,&_coor);CHKERRQ(ierr); 2420e2ec84fSDave May for (i=0; i<N; i++) { 2430e2ec84fSDave May for (b=0; b<bs; b++) { 2440e2ec84fSDave May gmin[b] = PetscMin(gmin[b],_coor[bs*i+b]); 2450e2ec84fSDave May gmax[b] = PetscMax(gmax[b],_coor[bs*i+b]); 2460e2ec84fSDave May } 2470e2ec84fSDave May } 2480e2ec84fSDave May ierr = VecRestoreArrayRead(coorlocal,&_coor);CHKERRQ(ierr); 2490e2ec84fSDave May 2500e2ec84fSDave May /* broadcast points from rank 0 if requested */ 2510e2ec84fSDave May if (redundant) { 2520e2ec84fSDave May my_npoints = npoints; 2530e2ec84fSDave May ierr = MPI_Bcast(&my_npoints,1,MPIU_INT,0,comm);CHKERRQ(ierr); 2540e2ec84fSDave May 2550e2ec84fSDave May if (rank > 0) { /* allocate space */ 256*71844bb8SDave May ierr = PetscMalloc1(bs*my_npoints,&my_coor);CHKERRQ(ierr); 2570e2ec84fSDave May } else { 2580e2ec84fSDave May my_coor = coor; 2590e2ec84fSDave May } 2600e2ec84fSDave May ierr = MPI_Bcast(my_coor,bs*my_npoints,MPIU_REAL,0,comm);CHKERRQ(ierr); 2610e2ec84fSDave May } else { 2620e2ec84fSDave May my_npoints = npoints; 2630e2ec84fSDave May my_coor = coor; 2640e2ec84fSDave May } 2650e2ec84fSDave May 2660e2ec84fSDave May /* determine the number of points living in the bounding box */ 2670e2ec84fSDave May n_estimate = 0; 2680e2ec84fSDave May for (i=0; i<my_npoints; i++) { 2690e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 2700e2ec84fSDave May 2710e2ec84fSDave May for (b=0; b<bs; b++) { 2720e2ec84fSDave May if (my_coor[bs*i+b] < gmin[b]) { point_inside = PETSC_FALSE; } 2730e2ec84fSDave May if (my_coor[bs*i+b] > gmax[b]) { point_inside = PETSC_FALSE; } 2740e2ec84fSDave May } 2750e2ec84fSDave May if (point_inside) { n_estimate++; } 2760e2ec84fSDave May } 2770e2ec84fSDave May 2780e2ec84fSDave May /* create candidate list */ 2790e2ec84fSDave May ierr = VecCreate(PetscObjectComm((PetscObject)dm),&pos);CHKERRQ(ierr); 2800e2ec84fSDave May ierr = VecSetSizes(pos,bs*n_estimate,PETSC_DECIDE);CHKERRQ(ierr); 2810e2ec84fSDave May ierr = VecSetBlockSize(pos,bs);CHKERRQ(ierr); 2820e2ec84fSDave May ierr = VecSetFromOptions(pos);CHKERRQ(ierr); 2830e2ec84fSDave May ierr = VecGetArray(pos,&_pos);CHKERRQ(ierr); 2840e2ec84fSDave May 2850e2ec84fSDave May n_estimate = 0; 2860e2ec84fSDave May for (i=0; i<my_npoints; i++) { 2870e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 2880e2ec84fSDave May 2890e2ec84fSDave May for (b=0; b<bs; b++) { 2900e2ec84fSDave May if (my_coor[bs*i+b] < gmin[b]) { point_inside = PETSC_FALSE; } 2910e2ec84fSDave May if (my_coor[bs*i+b] > gmax[b]) { point_inside = PETSC_FALSE; } 2920e2ec84fSDave May } 2930e2ec84fSDave May if (point_inside) { 2940e2ec84fSDave May for (b=0; b<bs; b++) { 2950e2ec84fSDave May _pos[bs*n_estimate+b] = my_coor[bs*i+b]; 2960e2ec84fSDave May } 2970e2ec84fSDave May n_estimate++; 2980e2ec84fSDave May } 2990e2ec84fSDave May } 3000e2ec84fSDave May ierr = VecRestoreArray(pos,&_pos);CHKERRQ(ierr); 3010e2ec84fSDave May 3020e2ec84fSDave May /* locate points */ 3030e2ec84fSDave May ierr = DMLocatePoints(celldm,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr); 3040e2ec84fSDave May 3050e2ec84fSDave May ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr); 3060e2ec84fSDave May n_found = 0; 3070e2ec84fSDave May for (p=0; p<n_estimate; p++) { 3080e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 3090e2ec84fSDave May n_found++; 3100e2ec84fSDave May } 3110e2ec84fSDave May } 3120e2ec84fSDave May 3130e2ec84fSDave May /* adjust size */ 3140e2ec84fSDave May if (mode == ADD_VALUES) { 3150e2ec84fSDave May ierr = DMSwarmGetLocalSize(dm,&n_curr);CHKERRQ(ierr); 3160e2ec84fSDave May n_new_est = n_curr + n_found; 3170e2ec84fSDave May ierr = DMSwarmSetLocalSizes(dm,n_new_est,-1);CHKERRQ(ierr); 3180e2ec84fSDave May } 3190e2ec84fSDave May if (mode == INSERT_VALUES) { 3200e2ec84fSDave May n_curr = 0; 3210e2ec84fSDave May n_new_est = n_found; 3220e2ec84fSDave May ierr = DMSwarmSetLocalSizes(dm,n_new_est,-1);CHKERRQ(ierr); 3230e2ec84fSDave May } 3240e2ec84fSDave May 3250e2ec84fSDave May /* initialize new coords, cell owners, pid */ 3260e2ec84fSDave May ierr = VecGetArrayRead(pos,&_coor);CHKERRQ(ierr); 3270e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 3280e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 3290e2ec84fSDave May n_found = 0; 3300e2ec84fSDave May for (p=0; p<n_estimate; p++) { 3310e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 3320e2ec84fSDave May for (b=0; b<bs; b++) { 3330e2ec84fSDave May swarm_coor[bs*(n_curr + n_found) + b] = _coor[bs*p+b]; 3340e2ec84fSDave May } 3350e2ec84fSDave May swarm_cellid[n_curr + n_found] = LA_sfcell[p].index; 3360e2ec84fSDave May n_found++; 3370e2ec84fSDave May } 3380e2ec84fSDave May } 3390e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 3400e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 3410e2ec84fSDave May ierr = VecRestoreArrayRead(pos,&_coor);CHKERRQ(ierr); 3420e2ec84fSDave May 3430e2ec84fSDave May if (redundant) { 3440e2ec84fSDave May if (rank > 0) { 3450e2ec84fSDave May ierr = PetscFree(my_coor);CHKERRQ(ierr); 3460e2ec84fSDave May } 3470e2ec84fSDave May } 3480e2ec84fSDave May ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr); 3490e2ec84fSDave May ierr = VecDestroy(&pos);CHKERRQ(ierr); 3500e2ec84fSDave May 3510e2ec84fSDave May PetscFunctionReturn(0); 3520e2ec84fSDave May } 3530e2ec84fSDave May 3540e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM,DM,DMSwarmPICLayoutType,PetscInt); 3550e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM,DM,DMSwarmPICLayoutType,PetscInt); 3560e2ec84fSDave May 3570e2ec84fSDave May /*@C 3580e2ec84fSDave May DMSwarmInsertPointsUsingCellDM - Insert point coordinates within each cell 3590e2ec84fSDave May 3600e2ec84fSDave May Not collective 3610e2ec84fSDave May 3620e2ec84fSDave May Input parameters: 3630e2ec84fSDave May + dm - the DMSwarm 3640e2ec84fSDave May . layout_type - method used to fill each cell with the cell DM 3650e2ec84fSDave May - fill_param - parameter controlling how many points per cell are added (the meaning of this parameter is dependent on the layout type) 3660e2ec84fSDave May 3670e2ec84fSDave May Level: beginner 3680e2ec84fSDave May 3690e2ec84fSDave May Notes: 3700e2ec84fSDave May The insert method will reset any previous defined points within the DMSwarm 3710e2ec84fSDave May 3720e2ec84fSDave May .seealso: DMSwarmPICLayoutType, DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType 3730e2ec84fSDave May @*/ 3740e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmInsertPointsUsingCellDM(DM dm,DMSwarmPICLayoutType layout_type,PetscInt fill_param) 3750e2ec84fSDave May { 3760e2ec84fSDave May PetscErrorCode ierr; 3770e2ec84fSDave May DM celldm; 3780e2ec84fSDave May PetscBool isDA,isPLEX; 3790e2ec84fSDave May 3800e2ec84fSDave May PetscFunctionBegin; 3810e2ec84fSDave May DMSWARMPICVALID(dm); 3820e2ec84fSDave May ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr); 3830e2ec84fSDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMDA,&isDA);CHKERRQ(ierr); 3840e2ec84fSDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMPLEX,&isPLEX);CHKERRQ(ierr); 3850e2ec84fSDave May if (isDA) { 3860e2ec84fSDave May ierr = private_DMSwarmInsertPointsUsingCellDM_DA(dm,celldm,layout_type,fill_param);CHKERRQ(ierr); 3870e2ec84fSDave May } else if (isPLEX) { 3880e2ec84fSDave May ierr = private_DMSwarmInsertPointsUsingCellDM_PLEX(dm,celldm,layout_type,fill_param);CHKERRQ(ierr); 3890e2ec84fSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only supported for cell DMs of type DMDA and DMPLEX"); 3900e2ec84fSDave May 3910e2ec84fSDave May PetscFunctionReturn(0); 3920e2ec84fSDave May } 3930e2ec84fSDave May 3940e2ec84fSDave May /* 3950e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmAddPointCoordinatesCellWise(DM dm,PetscInt cell,PetscInt npoints,PetscReal xi[],PetscBool proximity_initialization) 3960e2ec84fSDave May { 3970e2ec84fSDave May PetscFunctionBegin; 3980e2ec84fSDave May PetscFunctionReturn(0); 3990e2ec84fSDave May } 4000e2ec84fSDave May */ 4010e2ec84fSDave May 4020e2ec84fSDave May /* Field projection API */ 4030e2ec84fSDave May /* 4040e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM dm,PetscInt project_type,PetscInt nfields,const char *fieldnames[],Vec *fields) 4050e2ec84fSDave May { 4060e2ec84fSDave May PetscFunctionBegin; 4070e2ec84fSDave May PetscFunctionReturn(0); 4080e2ec84fSDave May } 4090e2ec84fSDave May */ 4100e2ec84fSDave May 4110e2ec84fSDave May /*@C 412b799feefSDave May DMSwarmCreatePointPerCellCount - Count the number of points within all cells in the cell DM 4130e2ec84fSDave May 4140e2ec84fSDave May Not collective 4150e2ec84fSDave May 4160e2ec84fSDave May Input parameter: 4170e2ec84fSDave May . dm - the DMSwarm 4180e2ec84fSDave May 4190e2ec84fSDave May Output parameters: 4200e2ec84fSDave May + ncells - the number of cells in the cell DM (optional argument, pass NULL to ignore) 421b799feefSDave May - count - array of length ncells containing the number of points per cell 4220e2ec84fSDave May 4230e2ec84fSDave May Level: beginner 4240e2ec84fSDave May 4250e2ec84fSDave May Notes: 4260e2ec84fSDave May The array count is allocated internally and must be free'd by the user. 4270e2ec84fSDave May 4280e2ec84fSDave May .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType 4290e2ec84fSDave May @*/ 4300e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmCreatePointPerCellCount(DM dm,PetscInt *ncells,PetscInt **count) 4310e2ec84fSDave May { 432b799feefSDave May PetscErrorCode ierr; 433b799feefSDave May PetscBool isvalid; 434b799feefSDave May PetscInt nel; 435b799feefSDave May PetscInt *sum; 436b799feefSDave May 4370e2ec84fSDave May PetscFunctionBegin; 438b799feefSDave May ierr = DMSwarmSortGetIsValid(dm,&isvalid);CHKERRQ(ierr); 439b799feefSDave May nel = 0; 440b799feefSDave May if (isvalid) { 441b799feefSDave May PetscInt e; 442b799feefSDave May 443b799feefSDave May ierr = DMSwarmSortGetSizes(dm,&nel,NULL);CHKERRQ(ierr); 444b799feefSDave May 445b799feefSDave May ierr = PetscMalloc1(nel,&sum);CHKERRQ(ierr); 446b799feefSDave May for (e=0; e<nel; e++) { 447b799feefSDave May ierr = DMSwarmSortGetNumberOfPointsPerCell(dm,e,&sum[e]);CHKERRQ(ierr); 448b799feefSDave May } 449b799feefSDave May } else { 450b799feefSDave May DM celldm; 451b799feefSDave May PetscBool isda,isplex,isshell; 452b799feefSDave May PetscInt p,npoints; 453b799feefSDave May PetscInt *swarm_cellid; 454b799feefSDave May 455b799feefSDave May /* get the number of cells */ 456b799feefSDave May ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr); 457b799feefSDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMDA,&isda);CHKERRQ(ierr); 458b799feefSDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMPLEX,&isplex);CHKERRQ(ierr); 459b799feefSDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMSHELL,&isshell);CHKERRQ(ierr); 460b799feefSDave May if (isda) { 461b799feefSDave May PetscInt _nel,_npe; 462b799feefSDave May const PetscInt *_element; 463b799feefSDave May 464b799feefSDave May ierr = DMDAGetElements(celldm,&_nel,&_npe,&_element);CHKERRQ(ierr); 465b799feefSDave May nel = _nel; 466b799feefSDave May ierr = DMDARestoreElements(celldm,&_nel,&_npe,&_element);CHKERRQ(ierr); 467b799feefSDave May } else if (isplex) { 468b799feefSDave May PetscInt ps,pe; 469b799feefSDave May 470b799feefSDave May ierr = DMPlexGetHeightStratum(celldm,0,&ps,&pe);CHKERRQ(ierr); 471b799feefSDave May nel = pe - ps; 472b799feefSDave May } else if (isshell) { 473b799feefSDave May PetscErrorCode (*method_DMShellGetNumberOfCells)(DM,PetscInt*); 474b799feefSDave May 475b799feefSDave May ierr = PetscObjectQueryFunction((PetscObject)celldm,"DMGetNumberOfCells_C",&method_DMShellGetNumberOfCells);CHKERRQ(ierr); 476b799feefSDave May if (method_DMShellGetNumberOfCells) { 477b799feefSDave May ierr = method_DMShellGetNumberOfCells(celldm,&nel);CHKERRQ(ierr); 478b799feefSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot determine the number of cells for the DMSHELL object. User must provide a method via PetscObjectComposeFunction( (PetscObject)shelldm, \"DMGetNumberOfCells_C\", your_function_to_compute_number_of_cells );"); 479b799feefSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot determine the number of cells for a DM not of type DA, PLEX or SHELL"); 480b799feefSDave May 481b799feefSDave May ierr = PetscMalloc1(nel,&sum);CHKERRQ(ierr); 482b799feefSDave May ierr = PetscMemzero(sum,sizeof(PetscInt)*nel);CHKERRQ(ierr); 483b799feefSDave May ierr = DMSwarmGetLocalSize(dm,&npoints);CHKERRQ(ierr); 484b799feefSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 485b799feefSDave May for (p=0; p<npoints; p++) { 486b799feefSDave May if (swarm_cellid[p] != DMLOCATEPOINT_POINT_NOT_FOUND) { 487b799feefSDave May sum[ swarm_cellid[p] ]++; 488b799feefSDave May } 489b799feefSDave May } 490b799feefSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 491b799feefSDave May } 492b799feefSDave May if (ncells) { *ncells = nel; } 493b799feefSDave May *count = sum; 4940e2ec84fSDave May PetscFunctionReturn(0); 4950e2ec84fSDave May } 496