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> 635b38c60SMatthew G. Knepley #include <petscdt.h> 7279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h" 80e2ec84fSDave May 935b38c60SMatthew G. Knepley #include <petsc/private/petscfeimpl.h> /* For CoordinatesRefToReal() */ 1035b38c60SMatthew G. Knepley 110e2ec84fSDave May /* Coordinate insertition/addition API */ 12cc4c1da9SBarry Smith /*@ 1320f4b53cSBarry Smith DMSwarmSetPointsUniformCoordinates - Set point coordinates in a `DMSWARM` on a regular (ijk) grid 140e2ec84fSDave May 1520f4b53cSBarry Smith Collective 160e2ec84fSDave May 1760225df5SJacob Faibussowitsch Input Parameters: 1820f4b53cSBarry Smith + dm - the `DMSWARM` 190e2ec84fSDave May . min - minimum coordinate values in the x, y, z directions (array of length dim) 200e2ec84fSDave May . max - maximum coordinate values in the x, y, z directions (array of length dim) 210e2ec84fSDave May . npoints - number of points in each spatial direction (array of length dim) 2220f4b53cSBarry Smith - mode - indicates whether to append points to the swarm (`ADD_VALUES`), or over-ride existing points (`INSERT_VALUES`) 230e2ec84fSDave May 240e2ec84fSDave May Level: beginner 250e2ec84fSDave May 260e2ec84fSDave May Notes: 2720f4b53cSBarry Smith When using mode = `INSERT_VALUES`, this method will reset the number of particles in the `DMSWARM` 28a3b724e8SBarry Smith to be `npoints[0]` x `npoints[1]` (2D) or `npoints[0]` x `npoints[1]` x `npoints[2]` (3D). When using mode = `ADD_VALUES`, 2920f4b53cSBarry Smith new points will be appended to any already existing in the `DMSWARM` 300e2ec84fSDave May 3120f4b53cSBarry Smith .seealso: `DM`, `DMSWARM`, `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType` 320e2ec84fSDave May @*/ 33cc4c1da9SBarry Smith PetscErrorCode DMSwarmSetPointsUniformCoordinates(DM dm, PetscReal min[], PetscReal max[], PetscInt npoints[], InsertMode mode) 34d71ae5a4SJacob Faibussowitsch { 35c448b993SMatthew G. Knepley PetscReal lmin[] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL}; 36c448b993SMatthew G. Knepley PetscReal lmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 37c448b993SMatthew G. Knepley PetscInt i, j, k, bs, b, n_estimate, n_curr, n_new_est, p, n_found; 380e2ec84fSDave May const PetscScalar *_coor; 390e2ec84fSDave May DM celldm; 400e2ec84fSDave May PetscReal dx[3]; 412e3d3761SDave May PetscInt _npoints[] = {0, 0, 1}; 420e2ec84fSDave May Vec pos; 430e2ec84fSDave May PetscScalar *_pos; 440e2ec84fSDave May PetscReal *swarm_coor; 450e2ec84fSDave May PetscInt *swarm_cellid; 460e2ec84fSDave May PetscSF sfcell = NULL; 470e2ec84fSDave May const PetscSFNode *LA_sfcell; 48*d52c2f21SMatthew G. Knepley const char *coordname; 490e2ec84fSDave May 500e2ec84fSDave May PetscFunctionBegin; 510e2ec84fSDave May DMSWARMPICVALID(dm); 529566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 53c448b993SMatthew G. Knepley PetscCall(DMGetLocalBoundingBox(celldm, lmin, lmax)); 54c448b993SMatthew G. Knepley PetscCall(DMGetCoordinateDim(celldm, &bs)); 55*d52c2f21SMatthew G. Knepley PetscCall(DMSwarmGetCoordinateField(dm, &coordname)); 560e2ec84fSDave May 570e2ec84fSDave May for (b = 0; b < bs; b++) { 5871844bb8SDave May if (npoints[b] > 1) { 590e2ec84fSDave May dx[b] = (max[b] - min[b]) / ((PetscReal)(npoints[b] - 1)); 6071844bb8SDave May } else { 6171844bb8SDave May dx[b] = 0.0; 6271844bb8SDave May } 632e3d3761SDave May _npoints[b] = npoints[b]; 640e2ec84fSDave May } 650e2ec84fSDave May 660e2ec84fSDave May /* determine number of points living in the bounding box */ 670e2ec84fSDave May n_estimate = 0; 682e3d3761SDave May for (k = 0; k < _npoints[2]; k++) { 692e3d3761SDave May for (j = 0; j < _npoints[1]; j++) { 702e3d3761SDave May for (i = 0; i < _npoints[0]; i++) { 710e2ec84fSDave May PetscReal xp[] = {0.0, 0.0, 0.0}; 720e2ec84fSDave May PetscInt ijk[3]; 730e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 740e2ec84fSDave May 750e2ec84fSDave May ijk[0] = i; 760e2ec84fSDave May ijk[1] = j; 770e2ec84fSDave May ijk[2] = k; 78ad540459SPierre Jolivet for (b = 0; b < bs; b++) xp[b] = min[b] + ijk[b] * dx[b]; 790e2ec84fSDave May for (b = 0; b < bs; b++) { 80c448b993SMatthew G. Knepley if (xp[b] < lmin[b]) point_inside = PETSC_FALSE; 81c448b993SMatthew G. Knepley if (xp[b] > lmax[b]) point_inside = PETSC_FALSE; 820e2ec84fSDave May } 83ad540459SPierre Jolivet if (point_inside) n_estimate++; 840e2ec84fSDave May } 850e2ec84fSDave May } 860e2ec84fSDave May } 870e2ec84fSDave May 880e2ec84fSDave May /* create candidate list */ 899566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &pos)); 909566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pos, bs * n_estimate, PETSC_DECIDE)); 919566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(pos, bs)); 929566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pos)); 939566063dSJacob Faibussowitsch PetscCall(VecGetArray(pos, &_pos)); 940e2ec84fSDave May 950e2ec84fSDave May n_estimate = 0; 962e3d3761SDave May for (k = 0; k < _npoints[2]; k++) { 972e3d3761SDave May for (j = 0; j < _npoints[1]; j++) { 982e3d3761SDave May for (i = 0; i < _npoints[0]; i++) { 990e2ec84fSDave May PetscReal xp[] = {0.0, 0.0, 0.0}; 1000e2ec84fSDave May PetscInt ijk[3]; 1010e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 1020e2ec84fSDave May 1030e2ec84fSDave May ijk[0] = i; 1040e2ec84fSDave May ijk[1] = j; 1050e2ec84fSDave May ijk[2] = k; 106ad540459SPierre Jolivet for (b = 0; b < bs; b++) xp[b] = min[b] + ijk[b] * dx[b]; 1070e2ec84fSDave May for (b = 0; b < bs; b++) { 108c448b993SMatthew G. Knepley if (xp[b] < lmin[b]) point_inside = PETSC_FALSE; 109c448b993SMatthew G. Knepley if (xp[b] > lmax[b]) point_inside = PETSC_FALSE; 1100e2ec84fSDave May } 1110e2ec84fSDave May if (point_inside) { 112ad540459SPierre Jolivet for (b = 0; b < bs; b++) _pos[bs * n_estimate + b] = xp[b]; 1130e2ec84fSDave May n_estimate++; 1140e2ec84fSDave May } 1150e2ec84fSDave May } 1160e2ec84fSDave May } 1170e2ec84fSDave May } 1189566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pos, &_pos)); 1190e2ec84fSDave May 1200e2ec84fSDave May /* locate points */ 1219566063dSJacob Faibussowitsch PetscCall(DMLocatePoints(celldm, pos, DM_POINTLOCATION_NONE, &sfcell)); 1229566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 1230e2ec84fSDave May n_found = 0; 1240e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 125ad540459SPierre Jolivet if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) n_found++; 1260e2ec84fSDave May } 1270e2ec84fSDave May 1280e2ec84fSDave May /* adjust size */ 1290e2ec84fSDave May if (mode == ADD_VALUES) { 1309566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dm, &n_curr)); 1310e2ec84fSDave May n_new_est = n_curr + n_found; 1329566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 1330e2ec84fSDave May } 1340e2ec84fSDave May if (mode == INSERT_VALUES) { 1350e2ec84fSDave May n_curr = 0; 1360e2ec84fSDave May n_new_est = n_found; 1379566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 1380e2ec84fSDave May } 1390e2ec84fSDave May 1400e2ec84fSDave May /* initialize new coords, cell owners, pid */ 1419566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pos, &_coor)); 142*d52c2f21SMatthew G. Knepley PetscCall(DMSwarmGetField(dm, coordname, NULL, NULL, (void **)&swarm_coor)); 1439566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 1440e2ec84fSDave May n_found = 0; 1450e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 1460e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 147ad540459SPierre Jolivet for (b = 0; b < bs; b++) swarm_coor[bs * (n_curr + n_found) + b] = PetscRealPart(_coor[bs * p + b]); 1480e2ec84fSDave May swarm_cellid[n_curr + n_found] = LA_sfcell[p].index; 1490e2ec84fSDave May n_found++; 1500e2ec84fSDave May } 1510e2ec84fSDave May } 1529566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 153*d52c2f21SMatthew G. Knepley PetscCall(DMSwarmRestoreField(dm, coordname, NULL, NULL, (void **)&swarm_coor)); 1549566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pos, &_coor)); 1550e2ec84fSDave May 1569566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfcell)); 1579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pos)); 1583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1590e2ec84fSDave May } 1600e2ec84fSDave May 161cc4c1da9SBarry Smith /*@ 16220f4b53cSBarry Smith DMSwarmSetPointCoordinates - Set point coordinates in a `DMSWARM` from a user defined list 1630e2ec84fSDave May 16420f4b53cSBarry Smith Collective 1650e2ec84fSDave May 16660225df5SJacob Faibussowitsch Input Parameters: 16720f4b53cSBarry Smith + dm - the `DMSWARM` 1680e2ec84fSDave May . npoints - the number of points to insert 1690e2ec84fSDave May . coor - the coordinate values 17020f4b53cSBarry Smith . 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 17120f4b53cSBarry Smith - mode - indicates whether to append points to the swarm (`ADD_VALUES`), or over-ride existing points (`INSERT_VALUES`) 1720e2ec84fSDave May 1730e2ec84fSDave May Level: beginner 1740e2ec84fSDave May 1750e2ec84fSDave May Notes: 17620f4b53cSBarry Smith If the user has specified `redundant` as `PETSC_FALSE`, the cell `DM` will attempt to locate the coordinates provided by `coor` within 17720f4b53cSBarry Smith its sub-domain. If they any values within `coor` are not located in the sub-domain, they will be ignored and will not get 17820f4b53cSBarry Smith added to the `DMSWARM`. 1790e2ec84fSDave May 18020f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType`, `DMSwarmSetPointsUniformCoordinates()` 1810e2ec84fSDave May @*/ 182cc4c1da9SBarry Smith PetscErrorCode DMSwarmSetPointCoordinates(DM dm, PetscInt npoints, PetscReal coor[], PetscBool redundant, InsertMode mode) 183d71ae5a4SJacob Faibussowitsch { 1840e2ec84fSDave May PetscReal gmin[] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL}; 1850e2ec84fSDave May PetscReal gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 1860e2ec84fSDave May PetscInt i, N, bs, b, n_estimate, n_curr, n_new_est, p, n_found; 1870e2ec84fSDave May Vec coorlocal; 1880e2ec84fSDave May const PetscScalar *_coor; 1890e2ec84fSDave May DM celldm; 1900e2ec84fSDave May Vec pos; 1910e2ec84fSDave May PetscScalar *_pos; 1920e2ec84fSDave May PetscReal *swarm_coor; 1930e2ec84fSDave May PetscInt *swarm_cellid; 1940e2ec84fSDave May PetscSF sfcell = NULL; 1950e2ec84fSDave May const PetscSFNode *LA_sfcell; 1960e2ec84fSDave May PetscReal *my_coor; 1970e2ec84fSDave May PetscInt my_npoints; 1980e2ec84fSDave May PetscMPIInt rank; 1990e2ec84fSDave May MPI_Comm comm; 200*d52c2f21SMatthew G. Knepley const char *coordname; 2010e2ec84fSDave May 2020e2ec84fSDave May PetscFunctionBegin; 2030e2ec84fSDave May DMSWARMPICVALID(dm); 2049566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 2059566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2060e2ec84fSDave May 2079566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 208*d52c2f21SMatthew G. Knepley PetscCall(DMSwarmGetCoordinateField(dm, &coordname)); 2099566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(celldm, &coorlocal)); 2109566063dSJacob Faibussowitsch PetscCall(VecGetSize(coorlocal, &N)); 2119566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(coorlocal, &bs)); 2120e2ec84fSDave May N = N / bs; 2139566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coorlocal, &_coor)); 2140e2ec84fSDave May for (i = 0; i < N; i++) { 2150e2ec84fSDave May for (b = 0; b < bs; b++) { 216a5f152d1SDave May gmin[b] = PetscMin(gmin[b], PetscRealPart(_coor[bs * i + b])); 217a5f152d1SDave May gmax[b] = PetscMax(gmax[b], PetscRealPart(_coor[bs * i + b])); 2180e2ec84fSDave May } 2190e2ec84fSDave May } 2209566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coorlocal, &_coor)); 2210e2ec84fSDave May 2220e2ec84fSDave May /* broadcast points from rank 0 if requested */ 2230e2ec84fSDave May if (redundant) { 2246497c311SBarry Smith PetscMPIInt imy; 2256497c311SBarry Smith 2260e2ec84fSDave May my_npoints = npoints; 2279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&my_npoints, 1, MPIU_INT, 0, comm)); 2280e2ec84fSDave May 2290e2ec84fSDave May if (rank > 0) { /* allocate space */ 2309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * my_npoints, &my_coor)); 2310e2ec84fSDave May } else { 2320e2ec84fSDave May my_coor = coor; 2330e2ec84fSDave May } 2346497c311SBarry Smith PetscCall(PetscMPIIntCast(bs * my_npoints, &imy)); 2356497c311SBarry Smith PetscCallMPI(MPI_Bcast(my_coor, imy, MPIU_REAL, 0, comm)); 2360e2ec84fSDave May } else { 2370e2ec84fSDave May my_npoints = npoints; 2380e2ec84fSDave May my_coor = coor; 2390e2ec84fSDave May } 2400e2ec84fSDave May 2410e2ec84fSDave May /* determine the number of points living in the bounding box */ 2420e2ec84fSDave May n_estimate = 0; 2430e2ec84fSDave May for (i = 0; i < my_npoints; i++) { 2440e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 2450e2ec84fSDave May 2460e2ec84fSDave May for (b = 0; b < bs; b++) { 247ad540459SPierre Jolivet if (my_coor[bs * i + b] < gmin[b]) point_inside = PETSC_FALSE; 248ad540459SPierre Jolivet if (my_coor[bs * i + b] > gmax[b]) point_inside = PETSC_FALSE; 2490e2ec84fSDave May } 250ad540459SPierre Jolivet if (point_inside) n_estimate++; 2510e2ec84fSDave May } 2520e2ec84fSDave May 2530e2ec84fSDave May /* create candidate list */ 2549566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &pos)); 2559566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pos, bs * n_estimate, PETSC_DECIDE)); 2569566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(pos, bs)); 2579566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pos)); 2589566063dSJacob Faibussowitsch PetscCall(VecGetArray(pos, &_pos)); 2590e2ec84fSDave May 2600e2ec84fSDave May n_estimate = 0; 2610e2ec84fSDave May for (i = 0; i < my_npoints; i++) { 2620e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 2630e2ec84fSDave May 2640e2ec84fSDave May for (b = 0; b < bs; b++) { 265ad540459SPierre Jolivet if (my_coor[bs * i + b] < gmin[b]) point_inside = PETSC_FALSE; 266ad540459SPierre Jolivet if (my_coor[bs * i + b] > gmax[b]) point_inside = PETSC_FALSE; 2670e2ec84fSDave May } 2680e2ec84fSDave May if (point_inside) { 269ad540459SPierre Jolivet for (b = 0; b < bs; b++) _pos[bs * n_estimate + b] = my_coor[bs * i + b]; 2700e2ec84fSDave May n_estimate++; 2710e2ec84fSDave May } 2720e2ec84fSDave May } 2739566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pos, &_pos)); 2740e2ec84fSDave May 2750e2ec84fSDave May /* locate points */ 2769566063dSJacob Faibussowitsch PetscCall(DMLocatePoints(celldm, pos, DM_POINTLOCATION_NONE, &sfcell)); 2770e2ec84fSDave May 2789566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 2790e2ec84fSDave May n_found = 0; 2800e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 281ad540459SPierre Jolivet if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) n_found++; 2820e2ec84fSDave May } 2830e2ec84fSDave May 2840e2ec84fSDave May /* adjust size */ 2850e2ec84fSDave May if (mode == ADD_VALUES) { 2869566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dm, &n_curr)); 2870e2ec84fSDave May n_new_est = n_curr + n_found; 2889566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 2890e2ec84fSDave May } 2900e2ec84fSDave May if (mode == INSERT_VALUES) { 2910e2ec84fSDave May n_curr = 0; 2920e2ec84fSDave May n_new_est = n_found; 2939566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 2940e2ec84fSDave May } 2950e2ec84fSDave May 2960e2ec84fSDave May /* initialize new coords, cell owners, pid */ 2979566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pos, &_coor)); 298*d52c2f21SMatthew G. Knepley PetscCall(DMSwarmGetField(dm, coordname, NULL, NULL, (void **)&swarm_coor)); 2999566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 3000e2ec84fSDave May n_found = 0; 3010e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 3020e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 303ad540459SPierre Jolivet for (b = 0; b < bs; b++) swarm_coor[bs * (n_curr + n_found) + b] = PetscRealPart(_coor[bs * p + b]); 3040e2ec84fSDave May swarm_cellid[n_curr + n_found] = LA_sfcell[p].index; 3050e2ec84fSDave May n_found++; 3060e2ec84fSDave May } 3070e2ec84fSDave May } 3089566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 309*d52c2f21SMatthew G. Knepley PetscCall(DMSwarmRestoreField(dm, coordname, NULL, NULL, (void **)&swarm_coor)); 3109566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pos, &_coor)); 3110e2ec84fSDave May 3120e2ec84fSDave May if (redundant) { 31348a46eb9SPierre Jolivet if (rank > 0) PetscCall(PetscFree(my_coor)); 3140e2ec84fSDave May } 3159566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfcell)); 3169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pos)); 3173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3180e2ec84fSDave May } 3190e2ec84fSDave May 3200e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM, DM, DMSwarmPICLayoutType, PetscInt); 3210e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM, DM, DMSwarmPICLayoutType, PetscInt); 3220e2ec84fSDave May 323cc4c1da9SBarry Smith /*@ 3240e2ec84fSDave May DMSwarmInsertPointsUsingCellDM - Insert point coordinates within each cell 3250e2ec84fSDave May 32620f4b53cSBarry Smith Not Collective 3270e2ec84fSDave May 32860225df5SJacob Faibussowitsch Input Parameters: 32920f4b53cSBarry Smith + dm - the `DMSWARM` 33020f4b53cSBarry Smith . layout_type - method used to fill each cell with the cell `DM` 3310e2ec84fSDave May - fill_param - parameter controlling how many points per cell are added (the meaning of this parameter is dependent on the layout type) 3320e2ec84fSDave May 3330e2ec84fSDave May Level: beginner 3340e2ec84fSDave May 3350e2ec84fSDave May Notes: 33620f4b53cSBarry Smith The insert method will reset any previous defined points within the `DMSWARM`. 337e7af74c8SDave May 33820f4b53cSBarry Smith When using a `DMDA` both 2D and 3D are supported for all layout types provided you are using `DMDA_ELEMENT_Q1`. 339e7af74c8SDave May 340a4e35b19SJacob Faibussowitsch When using a `DMPLEX` the following case are supported\: 34120f4b53cSBarry Smith .vb 342ea3d7275SDave May (i) DMSWARMPIC_LAYOUT_REGULAR: 2D (triangle), 343ea3d7275SDave May (ii) DMSWARMPIC_LAYOUT_GAUSS: 2D and 3D provided the cell is a tri/tet or a quad/hex, 344ea3d7275SDave May (iii) DMSWARMPIC_LAYOUT_SUBDIVISION: 2D and 3D for quad/hex and 2D tri. 34520f4b53cSBarry Smith .ve 3460e2ec84fSDave May 34720f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmPICLayoutType`, `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType` 3480e2ec84fSDave May @*/ 349cc4c1da9SBarry Smith PetscErrorCode DMSwarmInsertPointsUsingCellDM(DM dm, DMSwarmPICLayoutType layout_type, PetscInt fill_param) 350d71ae5a4SJacob Faibussowitsch { 3510e2ec84fSDave May DM celldm; 3520e2ec84fSDave May PetscBool isDA, isPLEX; 3530e2ec84fSDave May 3540e2ec84fSDave May PetscFunctionBegin; 3550e2ec84fSDave May DMSWARMPICVALID(dm); 3569566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 3579566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA)); 3589566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX)); 3590e2ec84fSDave May if (isDA) { 3609566063dSJacob Faibussowitsch PetscCall(private_DMSwarmInsertPointsUsingCellDM_DA(dm, celldm, layout_type, fill_param)); 3610e2ec84fSDave May } else if (isPLEX) { 3629566063dSJacob Faibussowitsch PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX(dm, celldm, layout_type, fill_param)); 3630e2ec84fSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX"); 3643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3650e2ec84fSDave May } 3660e2ec84fSDave May 367431382f2SDave May extern PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM, DM, PetscInt, PetscReal *); 368431382f2SDave May 369431382f2SDave May /*@C 370431382f2SDave May DMSwarmSetPointCoordinatesCellwise - Insert point coordinates (defined over the reference cell) within each cell 371431382f2SDave May 37220f4b53cSBarry Smith Not Collective 373431382f2SDave May 37460225df5SJacob Faibussowitsch Input Parameters: 37520f4b53cSBarry Smith + dm - the `DMSWARM` 376431382f2SDave May . npoints - the number of points to insert in each cell 377431382f2SDave May - xi - the coordinates (defined in the local coordinate system for each cell) to insert 378431382f2SDave May 379431382f2SDave May Level: beginner 380431382f2SDave May 381431382f2SDave May Notes: 38220f4b53cSBarry Smith The method will reset any previous defined points within the `DMSWARM`. 38320f4b53cSBarry Smith Only supported for `DMPLEX`. If you are using a `DMDA` it is recommended to either use 38420f4b53cSBarry Smith `DMSwarmInsertPointsUsingCellDM()`, or extract and set the coordinates yourself the following code 38520f4b53cSBarry Smith .vb 38620f4b53cSBarry Smith PetscReal *coor; 387*d52c2f21SMatthew G. Knepley const char *coordname; 388*d52c2f21SMatthew G. Knepley DMSwarmGetCoordinateField(dm, &coordname); 389*d52c2f21SMatthew G. Knepley DMSwarmGetField(dm,coordname,NULL,NULL,(void**)&coor); 39020f4b53cSBarry Smith // user code to define the coordinates here 391*d52c2f21SMatthew G. Knepley DMSwarmRestoreField(dm,coordname,NULL,NULL,(void**)&coor); 39220f4b53cSBarry Smith .ve 393e7af74c8SDave May 39420f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetCellDM()`, `DMSwarmInsertPointsUsingCellDM()` 395431382f2SDave May @*/ 396cc4c1da9SBarry Smith PetscErrorCode DMSwarmSetPointCoordinatesCellwise(DM dm, PetscInt npoints, PetscReal xi[]) 397d71ae5a4SJacob Faibussowitsch { 398431382f2SDave May DM celldm; 399431382f2SDave May PetscBool isDA, isPLEX; 400431382f2SDave May 4010e2ec84fSDave May PetscFunctionBegin; 402431382f2SDave May DMSWARMPICVALID(dm); 4039566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 4049566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA)); 4059566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX)); 40628b400f6SJacob Faibussowitsch PetscCheck(!isDA, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMPLEX. Recommended you use DMSwarmInsertPointsUsingCellDM()"); 407f7d195e4SLawrence Mitchell if (isPLEX) { 4089566063dSJacob Faibussowitsch PetscCall(private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm, celldm, npoints, xi)); 409431382f2SDave May } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX"); 4103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4110e2ec84fSDave May } 412431382f2SDave May 413cc4c1da9SBarry Smith /*@ 414b799feefSDave May DMSwarmCreatePointPerCellCount - Count the number of points within all cells in the cell DM 4150e2ec84fSDave May 41620f4b53cSBarry Smith Not Collective 4170e2ec84fSDave May 41860225df5SJacob Faibussowitsch Input Parameter: 41920f4b53cSBarry Smith . dm - the `DMSWARM` 4200e2ec84fSDave May 42160225df5SJacob Faibussowitsch Output Parameters: 42220f4b53cSBarry Smith + ncells - the number of cells in the cell `DM` (optional argument, pass `NULL` to ignore) 423b799feefSDave May - count - array of length ncells containing the number of points per cell 4240e2ec84fSDave May 4250e2ec84fSDave May Level: beginner 4260e2ec84fSDave May 4270e2ec84fSDave May Notes: 4280e2ec84fSDave May The array count is allocated internally and must be free'd by the user. 4290e2ec84fSDave May 43020f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType` 4310e2ec84fSDave May @*/ 432cc4c1da9SBarry Smith PetscErrorCode DMSwarmCreatePointPerCellCount(DM dm, PetscInt *ncells, PetscInt **count) 433d71ae5a4SJacob Faibussowitsch { 434b799feefSDave May PetscBool isvalid; 435b799feefSDave May PetscInt nel; 436b799feefSDave May PetscInt *sum; 437b799feefSDave May 4380e2ec84fSDave May PetscFunctionBegin; 4399566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetIsValid(dm, &isvalid)); 440b799feefSDave May nel = 0; 441b799feefSDave May if (isvalid) { 442b799feefSDave May PetscInt e; 443b799feefSDave May 4449566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetSizes(dm, &nel, NULL)); 445b799feefSDave May 4469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nel, &sum)); 44748a46eb9SPierre Jolivet for (e = 0; e < nel; e++) PetscCall(DMSwarmSortGetNumberOfPointsPerCell(dm, e, &sum[e])); 448b799feefSDave May } else { 449b799feefSDave May DM celldm; 450b799feefSDave May PetscBool isda, isplex, isshell; 451b799feefSDave May PetscInt p, npoints; 452b799feefSDave May PetscInt *swarm_cellid; 453b799feefSDave May 454b799feefSDave May /* get the number of cells */ 4559566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 4569566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isda)); 4579566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isplex)); 4589566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMSHELL, &isshell)); 459b799feefSDave May if (isda) { 460b799feefSDave May PetscInt _nel, _npe; 461b799feefSDave May const PetscInt *_element; 462b799feefSDave May 4639566063dSJacob Faibussowitsch PetscCall(DMDAGetElements(celldm, &_nel, &_npe, &_element)); 464b799feefSDave May nel = _nel; 4659566063dSJacob Faibussowitsch PetscCall(DMDARestoreElements(celldm, &_nel, &_npe, &_element)); 466b799feefSDave May } else if (isplex) { 467b799feefSDave May PetscInt ps, pe; 468b799feefSDave May 4699566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(celldm, 0, &ps, &pe)); 470b799feefSDave May nel = pe - ps; 471b799feefSDave May } else if (isshell) { 472b799feefSDave May PetscErrorCode (*method_DMShellGetNumberOfCells)(DM, PetscInt *); 473b799feefSDave May 4749566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)celldm, "DMGetNumberOfCells_C", &method_DMShellGetNumberOfCells)); 475b799feefSDave May if (method_DMShellGetNumberOfCells) { 4769566063dSJacob Faibussowitsch PetscCall(method_DMShellGetNumberOfCells(celldm, &nel)); 4779371c9d4SSatish Balay } else 4789371c9d4SSatish Balay 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 4819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nel, &sum)); 4829566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(sum, nel)); 4839566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dm, &npoints)); 4849566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 485b799feefSDave May for (p = 0; p < npoints; p++) { 486ad540459SPierre Jolivet if (swarm_cellid[p] != DMLOCATEPOINT_POINT_NOT_FOUND) sum[swarm_cellid[p]]++; 487b799feefSDave May } 4889566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 489b799feefSDave May } 490ad540459SPierre Jolivet if (ncells) *ncells = nel; 491b799feefSDave May *count = sum; 4923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4930e2ec84fSDave May } 49435b38c60SMatthew G. Knepley 49535b38c60SMatthew G. Knepley /*@ 49635b38c60SMatthew G. Knepley DMSwarmGetNumSpecies - Get the number of particle species 49735b38c60SMatthew G. Knepley 49820f4b53cSBarry Smith Not Collective 49935b38c60SMatthew G. Knepley 50060225df5SJacob Faibussowitsch Input Parameter: 50160225df5SJacob Faibussowitsch . sw - the `DMSWARM` 50235b38c60SMatthew G. Knepley 50360225df5SJacob Faibussowitsch Output Parameters: 50435b38c60SMatthew G. Knepley . Ns - the number of species 50535b38c60SMatthew G. Knepley 50635b38c60SMatthew G. Knepley Level: intermediate 50735b38c60SMatthew G. Knepley 50820f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetNumSpecies()`, `DMSwarmSetType()`, `DMSwarmType` 50935b38c60SMatthew G. Knepley @*/ 510d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetNumSpecies(DM sw, PetscInt *Ns) 511d71ae5a4SJacob Faibussowitsch { 51235b38c60SMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 51335b38c60SMatthew G. Knepley 51435b38c60SMatthew G. Knepley PetscFunctionBegin; 51535b38c60SMatthew G. Knepley *Ns = swarm->Ns; 5163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51735b38c60SMatthew G. Knepley } 51835b38c60SMatthew G. Knepley 51935b38c60SMatthew G. Knepley /*@ 52035b38c60SMatthew G. Knepley DMSwarmSetNumSpecies - Set the number of particle species 52135b38c60SMatthew G. Knepley 52220f4b53cSBarry Smith Not Collective 52335b38c60SMatthew G. Knepley 52460225df5SJacob Faibussowitsch Input Parameters: 52560225df5SJacob Faibussowitsch + sw - the `DMSWARM` 52635b38c60SMatthew G. Knepley - Ns - the number of species 52735b38c60SMatthew G. Knepley 52835b38c60SMatthew G. Knepley Level: intermediate 52935b38c60SMatthew G. Knepley 53020f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmGetNumSpecies()`, `DMSwarmSetType()`, `DMSwarmType` 53135b38c60SMatthew G. Knepley @*/ 532d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetNumSpecies(DM sw, PetscInt Ns) 533d71ae5a4SJacob Faibussowitsch { 53435b38c60SMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 53535b38c60SMatthew G. Knepley 53635b38c60SMatthew G. Knepley PetscFunctionBegin; 53735b38c60SMatthew G. Knepley swarm->Ns = Ns; 5383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 53935b38c60SMatthew G. Knepley } 54035b38c60SMatthew G. Knepley 54135b38c60SMatthew G. Knepley /*@C 5426c5a79ebSMatthew G. Knepley DMSwarmGetCoordinateFunction - Get the function setting initial particle positions, if it exists 5436c5a79ebSMatthew G. Knepley 54420f4b53cSBarry Smith Not Collective 5456c5a79ebSMatthew G. Knepley 54660225df5SJacob Faibussowitsch Input Parameter: 54760225df5SJacob Faibussowitsch . sw - the `DMSWARM` 5486c5a79ebSMatthew G. Knepley 5496c5a79ebSMatthew G. Knepley Output Parameter: 5508434afd1SBarry Smith . coordFunc - the function setting initial particle positions, or `NULL`, see `PetscSimplePointFn` for the calling sequence 5516c5a79ebSMatthew G. Knepley 5526c5a79ebSMatthew G. Knepley Level: intermediate 5536c5a79ebSMatthew G. Knepley 5548434afd1SBarry Smith .seealso: `DMSWARM`, `DMSwarmSetCoordinateFunction()`, `DMSwarmGetVelocityFunction()`, `DMSwarmInitializeCoordinates()`, `PetscSimplePointFn` 5556c5a79ebSMatthew G. Knepley @*/ 5568434afd1SBarry Smith PetscErrorCode DMSwarmGetCoordinateFunction(DM sw, PetscSimplePointFn **coordFunc) 557d71ae5a4SJacob Faibussowitsch { 5586c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 5596c5a79ebSMatthew G. Knepley 5606c5a79ebSMatthew G. Knepley PetscFunctionBegin; 5616c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 5624f572ea9SToby Isaac PetscAssertPointer(coordFunc, 2); 5636c5a79ebSMatthew G. Knepley *coordFunc = swarm->coordFunc; 5643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5656c5a79ebSMatthew G. Knepley } 5666c5a79ebSMatthew G. Knepley 5676c5a79ebSMatthew G. Knepley /*@C 5686c5a79ebSMatthew G. Knepley DMSwarmSetCoordinateFunction - Set the function setting initial particle positions 5696c5a79ebSMatthew G. Knepley 57020f4b53cSBarry Smith Not Collective 5716c5a79ebSMatthew G. Knepley 57260225df5SJacob Faibussowitsch Input Parameters: 57360225df5SJacob Faibussowitsch + sw - the `DMSWARM` 5748434afd1SBarry Smith - coordFunc - the function setting initial particle positions, see `PetscSimplePointFn` for the calling sequence 5756c5a79ebSMatthew G. Knepley 5766c5a79ebSMatthew G. Knepley Level: intermediate 5776c5a79ebSMatthew G. Knepley 5788434afd1SBarry Smith .seealso: `DMSWARM`, `DMSwarmGetCoordinateFunction()`, `DMSwarmSetVelocityFunction()`, `DMSwarmInitializeCoordinates()`, `PetscSimplePointFn` 5796c5a79ebSMatthew G. Knepley @*/ 5808434afd1SBarry Smith PetscErrorCode DMSwarmSetCoordinateFunction(DM sw, PetscSimplePointFn *coordFunc) 581d71ae5a4SJacob Faibussowitsch { 5826c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 5836c5a79ebSMatthew G. Knepley 5846c5a79ebSMatthew G. Knepley PetscFunctionBegin; 5856c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 5866c5a79ebSMatthew G. Knepley PetscValidFunction(coordFunc, 2); 5876c5a79ebSMatthew G. Knepley swarm->coordFunc = coordFunc; 5883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5896c5a79ebSMatthew G. Knepley } 5906c5a79ebSMatthew G. Knepley 5916c5a79ebSMatthew G. Knepley /*@C 59260225df5SJacob Faibussowitsch DMSwarmGetVelocityFunction - Get the function setting initial particle velocities, if it exists 5936c5a79ebSMatthew G. Knepley 59420f4b53cSBarry Smith Not Collective 5956c5a79ebSMatthew G. Knepley 59660225df5SJacob Faibussowitsch Input Parameter: 59760225df5SJacob Faibussowitsch . sw - the `DMSWARM` 5986c5a79ebSMatthew G. Knepley 5996c5a79ebSMatthew G. Knepley Output Parameter: 6008434afd1SBarry Smith . velFunc - the function setting initial particle velocities, or `NULL`, see `PetscSimplePointFn` for the calling sequence 6016c5a79ebSMatthew G. Knepley 6026c5a79ebSMatthew G. Knepley Level: intermediate 6036c5a79ebSMatthew G. Knepley 6048434afd1SBarry Smith .seealso: `DMSWARM`, `DMSwarmSetVelocityFunction()`, `DMSwarmGetCoordinateFunction()`, `DMSwarmInitializeVelocities()`, `PetscSimplePointFn` 6056c5a79ebSMatthew G. Knepley @*/ 6068434afd1SBarry Smith PetscErrorCode DMSwarmGetVelocityFunction(DM sw, PetscSimplePointFn **velFunc) 607d71ae5a4SJacob Faibussowitsch { 6086c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 6096c5a79ebSMatthew G. Knepley 6106c5a79ebSMatthew G. Knepley PetscFunctionBegin; 6116c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 6124f572ea9SToby Isaac PetscAssertPointer(velFunc, 2); 6136c5a79ebSMatthew G. Knepley *velFunc = swarm->velFunc; 6143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6156c5a79ebSMatthew G. Knepley } 6166c5a79ebSMatthew G. Knepley 6176c5a79ebSMatthew G. Knepley /*@C 6186c5a79ebSMatthew G. Knepley DMSwarmSetVelocityFunction - Set the function setting initial particle velocities 6196c5a79ebSMatthew G. Knepley 62020f4b53cSBarry Smith Not Collective 6216c5a79ebSMatthew G. Knepley 62260225df5SJacob Faibussowitsch Input Parameters: 623a4e35b19SJacob Faibussowitsch + sw - the `DMSWARM` 6248434afd1SBarry Smith - velFunc - the function setting initial particle velocities, see `PetscSimplePointFn` for the calling sequence 6256c5a79ebSMatthew G. Knepley 6266c5a79ebSMatthew G. Knepley Level: intermediate 6276c5a79ebSMatthew G. Knepley 6288434afd1SBarry Smith .seealso: `DMSWARM`, `DMSwarmGetVelocityFunction()`, `DMSwarmSetCoordinateFunction()`, `DMSwarmInitializeVelocities()`, `PetscSimplePointFn` 6296c5a79ebSMatthew G. Knepley @*/ 6308434afd1SBarry Smith PetscErrorCode DMSwarmSetVelocityFunction(DM sw, PetscSimplePointFn *velFunc) 631d71ae5a4SJacob Faibussowitsch { 6326c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 6336c5a79ebSMatthew G. Knepley 6346c5a79ebSMatthew G. Knepley PetscFunctionBegin; 6356c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 6366c5a79ebSMatthew G. Knepley PetscValidFunction(velFunc, 2); 6376c5a79ebSMatthew G. Knepley swarm->velFunc = velFunc; 6383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6396c5a79ebSMatthew G. Knepley } 6406c5a79ebSMatthew G. Knepley 6416c5a79ebSMatthew G. Knepley /*@C 64235b38c60SMatthew G. Knepley DMSwarmComputeLocalSize - Compute the local number and distribution of particles based upon a density function 64335b38c60SMatthew G. Knepley 64420f4b53cSBarry Smith Not Collective 64535b38c60SMatthew G. Knepley 64635b38c60SMatthew G. Knepley Input Parameters: 64720f4b53cSBarry Smith + sw - The `DMSWARM` 64835b38c60SMatthew G. Knepley . N - The target number of particles 64935b38c60SMatthew G. Knepley - density - The density field for the particle layout, normalized to unity 65035b38c60SMatthew G. Knepley 65135b38c60SMatthew G. Knepley Level: advanced 65235b38c60SMatthew G. Knepley 65320f4b53cSBarry Smith Note: 65420f4b53cSBarry Smith One particle will be created for each species. 65520f4b53cSBarry Smith 65620f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSizeFromOptions()` 65735b38c60SMatthew G. Knepley @*/ 658d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmComputeLocalSize(DM sw, PetscInt N, PetscProbFunc density) 659d71ae5a4SJacob Faibussowitsch { 66035b38c60SMatthew G. Knepley DM dm; 66135b38c60SMatthew G. Knepley PetscQuadrature quad; 66235b38c60SMatthew G. Knepley const PetscReal *xq, *wq; 663ea7032a0SMatthew G. Knepley PetscReal *n_int; 664ea7032a0SMatthew G. Knepley PetscInt *npc_s, *cellid, Ni; 665ea7032a0SMatthew G. Knepley PetscReal gmin[3], gmax[3], xi0[3]; 666ea7032a0SMatthew G. Knepley PetscInt Ns, cStart, cEnd, c, dim, d, Nq, q, Np = 0, p, s; 66735b38c60SMatthew G. Knepley PetscBool simplex; 66835b38c60SMatthew G. Knepley 66935b38c60SMatthew G. Knepley PetscFunctionBegin; 6709566063dSJacob Faibussowitsch PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 6719566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sw, &dm)); 6729566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 673ea7032a0SMatthew G. Knepley PetscCall(DMGetBoundingBox(dm, gmin, gmax)); 6749566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 6759566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 6766858538eSMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(dm)); 6779566063dSJacob Faibussowitsch if (simplex) PetscCall(PetscDTStroudConicalQuadrature(dim, 1, 5, -1.0, 1.0, &quad)); 6789566063dSJacob Faibussowitsch else PetscCall(PetscDTGaussTensorQuadrature(dim, 1, 5, -1.0, 1.0, &quad)); 6799566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &xq, &wq)); 680ea7032a0SMatthew G. Knepley PetscCall(PetscCalloc2(Ns, &n_int, (cEnd - cStart) * Ns, &npc_s)); 68135b38c60SMatthew G. Knepley /* Integrate the density function to get the number of particles in each cell */ 68235b38c60SMatthew G. Knepley for (d = 0; d < dim; ++d) xi0[d] = -1.0; 68335b38c60SMatthew G. Knepley for (c = 0; c < cEnd - cStart; ++c) { 68435b38c60SMatthew G. Knepley const PetscInt cell = c + cStart; 685ea7032a0SMatthew G. Knepley PetscReal v0[3], J[9], invJ[9], detJ, detJp = 2. / (gmax[0] - gmin[0]), xr[3], den; 68635b38c60SMatthew G. Knepley 687ea7032a0SMatthew G. Knepley /* Have to transform quadrature points/weights to cell domain */ 6889566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, NULL, v0, J, invJ, &detJ)); 689ea7032a0SMatthew G. Knepley PetscCall(PetscArrayzero(n_int, Ns)); 69035b38c60SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 69135b38c60SMatthew G. Knepley CoordinatesRefToReal(dim, dim, xi0, v0, J, &xq[q * dim], xr); 692ea7032a0SMatthew G. Knepley /* Have to transform mesh to domain of definition of PDF, [-1, 1], and weight PDF by |J|/2 */ 693ea7032a0SMatthew G. Knepley xr[0] = detJp * (xr[0] - gmin[0]) - 1.; 694ea7032a0SMatthew G. Knepley 695ea7032a0SMatthew G. Knepley for (s = 0; s < Ns; ++s) { 696ea7032a0SMatthew G. Knepley PetscCall(density(xr, NULL, &den)); 697ea7032a0SMatthew G. Knepley n_int[s] += (detJp * den) * (detJ * wq[q]) / (PetscReal)Ns; 69835b38c60SMatthew G. Knepley } 699ea7032a0SMatthew G. Knepley } 700ea7032a0SMatthew G. Knepley for (s = 0; s < Ns; ++s) { 701ea7032a0SMatthew G. Knepley Ni = N; 702ea7032a0SMatthew G. Knepley npc_s[c * Ns + s] += (PetscInt)(Ni * n_int[s]); 703ea7032a0SMatthew G. Knepley Np += npc_s[c * Ns + s]; 704ea7032a0SMatthew G. Knepley } 70535b38c60SMatthew G. Knepley } 7069566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quad)); 7079566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(sw, Np, 0)); 7089566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 70935b38c60SMatthew G. Knepley for (c = 0, p = 0; c < cEnd - cStart; ++c) { 710ea7032a0SMatthew G. Knepley for (s = 0; s < Ns; ++s) { 711ea7032a0SMatthew G. Knepley for (q = 0; q < npc_s[c * Ns + s]; ++q, ++p) cellid[p] = c; 712ea7032a0SMatthew G. Knepley } 71335b38c60SMatthew G. Knepley } 7149566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 715ea7032a0SMatthew G. Knepley PetscCall(PetscFree2(n_int, npc_s)); 7163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 71735b38c60SMatthew G. Knepley } 71835b38c60SMatthew G. Knepley 71935b38c60SMatthew G. Knepley /*@ 72035b38c60SMatthew G. Knepley DMSwarmComputeLocalSizeFromOptions - Compute the local number and distribution of particles based upon a density function determined by options 72135b38c60SMatthew G. Knepley 72220f4b53cSBarry Smith Not Collective 72335b38c60SMatthew G. Knepley 7242fe279fdSBarry Smith Input Parameter: 7252fe279fdSBarry Smith . sw - The `DMSWARM` 72635b38c60SMatthew G. Knepley 72735b38c60SMatthew G. Knepley Level: advanced 72835b38c60SMatthew G. Knepley 72920f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSize()` 73035b38c60SMatthew G. Knepley @*/ 731d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmComputeLocalSizeFromOptions(DM sw) 732d71ae5a4SJacob Faibussowitsch { 73335b38c60SMatthew G. Knepley PetscProbFunc pdf; 73435b38c60SMatthew G. Knepley const char *prefix; 7356c5a79ebSMatthew G. Knepley char funcname[PETSC_MAX_PATH_LEN]; 7366c5a79ebSMatthew G. Knepley PetscInt *N, Ns, dim, n; 7376c5a79ebSMatthew G. Knepley PetscBool flg; 7386c5a79ebSMatthew G. Knepley PetscMPIInt size, rank; 73935b38c60SMatthew G. Knepley 74035b38c60SMatthew G. Knepley PetscFunctionBegin; 7416c5a79ebSMatthew G. Knepley PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sw), &size)); 7426c5a79ebSMatthew G. Knepley PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sw), &rank)); 7436c5a79ebSMatthew G. Knepley PetscCall(PetscCalloc1(size, &N)); 744d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)sw), "", "DMSwarm Options", "DMSWARM"); 7456c5a79ebSMatthew G. Knepley n = size; 7466c5a79ebSMatthew G. Knepley PetscCall(PetscOptionsIntArray("-dm_swarm_num_particles", "The target number of particles", "", N, &n, NULL)); 7476c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 7489566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-dm_swarm_num_species", "The number of species", "DMSwarmSetNumSpecies", Ns, &Ns, &flg)); 7499566063dSJacob Faibussowitsch if (flg) PetscCall(DMSwarmSetNumSpecies(sw, Ns)); 7506c5a79ebSMatthew G. Knepley PetscCall(PetscOptionsString("-dm_swarm_coordinate_function", "Function to determine particle coordinates", "DMSwarmSetCoordinateFunction", funcname, funcname, sizeof(funcname), &flg)); 751d0609cedSBarry Smith PetscOptionsEnd(); 7526c5a79ebSMatthew G. Knepley if (flg) { 7538434afd1SBarry Smith PetscSimplePointFn *coordFunc; 75435b38c60SMatthew G. Knepley 7556c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 7566c5a79ebSMatthew G. Knepley PetscCall(PetscDLSym(NULL, funcname, (void **)&coordFunc)); 7576c5a79ebSMatthew G. Knepley PetscCheck(coordFunc, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname); 7586c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 7596c5a79ebSMatthew G. Knepley PetscCall(DMSwarmSetLocalSizes(sw, N[rank] * Ns, 0)); 7606c5a79ebSMatthew G. Knepley PetscCall(DMSwarmSetCoordinateFunction(sw, coordFunc)); 7616c5a79ebSMatthew G. Knepley } else { 7629566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 7639566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sw, &prefix)); 7649566063dSJacob Faibussowitsch PetscCall(PetscProbCreateFromOptions(dim, prefix, "-dm_swarm_coordinate_density", &pdf, NULL, NULL)); 7656c5a79ebSMatthew G. Knepley PetscCall(DMSwarmComputeLocalSize(sw, N[rank], pdf)); 7666c5a79ebSMatthew G. Knepley } 7676c5a79ebSMatthew G. Knepley PetscCall(PetscFree(N)); 7683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 76935b38c60SMatthew G. Knepley } 77035b38c60SMatthew G. Knepley 77135b38c60SMatthew G. Knepley /*@ 77235b38c60SMatthew G. Knepley DMSwarmInitializeCoordinates - Determine the initial coordinates of particles for a PIC method 77335b38c60SMatthew G. Knepley 77420f4b53cSBarry Smith Not Collective 77535b38c60SMatthew G. Knepley 77620f4b53cSBarry Smith Input Parameter: 77720f4b53cSBarry Smith . sw - The `DMSWARM` 77835b38c60SMatthew G. Knepley 77935b38c60SMatthew G. Knepley Level: advanced 78035b38c60SMatthew G. Knepley 78120f4b53cSBarry Smith Note: 78220f4b53cSBarry Smith Currently, we randomly place particles in their assigned cell 78320f4b53cSBarry Smith 78420f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSize()`, `DMSwarmInitializeVelocities()` 78535b38c60SMatthew G. Knepley @*/ 786d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmInitializeCoordinates(DM sw) 787d71ae5a4SJacob Faibussowitsch { 7888434afd1SBarry Smith PetscSimplePointFn *coordFunc; 78935b38c60SMatthew G. Knepley PetscScalar *weight; 7906c5a79ebSMatthew G. Knepley PetscReal *x; 79135b38c60SMatthew G. Knepley PetscInt *species; 7926c5a79ebSMatthew G. Knepley void *ctx; 79335b38c60SMatthew G. Knepley PetscBool removePoints = PETSC_TRUE; 79435b38c60SMatthew G. Knepley PetscDataType dtype; 7956c5a79ebSMatthew G. Knepley PetscInt Np, p, Ns, dim, d, bs; 796*d52c2f21SMatthew G. Knepley const char *coordname; 79735b38c60SMatthew G. Knepley 79835b38c60SMatthew G. Knepley PetscFunctionBeginUser; 7996c5a79ebSMatthew G. Knepley PetscCall(DMGetDimension(sw, &dim)); 8006c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetLocalSize(sw, &Np)); 8019566063dSJacob Faibussowitsch PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 802*d52c2f21SMatthew G. Knepley PetscCall(DMSwarmGetCoordinateField(sw, &coordname)); 8036c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetCoordinateFunction(sw, &coordFunc)); 8046c5a79ebSMatthew G. Knepley 805*d52c2f21SMatthew G. Knepley PetscCall(DMSwarmGetField(sw, coordname, &bs, &dtype, (void **)&x)); 8066c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&weight)); 8076c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 8086c5a79ebSMatthew G. Knepley if (coordFunc) { 8096c5a79ebSMatthew G. Knepley PetscCall(DMGetApplicationContext(sw, &ctx)); 8106c5a79ebSMatthew G. Knepley for (p = 0; p < Np; ++p) { 8116c5a79ebSMatthew G. Knepley PetscScalar X[3]; 8126c5a79ebSMatthew G. Knepley 8136c5a79ebSMatthew G. Knepley PetscCall((*coordFunc)(dim, 0., NULL, p, X, ctx)); 8146c5a79ebSMatthew G. Knepley for (d = 0; d < dim; ++d) x[p * dim + d] = PetscRealPart(X[d]); 8156c5a79ebSMatthew G. Knepley weight[p] = 1.0; 8166c5a79ebSMatthew G. Knepley species[p] = p % Ns; 8176c5a79ebSMatthew G. Knepley } 8186c5a79ebSMatthew G. Knepley } else { 8196c5a79ebSMatthew G. Knepley DM dm; 8206c5a79ebSMatthew G. Knepley PetscRandom rnd; 8216c5a79ebSMatthew G. Knepley PetscReal xi0[3]; 8226c5a79ebSMatthew G. Knepley PetscInt cStart, cEnd, c; 8236c5a79ebSMatthew G. Knepley 8249566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sw, &dm)); 8259566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 826ea7032a0SMatthew G. Knepley PetscCall(DMGetApplicationContext(sw, &ctx)); 82735b38c60SMatthew G. Knepley 82835b38c60SMatthew G. Knepley /* Set particle position randomly in cell, set weights to 1 */ 8299566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd)); 8309566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0)); 8319566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rnd)); 8329566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetAccess(sw)); 83335b38c60SMatthew G. Knepley for (d = 0; d < dim; ++d) xi0[d] = -1.0; 83435b38c60SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 83535b38c60SMatthew G. Knepley PetscReal v0[3], J[9], invJ[9], detJ; 83635b38c60SMatthew G. Knepley PetscInt *pidx, Npc, q; 83735b38c60SMatthew G. Knepley 8389566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx)); 8399566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); 84035b38c60SMatthew G. Knepley for (q = 0; q < Npc; ++q) { 84135b38c60SMatthew G. Knepley const PetscInt p = pidx[q]; 84235b38c60SMatthew G. Knepley PetscReal xref[3]; 84335b38c60SMatthew G. Knepley 8449566063dSJacob Faibussowitsch for (d = 0; d < dim; ++d) PetscCall(PetscRandomGetValueReal(rnd, &xref[d])); 84535b38c60SMatthew G. Knepley CoordinatesRefToReal(dim, dim, xi0, v0, J, xref, &x[p * dim]); 84635b38c60SMatthew G. Knepley 847ea7032a0SMatthew G. Knepley weight[p] = 1.0 / Np; 8486c5a79ebSMatthew G. Knepley species[p] = p % Ns; 84935b38c60SMatthew G. Knepley } 850fe11354eSMatthew G. Knepley PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Npc, &pidx)); 85135b38c60SMatthew G. Knepley } 8529566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rnd)); 8539566063dSJacob Faibussowitsch PetscCall(DMSwarmSortRestoreAccess(sw)); 8546c5a79ebSMatthew G. Knepley } 855*d52c2f21SMatthew G. Knepley PetscCall(DMSwarmRestoreField(sw, coordname, NULL, NULL, (void **)&x)); 8566c5a79ebSMatthew G. Knepley PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight)); 8579566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 8586c5a79ebSMatthew G. Knepley 8599566063dSJacob Faibussowitsch PetscCall(DMSwarmMigrate(sw, removePoints)); 8609566063dSJacob Faibussowitsch PetscCall(DMLocalizeCoordinates(sw)); 8613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 86235b38c60SMatthew G. Knepley } 86335b38c60SMatthew G. Knepley 86435b38c60SMatthew G. Knepley /*@C 86535b38c60SMatthew G. Knepley DMSwarmInitializeVelocities - Set the initial velocities of particles using a distribution. 86635b38c60SMatthew G. Knepley 86720f4b53cSBarry Smith Collective 86835b38c60SMatthew G. Knepley 86935b38c60SMatthew G. Knepley Input Parameters: 87020f4b53cSBarry Smith + sw - The `DMSWARM` object 87135b38c60SMatthew G. Knepley . sampler - A function which uniformly samples the velocity PDF 87235b38c60SMatthew G. Knepley - v0 - The velocity scale for nondimensionalization for each species 87335b38c60SMatthew G. Knepley 87435b38c60SMatthew G. Knepley Level: advanced 87535b38c60SMatthew G. Knepley 87620f4b53cSBarry Smith Note: 87720f4b53cSBarry Smith If `v0` is zero for the first species, all velocities are set to zero. If it is zero for any other species, the effect will be to give that species zero velocity. 87820f4b53cSBarry Smith 87920f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSize()`, `DMSwarmInitializeCoordinates()`, `DMSwarmInitializeVelocitiesFromOptions()` 88035b38c60SMatthew G. Knepley @*/ 881d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmInitializeVelocities(DM sw, PetscProbFunc sampler, const PetscReal v0[]) 882d71ae5a4SJacob Faibussowitsch { 8838434afd1SBarry Smith PetscSimplePointFn *velFunc; 88435b38c60SMatthew G. Knepley PetscReal *v; 88535b38c60SMatthew G. Knepley PetscInt *species; 8866c5a79ebSMatthew G. Knepley void *ctx; 88735b38c60SMatthew G. Knepley PetscInt dim, Np, p; 88835b38c60SMatthew G. Knepley 88935b38c60SMatthew G. Knepley PetscFunctionBegin; 8906c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetVelocityFunction(sw, &velFunc)); 89135b38c60SMatthew G. Knepley 8929566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 8939566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(sw, &Np)); 8949566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&v)); 8959566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 8966c5a79ebSMatthew G. Knepley if (v0[0] == 0.) { 8976c5a79ebSMatthew G. Knepley PetscCall(PetscArrayzero(v, Np * dim)); 8986c5a79ebSMatthew G. Knepley } else if (velFunc) { 8996c5a79ebSMatthew G. Knepley PetscCall(DMGetApplicationContext(sw, &ctx)); 9006c5a79ebSMatthew G. Knepley for (p = 0; p < Np; ++p) { 9016c5a79ebSMatthew G. Knepley PetscInt s = species[p], d; 9026c5a79ebSMatthew G. Knepley PetscScalar vel[3]; 9036c5a79ebSMatthew G. Knepley 9046c5a79ebSMatthew G. Knepley PetscCall((*velFunc)(dim, 0., NULL, p, vel, ctx)); 9056c5a79ebSMatthew G. Knepley for (d = 0; d < dim; ++d) v[p * dim + d] = (v0[s] / v0[0]) * PetscRealPart(vel[d]); 9066c5a79ebSMatthew G. Knepley } 9076c5a79ebSMatthew G. Knepley } else { 9086c5a79ebSMatthew G. Knepley PetscRandom rnd; 9096c5a79ebSMatthew G. Knepley 9106c5a79ebSMatthew G. Knepley PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rnd)); 9116c5a79ebSMatthew G. Knepley PetscCall(PetscRandomSetInterval(rnd, 0, 1.)); 9126c5a79ebSMatthew G. Knepley PetscCall(PetscRandomSetFromOptions(rnd)); 9136c5a79ebSMatthew G. Knepley 91435b38c60SMatthew G. Knepley for (p = 0; p < Np; ++p) { 91535b38c60SMatthew G. Knepley PetscInt s = species[p], d; 91635b38c60SMatthew G. Knepley PetscReal a[3], vel[3]; 91735b38c60SMatthew G. Knepley 9189566063dSJacob Faibussowitsch for (d = 0; d < dim; ++d) PetscCall(PetscRandomGetValueReal(rnd, &a[d])); 9199566063dSJacob Faibussowitsch PetscCall(sampler(a, NULL, vel)); 920ad540459SPierre Jolivet for (d = 0; d < dim; ++d) v[p * dim + d] = (v0[s] / v0[0]) * vel[d]; 92135b38c60SMatthew G. Knepley } 9226c5a79ebSMatthew G. Knepley PetscCall(PetscRandomDestroy(&rnd)); 9236c5a79ebSMatthew G. Knepley } 9249566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&v)); 9259566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 9263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 92735b38c60SMatthew G. Knepley } 92835b38c60SMatthew G. Knepley 92935b38c60SMatthew G. Knepley /*@ 93035b38c60SMatthew G. Knepley DMSwarmInitializeVelocitiesFromOptions - Set the initial velocities of particles using a distribution determined from options. 93135b38c60SMatthew G. Knepley 93220f4b53cSBarry Smith Collective 93335b38c60SMatthew G. Knepley 93435b38c60SMatthew G. Knepley Input Parameters: 93520f4b53cSBarry Smith + sw - The `DMSWARM` object 93635b38c60SMatthew G. Knepley - v0 - The velocity scale for nondimensionalization for each species 93735b38c60SMatthew G. Knepley 93835b38c60SMatthew G. Knepley Level: advanced 93935b38c60SMatthew G. Knepley 94020f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSize()`, `DMSwarmInitializeCoordinates()`, `DMSwarmInitializeVelocities()` 94135b38c60SMatthew G. Knepley @*/ 942d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmInitializeVelocitiesFromOptions(DM sw, const PetscReal v0[]) 943d71ae5a4SJacob Faibussowitsch { 94435b38c60SMatthew G. Knepley PetscProbFunc sampler; 94535b38c60SMatthew G. Knepley PetscInt dim; 94635b38c60SMatthew G. Knepley const char *prefix; 9476c5a79ebSMatthew G. Knepley char funcname[PETSC_MAX_PATH_LEN]; 9486c5a79ebSMatthew G. Knepley PetscBool flg; 94935b38c60SMatthew G. Knepley 95035b38c60SMatthew G. Knepley PetscFunctionBegin; 951d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)sw), "", "DMSwarm Options", "DMSWARM"); 9526c5a79ebSMatthew G. Knepley PetscCall(PetscOptionsString("-dm_swarm_velocity_function", "Function to determine particle velocities", "DMSwarmSetVelocityFunction", funcname, funcname, sizeof(funcname), &flg)); 953d0609cedSBarry Smith PetscOptionsEnd(); 9546c5a79ebSMatthew G. Knepley if (flg) { 9558434afd1SBarry Smith PetscSimplePointFn *velFunc; 9566c5a79ebSMatthew G. Knepley 9576c5a79ebSMatthew G. Knepley PetscCall(PetscDLSym(NULL, funcname, (void **)&velFunc)); 9586c5a79ebSMatthew G. Knepley PetscCheck(velFunc, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname); 9596c5a79ebSMatthew G. Knepley PetscCall(DMSwarmSetVelocityFunction(sw, velFunc)); 9606c5a79ebSMatthew G. Knepley } 9619566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 9629566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sw, &prefix)); 9639566063dSJacob Faibussowitsch PetscCall(PetscProbCreateFromOptions(dim, prefix, "-dm_swarm_velocity_density", NULL, NULL, &sampler)); 9649566063dSJacob Faibussowitsch PetscCall(DMSwarmInitializeVelocities(sw, sampler, v0)); 9653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 96635b38c60SMatthew G. Knepley } 967