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 */ 120e2ec84fSDave May /*@C 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` 28*a3b724e8SBarry 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 @*/ 33d71ae5a4SJacob Faibussowitsch PETSC_EXTERN 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; 480e2ec84fSDave May 490e2ec84fSDave May PetscFunctionBegin; 500e2ec84fSDave May DMSWARMPICVALID(dm); 519566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 52c448b993SMatthew G. Knepley PetscCall(DMGetLocalBoundingBox(celldm, lmin, lmax)); 53c448b993SMatthew G. Knepley PetscCall(DMGetCoordinateDim(celldm, &bs)); 540e2ec84fSDave May 550e2ec84fSDave May for (b = 0; b < bs; b++) { 5671844bb8SDave May if (npoints[b] > 1) { 570e2ec84fSDave May dx[b] = (max[b] - min[b]) / ((PetscReal)(npoints[b] - 1)); 5871844bb8SDave May } else { 5971844bb8SDave May dx[b] = 0.0; 6071844bb8SDave May } 612e3d3761SDave May _npoints[b] = npoints[b]; 620e2ec84fSDave May } 630e2ec84fSDave May 640e2ec84fSDave May /* determine number of points living in the bounding box */ 650e2ec84fSDave May n_estimate = 0; 662e3d3761SDave May for (k = 0; k < _npoints[2]; k++) { 672e3d3761SDave May for (j = 0; j < _npoints[1]; j++) { 682e3d3761SDave May for (i = 0; i < _npoints[0]; i++) { 690e2ec84fSDave May PetscReal xp[] = {0.0, 0.0, 0.0}; 700e2ec84fSDave May PetscInt ijk[3]; 710e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 720e2ec84fSDave May 730e2ec84fSDave May ijk[0] = i; 740e2ec84fSDave May ijk[1] = j; 750e2ec84fSDave May ijk[2] = k; 76ad540459SPierre Jolivet for (b = 0; b < bs; b++) xp[b] = min[b] + ijk[b] * dx[b]; 770e2ec84fSDave May for (b = 0; b < bs; b++) { 78c448b993SMatthew G. Knepley if (xp[b] < lmin[b]) point_inside = PETSC_FALSE; 79c448b993SMatthew G. Knepley if (xp[b] > lmax[b]) point_inside = PETSC_FALSE; 800e2ec84fSDave May } 81ad540459SPierre Jolivet if (point_inside) n_estimate++; 820e2ec84fSDave May } 830e2ec84fSDave May } 840e2ec84fSDave May } 850e2ec84fSDave May 860e2ec84fSDave May /* create candidate list */ 879566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &pos)); 889566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pos, bs * n_estimate, PETSC_DECIDE)); 899566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(pos, bs)); 909566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pos)); 919566063dSJacob Faibussowitsch PetscCall(VecGetArray(pos, &_pos)); 920e2ec84fSDave May 930e2ec84fSDave May n_estimate = 0; 942e3d3761SDave May for (k = 0; k < _npoints[2]; k++) { 952e3d3761SDave May for (j = 0; j < _npoints[1]; j++) { 962e3d3761SDave May for (i = 0; i < _npoints[0]; i++) { 970e2ec84fSDave May PetscReal xp[] = {0.0, 0.0, 0.0}; 980e2ec84fSDave May PetscInt ijk[3]; 990e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 1000e2ec84fSDave May 1010e2ec84fSDave May ijk[0] = i; 1020e2ec84fSDave May ijk[1] = j; 1030e2ec84fSDave May ijk[2] = k; 104ad540459SPierre Jolivet for (b = 0; b < bs; b++) xp[b] = min[b] + ijk[b] * dx[b]; 1050e2ec84fSDave May for (b = 0; b < bs; b++) { 106c448b993SMatthew G. Knepley if (xp[b] < lmin[b]) point_inside = PETSC_FALSE; 107c448b993SMatthew G. Knepley if (xp[b] > lmax[b]) point_inside = PETSC_FALSE; 1080e2ec84fSDave May } 1090e2ec84fSDave May if (point_inside) { 110ad540459SPierre Jolivet for (b = 0; b < bs; b++) _pos[bs * n_estimate + b] = xp[b]; 1110e2ec84fSDave May n_estimate++; 1120e2ec84fSDave May } 1130e2ec84fSDave May } 1140e2ec84fSDave May } 1150e2ec84fSDave May } 1169566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pos, &_pos)); 1170e2ec84fSDave May 1180e2ec84fSDave May /* locate points */ 1199566063dSJacob Faibussowitsch PetscCall(DMLocatePoints(celldm, pos, DM_POINTLOCATION_NONE, &sfcell)); 1209566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 1210e2ec84fSDave May n_found = 0; 1220e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 123ad540459SPierre Jolivet if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) n_found++; 1240e2ec84fSDave May } 1250e2ec84fSDave May 1260e2ec84fSDave May /* adjust size */ 1270e2ec84fSDave May if (mode == ADD_VALUES) { 1289566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dm, &n_curr)); 1290e2ec84fSDave May n_new_est = n_curr + n_found; 1309566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 1310e2ec84fSDave May } 1320e2ec84fSDave May if (mode == INSERT_VALUES) { 1330e2ec84fSDave May n_curr = 0; 1340e2ec84fSDave May n_new_est = n_found; 1359566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 1360e2ec84fSDave May } 1370e2ec84fSDave May 1380e2ec84fSDave May /* initialize new coords, cell owners, pid */ 1399566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pos, &_coor)); 1409566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 1419566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 1420e2ec84fSDave May n_found = 0; 1430e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 1440e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 145ad540459SPierre Jolivet for (b = 0; b < bs; b++) swarm_coor[bs * (n_curr + n_found) + b] = PetscRealPart(_coor[bs * p + b]); 1460e2ec84fSDave May swarm_cellid[n_curr + n_found] = LA_sfcell[p].index; 1470e2ec84fSDave May n_found++; 1480e2ec84fSDave May } 1490e2ec84fSDave May } 1509566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 1519566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 1529566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pos, &_coor)); 1530e2ec84fSDave May 1549566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfcell)); 1559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pos)); 1563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1570e2ec84fSDave May } 1580e2ec84fSDave May 1590e2ec84fSDave May /*@C 16020f4b53cSBarry Smith DMSwarmSetPointCoordinates - Set point coordinates in a `DMSWARM` from a user defined list 1610e2ec84fSDave May 16220f4b53cSBarry Smith Collective 1630e2ec84fSDave May 16460225df5SJacob Faibussowitsch Input Parameters: 16520f4b53cSBarry Smith + dm - the `DMSWARM` 1660e2ec84fSDave May . npoints - the number of points to insert 1670e2ec84fSDave May . coor - the coordinate values 16820f4b53cSBarry 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 16920f4b53cSBarry Smith - mode - indicates whether to append points to the swarm (`ADD_VALUES`), or over-ride existing points (`INSERT_VALUES`) 1700e2ec84fSDave May 1710e2ec84fSDave May Level: beginner 1720e2ec84fSDave May 1730e2ec84fSDave May Notes: 17420f4b53cSBarry Smith If the user has specified `redundant` as `PETSC_FALSE`, the cell `DM` will attempt to locate the coordinates provided by `coor` within 17520f4b53cSBarry 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 17620f4b53cSBarry Smith added to the `DMSWARM`. 1770e2ec84fSDave May 17820f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType`, `DMSwarmSetPointsUniformCoordinates()` 1790e2ec84fSDave May @*/ 180d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinates(DM dm, PetscInt npoints, PetscReal coor[], PetscBool redundant, InsertMode mode) 181d71ae5a4SJacob Faibussowitsch { 1820e2ec84fSDave May PetscReal gmin[] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL}; 1830e2ec84fSDave May PetscReal gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 1840e2ec84fSDave May PetscInt i, N, bs, b, n_estimate, n_curr, n_new_est, p, n_found; 1850e2ec84fSDave May Vec coorlocal; 1860e2ec84fSDave May const PetscScalar *_coor; 1870e2ec84fSDave May DM celldm; 1880e2ec84fSDave May Vec pos; 1890e2ec84fSDave May PetscScalar *_pos; 1900e2ec84fSDave May PetscReal *swarm_coor; 1910e2ec84fSDave May PetscInt *swarm_cellid; 1920e2ec84fSDave May PetscSF sfcell = NULL; 1930e2ec84fSDave May const PetscSFNode *LA_sfcell; 1940e2ec84fSDave May PetscReal *my_coor; 1950e2ec84fSDave May PetscInt my_npoints; 1960e2ec84fSDave May PetscMPIInt rank; 1970e2ec84fSDave May MPI_Comm comm; 1980e2ec84fSDave May 1990e2ec84fSDave May PetscFunctionBegin; 2000e2ec84fSDave May DMSWARMPICVALID(dm); 2019566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 2029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2030e2ec84fSDave May 2049566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 2059566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(celldm, &coorlocal)); 2069566063dSJacob Faibussowitsch PetscCall(VecGetSize(coorlocal, &N)); 2079566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(coorlocal, &bs)); 2080e2ec84fSDave May N = N / bs; 2099566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coorlocal, &_coor)); 2100e2ec84fSDave May for (i = 0; i < N; i++) { 2110e2ec84fSDave May for (b = 0; b < bs; b++) { 212a5f152d1SDave May gmin[b] = PetscMin(gmin[b], PetscRealPart(_coor[bs * i + b])); 213a5f152d1SDave May gmax[b] = PetscMax(gmax[b], PetscRealPart(_coor[bs * i + b])); 2140e2ec84fSDave May } 2150e2ec84fSDave May } 2169566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coorlocal, &_coor)); 2170e2ec84fSDave May 2180e2ec84fSDave May /* broadcast points from rank 0 if requested */ 2190e2ec84fSDave May if (redundant) { 2200e2ec84fSDave May my_npoints = npoints; 2219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&my_npoints, 1, MPIU_INT, 0, comm)); 2220e2ec84fSDave May 2230e2ec84fSDave May if (rank > 0) { /* allocate space */ 2249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * my_npoints, &my_coor)); 2250e2ec84fSDave May } else { 2260e2ec84fSDave May my_coor = coor; 2270e2ec84fSDave May } 2289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(my_coor, bs * my_npoints, MPIU_REAL, 0, comm)); 2290e2ec84fSDave May } else { 2300e2ec84fSDave May my_npoints = npoints; 2310e2ec84fSDave May my_coor = coor; 2320e2ec84fSDave May } 2330e2ec84fSDave May 2340e2ec84fSDave May /* determine the number of points living in the bounding box */ 2350e2ec84fSDave May n_estimate = 0; 2360e2ec84fSDave May for (i = 0; i < my_npoints; i++) { 2370e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 2380e2ec84fSDave May 2390e2ec84fSDave May for (b = 0; b < bs; b++) { 240ad540459SPierre Jolivet if (my_coor[bs * i + b] < gmin[b]) point_inside = PETSC_FALSE; 241ad540459SPierre Jolivet if (my_coor[bs * i + b] > gmax[b]) point_inside = PETSC_FALSE; 2420e2ec84fSDave May } 243ad540459SPierre Jolivet if (point_inside) n_estimate++; 2440e2ec84fSDave May } 2450e2ec84fSDave May 2460e2ec84fSDave May /* create candidate list */ 2479566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &pos)); 2489566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pos, bs * n_estimate, PETSC_DECIDE)); 2499566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(pos, bs)); 2509566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pos)); 2519566063dSJacob Faibussowitsch PetscCall(VecGetArray(pos, &_pos)); 2520e2ec84fSDave May 2530e2ec84fSDave May n_estimate = 0; 2540e2ec84fSDave May for (i = 0; i < my_npoints; i++) { 2550e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 2560e2ec84fSDave May 2570e2ec84fSDave May for (b = 0; b < bs; b++) { 258ad540459SPierre Jolivet if (my_coor[bs * i + b] < gmin[b]) point_inside = PETSC_FALSE; 259ad540459SPierre Jolivet if (my_coor[bs * i + b] > gmax[b]) point_inside = PETSC_FALSE; 2600e2ec84fSDave May } 2610e2ec84fSDave May if (point_inside) { 262ad540459SPierre Jolivet for (b = 0; b < bs; b++) _pos[bs * n_estimate + b] = my_coor[bs * i + b]; 2630e2ec84fSDave May n_estimate++; 2640e2ec84fSDave May } 2650e2ec84fSDave May } 2669566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pos, &_pos)); 2670e2ec84fSDave May 2680e2ec84fSDave May /* locate points */ 2699566063dSJacob Faibussowitsch PetscCall(DMLocatePoints(celldm, pos, DM_POINTLOCATION_NONE, &sfcell)); 2700e2ec84fSDave May 2719566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 2720e2ec84fSDave May n_found = 0; 2730e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 274ad540459SPierre Jolivet if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) n_found++; 2750e2ec84fSDave May } 2760e2ec84fSDave May 2770e2ec84fSDave May /* adjust size */ 2780e2ec84fSDave May if (mode == ADD_VALUES) { 2799566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dm, &n_curr)); 2800e2ec84fSDave May n_new_est = n_curr + n_found; 2819566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 2820e2ec84fSDave May } 2830e2ec84fSDave May if (mode == INSERT_VALUES) { 2840e2ec84fSDave May n_curr = 0; 2850e2ec84fSDave May n_new_est = n_found; 2869566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 2870e2ec84fSDave May } 2880e2ec84fSDave May 2890e2ec84fSDave May /* initialize new coords, cell owners, pid */ 2909566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pos, &_coor)); 2919566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 2929566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 2930e2ec84fSDave May n_found = 0; 2940e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 2950e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 296ad540459SPierre Jolivet for (b = 0; b < bs; b++) swarm_coor[bs * (n_curr + n_found) + b] = PetscRealPart(_coor[bs * p + b]); 2970e2ec84fSDave May swarm_cellid[n_curr + n_found] = LA_sfcell[p].index; 2980e2ec84fSDave May n_found++; 2990e2ec84fSDave May } 3000e2ec84fSDave May } 3019566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 3029566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 3039566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pos, &_coor)); 3040e2ec84fSDave May 3050e2ec84fSDave May if (redundant) { 30648a46eb9SPierre Jolivet if (rank > 0) PetscCall(PetscFree(my_coor)); 3070e2ec84fSDave May } 3089566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfcell)); 3099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pos)); 3103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3110e2ec84fSDave May } 3120e2ec84fSDave May 3130e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM, DM, DMSwarmPICLayoutType, PetscInt); 3140e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM, DM, DMSwarmPICLayoutType, PetscInt); 3150e2ec84fSDave May 3160e2ec84fSDave May /*@C 3170e2ec84fSDave May DMSwarmInsertPointsUsingCellDM - Insert point coordinates within each cell 3180e2ec84fSDave May 31920f4b53cSBarry Smith Not Collective 3200e2ec84fSDave May 32160225df5SJacob Faibussowitsch Input Parameters: 32220f4b53cSBarry Smith + dm - the `DMSWARM` 32320f4b53cSBarry Smith . layout_type - method used to fill each cell with the cell `DM` 3240e2ec84fSDave May - fill_param - parameter controlling how many points per cell are added (the meaning of this parameter is dependent on the layout type) 3250e2ec84fSDave May 3260e2ec84fSDave May Level: beginner 3270e2ec84fSDave May 3280e2ec84fSDave May Notes: 32920f4b53cSBarry Smith The insert method will reset any previous defined points within the `DMSWARM`. 330e7af74c8SDave May 33120f4b53cSBarry Smith When using a `DMDA` both 2D and 3D are supported for all layout types provided you are using `DMDA_ELEMENT_Q1`. 332e7af74c8SDave May 333a4e35b19SJacob Faibussowitsch When using a `DMPLEX` the following case are supported\: 33420f4b53cSBarry Smith .vb 335ea3d7275SDave May (i) DMSWARMPIC_LAYOUT_REGULAR: 2D (triangle), 336ea3d7275SDave May (ii) DMSWARMPIC_LAYOUT_GAUSS: 2D and 3D provided the cell is a tri/tet or a quad/hex, 337ea3d7275SDave May (iii) DMSWARMPIC_LAYOUT_SUBDIVISION: 2D and 3D for quad/hex and 2D tri. 33820f4b53cSBarry Smith .ve 3390e2ec84fSDave May 34020f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmPICLayoutType`, `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType` 3410e2ec84fSDave May @*/ 342d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSwarmInsertPointsUsingCellDM(DM dm, DMSwarmPICLayoutType layout_type, PetscInt fill_param) 343d71ae5a4SJacob Faibussowitsch { 3440e2ec84fSDave May DM celldm; 3450e2ec84fSDave May PetscBool isDA, isPLEX; 3460e2ec84fSDave May 3470e2ec84fSDave May PetscFunctionBegin; 3480e2ec84fSDave May DMSWARMPICVALID(dm); 3499566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 3509566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA)); 3519566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX)); 3520e2ec84fSDave May if (isDA) { 3539566063dSJacob Faibussowitsch PetscCall(private_DMSwarmInsertPointsUsingCellDM_DA(dm, celldm, layout_type, fill_param)); 3540e2ec84fSDave May } else if (isPLEX) { 3559566063dSJacob Faibussowitsch PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX(dm, celldm, layout_type, fill_param)); 3560e2ec84fSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX"); 3573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3580e2ec84fSDave May } 3590e2ec84fSDave May 360431382f2SDave May extern PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM, DM, PetscInt, PetscReal *); 361431382f2SDave May 362431382f2SDave May /*@C 363431382f2SDave May DMSwarmSetPointCoordinatesCellwise - Insert point coordinates (defined over the reference cell) within each cell 364431382f2SDave May 36520f4b53cSBarry Smith Not Collective 366431382f2SDave May 36760225df5SJacob Faibussowitsch Input Parameters: 36820f4b53cSBarry Smith + dm - the `DMSWARM` 369431382f2SDave May . npoints - the number of points to insert in each cell 370431382f2SDave May - xi - the coordinates (defined in the local coordinate system for each cell) to insert 371431382f2SDave May 372431382f2SDave May Level: beginner 373431382f2SDave May 374431382f2SDave May Notes: 37520f4b53cSBarry Smith The method will reset any previous defined points within the `DMSWARM`. 37620f4b53cSBarry Smith Only supported for `DMPLEX`. If you are using a `DMDA` it is recommended to either use 37720f4b53cSBarry Smith `DMSwarmInsertPointsUsingCellDM()`, or extract and set the coordinates yourself the following code 37820f4b53cSBarry Smith .vb 37920f4b53cSBarry Smith PetscReal *coor; 38020f4b53cSBarry Smith DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&coor); 38120f4b53cSBarry Smith // user code to define the coordinates here 38220f4b53cSBarry Smith DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&coor); 38320f4b53cSBarry Smith .ve 384e7af74c8SDave May 38520f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetCellDM()`, `DMSwarmInsertPointsUsingCellDM()` 386431382f2SDave May @*/ 387d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinatesCellwise(DM dm, PetscInt npoints, PetscReal xi[]) 388d71ae5a4SJacob Faibussowitsch { 389431382f2SDave May DM celldm; 390431382f2SDave May PetscBool isDA, isPLEX; 391431382f2SDave May 3920e2ec84fSDave May PetscFunctionBegin; 393431382f2SDave May DMSWARMPICVALID(dm); 3949566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 3959566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA)); 3969566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX)); 39728b400f6SJacob Faibussowitsch PetscCheck(!isDA, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMPLEX. Recommended you use DMSwarmInsertPointsUsingCellDM()"); 398f7d195e4SLawrence Mitchell if (isPLEX) { 3999566063dSJacob Faibussowitsch PetscCall(private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm, celldm, npoints, xi)); 400431382f2SDave May } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX"); 4013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4020e2ec84fSDave May } 403431382f2SDave May 4040e2ec84fSDave May /*@C 405b799feefSDave May DMSwarmCreatePointPerCellCount - Count the number of points within all cells in the cell DM 4060e2ec84fSDave May 40720f4b53cSBarry Smith Not Collective 4080e2ec84fSDave May 40960225df5SJacob Faibussowitsch Input Parameter: 41020f4b53cSBarry Smith . dm - the `DMSWARM` 4110e2ec84fSDave May 41260225df5SJacob Faibussowitsch Output Parameters: 41320f4b53cSBarry Smith + ncells - the number of cells in the cell `DM` (optional argument, pass `NULL` to ignore) 414b799feefSDave May - count - array of length ncells containing the number of points per cell 4150e2ec84fSDave May 4160e2ec84fSDave May Level: beginner 4170e2ec84fSDave May 4180e2ec84fSDave May Notes: 4190e2ec84fSDave May The array count is allocated internally and must be free'd by the user. 4200e2ec84fSDave May 42120f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType` 4220e2ec84fSDave May @*/ 423d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSwarmCreatePointPerCellCount(DM dm, PetscInt *ncells, PetscInt **count) 424d71ae5a4SJacob Faibussowitsch { 425b799feefSDave May PetscBool isvalid; 426b799feefSDave May PetscInt nel; 427b799feefSDave May PetscInt *sum; 428b799feefSDave May 4290e2ec84fSDave May PetscFunctionBegin; 4309566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetIsValid(dm, &isvalid)); 431b799feefSDave May nel = 0; 432b799feefSDave May if (isvalid) { 433b799feefSDave May PetscInt e; 434b799feefSDave May 4359566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetSizes(dm, &nel, NULL)); 436b799feefSDave May 4379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nel, &sum)); 43848a46eb9SPierre Jolivet for (e = 0; e < nel; e++) PetscCall(DMSwarmSortGetNumberOfPointsPerCell(dm, e, &sum[e])); 439b799feefSDave May } else { 440b799feefSDave May DM celldm; 441b799feefSDave May PetscBool isda, isplex, isshell; 442b799feefSDave May PetscInt p, npoints; 443b799feefSDave May PetscInt *swarm_cellid; 444b799feefSDave May 445b799feefSDave May /* get the number of cells */ 4469566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 4479566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isda)); 4489566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isplex)); 4499566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMSHELL, &isshell)); 450b799feefSDave May if (isda) { 451b799feefSDave May PetscInt _nel, _npe; 452b799feefSDave May const PetscInt *_element; 453b799feefSDave May 4549566063dSJacob Faibussowitsch PetscCall(DMDAGetElements(celldm, &_nel, &_npe, &_element)); 455b799feefSDave May nel = _nel; 4569566063dSJacob Faibussowitsch PetscCall(DMDARestoreElements(celldm, &_nel, &_npe, &_element)); 457b799feefSDave May } else if (isplex) { 458b799feefSDave May PetscInt ps, pe; 459b799feefSDave May 4609566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(celldm, 0, &ps, &pe)); 461b799feefSDave May nel = pe - ps; 462b799feefSDave May } else if (isshell) { 463b799feefSDave May PetscErrorCode (*method_DMShellGetNumberOfCells)(DM, PetscInt *); 464b799feefSDave May 4659566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)celldm, "DMGetNumberOfCells_C", &method_DMShellGetNumberOfCells)); 466b799feefSDave May if (method_DMShellGetNumberOfCells) { 4679566063dSJacob Faibussowitsch PetscCall(method_DMShellGetNumberOfCells(celldm, &nel)); 4689371c9d4SSatish Balay } else 4699371c9d4SSatish 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);"); 470b799feefSDave 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"); 471b799feefSDave May 4729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nel, &sum)); 4739566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(sum, nel)); 4749566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dm, &npoints)); 4759566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 476b799feefSDave May for (p = 0; p < npoints; p++) { 477ad540459SPierre Jolivet if (swarm_cellid[p] != DMLOCATEPOINT_POINT_NOT_FOUND) sum[swarm_cellid[p]]++; 478b799feefSDave May } 4799566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 480b799feefSDave May } 481ad540459SPierre Jolivet if (ncells) *ncells = nel; 482b799feefSDave May *count = sum; 4833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4840e2ec84fSDave May } 48535b38c60SMatthew G. Knepley 48635b38c60SMatthew G. Knepley /*@ 48735b38c60SMatthew G. Knepley DMSwarmGetNumSpecies - Get the number of particle species 48835b38c60SMatthew G. Knepley 48920f4b53cSBarry Smith Not Collective 49035b38c60SMatthew G. Knepley 49160225df5SJacob Faibussowitsch Input Parameter: 49260225df5SJacob Faibussowitsch . sw - the `DMSWARM` 49335b38c60SMatthew G. Knepley 49460225df5SJacob Faibussowitsch Output Parameters: 49535b38c60SMatthew G. Knepley . Ns - the number of species 49635b38c60SMatthew G. Knepley 49735b38c60SMatthew G. Knepley Level: intermediate 49835b38c60SMatthew G. Knepley 49920f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmSetNumSpecies()`, `DMSwarmSetType()`, `DMSwarmType` 50035b38c60SMatthew G. Knepley @*/ 501d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetNumSpecies(DM sw, PetscInt *Ns) 502d71ae5a4SJacob Faibussowitsch { 50335b38c60SMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 50435b38c60SMatthew G. Knepley 50535b38c60SMatthew G. Knepley PetscFunctionBegin; 50635b38c60SMatthew G. Knepley *Ns = swarm->Ns; 5073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 50835b38c60SMatthew G. Knepley } 50935b38c60SMatthew G. Knepley 51035b38c60SMatthew G. Knepley /*@ 51135b38c60SMatthew G. Knepley DMSwarmSetNumSpecies - Set the number of particle species 51235b38c60SMatthew G. Knepley 51320f4b53cSBarry Smith Not Collective 51435b38c60SMatthew G. Knepley 51560225df5SJacob Faibussowitsch Input Parameters: 51660225df5SJacob Faibussowitsch + sw - the `DMSWARM` 51735b38c60SMatthew G. Knepley - Ns - the number of species 51835b38c60SMatthew G. Knepley 51935b38c60SMatthew G. Knepley Level: intermediate 52035b38c60SMatthew G. Knepley 52120f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmGetNumSpecies()`, `DMSwarmSetType()`, `DMSwarmType` 52235b38c60SMatthew G. Knepley @*/ 523d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetNumSpecies(DM sw, PetscInt Ns) 524d71ae5a4SJacob Faibussowitsch { 52535b38c60SMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 52635b38c60SMatthew G. Knepley 52735b38c60SMatthew G. Knepley PetscFunctionBegin; 52835b38c60SMatthew G. Knepley swarm->Ns = Ns; 5293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 53035b38c60SMatthew G. Knepley } 53135b38c60SMatthew G. Knepley 53235b38c60SMatthew G. Knepley /*@C 5336c5a79ebSMatthew G. Knepley DMSwarmGetCoordinateFunction - Get the function setting initial particle positions, if it exists 5346c5a79ebSMatthew G. Knepley 53520f4b53cSBarry Smith Not Collective 5366c5a79ebSMatthew G. Knepley 53760225df5SJacob Faibussowitsch Input Parameter: 53860225df5SJacob Faibussowitsch . sw - the `DMSWARM` 5396c5a79ebSMatthew G. Knepley 5406c5a79ebSMatthew G. Knepley Output Parameter: 5418434afd1SBarry Smith . coordFunc - the function setting initial particle positions, or `NULL`, see `PetscSimplePointFn` for the calling sequence 5426c5a79ebSMatthew G. Knepley 5436c5a79ebSMatthew G. Knepley Level: intermediate 5446c5a79ebSMatthew G. Knepley 5458434afd1SBarry Smith .seealso: `DMSWARM`, `DMSwarmSetCoordinateFunction()`, `DMSwarmGetVelocityFunction()`, `DMSwarmInitializeCoordinates()`, `PetscSimplePointFn` 5466c5a79ebSMatthew G. Knepley @*/ 5478434afd1SBarry Smith PetscErrorCode DMSwarmGetCoordinateFunction(DM sw, PetscSimplePointFn **coordFunc) 548d71ae5a4SJacob Faibussowitsch { 5496c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 5506c5a79ebSMatthew G. Knepley 5516c5a79ebSMatthew G. Knepley PetscFunctionBegin; 5526c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 5534f572ea9SToby Isaac PetscAssertPointer(coordFunc, 2); 5546c5a79ebSMatthew G. Knepley *coordFunc = swarm->coordFunc; 5553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5566c5a79ebSMatthew G. Knepley } 5576c5a79ebSMatthew G. Knepley 5586c5a79ebSMatthew G. Knepley /*@C 5596c5a79ebSMatthew G. Knepley DMSwarmSetCoordinateFunction - Set the function setting initial particle positions 5606c5a79ebSMatthew G. Knepley 56120f4b53cSBarry Smith Not Collective 5626c5a79ebSMatthew G. Knepley 56360225df5SJacob Faibussowitsch Input Parameters: 56460225df5SJacob Faibussowitsch + sw - the `DMSWARM` 5658434afd1SBarry Smith - coordFunc - the function setting initial particle positions, see `PetscSimplePointFn` for the calling sequence 5666c5a79ebSMatthew G. Knepley 5676c5a79ebSMatthew G. Knepley Level: intermediate 5686c5a79ebSMatthew G. Knepley 5698434afd1SBarry Smith .seealso: `DMSWARM`, `DMSwarmGetCoordinateFunction()`, `DMSwarmSetVelocityFunction()`, `DMSwarmInitializeCoordinates()`, `PetscSimplePointFn` 5706c5a79ebSMatthew G. Knepley @*/ 5718434afd1SBarry Smith PetscErrorCode DMSwarmSetCoordinateFunction(DM sw, PetscSimplePointFn *coordFunc) 572d71ae5a4SJacob Faibussowitsch { 5736c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 5746c5a79ebSMatthew G. Knepley 5756c5a79ebSMatthew G. Knepley PetscFunctionBegin; 5766c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 5776c5a79ebSMatthew G. Knepley PetscValidFunction(coordFunc, 2); 5786c5a79ebSMatthew G. Knepley swarm->coordFunc = coordFunc; 5793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5806c5a79ebSMatthew G. Knepley } 5816c5a79ebSMatthew G. Knepley 5826c5a79ebSMatthew G. Knepley /*@C 58360225df5SJacob Faibussowitsch DMSwarmGetVelocityFunction - Get the function setting initial particle velocities, if it exists 5846c5a79ebSMatthew G. Knepley 58520f4b53cSBarry Smith Not Collective 5866c5a79ebSMatthew G. Knepley 58760225df5SJacob Faibussowitsch Input Parameter: 58860225df5SJacob Faibussowitsch . sw - the `DMSWARM` 5896c5a79ebSMatthew G. Knepley 5906c5a79ebSMatthew G. Knepley Output Parameter: 5918434afd1SBarry Smith . velFunc - the function setting initial particle velocities, or `NULL`, see `PetscSimplePointFn` for the calling sequence 5926c5a79ebSMatthew G. Knepley 5936c5a79ebSMatthew G. Knepley Level: intermediate 5946c5a79ebSMatthew G. Knepley 5958434afd1SBarry Smith .seealso: `DMSWARM`, `DMSwarmSetVelocityFunction()`, `DMSwarmGetCoordinateFunction()`, `DMSwarmInitializeVelocities()`, `PetscSimplePointFn` 5966c5a79ebSMatthew G. Knepley @*/ 5978434afd1SBarry Smith PetscErrorCode DMSwarmGetVelocityFunction(DM sw, PetscSimplePointFn **velFunc) 598d71ae5a4SJacob Faibussowitsch { 5996c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 6006c5a79ebSMatthew G. Knepley 6016c5a79ebSMatthew G. Knepley PetscFunctionBegin; 6026c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 6034f572ea9SToby Isaac PetscAssertPointer(velFunc, 2); 6046c5a79ebSMatthew G. Knepley *velFunc = swarm->velFunc; 6053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6066c5a79ebSMatthew G. Knepley } 6076c5a79ebSMatthew G. Knepley 6086c5a79ebSMatthew G. Knepley /*@C 6096c5a79ebSMatthew G. Knepley DMSwarmSetVelocityFunction - Set the function setting initial particle velocities 6106c5a79ebSMatthew G. Knepley 61120f4b53cSBarry Smith Not Collective 6126c5a79ebSMatthew G. Knepley 61360225df5SJacob Faibussowitsch Input Parameters: 614a4e35b19SJacob Faibussowitsch + sw - the `DMSWARM` 6158434afd1SBarry Smith - velFunc - the function setting initial particle velocities, see `PetscSimplePointFn` for the calling sequence 6166c5a79ebSMatthew G. Knepley 6176c5a79ebSMatthew G. Knepley Level: intermediate 6186c5a79ebSMatthew G. Knepley 6198434afd1SBarry Smith .seealso: `DMSWARM`, `DMSwarmGetVelocityFunction()`, `DMSwarmSetCoordinateFunction()`, `DMSwarmInitializeVelocities()`, `PetscSimplePointFn` 6206c5a79ebSMatthew G. Knepley @*/ 6218434afd1SBarry Smith PetscErrorCode DMSwarmSetVelocityFunction(DM sw, PetscSimplePointFn *velFunc) 622d71ae5a4SJacob Faibussowitsch { 6236c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 6246c5a79ebSMatthew G. Knepley 6256c5a79ebSMatthew G. Knepley PetscFunctionBegin; 6266c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 6276c5a79ebSMatthew G. Knepley PetscValidFunction(velFunc, 2); 6286c5a79ebSMatthew G. Knepley swarm->velFunc = velFunc; 6293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6306c5a79ebSMatthew G. Knepley } 6316c5a79ebSMatthew G. Knepley 6326c5a79ebSMatthew G. Knepley /*@C 63335b38c60SMatthew G. Knepley DMSwarmComputeLocalSize - Compute the local number and distribution of particles based upon a density function 63435b38c60SMatthew G. Knepley 63520f4b53cSBarry Smith Not Collective 63635b38c60SMatthew G. Knepley 63735b38c60SMatthew G. Knepley Input Parameters: 63820f4b53cSBarry Smith + sw - The `DMSWARM` 63935b38c60SMatthew G. Knepley . N - The target number of particles 64035b38c60SMatthew G. Knepley - density - The density field for the particle layout, normalized to unity 64135b38c60SMatthew G. Knepley 64235b38c60SMatthew G. Knepley Level: advanced 64335b38c60SMatthew G. Knepley 64420f4b53cSBarry Smith Note: 64520f4b53cSBarry Smith One particle will be created for each species. 64620f4b53cSBarry Smith 64720f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSizeFromOptions()` 64835b38c60SMatthew G. Knepley @*/ 649d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmComputeLocalSize(DM sw, PetscInt N, PetscProbFunc density) 650d71ae5a4SJacob Faibussowitsch { 65135b38c60SMatthew G. Knepley DM dm; 65235b38c60SMatthew G. Knepley PetscQuadrature quad; 65335b38c60SMatthew G. Knepley const PetscReal *xq, *wq; 654ea7032a0SMatthew G. Knepley PetscReal *n_int; 655ea7032a0SMatthew G. Knepley PetscInt *npc_s, *cellid, Ni; 656ea7032a0SMatthew G. Knepley PetscReal gmin[3], gmax[3], xi0[3]; 657ea7032a0SMatthew G. Knepley PetscInt Ns, cStart, cEnd, c, dim, d, Nq, q, Np = 0, p, s; 65835b38c60SMatthew G. Knepley PetscBool simplex; 65935b38c60SMatthew G. Knepley 66035b38c60SMatthew G. Knepley PetscFunctionBegin; 6619566063dSJacob Faibussowitsch PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 6629566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sw, &dm)); 6639566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 664ea7032a0SMatthew G. Knepley PetscCall(DMGetBoundingBox(dm, gmin, gmax)); 6659566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 6669566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 6676858538eSMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(dm)); 6689566063dSJacob Faibussowitsch if (simplex) PetscCall(PetscDTStroudConicalQuadrature(dim, 1, 5, -1.0, 1.0, &quad)); 6699566063dSJacob Faibussowitsch else PetscCall(PetscDTGaussTensorQuadrature(dim, 1, 5, -1.0, 1.0, &quad)); 6709566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &xq, &wq)); 671ea7032a0SMatthew G. Knepley PetscCall(PetscCalloc2(Ns, &n_int, (cEnd - cStart) * Ns, &npc_s)); 67235b38c60SMatthew G. Knepley /* Integrate the density function to get the number of particles in each cell */ 67335b38c60SMatthew G. Knepley for (d = 0; d < dim; ++d) xi0[d] = -1.0; 67435b38c60SMatthew G. Knepley for (c = 0; c < cEnd - cStart; ++c) { 67535b38c60SMatthew G. Knepley const PetscInt cell = c + cStart; 676ea7032a0SMatthew G. Knepley PetscReal v0[3], J[9], invJ[9], detJ, detJp = 2. / (gmax[0] - gmin[0]), xr[3], den; 67735b38c60SMatthew G. Knepley 678ea7032a0SMatthew G. Knepley /* Have to transform quadrature points/weights to cell domain */ 6799566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, NULL, v0, J, invJ, &detJ)); 680ea7032a0SMatthew G. Knepley PetscCall(PetscArrayzero(n_int, Ns)); 68135b38c60SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 68235b38c60SMatthew G. Knepley CoordinatesRefToReal(dim, dim, xi0, v0, J, &xq[q * dim], xr); 683ea7032a0SMatthew G. Knepley /* Have to transform mesh to domain of definition of PDF, [-1, 1], and weight PDF by |J|/2 */ 684ea7032a0SMatthew G. Knepley xr[0] = detJp * (xr[0] - gmin[0]) - 1.; 685ea7032a0SMatthew G. Knepley 686ea7032a0SMatthew G. Knepley for (s = 0; s < Ns; ++s) { 687ea7032a0SMatthew G. Knepley PetscCall(density(xr, NULL, &den)); 688ea7032a0SMatthew G. Knepley n_int[s] += (detJp * den) * (detJ * wq[q]) / (PetscReal)Ns; 68935b38c60SMatthew G. Knepley } 690ea7032a0SMatthew G. Knepley } 691ea7032a0SMatthew G. Knepley for (s = 0; s < Ns; ++s) { 692ea7032a0SMatthew G. Knepley Ni = N; 693ea7032a0SMatthew G. Knepley npc_s[c * Ns + s] += (PetscInt)(Ni * n_int[s]); 694ea7032a0SMatthew G. Knepley Np += npc_s[c * Ns + s]; 695ea7032a0SMatthew G. Knepley } 69635b38c60SMatthew G. Knepley } 6979566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quad)); 6989566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(sw, Np, 0)); 6999566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 70035b38c60SMatthew G. Knepley for (c = 0, p = 0; c < cEnd - cStart; ++c) { 701ea7032a0SMatthew G. Knepley for (s = 0; s < Ns; ++s) { 702ea7032a0SMatthew G. Knepley for (q = 0; q < npc_s[c * Ns + s]; ++q, ++p) cellid[p] = c; 703ea7032a0SMatthew G. Knepley } 70435b38c60SMatthew G. Knepley } 7059566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 706ea7032a0SMatthew G. Knepley PetscCall(PetscFree2(n_int, npc_s)); 7073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 70835b38c60SMatthew G. Knepley } 70935b38c60SMatthew G. Knepley 71035b38c60SMatthew G. Knepley /*@ 71135b38c60SMatthew G. Knepley DMSwarmComputeLocalSizeFromOptions - Compute the local number and distribution of particles based upon a density function determined by options 71235b38c60SMatthew G. Knepley 71320f4b53cSBarry Smith Not Collective 71435b38c60SMatthew G. Knepley 7152fe279fdSBarry Smith Input Parameter: 7162fe279fdSBarry Smith . sw - The `DMSWARM` 71735b38c60SMatthew G. Knepley 71835b38c60SMatthew G. Knepley Level: advanced 71935b38c60SMatthew G. Knepley 72020f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSize()` 72135b38c60SMatthew G. Knepley @*/ 722d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmComputeLocalSizeFromOptions(DM sw) 723d71ae5a4SJacob Faibussowitsch { 72435b38c60SMatthew G. Knepley PetscProbFunc pdf; 72535b38c60SMatthew G. Knepley const char *prefix; 7266c5a79ebSMatthew G. Knepley char funcname[PETSC_MAX_PATH_LEN]; 7276c5a79ebSMatthew G. Knepley PetscInt *N, Ns, dim, n; 7286c5a79ebSMatthew G. Knepley PetscBool flg; 7296c5a79ebSMatthew G. Knepley PetscMPIInt size, rank; 73035b38c60SMatthew G. Knepley 73135b38c60SMatthew G. Knepley PetscFunctionBegin; 7326c5a79ebSMatthew G. Knepley PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sw), &size)); 7336c5a79ebSMatthew G. Knepley PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sw), &rank)); 7346c5a79ebSMatthew G. Knepley PetscCall(PetscCalloc1(size, &N)); 735d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)sw), "", "DMSwarm Options", "DMSWARM"); 7366c5a79ebSMatthew G. Knepley n = size; 7376c5a79ebSMatthew G. Knepley PetscCall(PetscOptionsIntArray("-dm_swarm_num_particles", "The target number of particles", "", N, &n, NULL)); 7386c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 7399566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-dm_swarm_num_species", "The number of species", "DMSwarmSetNumSpecies", Ns, &Ns, &flg)); 7409566063dSJacob Faibussowitsch if (flg) PetscCall(DMSwarmSetNumSpecies(sw, Ns)); 7416c5a79ebSMatthew G. Knepley PetscCall(PetscOptionsString("-dm_swarm_coordinate_function", "Function to determine particle coordinates", "DMSwarmSetCoordinateFunction", funcname, funcname, sizeof(funcname), &flg)); 742d0609cedSBarry Smith PetscOptionsEnd(); 7436c5a79ebSMatthew G. Knepley if (flg) { 7448434afd1SBarry Smith PetscSimplePointFn *coordFunc; 74535b38c60SMatthew G. Knepley 7466c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 7476c5a79ebSMatthew G. Knepley PetscCall(PetscDLSym(NULL, funcname, (void **)&coordFunc)); 7486c5a79ebSMatthew G. Knepley PetscCheck(coordFunc, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname); 7496c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 7506c5a79ebSMatthew G. Knepley PetscCall(DMSwarmSetLocalSizes(sw, N[rank] * Ns, 0)); 7516c5a79ebSMatthew G. Knepley PetscCall(DMSwarmSetCoordinateFunction(sw, coordFunc)); 7526c5a79ebSMatthew G. Knepley } else { 7539566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 7549566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sw, &prefix)); 7559566063dSJacob Faibussowitsch PetscCall(PetscProbCreateFromOptions(dim, prefix, "-dm_swarm_coordinate_density", &pdf, NULL, NULL)); 7566c5a79ebSMatthew G. Knepley PetscCall(DMSwarmComputeLocalSize(sw, N[rank], pdf)); 7576c5a79ebSMatthew G. Knepley } 7586c5a79ebSMatthew G. Knepley PetscCall(PetscFree(N)); 7593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 76035b38c60SMatthew G. Knepley } 76135b38c60SMatthew G. Knepley 76235b38c60SMatthew G. Knepley /*@ 76335b38c60SMatthew G. Knepley DMSwarmInitializeCoordinates - Determine the initial coordinates of particles for a PIC method 76435b38c60SMatthew G. Knepley 76520f4b53cSBarry Smith Not Collective 76635b38c60SMatthew G. Knepley 76720f4b53cSBarry Smith Input Parameter: 76820f4b53cSBarry Smith . sw - The `DMSWARM` 76935b38c60SMatthew G. Knepley 77035b38c60SMatthew G. Knepley Level: advanced 77135b38c60SMatthew G. Knepley 77220f4b53cSBarry Smith Note: 77320f4b53cSBarry Smith Currently, we randomly place particles in their assigned cell 77420f4b53cSBarry Smith 77520f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSize()`, `DMSwarmInitializeVelocities()` 77635b38c60SMatthew G. Knepley @*/ 777d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmInitializeCoordinates(DM sw) 778d71ae5a4SJacob Faibussowitsch { 7798434afd1SBarry Smith PetscSimplePointFn *coordFunc; 78035b38c60SMatthew G. Knepley PetscScalar *weight; 7816c5a79ebSMatthew G. Knepley PetscReal *x; 78235b38c60SMatthew G. Knepley PetscInt *species; 7836c5a79ebSMatthew G. Knepley void *ctx; 78435b38c60SMatthew G. Knepley PetscBool removePoints = PETSC_TRUE; 78535b38c60SMatthew G. Knepley PetscDataType dtype; 7866c5a79ebSMatthew G. Knepley PetscInt Np, p, Ns, dim, d, bs; 78735b38c60SMatthew G. Knepley 78835b38c60SMatthew G. Knepley PetscFunctionBeginUser; 7896c5a79ebSMatthew G. Knepley PetscCall(DMGetDimension(sw, &dim)); 7906c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetLocalSize(sw, &Np)); 7919566063dSJacob Faibussowitsch PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 7926c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetCoordinateFunction(sw, &coordFunc)); 7936c5a79ebSMatthew G. Knepley 7946c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, &bs, &dtype, (void **)&x)); 7956c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&weight)); 7966c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 7976c5a79ebSMatthew G. Knepley if (coordFunc) { 7986c5a79ebSMatthew G. Knepley PetscCall(DMGetApplicationContext(sw, &ctx)); 7996c5a79ebSMatthew G. Knepley for (p = 0; p < Np; ++p) { 8006c5a79ebSMatthew G. Knepley PetscScalar X[3]; 8016c5a79ebSMatthew G. Knepley 8026c5a79ebSMatthew G. Knepley PetscCall((*coordFunc)(dim, 0., NULL, p, X, ctx)); 8036c5a79ebSMatthew G. Knepley for (d = 0; d < dim; ++d) x[p * dim + d] = PetscRealPart(X[d]); 8046c5a79ebSMatthew G. Knepley weight[p] = 1.0; 8056c5a79ebSMatthew G. Knepley species[p] = p % Ns; 8066c5a79ebSMatthew G. Knepley } 8076c5a79ebSMatthew G. Knepley } else { 8086c5a79ebSMatthew G. Knepley DM dm; 8096c5a79ebSMatthew G. Knepley PetscRandom rnd; 8106c5a79ebSMatthew G. Knepley PetscReal xi0[3]; 8116c5a79ebSMatthew G. Knepley PetscInt cStart, cEnd, c; 8126c5a79ebSMatthew G. Knepley 8139566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sw, &dm)); 8149566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 815ea7032a0SMatthew G. Knepley PetscCall(DMGetApplicationContext(sw, &ctx)); 81635b38c60SMatthew G. Knepley 81735b38c60SMatthew G. Knepley /* Set particle position randomly in cell, set weights to 1 */ 8189566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd)); 8199566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0)); 8209566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rnd)); 8219566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetAccess(sw)); 82235b38c60SMatthew G. Knepley for (d = 0; d < dim; ++d) xi0[d] = -1.0; 82335b38c60SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 82435b38c60SMatthew G. Knepley PetscReal v0[3], J[9], invJ[9], detJ; 82535b38c60SMatthew G. Knepley PetscInt *pidx, Npc, q; 82635b38c60SMatthew G. Knepley 8279566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx)); 8289566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); 82935b38c60SMatthew G. Knepley for (q = 0; q < Npc; ++q) { 83035b38c60SMatthew G. Knepley const PetscInt p = pidx[q]; 83135b38c60SMatthew G. Knepley PetscReal xref[3]; 83235b38c60SMatthew G. Knepley 8339566063dSJacob Faibussowitsch for (d = 0; d < dim; ++d) PetscCall(PetscRandomGetValueReal(rnd, &xref[d])); 83435b38c60SMatthew G. Knepley CoordinatesRefToReal(dim, dim, xi0, v0, J, xref, &x[p * dim]); 83535b38c60SMatthew G. Knepley 836ea7032a0SMatthew G. Knepley weight[p] = 1.0 / Np; 8376c5a79ebSMatthew G. Knepley species[p] = p % Ns; 83835b38c60SMatthew G. Knepley } 8399566063dSJacob Faibussowitsch PetscCall(PetscFree(pidx)); 84035b38c60SMatthew G. Knepley } 8419566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rnd)); 8429566063dSJacob Faibussowitsch PetscCall(DMSwarmSortRestoreAccess(sw)); 8436c5a79ebSMatthew G. Knepley } 8449566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x)); 8456c5a79ebSMatthew G. Knepley PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight)); 8469566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 8476c5a79ebSMatthew G. Knepley 8489566063dSJacob Faibussowitsch PetscCall(DMSwarmMigrate(sw, removePoints)); 8499566063dSJacob Faibussowitsch PetscCall(DMLocalizeCoordinates(sw)); 8503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 85135b38c60SMatthew G. Knepley } 85235b38c60SMatthew G. Knepley 85335b38c60SMatthew G. Knepley /*@C 85435b38c60SMatthew G. Knepley DMSwarmInitializeVelocities - Set the initial velocities of particles using a distribution. 85535b38c60SMatthew G. Knepley 85620f4b53cSBarry Smith Collective 85735b38c60SMatthew G. Knepley 85835b38c60SMatthew G. Knepley Input Parameters: 85920f4b53cSBarry Smith + sw - The `DMSWARM` object 86035b38c60SMatthew G. Knepley . sampler - A function which uniformly samples the velocity PDF 86135b38c60SMatthew G. Knepley - v0 - The velocity scale for nondimensionalization for each species 86235b38c60SMatthew G. Knepley 86335b38c60SMatthew G. Knepley Level: advanced 86435b38c60SMatthew G. Knepley 86520f4b53cSBarry Smith Note: 86620f4b53cSBarry 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. 86720f4b53cSBarry Smith 86820f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSize()`, `DMSwarmInitializeCoordinates()`, `DMSwarmInitializeVelocitiesFromOptions()` 86935b38c60SMatthew G. Knepley @*/ 870d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmInitializeVelocities(DM sw, PetscProbFunc sampler, const PetscReal v0[]) 871d71ae5a4SJacob Faibussowitsch { 8728434afd1SBarry Smith PetscSimplePointFn *velFunc; 87335b38c60SMatthew G. Knepley PetscReal *v; 87435b38c60SMatthew G. Knepley PetscInt *species; 8756c5a79ebSMatthew G. Knepley void *ctx; 87635b38c60SMatthew G. Knepley PetscInt dim, Np, p; 87735b38c60SMatthew G. Knepley 87835b38c60SMatthew G. Knepley PetscFunctionBegin; 8796c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetVelocityFunction(sw, &velFunc)); 88035b38c60SMatthew G. Knepley 8819566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 8829566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(sw, &Np)); 8839566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&v)); 8849566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 8856c5a79ebSMatthew G. Knepley if (v0[0] == 0.) { 8866c5a79ebSMatthew G. Knepley PetscCall(PetscArrayzero(v, Np * dim)); 8876c5a79ebSMatthew G. Knepley } else if (velFunc) { 8886c5a79ebSMatthew G. Knepley PetscCall(DMGetApplicationContext(sw, &ctx)); 8896c5a79ebSMatthew G. Knepley for (p = 0; p < Np; ++p) { 8906c5a79ebSMatthew G. Knepley PetscInt s = species[p], d; 8916c5a79ebSMatthew G. Knepley PetscScalar vel[3]; 8926c5a79ebSMatthew G. Knepley 8936c5a79ebSMatthew G. Knepley PetscCall((*velFunc)(dim, 0., NULL, p, vel, ctx)); 8946c5a79ebSMatthew G. Knepley for (d = 0; d < dim; ++d) v[p * dim + d] = (v0[s] / v0[0]) * PetscRealPart(vel[d]); 8956c5a79ebSMatthew G. Knepley } 8966c5a79ebSMatthew G. Knepley } else { 8976c5a79ebSMatthew G. Knepley PetscRandom rnd; 8986c5a79ebSMatthew G. Knepley 8996c5a79ebSMatthew G. Knepley PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rnd)); 9006c5a79ebSMatthew G. Knepley PetscCall(PetscRandomSetInterval(rnd, 0, 1.)); 9016c5a79ebSMatthew G. Knepley PetscCall(PetscRandomSetFromOptions(rnd)); 9026c5a79ebSMatthew G. Knepley 90335b38c60SMatthew G. Knepley for (p = 0; p < Np; ++p) { 90435b38c60SMatthew G. Knepley PetscInt s = species[p], d; 90535b38c60SMatthew G. Knepley PetscReal a[3], vel[3]; 90635b38c60SMatthew G. Knepley 9079566063dSJacob Faibussowitsch for (d = 0; d < dim; ++d) PetscCall(PetscRandomGetValueReal(rnd, &a[d])); 9089566063dSJacob Faibussowitsch PetscCall(sampler(a, NULL, vel)); 909ad540459SPierre Jolivet for (d = 0; d < dim; ++d) v[p * dim + d] = (v0[s] / v0[0]) * vel[d]; 91035b38c60SMatthew G. Knepley } 9116c5a79ebSMatthew G. Knepley PetscCall(PetscRandomDestroy(&rnd)); 9126c5a79ebSMatthew G. Knepley } 9139566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&v)); 9149566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 9153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 91635b38c60SMatthew G. Knepley } 91735b38c60SMatthew G. Knepley 91835b38c60SMatthew G. Knepley /*@ 91935b38c60SMatthew G. Knepley DMSwarmInitializeVelocitiesFromOptions - Set the initial velocities of particles using a distribution determined from options. 92035b38c60SMatthew G. Knepley 92120f4b53cSBarry Smith Collective 92235b38c60SMatthew G. Knepley 92335b38c60SMatthew G. Knepley Input Parameters: 92420f4b53cSBarry Smith + sw - The `DMSWARM` object 92535b38c60SMatthew G. Knepley - v0 - The velocity scale for nondimensionalization for each species 92635b38c60SMatthew G. Knepley 92735b38c60SMatthew G. Knepley Level: advanced 92835b38c60SMatthew G. Knepley 92920f4b53cSBarry Smith .seealso: `DMSWARM`, `DMSwarmComputeLocalSize()`, `DMSwarmInitializeCoordinates()`, `DMSwarmInitializeVelocities()` 93035b38c60SMatthew G. Knepley @*/ 931d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmInitializeVelocitiesFromOptions(DM sw, const PetscReal v0[]) 932d71ae5a4SJacob Faibussowitsch { 93335b38c60SMatthew G. Knepley PetscProbFunc sampler; 93435b38c60SMatthew G. Knepley PetscInt dim; 93535b38c60SMatthew G. Knepley const char *prefix; 9366c5a79ebSMatthew G. Knepley char funcname[PETSC_MAX_PATH_LEN]; 9376c5a79ebSMatthew G. Knepley PetscBool flg; 93835b38c60SMatthew G. Knepley 93935b38c60SMatthew G. Knepley PetscFunctionBegin; 940d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)sw), "", "DMSwarm Options", "DMSWARM"); 9416c5a79ebSMatthew G. Knepley PetscCall(PetscOptionsString("-dm_swarm_velocity_function", "Function to determine particle velocities", "DMSwarmSetVelocityFunction", funcname, funcname, sizeof(funcname), &flg)); 942d0609cedSBarry Smith PetscOptionsEnd(); 9436c5a79ebSMatthew G. Knepley if (flg) { 9448434afd1SBarry Smith PetscSimplePointFn *velFunc; 9456c5a79ebSMatthew G. Knepley 9466c5a79ebSMatthew G. Knepley PetscCall(PetscDLSym(NULL, funcname, (void **)&velFunc)); 9476c5a79ebSMatthew G. Knepley PetscCheck(velFunc, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname); 9486c5a79ebSMatthew G. Knepley PetscCall(DMSwarmSetVelocityFunction(sw, velFunc)); 9496c5a79ebSMatthew G. Knepley } 9509566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 9519566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sw, &prefix)); 9529566063dSJacob Faibussowitsch PetscCall(PetscProbCreateFromOptions(dim, prefix, "-dm_swarm_velocity_density", NULL, NULL, &sampler)); 9539566063dSJacob Faibussowitsch PetscCall(DMSwarmInitializeVelocities(sw, sampler, v0)); 9543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 95535b38c60SMatthew G. Knepley } 956