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> 7*94f7d2dcSDave May #include "data_bucket.h" 80e2ec84fSDave May 90e2ec84fSDave May /* 100e2ec84fSDave May Error chceking macto to ensure the swarm type is correct and that a cell DM has been set 110e2ec84fSDave May */ 120e2ec84fSDave May #define DMSWARMPICVALID(dm) \ 130e2ec84fSDave May { \ 140e2ec84fSDave May DM_Swarm *_swarm = (DM_Swarm*)(dm)->data; \ 150e2ec84fSDave 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)"); \ 160e2ec84fSDave May else \ 170e2ec84fSDave 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)"); \ 180e2ec84fSDave May } 190e2ec84fSDave May 200e2ec84fSDave May /* Coordinate insertition/addition API */ 210e2ec84fSDave May /*@C 220e2ec84fSDave May DMSwarmSetPointsUniformCoordinates - Set point coordinates in a DMSwarm on a regular (ijk) grid 230e2ec84fSDave May 240e2ec84fSDave May Collective on DM 250e2ec84fSDave May 260e2ec84fSDave May Input parameters: 270e2ec84fSDave May + dm - the DMSwarm 280e2ec84fSDave May . min - minimum coordinate values in the x, y, z directions (array of length dim) 290e2ec84fSDave May . max - maximum coordinate values in the x, y, z directions (array of length dim) 300e2ec84fSDave May . npoints - number of points in each spatial direction (array of length dim) 310e2ec84fSDave May - mode - indicates whether to append points to the swarm (ADD_VALUES), or over-ride existing points (INSERT_VALUES) 320e2ec84fSDave May 330e2ec84fSDave May Level: beginner 340e2ec84fSDave May 350e2ec84fSDave May Notes: 360e2ec84fSDave May When using mode = INSERT_VALUES, this method will reset the number of particles in the DMSwarm 370e2ec84fSDave May to be npoints[0]*npoints[1] (2D) or npoints[0]*npoints[1]*npoints[2] (3D). When using mode = ADD_VALUES, 380e2ec84fSDave May new points will be appended to any already existing in the DMSwarm 390e2ec84fSDave May 400e2ec84fSDave May .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType 410e2ec84fSDave May @*/ 420e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmSetPointsUniformCoordinates(DM dm,PetscReal min[],PetscReal max[],PetscInt npoints[],InsertMode mode) 430e2ec84fSDave May { 440e2ec84fSDave May PetscErrorCode ierr; 450e2ec84fSDave May PetscReal gmin[] = {PETSC_MAX_REAL ,PETSC_MAX_REAL, PETSC_MAX_REAL}; 460e2ec84fSDave May PetscReal gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 470e2ec84fSDave May PetscInt i,j,k,N,bs,b,n_estimate,n_curr,n_new_est,p,n_found; 480e2ec84fSDave May Vec coorlocal; 490e2ec84fSDave May const PetscScalar *_coor; 500e2ec84fSDave May DM celldm; 510e2ec84fSDave May PetscReal dx[3]; 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++) { 690e2ec84fSDave May gmin[b] = PetscMin(gmin[b],_coor[bs*i+b]); 700e2ec84fSDave May gmax[b] = PetscMax(gmax[b],_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 } 810e2ec84fSDave May } 820e2ec84fSDave May 830e2ec84fSDave May /* determine number of points living in the bounding box */ 840e2ec84fSDave May n_estimate = 0; 850e2ec84fSDave May if (bs == 2) { npoints[2] = 1; } 860e2ec84fSDave May for (k=0; k<npoints[2]; k++) { 870e2ec84fSDave May for (j=0; j<npoints[1]; j++) { 880e2ec84fSDave 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 */ 1090e2ec84fSDave May ierr = VecCreate(PetscObjectComm((PetscObject)dm),&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; 1160e2ec84fSDave May for (k=0; k<npoints[2]; k++) { 1170e2ec84fSDave May for (j=0; j<npoints[1]; j++) { 1180e2ec84fSDave 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 1470e2ec84fSDave May ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr); 1480e2ec84fSDave May n_found = 0; 1490e2ec84fSDave May for (p=0; p<n_estimate; p++) { 1500e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 1510e2ec84fSDave May n_found++; 1520e2ec84fSDave May } 1530e2ec84fSDave May } 1540e2ec84fSDave May 1550e2ec84fSDave May /* adjust size */ 1560e2ec84fSDave May if (mode == ADD_VALUES) { 1570e2ec84fSDave May ierr = DMSwarmGetLocalSize(dm,&n_curr);CHKERRQ(ierr); 1580e2ec84fSDave May n_new_est = n_curr + n_found; 1590e2ec84fSDave May ierr = DMSwarmSetLocalSizes(dm,n_new_est,-1);CHKERRQ(ierr); 1600e2ec84fSDave May } 1610e2ec84fSDave May if (mode == INSERT_VALUES) { 1620e2ec84fSDave May n_curr = 0; 1630e2ec84fSDave May n_new_est = n_found; 1640e2ec84fSDave May ierr = DMSwarmSetLocalSizes(dm,n_new_est,-1);CHKERRQ(ierr); 1650e2ec84fSDave May } 1660e2ec84fSDave May 1670e2ec84fSDave May /* initialize new coords, cell owners, pid */ 1680e2ec84fSDave May ierr = VecGetArrayRead(pos,&_coor);CHKERRQ(ierr); 1690e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 1700e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 1710e2ec84fSDave May n_found = 0; 1720e2ec84fSDave May for (p=0; p<n_estimate; p++) { 1730e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 1740e2ec84fSDave May for (b=0; b<bs; b++) { 1750e2ec84fSDave May swarm_coor[bs*(n_curr + n_found) + b] = _coor[bs*p+b]; 1760e2ec84fSDave May } 1770e2ec84fSDave May swarm_cellid[n_curr + n_found] = LA_sfcell[p].index; 1780e2ec84fSDave May n_found++; 1790e2ec84fSDave May } 1800e2ec84fSDave May } 1810e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 1820e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 1830e2ec84fSDave May ierr = VecRestoreArrayRead(pos,&_coor);CHKERRQ(ierr); 1840e2ec84fSDave May 1850e2ec84fSDave May ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr); 1860e2ec84fSDave May ierr = VecDestroy(&pos);CHKERRQ(ierr); 1870e2ec84fSDave May 1880e2ec84fSDave May PetscFunctionReturn(0); 1890e2ec84fSDave May } 1900e2ec84fSDave May 1910e2ec84fSDave May /*@C 1920e2ec84fSDave May DMSwarmSetPointCoordinates - Set point coordinates in a DMSwarm from a user defined list 1930e2ec84fSDave May 1940e2ec84fSDave May Collective on DM 1950e2ec84fSDave May 1960e2ec84fSDave May Input parameters: 1970e2ec84fSDave May + dm - the DMSwarm 1980e2ec84fSDave May . npoints - the number of points to insert 1990e2ec84fSDave May . coor - the coordinate values 2000e2ec84fSDave 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 2010e2ec84fSDave May - mode - indicates whether to append points to the swarm (ADD_VALUES), or over-ride existing points (INSERT_VALUES) 2020e2ec84fSDave May 2030e2ec84fSDave May Level: beginner 2040e2ec84fSDave May 2050e2ec84fSDave May Notes: 2060e2ec84fSDave May If the user has specified redundant = PETSC_FALSE, the cell DM will attempt to locate the coordinates provided by coor[] within 2070e2ec84fSDave 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 2080e2ec84fSDave May added to the DMSwarm. 2090e2ec84fSDave May 2100e2ec84fSDave May .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType, DMSwarmSetPointsUniformCoordinates() 2110e2ec84fSDave May @*/ 2120e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinates(DM dm,PetscInt npoints,PetscReal coor[],PetscBool redundant,InsertMode mode) 2130e2ec84fSDave May { 2140e2ec84fSDave May PetscErrorCode ierr; 2150e2ec84fSDave May PetscReal gmin[] = {PETSC_MAX_REAL ,PETSC_MAX_REAL, PETSC_MAX_REAL}; 2160e2ec84fSDave May PetscReal gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 2170e2ec84fSDave May PetscInt i,N,bs,b,n_estimate,n_curr,n_new_est,p,n_found; 2180e2ec84fSDave May Vec coorlocal; 2190e2ec84fSDave May const PetscScalar *_coor; 2200e2ec84fSDave May DM celldm; 2210e2ec84fSDave May Vec pos; 2220e2ec84fSDave May PetscScalar *_pos; 2230e2ec84fSDave May PetscReal *swarm_coor; 2240e2ec84fSDave May PetscInt *swarm_cellid; 2250e2ec84fSDave May PetscSF sfcell = NULL; 2260e2ec84fSDave May const PetscSFNode *LA_sfcell; 2270e2ec84fSDave May PetscReal *my_coor; 2280e2ec84fSDave May PetscInt my_npoints; 2290e2ec84fSDave May PetscMPIInt rank; 2300e2ec84fSDave May MPI_Comm comm; 2310e2ec84fSDave May 2320e2ec84fSDave May PetscFunctionBegin; 2330e2ec84fSDave May DMSWARMPICVALID(dm); 2340e2ec84fSDave May ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 2350e2ec84fSDave May ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 2360e2ec84fSDave May 2370e2ec84fSDave May ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr); 2380e2ec84fSDave May ierr = DMGetCoordinatesLocal(celldm,&coorlocal);CHKERRQ(ierr); 2390e2ec84fSDave May ierr = VecGetSize(coorlocal,&N);CHKERRQ(ierr); 2400e2ec84fSDave May ierr = VecGetBlockSize(coorlocal,&bs);CHKERRQ(ierr); 2410e2ec84fSDave May N = N / bs; 2420e2ec84fSDave May ierr = VecGetArrayRead(coorlocal,&_coor);CHKERRQ(ierr); 2430e2ec84fSDave May for (i=0; i<N; i++) { 2440e2ec84fSDave May for (b=0; b<bs; b++) { 2450e2ec84fSDave May gmin[b] = PetscMin(gmin[b],_coor[bs*i+b]); 2460e2ec84fSDave May gmax[b] = PetscMax(gmax[b],_coor[bs*i+b]); 2470e2ec84fSDave May } 2480e2ec84fSDave May } 2490e2ec84fSDave May ierr = VecRestoreArrayRead(coorlocal,&_coor);CHKERRQ(ierr); 2500e2ec84fSDave May 2510e2ec84fSDave May /* broadcast points from rank 0 if requested */ 2520e2ec84fSDave May if (redundant) { 2530e2ec84fSDave May my_npoints = npoints; 2540e2ec84fSDave May ierr = MPI_Bcast(&my_npoints,1,MPIU_INT,0,comm);CHKERRQ(ierr); 2550e2ec84fSDave May 2560e2ec84fSDave May if (rank > 0) { /* allocate space */ 25771844bb8SDave May ierr = PetscMalloc1(bs*my_npoints,&my_coor);CHKERRQ(ierr); 2580e2ec84fSDave May } else { 2590e2ec84fSDave May my_coor = coor; 2600e2ec84fSDave May } 2610e2ec84fSDave May ierr = MPI_Bcast(my_coor,bs*my_npoints,MPIU_REAL,0,comm);CHKERRQ(ierr); 2620e2ec84fSDave May } else { 2630e2ec84fSDave May my_npoints = npoints; 2640e2ec84fSDave May my_coor = coor; 2650e2ec84fSDave May } 2660e2ec84fSDave May 2670e2ec84fSDave May /* determine the number of points living in the bounding box */ 2680e2ec84fSDave May n_estimate = 0; 2690e2ec84fSDave May for (i=0; i<my_npoints; i++) { 2700e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 2710e2ec84fSDave May 2720e2ec84fSDave May for (b=0; b<bs; b++) { 2730e2ec84fSDave May if (my_coor[bs*i+b] < gmin[b]) { point_inside = PETSC_FALSE; } 2740e2ec84fSDave May if (my_coor[bs*i+b] > gmax[b]) { point_inside = PETSC_FALSE; } 2750e2ec84fSDave May } 2760e2ec84fSDave May if (point_inside) { n_estimate++; } 2770e2ec84fSDave May } 2780e2ec84fSDave May 2790e2ec84fSDave May /* create candidate list */ 2800e2ec84fSDave May ierr = VecCreate(PetscObjectComm((PetscObject)dm),&pos);CHKERRQ(ierr); 2810e2ec84fSDave May ierr = VecSetSizes(pos,bs*n_estimate,PETSC_DECIDE);CHKERRQ(ierr); 2820e2ec84fSDave May ierr = VecSetBlockSize(pos,bs);CHKERRQ(ierr); 2830e2ec84fSDave May ierr = VecSetFromOptions(pos);CHKERRQ(ierr); 2840e2ec84fSDave May ierr = VecGetArray(pos,&_pos);CHKERRQ(ierr); 2850e2ec84fSDave May 2860e2ec84fSDave May n_estimate = 0; 2870e2ec84fSDave May for (i=0; i<my_npoints; i++) { 2880e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 2890e2ec84fSDave May 2900e2ec84fSDave May for (b=0; b<bs; b++) { 2910e2ec84fSDave May if (my_coor[bs*i+b] < gmin[b]) { point_inside = PETSC_FALSE; } 2920e2ec84fSDave May if (my_coor[bs*i+b] > gmax[b]) { point_inside = PETSC_FALSE; } 2930e2ec84fSDave May } 2940e2ec84fSDave May if (point_inside) { 2950e2ec84fSDave May for (b=0; b<bs; b++) { 2960e2ec84fSDave May _pos[bs*n_estimate+b] = my_coor[bs*i+b]; 2970e2ec84fSDave May } 2980e2ec84fSDave May n_estimate++; 2990e2ec84fSDave May } 3000e2ec84fSDave May } 3010e2ec84fSDave May ierr = VecRestoreArray(pos,&_pos);CHKERRQ(ierr); 3020e2ec84fSDave May 3030e2ec84fSDave May /* locate points */ 3040e2ec84fSDave May ierr = DMLocatePoints(celldm,pos,DM_POINTLOCATION_NONE,&sfcell);CHKERRQ(ierr); 3050e2ec84fSDave May 3060e2ec84fSDave May ierr = PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell);CHKERRQ(ierr); 3070e2ec84fSDave May n_found = 0; 3080e2ec84fSDave May for (p=0; p<n_estimate; p++) { 3090e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 3100e2ec84fSDave May n_found++; 3110e2ec84fSDave May } 3120e2ec84fSDave May } 3130e2ec84fSDave May 3140e2ec84fSDave May /* adjust size */ 3150e2ec84fSDave May if (mode == ADD_VALUES) { 3160e2ec84fSDave May ierr = DMSwarmGetLocalSize(dm,&n_curr);CHKERRQ(ierr); 3170e2ec84fSDave May n_new_est = n_curr + n_found; 3180e2ec84fSDave May ierr = DMSwarmSetLocalSizes(dm,n_new_est,-1);CHKERRQ(ierr); 3190e2ec84fSDave May } 3200e2ec84fSDave May if (mode == INSERT_VALUES) { 3210e2ec84fSDave May n_curr = 0; 3220e2ec84fSDave May n_new_est = n_found; 3230e2ec84fSDave May ierr = DMSwarmSetLocalSizes(dm,n_new_est,-1);CHKERRQ(ierr); 3240e2ec84fSDave May } 3250e2ec84fSDave May 3260e2ec84fSDave May /* initialize new coords, cell owners, pid */ 3270e2ec84fSDave May ierr = VecGetArrayRead(pos,&_coor);CHKERRQ(ierr); 3280e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 3290e2ec84fSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 3300e2ec84fSDave May n_found = 0; 3310e2ec84fSDave May for (p=0; p<n_estimate; p++) { 3320e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 3330e2ec84fSDave May for (b=0; b<bs; b++) { 3340e2ec84fSDave May swarm_coor[bs*(n_curr + n_found) + b] = _coor[bs*p+b]; 3350e2ec84fSDave May } 3360e2ec84fSDave May swarm_cellid[n_curr + n_found] = LA_sfcell[p].index; 3370e2ec84fSDave May n_found++; 3380e2ec84fSDave May } 3390e2ec84fSDave May } 3400e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 3410e2ec84fSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&swarm_coor);CHKERRQ(ierr); 3420e2ec84fSDave May ierr = VecRestoreArrayRead(pos,&_coor);CHKERRQ(ierr); 3430e2ec84fSDave May 3440e2ec84fSDave May if (redundant) { 3450e2ec84fSDave May if (rank > 0) { 3460e2ec84fSDave May ierr = PetscFree(my_coor);CHKERRQ(ierr); 3470e2ec84fSDave May } 3480e2ec84fSDave May } 3490e2ec84fSDave May ierr = PetscSFDestroy(&sfcell);CHKERRQ(ierr); 3500e2ec84fSDave May ierr = VecDestroy(&pos);CHKERRQ(ierr); 3510e2ec84fSDave May 3520e2ec84fSDave May PetscFunctionReturn(0); 3530e2ec84fSDave May } 3540e2ec84fSDave May 3550e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM,DM,DMSwarmPICLayoutType,PetscInt); 3560e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM,DM,DMSwarmPICLayoutType,PetscInt); 3570e2ec84fSDave May 3580e2ec84fSDave May /*@C 3590e2ec84fSDave May DMSwarmInsertPointsUsingCellDM - Insert point coordinates within each cell 3600e2ec84fSDave May 3610e2ec84fSDave May Not collective 3620e2ec84fSDave May 3630e2ec84fSDave May Input parameters: 3640e2ec84fSDave May + dm - the DMSwarm 3650e2ec84fSDave May . layout_type - method used to fill each cell with the cell DM 3660e2ec84fSDave May - fill_param - parameter controlling how many points per cell are added (the meaning of this parameter is dependent on the layout type) 3670e2ec84fSDave May 3680e2ec84fSDave May Level: beginner 3690e2ec84fSDave May 3700e2ec84fSDave May Notes: 3710e2ec84fSDave May The insert method will reset any previous defined points within the DMSwarm 3720e2ec84fSDave May 3730e2ec84fSDave May .seealso: DMSwarmPICLayoutType, DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType 3740e2ec84fSDave May @*/ 3750e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmInsertPointsUsingCellDM(DM dm,DMSwarmPICLayoutType layout_type,PetscInt fill_param) 3760e2ec84fSDave May { 3770e2ec84fSDave May PetscErrorCode ierr; 3780e2ec84fSDave May DM celldm; 3790e2ec84fSDave May PetscBool isDA,isPLEX; 3800e2ec84fSDave May 3810e2ec84fSDave May PetscFunctionBegin; 3820e2ec84fSDave May DMSWARMPICVALID(dm); 3830e2ec84fSDave May ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr); 3840e2ec84fSDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMDA,&isDA);CHKERRQ(ierr); 3850e2ec84fSDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMPLEX,&isPLEX);CHKERRQ(ierr); 3860e2ec84fSDave May if (isDA) { 3870e2ec84fSDave May ierr = private_DMSwarmInsertPointsUsingCellDM_DA(dm,celldm,layout_type,fill_param);CHKERRQ(ierr); 3880e2ec84fSDave May } else if (isPLEX) { 3890e2ec84fSDave May ierr = private_DMSwarmInsertPointsUsingCellDM_PLEX(dm,celldm,layout_type,fill_param);CHKERRQ(ierr); 3900e2ec84fSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only supported for cell DMs of type DMDA and DMPLEX"); 3910e2ec84fSDave May 3920e2ec84fSDave May PetscFunctionReturn(0); 3930e2ec84fSDave May } 3940e2ec84fSDave May 3950e2ec84fSDave May /* 3960e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmAddPointCoordinatesCellWise(DM dm,PetscInt cell,PetscInt npoints,PetscReal xi[],PetscBool proximity_initialization) 3970e2ec84fSDave May { 3980e2ec84fSDave May PetscFunctionBegin; 3990e2ec84fSDave May PetscFunctionReturn(0); 4000e2ec84fSDave May } 4010e2ec84fSDave May */ 4020e2ec84fSDave May 4030e2ec84fSDave May /* Field projection API */ 404*94f7d2dcSDave May extern PetscErrorCode private_DMSwarmProjectFields_DA(DM swarm,DM celldm,PetscInt project_type,PetscInt nfields,DataField dfield[],Vec vecs[]); 405*94f7d2dcSDave May extern PetscErrorCode private_DMSwarmProjectFields_PLEX(DM swarm,DM celldm,PetscInt project_type,PetscInt nfields,DataField dfield[],Vec vecs[]); 406*94f7d2dcSDave May 407*94f7d2dcSDave May /*@C 408*94f7d2dcSDave May DMSwarmProjectFields - Project a set of swarm fields onto the cell DM 409*94f7d2dcSDave May 410*94f7d2dcSDave May Collective on Vec 411*94f7d2dcSDave May 412*94f7d2dcSDave May Input parameters: 413*94f7d2dcSDave May + dm - the DMSwarm 414*94f7d2dcSDave May . nfields - the number of swarm fields to project 415*94f7d2dcSDave May . fieldnames - the textual names of the swarm fields to project 416*94f7d2dcSDave May . fields - an array of Vec's of length nfields 417*94f7d2dcSDave May - reuse - flag indicating whether the array and contents of fields should be re-used or internally allocated 418*94f7d2dcSDave May 419*94f7d2dcSDave May Currently, the only available projection method consists of 420*94f7d2dcSDave May phi_i = \sum_{p=0}^{np} N_i(x_p) phi_p dJ / \sum_{p=0}^{np} N_i(x_p) dJ 421*94f7d2dcSDave May where phi_p is the swarm field at point p, 422*94f7d2dcSDave May N_i() is the cell DM basis function at vertex i, 423*94f7d2dcSDave May dJ is the determinant of the cell Jacobian and 424*94f7d2dcSDave May phi_i is the projected vertex value of the field phi. 425*94f7d2dcSDave May 426*94f7d2dcSDave May Level: beginner 427*94f7d2dcSDave May 428*94f7d2dcSDave May Notes: 429*94f7d2dcSDave May - If reuse = PETSC_FALSE, this function will allocate the array of Vec's, and each individual Vec. 430*94f7d2dcSDave May The user is responsible for destroying both the array and the individual Vec objects. 431*94f7d2dcSDave May - Only swarm fields registered with data type = PETSC_REAL can be projected onto the cell DM. 432*94f7d2dcSDave May - Only swarm fields of block size = 1 can currently be projected. 433*94f7d2dcSDave May - The only projection methods currently only support the DA (2D) and PLEX (triangles 2D). 434*94f7d2dcSDave May 435*94f7d2dcSDave May .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType 436*94f7d2dcSDave May @*/ 437*94f7d2dcSDave May PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM dm,PetscInt nfields,const char *fieldnames[],Vec **fields,PetscBool reuse) 4380e2ec84fSDave May { 439*94f7d2dcSDave May DM_Swarm *swarm = (DM_Swarm*)dm->data; 440*94f7d2dcSDave May DataField *gfield; 441*94f7d2dcSDave May DM celldm; 442*94f7d2dcSDave May PetscBool isDA,isPLEX; 443*94f7d2dcSDave May Vec *vecs; 444*94f7d2dcSDave May PetscInt f,nvecs; 445*94f7d2dcSDave May PetscInt project_type = 0; 446*94f7d2dcSDave May PetscErrorCode ierr; 447*94f7d2dcSDave May 4480e2ec84fSDave May PetscFunctionBegin; 449*94f7d2dcSDave May DMSWARMPICVALID(dm); 450*94f7d2dcSDave May ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr); 451*94f7d2dcSDave May ierr = PetscMalloc1(nfields,&gfield);CHKERRQ(ierr); 452*94f7d2dcSDave May nvecs = 0; 453*94f7d2dcSDave May for (f=0; f<nfields; f++) { 454*94f7d2dcSDave May ierr = DataBucketGetDataFieldByName(swarm->db,fieldnames[f],&gfield[f]);CHKERRQ(ierr); 455*94f7d2dcSDave 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"); 456*94f7d2dcSDave May if (gfield[f]->bs != 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Projection only valid for fields with block size = 1"); 457*94f7d2dcSDave May nvecs += gfield[f]->bs; 458*94f7d2dcSDave May } 459*94f7d2dcSDave May if (!reuse) { 460*94f7d2dcSDave May ierr = PetscMalloc1(nvecs,&vecs);CHKERRQ(ierr); 461*94f7d2dcSDave May for (f=0; f<nvecs; f++) { 462*94f7d2dcSDave May ierr = DMCreateGlobalVector(celldm,&vecs[f]);CHKERRQ(ierr); 463*94f7d2dcSDave May ierr = PetscObjectSetName((PetscObject)vecs[f],gfield[f]->name);CHKERRQ(ierr); 464*94f7d2dcSDave May } 465*94f7d2dcSDave May } else { 466*94f7d2dcSDave May vecs = *fields; 467*94f7d2dcSDave May } 468*94f7d2dcSDave May 469*94f7d2dcSDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMDA,&isDA);CHKERRQ(ierr); 470*94f7d2dcSDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMPLEX,&isPLEX);CHKERRQ(ierr); 471*94f7d2dcSDave May if (isDA) { 472*94f7d2dcSDave May ierr = private_DMSwarmProjectFields_DA(dm,celldm,project_type,nfields,gfield,vecs);CHKERRQ(ierr); 473*94f7d2dcSDave May } else if (isPLEX) { 474*94f7d2dcSDave May ierr = private_DMSwarmProjectFields_PLEX(dm,celldm,project_type,nfields,gfield,vecs);CHKERRQ(ierr); 475*94f7d2dcSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only supported for cell DMs of type DMDA and DMPLEX"); 476*94f7d2dcSDave May 477*94f7d2dcSDave May ierr = PetscFree(gfield);CHKERRQ(ierr); 478*94f7d2dcSDave May if (!reuse) { 479*94f7d2dcSDave May *fields = vecs; 480*94f7d2dcSDave May } 481*94f7d2dcSDave May 4820e2ec84fSDave May PetscFunctionReturn(0); 4830e2ec84fSDave May } 4840e2ec84fSDave May 4850e2ec84fSDave May /*@C 486b799feefSDave May DMSwarmCreatePointPerCellCount - Count the number of points within all cells in the cell DM 4870e2ec84fSDave May 4880e2ec84fSDave May Not collective 4890e2ec84fSDave May 4900e2ec84fSDave May Input parameter: 4910e2ec84fSDave May . dm - the DMSwarm 4920e2ec84fSDave May 4930e2ec84fSDave May Output parameters: 4940e2ec84fSDave May + ncells - the number of cells in the cell DM (optional argument, pass NULL to ignore) 495b799feefSDave May - count - array of length ncells containing the number of points per cell 4960e2ec84fSDave May 4970e2ec84fSDave May Level: beginner 4980e2ec84fSDave May 4990e2ec84fSDave May Notes: 5000e2ec84fSDave May The array count is allocated internally and must be free'd by the user. 5010e2ec84fSDave May 5020e2ec84fSDave May .seealso: DMSwarmSetType(), DMSwarmSetCellDM(), DMSwarmType 5030e2ec84fSDave May @*/ 5040e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmCreatePointPerCellCount(DM dm,PetscInt *ncells,PetscInt **count) 5050e2ec84fSDave May { 506b799feefSDave May PetscErrorCode ierr; 507b799feefSDave May PetscBool isvalid; 508b799feefSDave May PetscInt nel; 509b799feefSDave May PetscInt *sum; 510b799feefSDave May 5110e2ec84fSDave May PetscFunctionBegin; 512b799feefSDave May ierr = DMSwarmSortGetIsValid(dm,&isvalid);CHKERRQ(ierr); 513b799feefSDave May nel = 0; 514b799feefSDave May if (isvalid) { 515b799feefSDave May PetscInt e; 516b799feefSDave May 517b799feefSDave May ierr = DMSwarmSortGetSizes(dm,&nel,NULL);CHKERRQ(ierr); 518b799feefSDave May 519b799feefSDave May ierr = PetscMalloc1(nel,&sum);CHKERRQ(ierr); 520b799feefSDave May for (e=0; e<nel; e++) { 521b799feefSDave May ierr = DMSwarmSortGetNumberOfPointsPerCell(dm,e,&sum[e]);CHKERRQ(ierr); 522b799feefSDave May } 523b799feefSDave May } else { 524b799feefSDave May DM celldm; 525b799feefSDave May PetscBool isda,isplex,isshell; 526b799feefSDave May PetscInt p,npoints; 527b799feefSDave May PetscInt *swarm_cellid; 528b799feefSDave May 529b799feefSDave May /* get the number of cells */ 530b799feefSDave May ierr = DMSwarmGetCellDM(dm,&celldm);CHKERRQ(ierr); 531b799feefSDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMDA,&isda);CHKERRQ(ierr); 532b799feefSDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMPLEX,&isplex);CHKERRQ(ierr); 533b799feefSDave May ierr = PetscObjectTypeCompare((PetscObject)celldm,DMSHELL,&isshell);CHKERRQ(ierr); 534b799feefSDave May if (isda) { 535b799feefSDave May PetscInt _nel,_npe; 536b799feefSDave May const PetscInt *_element; 537b799feefSDave May 538b799feefSDave May ierr = DMDAGetElements(celldm,&_nel,&_npe,&_element);CHKERRQ(ierr); 539b799feefSDave May nel = _nel; 540b799feefSDave May ierr = DMDARestoreElements(celldm,&_nel,&_npe,&_element);CHKERRQ(ierr); 541b799feefSDave May } else if (isplex) { 542b799feefSDave May PetscInt ps,pe; 543b799feefSDave May 544b799feefSDave May ierr = DMPlexGetHeightStratum(celldm,0,&ps,&pe);CHKERRQ(ierr); 545b799feefSDave May nel = pe - ps; 546b799feefSDave May } else if (isshell) { 547b799feefSDave May PetscErrorCode (*method_DMShellGetNumberOfCells)(DM,PetscInt*); 548b799feefSDave May 549b799feefSDave May ierr = PetscObjectQueryFunction((PetscObject)celldm,"DMGetNumberOfCells_C",&method_DMShellGetNumberOfCells);CHKERRQ(ierr); 550b799feefSDave May if (method_DMShellGetNumberOfCells) { 551b799feefSDave May ierr = method_DMShellGetNumberOfCells(celldm,&nel);CHKERRQ(ierr); 552b799feefSDave 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 );"); 553b799feefSDave 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"); 554b799feefSDave May 555b799feefSDave May ierr = PetscMalloc1(nel,&sum);CHKERRQ(ierr); 556b799feefSDave May ierr = PetscMemzero(sum,sizeof(PetscInt)*nel);CHKERRQ(ierr); 557b799feefSDave May ierr = DMSwarmGetLocalSize(dm,&npoints);CHKERRQ(ierr); 558b799feefSDave May ierr = DMSwarmGetField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 559b799feefSDave May for (p=0; p<npoints; p++) { 560b799feefSDave May if (swarm_cellid[p] != DMLOCATEPOINT_POINT_NOT_FOUND) { 561b799feefSDave May sum[ swarm_cellid[p] ]++; 562b799feefSDave May } 563b799feefSDave May } 564b799feefSDave May ierr = DMSwarmRestoreField(dm,DMSwarmPICField_cellid,NULL,NULL,(void**)&swarm_cellid);CHKERRQ(ierr); 565b799feefSDave May } 566b799feefSDave May if (ncells) { *ncells = nel; } 567b799feefSDave May *count = sum; 5680e2ec84fSDave May PetscFunctionReturn(0); 5690e2ec84fSDave May } 570