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