10e2ec84fSDave May #define PETSCDM_DLL 20e2ec84fSDave May #include <petsc/private/dmswarmimpl.h> /*I "petscdmswarm.h" I*/ 30e2ec84fSDave May #include <petscsf.h> 4b799feefSDave May #include <petscdmda.h> 5b799feefSDave May #include <petscdmplex.h> 6279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h" 70e2ec84fSDave May 80e2ec84fSDave May /* 9*5627991aSBarry Smith Error checking 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; \ 14*5627991aSBarry Smith if (_swarm->swarm_type != DMSWARM_PIC) SETERRQ(PetscObjectComm((PetscObject)(dm)),PETSC_ERR_SUP,"Valid only for DMSwarm-PIC. You must call DMSwarmSetType(dm,DMSWARM_PIC)"); \ 150e2ec84fSDave May else \ 16*5627991aSBarry Smith if (!_swarm->dmcell) SETERRQ(PetscObjectComm((PetscObject)(dm)),PETSC_ERR_SUP,"Valid only 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 23d083f849SBarry Smith 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]; 512e3d3761SDave May PetscInt _npoints[] = { 0, 0, 1 }; 520e2ec84fSDave May Vec pos; 530e2ec84fSDave May PetscScalar *_pos; 540e2ec84fSDave May PetscReal *swarm_coor; 550e2ec84fSDave May PetscInt *swarm_cellid; 560e2ec84fSDave May PetscSF sfcell = NULL; 570e2ec84fSDave May const PetscSFNode *LA_sfcell; 580e2ec84fSDave May 590e2ec84fSDave May PetscFunctionBegin; 600e2ec84fSDave May DMSWARMPICVALID(dm); 610e2ec84fSDave May ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr); 620e2ec84fSDave May ierr = DMGetCoordinatesLocal(celldm,&coorlocal);CHKERRQ(ierr); 630e2ec84fSDave May ierr = VecGetSize(coorlocal,&N);CHKERRQ(ierr); 640e2ec84fSDave May ierr = VecGetBlockSize(coorlocal,&bs);CHKERRQ(ierr); 650e2ec84fSDave May N = N / bs; 660e2ec84fSDave May ierr = VecGetArrayRead(coorlocal,&_coor);CHKERRQ(ierr); 670e2ec84fSDave May for (i=0; i<N; i++) { 680e2ec84fSDave May for (b=0; b<bs; b++) { 69a5f152d1SDave May gmin[b] = PetscMin(gmin[b],PetscRealPart(_coor[bs*i+b])); 70a5f152d1SDave May gmax[b] = PetscMax(gmax[b],PetscRealPart(_coor[bs*i+b])); 710e2ec84fSDave May } 720e2ec84fSDave May } 730e2ec84fSDave May ierr = VecRestoreArrayRead(coorlocal,&_coor);CHKERRQ(ierr); 740e2ec84fSDave May 750e2ec84fSDave May for (b=0; b<bs; b++) { 7671844bb8SDave May if (npoints[b] > 1) { 770e2ec84fSDave May dx[b] = (max[b] - min[b])/((PetscReal)(npoints[b]-1)); 7871844bb8SDave May } else { 7971844bb8SDave May dx[b] = 0.0; 8071844bb8SDave May } 812e3d3761SDave May _npoints[b] = npoints[b]; 820e2ec84fSDave May } 830e2ec84fSDave May 840e2ec84fSDave May /* determine number of points living in the bounding box */ 850e2ec84fSDave May n_estimate = 0; 862e3d3761SDave May for (k=0; k<_npoints[2]; k++) { 872e3d3761SDave May for (j=0; j<_npoints[1]; j++) { 882e3d3761SDave May for (i=0; i<_npoints[0]; i++) { 890e2ec84fSDave May PetscReal xp[] = {0.0,0.0,0.0}; 900e2ec84fSDave May PetscInt ijk[3]; 910e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 920e2ec84fSDave May 930e2ec84fSDave May ijk[0] = i; 940e2ec84fSDave May ijk[1] = j; 950e2ec84fSDave May ijk[2] = k; 960e2ec84fSDave May for (b=0; b<bs; b++) { 970e2ec84fSDave May xp[b] = min[b] + ijk[b] * dx[b]; 980e2ec84fSDave May } 990e2ec84fSDave May for (b=0; b<bs; b++) { 1000e2ec84fSDave May if (xp[b] < gmin[b]) { point_inside = PETSC_FALSE; } 1010e2ec84fSDave May if (xp[b] > gmax[b]) { point_inside = PETSC_FALSE; } 1020e2ec84fSDave May } 1030e2ec84fSDave May if (point_inside) { n_estimate++; } 1040e2ec84fSDave May } 1050e2ec84fSDave May } 1060e2ec84fSDave May } 1070e2ec84fSDave May 1080e2ec84fSDave May /* create candidate list */ 109c28a3c15SDave May ierr = VecCreate(PETSC_COMM_SELF,&pos);CHKERRQ(ierr); 1100e2ec84fSDave May ierr = VecSetSizes(pos,bs*n_estimate,PETSC_DECIDE);CHKERRQ(ierr); 1110e2ec84fSDave May ierr = VecSetBlockSize(pos,bs);CHKERRQ(ierr); 1120e2ec84fSDave May ierr = VecSetFromOptions(pos);CHKERRQ(ierr); 1130e2ec84fSDave May ierr = VecGetArray(pos,&_pos);CHKERRQ(ierr); 1140e2ec84fSDave May 1150e2ec84fSDave May n_estimate = 0; 1162e3d3761SDave May for (k=0; k<_npoints[2]; k++) { 1172e3d3761SDave May for (j=0; j<_npoints[1]; j++) { 1182e3d3761SDave May for (i=0; i<_npoints[0]; i++) { 1190e2ec84fSDave May PetscReal xp[] = {0.0,0.0,0.0}; 1200e2ec84fSDave May PetscInt ijk[3]; 1210e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 1220e2ec84fSDave May 1230e2ec84fSDave May ijk[0] = i; 1240e2ec84fSDave May ijk[1] = j; 1250e2ec84fSDave May ijk[2] = k; 1260e2ec84fSDave May for (b=0; b<bs; b++) { 1270e2ec84fSDave May xp[b] = min[b] + ijk[b] * dx[b]; 1280e2ec84fSDave May } 1290e2ec84fSDave May for (b=0; b<bs; b++) { 1300e2ec84fSDave May if (xp[b] < gmin[b]) { point_inside = PETSC_FALSE; } 1310e2ec84fSDave May if (xp[b] > gmax[b]) { point_inside = PETSC_FALSE; } 1320e2ec84fSDave May } 1330e2ec84fSDave May if (point_inside) { 1340e2ec84fSDave May for (b=0; b<bs; b++) { 1350e2ec84fSDave May _pos[bs*n_estimate+b] = xp[b]; 1360e2ec84fSDave May } 1370e2ec84fSDave May n_estimate++; 1380e2ec84fSDave May } 1390e2ec84fSDave May } 1400e2ec84fSDave May } 1410e2ec84fSDave May } 1420e2ec84fSDave May ierr = VecRestoreArray(pos,&_pos);CHKERRQ(ierr); 1430e2ec84fSDave May 1440e2ec84fSDave May /* locate points */ 1450e2ec84fSDave May ierr = DMLocatePoints(celldm,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr); 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++) { 1740ca99263SDave May swarm_coor[bs*(n_curr + n_found) + b] = PetscRealPart(_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 PetscFunctionReturn(0); 1870e2ec84fSDave May } 1880e2ec84fSDave May 1890e2ec84fSDave May /*@C 1900e2ec84fSDave May DMSwarmSetPointCoordinates - Set point coordinates in a DMSwarm from a user defined list 1910e2ec84fSDave May 192d083f849SBarry Smith Collective on dm 1930e2ec84fSDave May 1940e2ec84fSDave May Input parameters: 1950e2ec84fSDave May + dm - the DMSwarm 1960e2ec84fSDave May . npoints - the number of points to insert 1970e2ec84fSDave May . coor - the coordinate values 1980e2ec84fSDave 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 1990e2ec84fSDave May - mode - indicates whether to append points to the swarm (ADD_VALUES), or over-ride existing points (INSERT_VALUES) 2000e2ec84fSDave May 2010e2ec84fSDave May Level: beginner 2020e2ec84fSDave May 2030e2ec84fSDave May Notes: 2040e2ec84fSDave May If the user has specified redundant = PETSC_FALSE, the cell DM will attempt to locate the coordinates provided by coor[] within 2050e2ec84fSDave 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 2060e2ec84fSDave May added to the DMSwarm. 2070e2ec84fSDave May 2080e2ec84fSDave May .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType, DMSwarmSetPointsUniformCoordinates() 2090e2ec84fSDave May @*/ 2100e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinates(DM dm,PetscInt npoints,PetscReal coor[],PetscBool redundant,InsertMode mode) 2110e2ec84fSDave May { 2120e2ec84fSDave May PetscErrorCode ierr; 2130e2ec84fSDave May PetscReal gmin[] = {PETSC_MAX_REAL ,PETSC_MAX_REAL, PETSC_MAX_REAL}; 2140e2ec84fSDave May PetscReal gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 2150e2ec84fSDave May PetscInt i,N,bs,b,n_estimate,n_curr,n_new_est,p,n_found; 2160e2ec84fSDave May Vec coorlocal; 2170e2ec84fSDave May const PetscScalar *_coor; 2180e2ec84fSDave May DM celldm; 2190e2ec84fSDave May Vec pos; 2200e2ec84fSDave May PetscScalar *_pos; 2210e2ec84fSDave May PetscReal *swarm_coor; 2220e2ec84fSDave May PetscInt *swarm_cellid; 2230e2ec84fSDave May PetscSF sfcell = NULL; 2240e2ec84fSDave May const PetscSFNode *LA_sfcell; 2250e2ec84fSDave May PetscReal *my_coor; 2260e2ec84fSDave May PetscInt my_npoints; 2270e2ec84fSDave May PetscMPIInt rank; 2280e2ec84fSDave May MPI_Comm comm; 2290e2ec84fSDave May 2300e2ec84fSDave May PetscFunctionBegin; 2310e2ec84fSDave May DMSWARMPICVALID(dm); 2320e2ec84fSDave May ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 233ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 2340e2ec84fSDave May 2350e2ec84fSDave May ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr); 2360e2ec84fSDave May ierr = DMGetCoordinatesLocal(celldm,&coorlocal);CHKERRQ(ierr); 2370e2ec84fSDave May ierr = VecGetSize(coorlocal,&N);CHKERRQ(ierr); 2380e2ec84fSDave May ierr = VecGetBlockSize(coorlocal,&bs);CHKERRQ(ierr); 2390e2ec84fSDave May N = N / bs; 2400e2ec84fSDave May ierr = VecGetArrayRead(coorlocal,&_coor);CHKERRQ(ierr); 2410e2ec84fSDave May for (i=0; i<N; i++) { 2420e2ec84fSDave May for (b=0; b<bs; b++) { 243a5f152d1SDave May gmin[b] = PetscMin(gmin[b],PetscRealPart(_coor[bs*i+b])); 244a5f152d1SDave May gmax[b] = PetscMax(gmax[b],PetscRealPart(_coor[bs*i+b])); 2450e2ec84fSDave May } 2460e2ec84fSDave May } 2470e2ec84fSDave May ierr = VecRestoreArrayRead(coorlocal,&_coor);CHKERRQ(ierr); 2480e2ec84fSDave May 2490e2ec84fSDave May /* broadcast points from rank 0 if requested */ 2500e2ec84fSDave May if (redundant) { 2510e2ec84fSDave May my_npoints = npoints; 252ffc4695bSBarry Smith ierr = MPI_Bcast(&my_npoints,1,MPIU_INT,0,comm);CHKERRMPI(ierr); 2530e2ec84fSDave May 2540e2ec84fSDave May if (rank > 0) { /* allocate space */ 25571844bb8SDave May ierr = PetscMalloc1(bs*my_npoints,&my_coor);CHKERRQ(ierr); 2560e2ec84fSDave May } else { 2570e2ec84fSDave May my_coor = coor; 2580e2ec84fSDave May } 259ffc4695bSBarry Smith ierr = MPI_Bcast(my_coor,bs*my_npoints,MPIU_REAL,0,comm);CHKERRMPI(ierr); 2600e2ec84fSDave May } else { 2610e2ec84fSDave May my_npoints = npoints; 2620e2ec84fSDave May my_coor = coor; 2630e2ec84fSDave May } 2640e2ec84fSDave May 2650e2ec84fSDave May /* determine the number of points living in the bounding box */ 2660e2ec84fSDave May n_estimate = 0; 2670e2ec84fSDave May for (i=0; i<my_npoints; i++) { 2680e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 2690e2ec84fSDave May 2700e2ec84fSDave May for (b=0; b<bs; b++) { 2710e2ec84fSDave May if (my_coor[bs*i+b] < gmin[b]) { point_inside = PETSC_FALSE; } 2720e2ec84fSDave May if (my_coor[bs*i+b] > gmax[b]) { point_inside = PETSC_FALSE; } 2730e2ec84fSDave May } 2740e2ec84fSDave May if (point_inside) { n_estimate++; } 2750e2ec84fSDave May } 2760e2ec84fSDave May 2770e2ec84fSDave May /* create candidate list */ 278c28a3c15SDave May ierr = VecCreate(PETSC_COMM_SELF,&pos);CHKERRQ(ierr); 2790e2ec84fSDave May ierr = VecSetSizes(pos,bs*n_estimate,PETSC_DECIDE);CHKERRQ(ierr); 2800e2ec84fSDave May ierr = VecSetBlockSize(pos,bs);CHKERRQ(ierr); 2810e2ec84fSDave May ierr = VecSetFromOptions(pos);CHKERRQ(ierr); 2820e2ec84fSDave May ierr = VecGetArray(pos,&_pos);CHKERRQ(ierr); 2830e2ec84fSDave May 2840e2ec84fSDave May n_estimate = 0; 2850e2ec84fSDave May for (i=0; i<my_npoints; i++) { 2860e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 2870e2ec84fSDave May 2880e2ec84fSDave May for (b=0; b<bs; b++) { 2890e2ec84fSDave May if (my_coor[bs*i+b] < gmin[b]) { point_inside = PETSC_FALSE; } 2900e2ec84fSDave May if (my_coor[bs*i+b] > gmax[b]) { point_inside = PETSC_FALSE; } 2910e2ec84fSDave May } 2920e2ec84fSDave May if (point_inside) { 2930e2ec84fSDave May for (b=0; b<bs; b++) { 2940e2ec84fSDave May _pos[bs*n_estimate+b] = my_coor[bs*i+b]; 2950e2ec84fSDave May } 2960e2ec84fSDave May n_estimate++; 2970e2ec84fSDave May } 2980e2ec84fSDave May } 2990e2ec84fSDave May ierr = VecRestoreArray(pos,&_pos);CHKERRQ(ierr); 3000e2ec84fSDave May 3010e2ec84fSDave May /* locate points */ 3020e2ec84fSDave May ierr = DMLocatePoints(celldm,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr); 3030e2ec84fSDave May 3040e2ec84fSDave May ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr); 3050e2ec84fSDave May n_found = 0; 3060e2ec84fSDave May for (p=0; p<n_estimate; p++) { 3070e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 3080e2ec84fSDave May n_found++; 3090e2ec84fSDave May } 3100e2ec84fSDave May } 3110e2ec84fSDave May 3120e2ec84fSDave May /* adjust size */ 3130e2ec84fSDave May if (mode == ADD_VALUES) { 3140e2ec84fSDave May ierr = DMSwarmGetLocalSize(dm,&n_curr);CHKERRQ(ierr); 3150e2ec84fSDave May n_new_est = n_curr + n_found; 3160e2ec84fSDave May ierr = DMSwarmSetLocalSizes(dm,n_new_est,-1);CHKERRQ(ierr); 3170e2ec84fSDave May } 3180e2ec84fSDave May if (mode == INSERT_VALUES) { 3190e2ec84fSDave May n_curr = 0; 3200e2ec84fSDave May n_new_est = n_found; 3210e2ec84fSDave May ierr = DMSwarmSetLocalSizes(dm,n_new_est,-1);CHKERRQ(ierr); 3220e2ec84fSDave May } 3230e2ec84fSDave May 3240e2ec84fSDave May /* initialize new coords, cell owners, pid */ 3250e2ec84fSDave May ierr = VecGetArrayRead(pos,&_coor);CHKERRQ(ierr); 3260e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 3270e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 3280e2ec84fSDave May n_found = 0; 3290e2ec84fSDave May for (p=0; p<n_estimate; p++) { 3300e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 3310e2ec84fSDave May for (b=0; b<bs; b++) { 3320ca99263SDave May swarm_coor[bs*(n_curr + n_found) + b] = PetscRealPart(_coor[bs*p+b]); 3330e2ec84fSDave May } 3340e2ec84fSDave May swarm_cellid[n_curr + n_found] = LA_sfcell[p].index; 3350e2ec84fSDave May n_found++; 3360e2ec84fSDave May } 3370e2ec84fSDave May } 3380e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 3390e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 3400e2ec84fSDave May ierr = VecRestoreArrayRead(pos,&_coor);CHKERRQ(ierr); 3410e2ec84fSDave May 3420e2ec84fSDave May if (redundant) { 3430e2ec84fSDave May if (rank > 0) { 3440e2ec84fSDave May ierr = PetscFree(my_coor);CHKERRQ(ierr); 3450e2ec84fSDave May } 3460e2ec84fSDave May } 3470e2ec84fSDave May ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr); 3480e2ec84fSDave May ierr = VecDestroy(&pos);CHKERRQ(ierr); 3490e2ec84fSDave May PetscFunctionReturn(0); 3500e2ec84fSDave May } 3510e2ec84fSDave May 3520e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM,DM,DMSwarmPICLayoutType,PetscInt); 3530e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM,DM,DMSwarmPICLayoutType,PetscInt); 3540e2ec84fSDave May 3550e2ec84fSDave May /*@C 3560e2ec84fSDave May DMSwarmInsertPointsUsingCellDM - Insert point coordinates within each cell 3570e2ec84fSDave May 3580e2ec84fSDave May Not collective 3590e2ec84fSDave May 3600e2ec84fSDave May Input parameters: 3610e2ec84fSDave May + dm - the DMSwarm 3620e2ec84fSDave May . layout_type - method used to fill each cell with the cell DM 3630e2ec84fSDave May - fill_param - parameter controlling how many points per cell are added (the meaning of this parameter is dependent on the layout type) 3640e2ec84fSDave May 3650e2ec84fSDave May Level: beginner 3660e2ec84fSDave May 3670e2ec84fSDave May Notes: 368e7af74c8SDave May 369e7af74c8SDave May The insert method will reset any previous defined points within the DMSwarm. 370e7af74c8SDave May 371e7af74c8SDave May When using a DMDA both 2D and 3D are supported for all layout types provided you are using DMDA_ELEMENT_Q1. 372e7af74c8SDave May 373e7af74c8SDave May When using a DMPLEX the following case are supported: 374ea3d7275SDave May (i) DMSWARMPIC_LAYOUT_REGULAR: 2D (triangle), 375ea3d7275SDave May (ii) DMSWARMPIC_LAYOUT_GAUSS: 2D and 3D provided the cell is a tri/tet or a quad/hex, 376ea3d7275SDave May (iii) DMSWARMPIC_LAYOUT_SUBDIVISION: 2D and 3D for quad/hex and 2D tri. 3770e2ec84fSDave May 3780e2ec84fSDave May .seealso: DMSwarmPICLayoutType, DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType 3790e2ec84fSDave May @*/ 3800e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmInsertPointsUsingCellDM(DM dm,DMSwarmPICLayoutType layout_type,PetscInt fill_param) 3810e2ec84fSDave May { 3820e2ec84fSDave May PetscErrorCode ierr; 3830e2ec84fSDave May DM celldm; 3840e2ec84fSDave May PetscBool isDA,isPLEX; 3850e2ec84fSDave May 3860e2ec84fSDave May PetscFunctionBegin; 3870e2ec84fSDave May DMSWARMPICVALID(dm); 3880e2ec84fSDave May ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr); 3890e2ec84fSDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMDA,&isDA);CHKERRQ(ierr); 3900e2ec84fSDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMPLEX,&isPLEX);CHKERRQ(ierr); 3910e2ec84fSDave May if (isDA) { 3920e2ec84fSDave May ierr = private_DMSwarmInsertPointsUsingCellDM_DA(dm,celldm,layout_type,fill_param);CHKERRQ(ierr); 3930e2ec84fSDave May } else if (isPLEX) { 3940e2ec84fSDave May ierr = private_DMSwarmInsertPointsUsingCellDM_PLEX(dm,celldm,layout_type,fill_param);CHKERRQ(ierr); 3950e2ec84fSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only supported for cell DMs of type DMDA and DMPLEX"); 3960e2ec84fSDave May PetscFunctionReturn(0); 3970e2ec84fSDave May } 3980e2ec84fSDave May 399431382f2SDave May extern PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM,DM,PetscInt,PetscReal*); 400431382f2SDave May 401431382f2SDave May /*@C 402431382f2SDave May DMSwarmSetPointCoordinatesCellwise - Insert point coordinates (defined over the reference cell) within each cell 403431382f2SDave May 404431382f2SDave May Not collective 405431382f2SDave May 406431382f2SDave May Input parameters: 407431382f2SDave May + dm - the DMSwarm 408431382f2SDave May . celldm - the cell DM 409431382f2SDave May . npoints - the number of points to insert in each cell 410431382f2SDave May - xi - the coordinates (defined in the local coordinate system for each cell) to insert 411431382f2SDave May 412431382f2SDave May Level: beginner 413431382f2SDave May 414431382f2SDave May Notes: 415431382f2SDave May The method will reset any previous defined points within the DMSwarm. 416ea3d7275SDave May Only supported for DMPLEX. If you are using a DMDA it is recommended to either use 417e7af74c8SDave May DMSwarmInsertPointsUsingCellDM(), or extract and set the coordinates yourself the following code 418e7af74c8SDave May 419e7af74c8SDave May $ PetscReal *coor; 420e7af74c8SDave May $ DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&coor); 421e7af74c8SDave May $ // user code to define the coordinates here 422e7af74c8SDave May $ DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&coor); 423431382f2SDave May 424431382f2SDave May .seealso: DMSwarmSetCellDM(), DMSwarmInsertPointsUsingCellDM() 425431382f2SDave May @*/ 426431382f2SDave May PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinatesCellwise(DM dm,PetscInt npoints,PetscReal xi[]) 4270e2ec84fSDave May { 428431382f2SDave May PetscErrorCode ierr; 429431382f2SDave May DM celldm; 430431382f2SDave May PetscBool isDA,isPLEX; 431431382f2SDave May 4320e2ec84fSDave May PetscFunctionBegin; 433431382f2SDave May DMSWARMPICVALID(dm); 434431382f2SDave May ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr); 435431382f2SDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMDA,&isDA);CHKERRQ(ierr); 436431382f2SDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMPLEX,&isPLEX);CHKERRQ(ierr); 437431382f2SDave May if (isDA) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only supported for cell DMs of type DMPLEX. Recommended you use DMSwarmInsertPointsUsingCellDM()"); 438431382f2SDave May else if (isPLEX) { 439431382f2SDave May ierr = private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm,celldm,npoints,xi);CHKERRQ(ierr); 440431382f2SDave May } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only supported for cell DMs of type DMDA and DMPLEX"); 4410e2ec84fSDave May PetscFunctionReturn(0); 4420e2ec84fSDave May } 443431382f2SDave May 4440e2ec84fSDave May /* Field projection API */ 44577048351SPatrick Sanan extern PetscErrorCode private_DMSwarmProjectFields_DA(DM swarm,DM celldm,PetscInt project_type,PetscInt nfields,DMSwarmDataField dfield[],Vec vecs[]); 44677048351SPatrick Sanan extern PetscErrorCode private_DMSwarmProjectFields_PLEX(DM swarm,DM celldm,PetscInt project_type,PetscInt nfields,DMSwarmDataField dfield[],Vec vecs[]); 44794f7d2dcSDave May 44894f7d2dcSDave May /*@C 44994f7d2dcSDave May DMSwarmProjectFields - Project a set of swarm fields onto the cell DM 45094f7d2dcSDave May 451d083f849SBarry Smith Collective on dm 45294f7d2dcSDave May 45394f7d2dcSDave May Input parameters: 45494f7d2dcSDave May + dm - the DMSwarm 45594f7d2dcSDave May . nfields - the number of swarm fields to project 45694f7d2dcSDave May . fieldnames - the textual names of the swarm fields to project 45794f7d2dcSDave May . fields - an array of Vec's of length nfields 45894f7d2dcSDave May - reuse - flag indicating whether the array and contents of fields should be re-used or internally allocated 45994f7d2dcSDave May 46094f7d2dcSDave May Currently, the only available projection method consists of 46194f7d2dcSDave May phi_i = \sum_{p=0}^{np} N_i(x_p) phi_p dJ / \sum_{p=0}^{np} N_i(x_p) dJ 46294f7d2dcSDave May where phi_p is the swarm field at point p, 46394f7d2dcSDave May N_i() is the cell DM basis function at vertex i, 46494f7d2dcSDave May dJ is the determinant of the cell Jacobian and 46594f7d2dcSDave May phi_i is the projected vertex value of the field phi. 46694f7d2dcSDave May 46794f7d2dcSDave May Level: beginner 46894f7d2dcSDave May 46994f7d2dcSDave May Notes: 470e7af74c8SDave May 471e7af74c8SDave May If reuse = PETSC_FALSE, this function will allocate the array of Vec's, and each individual Vec. 47294f7d2dcSDave May The user is responsible for destroying both the array and the individual Vec objects. 473e7af74c8SDave May 474e7af74c8SDave May Only swarm fields registered with data type = PETSC_REAL can be projected onto the cell DM. 475e7af74c8SDave May 476e7af74c8SDave May Only swarm fields of block size = 1 can currently be projected. 477e7af74c8SDave May 478e7af74c8SDave May The only projection methods currently only support the DA (2D) and PLEX (triangles 2D). 47994f7d2dcSDave May 48094f7d2dcSDave May .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType 48194f7d2dcSDave May @*/ 48294f7d2dcSDave May PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM dm,PetscInt nfields,const char *fieldnames[],Vec **fields,PetscBool reuse) 4830e2ec84fSDave May { 48494f7d2dcSDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 48577048351SPatrick Sanan DMSwarmDataField *gfield; 48694f7d2dcSDave May DM celldm; 48794f7d2dcSDave May PetscBool isDA,isPLEX; 48894f7d2dcSDave May Vec *vecs; 48994f7d2dcSDave May PetscInt f,nvecs; 49094f7d2dcSDave May PetscInt project_type = 0; 49194f7d2dcSDave May PetscErrorCode ierr; 49294f7d2dcSDave May 4930e2ec84fSDave May PetscFunctionBegin; 49494f7d2dcSDave May DMSWARMPICVALID(dm); 49594f7d2dcSDave May ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr); 49694f7d2dcSDave May ierr = PetscMalloc1(nfields,&gfield);CHKERRQ(ierr); 49794f7d2dcSDave May nvecs = 0; 49894f7d2dcSDave May for (f=0; f<nfields; f++) { 49977048351SPatrick Sanan ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldnames[f],&gfield[f]);CHKERRQ(ierr); 50094f7d2dcSDave May if (gfield[f]->petsc_type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Projection only valid for fields using a data type = PETSC_REAL"); 50194f7d2dcSDave May if (gfield[f]->bs != 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Projection only valid for fields with block size = 1"); 50294f7d2dcSDave May nvecs += gfield[f]->bs; 50394f7d2dcSDave May } 50494f7d2dcSDave May if (!reuse) { 50594f7d2dcSDave May ierr = PetscMalloc1(nvecs,&vecs);CHKERRQ(ierr); 50694f7d2dcSDave May for (f=0; f<nvecs; f++) { 50794f7d2dcSDave May ierr = DMCreateGlobalVector(celldm,&vecs[f]);CHKERRQ(ierr); 50894f7d2dcSDave May ierr = PetscObjectSetName((PetscObject)vecs[f],gfield[f]->name);CHKERRQ(ierr); 50994f7d2dcSDave May } 51094f7d2dcSDave May } else { 51194f7d2dcSDave May vecs = *fields; 51294f7d2dcSDave May } 51394f7d2dcSDave May 51494f7d2dcSDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMDA,&isDA);CHKERRQ(ierr); 51594f7d2dcSDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMPLEX,&isPLEX);CHKERRQ(ierr); 51694f7d2dcSDave May if (isDA) { 51794f7d2dcSDave May ierr = private_DMSwarmProjectFields_DA(dm,celldm,project_type,nfields,gfield,vecs);CHKERRQ(ierr); 51894f7d2dcSDave May } else if (isPLEX) { 51994f7d2dcSDave May ierr = private_DMSwarmProjectFields_PLEX(dm,celldm,project_type,nfields,gfield,vecs);CHKERRQ(ierr); 52094f7d2dcSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only supported for cell DMs of type DMDA and DMPLEX"); 52194f7d2dcSDave May 52294f7d2dcSDave May ierr = PetscFree(gfield);CHKERRQ(ierr); 52394f7d2dcSDave May if (!reuse) { 52494f7d2dcSDave May *fields = vecs; 52594f7d2dcSDave May } 5260e2ec84fSDave May PetscFunctionReturn(0); 5270e2ec84fSDave May } 5280e2ec84fSDave May 5290e2ec84fSDave May /*@C 530b799feefSDave May DMSwarmCreatePointPerCellCount - Count the number of points within all cells in the cell DM 5310e2ec84fSDave May 5320e2ec84fSDave May Not collective 5330e2ec84fSDave May 5340e2ec84fSDave May Input parameter: 5350e2ec84fSDave May . dm - the DMSwarm 5360e2ec84fSDave May 5370e2ec84fSDave May Output parameters: 5380e2ec84fSDave May + ncells - the number of cells in the cell DM (optional argument, pass NULL to ignore) 539b799feefSDave May - count - array of length ncells containing the number of points per cell 5400e2ec84fSDave May 5410e2ec84fSDave May Level: beginner 5420e2ec84fSDave May 5430e2ec84fSDave May Notes: 5440e2ec84fSDave May The array count is allocated internally and must be free'd by the user. 5450e2ec84fSDave May 5460e2ec84fSDave May .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType 5470e2ec84fSDave May @*/ 5480e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmCreatePointPerCellCount(DM dm,PetscInt *ncells,PetscInt **count) 5490e2ec84fSDave May { 550b799feefSDave May PetscErrorCode ierr; 551b799feefSDave May PetscBool isvalid; 552b799feefSDave May PetscInt nel; 553b799feefSDave May PetscInt *sum; 554b799feefSDave May 5550e2ec84fSDave May PetscFunctionBegin; 556b799feefSDave May ierr = DMSwarmSortGetIsValid(dm,&isvalid);CHKERRQ(ierr); 557b799feefSDave May nel = 0; 558b799feefSDave May if (isvalid) { 559b799feefSDave May PetscInt e; 560b799feefSDave May 561b799feefSDave May ierr = DMSwarmSortGetSizes(dm,&nel,NULL);CHKERRQ(ierr); 562b799feefSDave May 563b799feefSDave May ierr = PetscMalloc1(nel,&sum);CHKERRQ(ierr); 564b799feefSDave May for (e=0; e<nel; e++) { 565b799feefSDave May ierr = DMSwarmSortGetNumberOfPointsPerCell(dm,e,&sum[e]);CHKERRQ(ierr); 566b799feefSDave May } 567b799feefSDave May } else { 568b799feefSDave May DM celldm; 569b799feefSDave May PetscBool isda,isplex,isshell; 570b799feefSDave May PetscInt p,npoints; 571b799feefSDave May PetscInt *swarm_cellid; 572b799feefSDave May 573b799feefSDave May /* get the number of cells */ 574b799feefSDave May ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr); 575b799feefSDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMDA,&isda);CHKERRQ(ierr); 576b799feefSDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMPLEX,&isplex);CHKERRQ(ierr); 577b799feefSDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMSHELL,&isshell);CHKERRQ(ierr); 578b799feefSDave May if (isda) { 579b799feefSDave May PetscInt _nel,_npe; 580b799feefSDave May const PetscInt *_element; 581b799feefSDave May 582b799feefSDave May ierr = DMDAGetElements(celldm,&_nel,&_npe,&_element);CHKERRQ(ierr); 583b799feefSDave May nel = _nel; 584b799feefSDave May ierr = DMDARestoreElements(celldm,&_nel,&_npe,&_element);CHKERRQ(ierr); 585b799feefSDave May } else if (isplex) { 586b799feefSDave May PetscInt ps,pe; 587b799feefSDave May 588b799feefSDave May ierr = DMPlexGetHeightStratum(celldm,0,&ps,&pe);CHKERRQ(ierr); 589b799feefSDave May nel = pe - ps; 590b799feefSDave May } else if (isshell) { 591b799feefSDave May PetscErrorCode (*method_DMShellGetNumberOfCells)(DM,PetscInt*); 592b799feefSDave May 593b799feefSDave May ierr = PetscObjectQueryFunction((PetscObject)celldm,"DMGetNumberOfCells_C",&method_DMShellGetNumberOfCells);CHKERRQ(ierr); 594b799feefSDave May if (method_DMShellGetNumberOfCells) { 595b799feefSDave May ierr = method_DMShellGetNumberOfCells(celldm,&nel);CHKERRQ(ierr); 596b799feefSDave 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);"); 597b799feefSDave 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"); 598b799feefSDave May 599b799feefSDave May ierr = PetscMalloc1(nel,&sum);CHKERRQ(ierr); 600580bdb30SBarry Smith ierr = PetscArrayzero(sum,nel);CHKERRQ(ierr); 601b799feefSDave May ierr = DMSwarmGetLocalSize(dm,&npoints);CHKERRQ(ierr); 602b799feefSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 603b799feefSDave May for (p=0; p<npoints; p++) { 604b799feefSDave May if (swarm_cellid[p] != DMLOCATEPOINT_POINT_NOT_FOUND) { 605b799feefSDave May sum[ swarm_cellid[p] ]++; 606b799feefSDave May } 607b799feefSDave May } 608b799feefSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 609b799feefSDave May } 610b799feefSDave May if (ncells) { *ncells = nel; } 611b799feefSDave May *count = sum; 6120e2ec84fSDave May PetscFunctionReturn(0); 6130e2ec84fSDave May } 614