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 /* 125627991aSBarry Smith Error checking to ensure the swarm type is correct and that a cell DM has been set 130e2ec84fSDave May */ 140e2ec84fSDave May #define DMSWARMPICVALID(dm) \ 15f7d195e4SLawrence Mitchell do { \ 160e2ec84fSDave May DM_Swarm *_swarm = (DM_Swarm *)(dm)->data; \ 175f80ce2aSJacob Faibussowitsch PetscCheck(_swarm->swarm_type == DMSWARM_PIC, PetscObjectComm((PetscObject)(dm)), PETSC_ERR_SUP, "Valid only for DMSwarm-PIC. You must call DMSwarmSetType(dm,DMSWARM_PIC)"); \ 1828b400f6SJacob Faibussowitsch PetscCheck(_swarm->dmcell, PetscObjectComm((PetscObject)(dm)), PETSC_ERR_SUP, "Valid only for DMSwarmPIC if the cell DM is set. You must call DMSwarmSetCellDM(dm,celldm)"); \ 19f7d195e4SLawrence Mitchell } while (0) 200e2ec84fSDave May 210e2ec84fSDave May /* Coordinate insertition/addition API */ 220e2ec84fSDave May /*@C 230e2ec84fSDave May DMSwarmSetPointsUniformCoordinates - Set point coordinates in a DMSwarm on a regular (ijk) grid 240e2ec84fSDave May 25d083f849SBarry Smith Collective on dm 260e2ec84fSDave May 270e2ec84fSDave May Input parameters: 280e2ec84fSDave May + dm - the DMSwarm 290e2ec84fSDave May . min - minimum coordinate values in the x, y, z directions (array of length dim) 300e2ec84fSDave May . max - maximum coordinate values in the x, y, z directions (array of length dim) 310e2ec84fSDave May . npoints - number of points in each spatial direction (array of length dim) 320e2ec84fSDave May - mode - indicates whether to append points to the swarm (ADD_VALUES), or over-ride existing points (INSERT_VALUES) 330e2ec84fSDave May 340e2ec84fSDave May Level: beginner 350e2ec84fSDave May 360e2ec84fSDave May Notes: 370e2ec84fSDave May When using mode = INSERT_VALUES, this method will reset the number of particles in the DMSwarm 380e2ec84fSDave May to be npoints[0]*npoints[1] (2D) or npoints[0]*npoints[1]*npoints[2] (3D). When using mode = ADD_VALUES, 390e2ec84fSDave May new points will be appended to any already existing in the DMSwarm 400e2ec84fSDave May 41db781477SPatrick Sanan .seealso: `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType` 420e2ec84fSDave May @*/ 43*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSwarmSetPointsUniformCoordinates(DM dm, PetscReal min[], PetscReal max[], PetscInt npoints[], InsertMode mode) 44*d71ae5a4SJacob Faibussowitsch { 450e2ec84fSDave May PetscReal gmin[] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL}; 460e2ec84fSDave May PetscReal gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 470e2ec84fSDave May PetscInt i, j, k, N, bs, b, n_estimate, n_curr, n_new_est, p, n_found; 480e2ec84fSDave May Vec coorlocal; 490e2ec84fSDave May const PetscScalar *_coor; 500e2ec84fSDave May DM celldm; 510e2ec84fSDave May PetscReal dx[3]; 522e3d3761SDave May PetscInt _npoints[] = {0, 0, 1}; 530e2ec84fSDave May Vec pos; 540e2ec84fSDave May PetscScalar *_pos; 550e2ec84fSDave May PetscReal *swarm_coor; 560e2ec84fSDave May PetscInt *swarm_cellid; 570e2ec84fSDave May PetscSF sfcell = NULL; 580e2ec84fSDave May const PetscSFNode *LA_sfcell; 590e2ec84fSDave May 600e2ec84fSDave May PetscFunctionBegin; 610e2ec84fSDave May DMSWARMPICVALID(dm); 629566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 639566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(celldm, &coorlocal)); 649566063dSJacob Faibussowitsch PetscCall(VecGetSize(coorlocal, &N)); 659566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(coorlocal, &bs)); 660e2ec84fSDave May N = N / bs; 679566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coorlocal, &_coor)); 680e2ec84fSDave May for (i = 0; i < N; i++) { 690e2ec84fSDave May for (b = 0; b < bs; b++) { 70a5f152d1SDave May gmin[b] = PetscMin(gmin[b], PetscRealPart(_coor[bs * i + b])); 71a5f152d1SDave May gmax[b] = PetscMax(gmax[b], PetscRealPart(_coor[bs * i + b])); 720e2ec84fSDave May } 730e2ec84fSDave May } 749566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coorlocal, &_coor)); 750e2ec84fSDave May 760e2ec84fSDave May for (b = 0; b < bs; b++) { 7771844bb8SDave May if (npoints[b] > 1) { 780e2ec84fSDave May dx[b] = (max[b] - min[b]) / ((PetscReal)(npoints[b] - 1)); 7971844bb8SDave May } else { 8071844bb8SDave May dx[b] = 0.0; 8171844bb8SDave May } 822e3d3761SDave May _npoints[b] = npoints[b]; 830e2ec84fSDave May } 840e2ec84fSDave May 850e2ec84fSDave May /* determine number of points living in the bounding box */ 860e2ec84fSDave May n_estimate = 0; 872e3d3761SDave May for (k = 0; k < _npoints[2]; k++) { 882e3d3761SDave May for (j = 0; j < _npoints[1]; j++) { 892e3d3761SDave May for (i = 0; i < _npoints[0]; i++) { 900e2ec84fSDave May PetscReal xp[] = {0.0, 0.0, 0.0}; 910e2ec84fSDave May PetscInt ijk[3]; 920e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 930e2ec84fSDave May 940e2ec84fSDave May ijk[0] = i; 950e2ec84fSDave May ijk[1] = j; 960e2ec84fSDave May ijk[2] = k; 97ad540459SPierre Jolivet for (b = 0; b < bs; b++) xp[b] = min[b] + ijk[b] * dx[b]; 980e2ec84fSDave May for (b = 0; b < bs; b++) { 99ad540459SPierre Jolivet if (xp[b] < gmin[b]) point_inside = PETSC_FALSE; 100ad540459SPierre Jolivet if (xp[b] > gmax[b]) point_inside = PETSC_FALSE; 1010e2ec84fSDave May } 102ad540459SPierre Jolivet if (point_inside) n_estimate++; 1030e2ec84fSDave May } 1040e2ec84fSDave May } 1050e2ec84fSDave May } 1060e2ec84fSDave May 1070e2ec84fSDave May /* create candidate list */ 1089566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &pos)); 1099566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pos, bs * n_estimate, PETSC_DECIDE)); 1109566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(pos, bs)); 1119566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pos)); 1129566063dSJacob Faibussowitsch PetscCall(VecGetArray(pos, &_pos)); 1130e2ec84fSDave May 1140e2ec84fSDave May n_estimate = 0; 1152e3d3761SDave May for (k = 0; k < _npoints[2]; k++) { 1162e3d3761SDave May for (j = 0; j < _npoints[1]; j++) { 1172e3d3761SDave May for (i = 0; i < _npoints[0]; i++) { 1180e2ec84fSDave May PetscReal xp[] = {0.0, 0.0, 0.0}; 1190e2ec84fSDave May PetscInt ijk[3]; 1200e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 1210e2ec84fSDave May 1220e2ec84fSDave May ijk[0] = i; 1230e2ec84fSDave May ijk[1] = j; 1240e2ec84fSDave May ijk[2] = k; 125ad540459SPierre Jolivet for (b = 0; b < bs; b++) xp[b] = min[b] + ijk[b] * dx[b]; 1260e2ec84fSDave May for (b = 0; b < bs; b++) { 127ad540459SPierre Jolivet if (xp[b] < gmin[b]) point_inside = PETSC_FALSE; 128ad540459SPierre Jolivet if (xp[b] > gmax[b]) point_inside = PETSC_FALSE; 1290e2ec84fSDave May } 1300e2ec84fSDave May if (point_inside) { 131ad540459SPierre Jolivet for (b = 0; b < bs; b++) _pos[bs * n_estimate + b] = xp[b]; 1320e2ec84fSDave May n_estimate++; 1330e2ec84fSDave May } 1340e2ec84fSDave May } 1350e2ec84fSDave May } 1360e2ec84fSDave May } 1379566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pos, &_pos)); 1380e2ec84fSDave May 1390e2ec84fSDave May /* locate points */ 1409566063dSJacob Faibussowitsch PetscCall(DMLocatePoints(celldm, pos, DM_POINTLOCATION_NONE, &sfcell)); 1419566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 1420e2ec84fSDave May n_found = 0; 1430e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 144ad540459SPierre Jolivet if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) n_found++; 1450e2ec84fSDave May } 1460e2ec84fSDave May 1470e2ec84fSDave May /* adjust size */ 1480e2ec84fSDave May if (mode == ADD_VALUES) { 1499566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dm, &n_curr)); 1500e2ec84fSDave May n_new_est = n_curr + n_found; 1519566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 1520e2ec84fSDave May } 1530e2ec84fSDave May if (mode == INSERT_VALUES) { 1540e2ec84fSDave May n_curr = 0; 1550e2ec84fSDave May n_new_est = n_found; 1569566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 1570e2ec84fSDave May } 1580e2ec84fSDave May 1590e2ec84fSDave May /* initialize new coords, cell owners, pid */ 1609566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pos, &_coor)); 1619566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 1629566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 1630e2ec84fSDave May n_found = 0; 1640e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 1650e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 166ad540459SPierre Jolivet for (b = 0; b < bs; b++) swarm_coor[bs * (n_curr + n_found) + b] = PetscRealPart(_coor[bs * p + b]); 1670e2ec84fSDave May swarm_cellid[n_curr + n_found] = LA_sfcell[p].index; 1680e2ec84fSDave May n_found++; 1690e2ec84fSDave May } 1700e2ec84fSDave May } 1719566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 1729566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 1739566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pos, &_coor)); 1740e2ec84fSDave May 1759566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfcell)); 1769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pos)); 1770e2ec84fSDave May PetscFunctionReturn(0); 1780e2ec84fSDave May } 1790e2ec84fSDave May 1800e2ec84fSDave May /*@C 1810e2ec84fSDave May DMSwarmSetPointCoordinates - Set point coordinates in a DMSwarm from a user defined list 1820e2ec84fSDave May 183d083f849SBarry Smith Collective on dm 1840e2ec84fSDave May 1850e2ec84fSDave May Input parameters: 1860e2ec84fSDave May + dm - the DMSwarm 1870e2ec84fSDave May . npoints - the number of points to insert 1880e2ec84fSDave May . coor - the coordinate values 1890e2ec84fSDave May . redundant - if set to PETSC_TRUE, it is assumed that npoints and coor[] are only valid on rank 0 and should be broadcast to other ranks 1900e2ec84fSDave May - mode - indicates whether to append points to the swarm (ADD_VALUES), or over-ride existing points (INSERT_VALUES) 1910e2ec84fSDave May 1920e2ec84fSDave May Level: beginner 1930e2ec84fSDave May 1940e2ec84fSDave May Notes: 1950e2ec84fSDave May If the user has specified redundant = PETSC_FALSE, the cell DM will attempt to locate the coordinates provided by coor[] within 1960e2ec84fSDave May its sub-domain. If they any values within coor[] are not located in the sub-domain, they will be ignored and will not get 1970e2ec84fSDave May added to the DMSwarm. 1980e2ec84fSDave May 199db781477SPatrick Sanan .seealso: `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType`, `DMSwarmSetPointsUniformCoordinates()` 2000e2ec84fSDave May @*/ 201*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinates(DM dm, PetscInt npoints, PetscReal coor[], PetscBool redundant, InsertMode mode) 202*d71ae5a4SJacob Faibussowitsch { 2030e2ec84fSDave May PetscReal gmin[] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL}; 2040e2ec84fSDave May PetscReal gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 2050e2ec84fSDave May PetscInt i, N, bs, b, n_estimate, n_curr, n_new_est, p, n_found; 2060e2ec84fSDave May Vec coorlocal; 2070e2ec84fSDave May const PetscScalar *_coor; 2080e2ec84fSDave May DM celldm; 2090e2ec84fSDave May Vec pos; 2100e2ec84fSDave May PetscScalar *_pos; 2110e2ec84fSDave May PetscReal *swarm_coor; 2120e2ec84fSDave May PetscInt *swarm_cellid; 2130e2ec84fSDave May PetscSF sfcell = NULL; 2140e2ec84fSDave May const PetscSFNode *LA_sfcell; 2150e2ec84fSDave May PetscReal *my_coor; 2160e2ec84fSDave May PetscInt my_npoints; 2170e2ec84fSDave May PetscMPIInt rank; 2180e2ec84fSDave May MPI_Comm comm; 2190e2ec84fSDave May 2200e2ec84fSDave May PetscFunctionBegin; 2210e2ec84fSDave May DMSWARMPICVALID(dm); 2229566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 2239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2240e2ec84fSDave May 2259566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 2269566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(celldm, &coorlocal)); 2279566063dSJacob Faibussowitsch PetscCall(VecGetSize(coorlocal, &N)); 2289566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(coorlocal, &bs)); 2290e2ec84fSDave May N = N / bs; 2309566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coorlocal, &_coor)); 2310e2ec84fSDave May for (i = 0; i < N; i++) { 2320e2ec84fSDave May for (b = 0; b < bs; b++) { 233a5f152d1SDave May gmin[b] = PetscMin(gmin[b], PetscRealPart(_coor[bs * i + b])); 234a5f152d1SDave May gmax[b] = PetscMax(gmax[b], PetscRealPart(_coor[bs * i + b])); 2350e2ec84fSDave May } 2360e2ec84fSDave May } 2379566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coorlocal, &_coor)); 2380e2ec84fSDave May 2390e2ec84fSDave May /* broadcast points from rank 0 if requested */ 2400e2ec84fSDave May if (redundant) { 2410e2ec84fSDave May my_npoints = npoints; 2429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&my_npoints, 1, MPIU_INT, 0, comm)); 2430e2ec84fSDave May 2440e2ec84fSDave May if (rank > 0) { /* allocate space */ 2459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * my_npoints, &my_coor)); 2460e2ec84fSDave May } else { 2470e2ec84fSDave May my_coor = coor; 2480e2ec84fSDave May } 2499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(my_coor, bs * my_npoints, MPIU_REAL, 0, comm)); 2500e2ec84fSDave May } else { 2510e2ec84fSDave May my_npoints = npoints; 2520e2ec84fSDave May my_coor = coor; 2530e2ec84fSDave May } 2540e2ec84fSDave May 2550e2ec84fSDave May /* determine the number of points living in the bounding box */ 2560e2ec84fSDave May n_estimate = 0; 2570e2ec84fSDave May for (i = 0; i < my_npoints; i++) { 2580e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 2590e2ec84fSDave May 2600e2ec84fSDave May for (b = 0; b < bs; b++) { 261ad540459SPierre Jolivet if (my_coor[bs * i + b] < gmin[b]) point_inside = PETSC_FALSE; 262ad540459SPierre Jolivet if (my_coor[bs * i + b] > gmax[b]) point_inside = PETSC_FALSE; 2630e2ec84fSDave May } 264ad540459SPierre Jolivet if (point_inside) n_estimate++; 2650e2ec84fSDave May } 2660e2ec84fSDave May 2670e2ec84fSDave May /* create candidate list */ 2689566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &pos)); 2699566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pos, bs * n_estimate, PETSC_DECIDE)); 2709566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(pos, bs)); 2719566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pos)); 2729566063dSJacob Faibussowitsch PetscCall(VecGetArray(pos, &_pos)); 2730e2ec84fSDave May 2740e2ec84fSDave May n_estimate = 0; 2750e2ec84fSDave May for (i = 0; i < my_npoints; i++) { 2760e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 2770e2ec84fSDave May 2780e2ec84fSDave May for (b = 0; b < bs; b++) { 279ad540459SPierre Jolivet if (my_coor[bs * i + b] < gmin[b]) point_inside = PETSC_FALSE; 280ad540459SPierre Jolivet if (my_coor[bs * i + b] > gmax[b]) point_inside = PETSC_FALSE; 2810e2ec84fSDave May } 2820e2ec84fSDave May if (point_inside) { 283ad540459SPierre Jolivet for (b = 0; b < bs; b++) _pos[bs * n_estimate + b] = my_coor[bs * i + b]; 2840e2ec84fSDave May n_estimate++; 2850e2ec84fSDave May } 2860e2ec84fSDave May } 2879566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pos, &_pos)); 2880e2ec84fSDave May 2890e2ec84fSDave May /* locate points */ 2909566063dSJacob Faibussowitsch PetscCall(DMLocatePoints(celldm, pos, DM_POINTLOCATION_NONE, &sfcell)); 2910e2ec84fSDave May 2929566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 2930e2ec84fSDave May n_found = 0; 2940e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 295ad540459SPierre Jolivet if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) n_found++; 2960e2ec84fSDave May } 2970e2ec84fSDave May 2980e2ec84fSDave May /* adjust size */ 2990e2ec84fSDave May if (mode == ADD_VALUES) { 3009566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dm, &n_curr)); 3010e2ec84fSDave May n_new_est = n_curr + n_found; 3029566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 3030e2ec84fSDave May } 3040e2ec84fSDave May if (mode == INSERT_VALUES) { 3050e2ec84fSDave May n_curr = 0; 3060e2ec84fSDave May n_new_est = n_found; 3079566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 3080e2ec84fSDave May } 3090e2ec84fSDave May 3100e2ec84fSDave May /* initialize new coords, cell owners, pid */ 3119566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pos, &_coor)); 3129566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 3139566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 3140e2ec84fSDave May n_found = 0; 3150e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 3160e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 317ad540459SPierre Jolivet for (b = 0; b < bs; b++) swarm_coor[bs * (n_curr + n_found) + b] = PetscRealPart(_coor[bs * p + b]); 3180e2ec84fSDave May swarm_cellid[n_curr + n_found] = LA_sfcell[p].index; 3190e2ec84fSDave May n_found++; 3200e2ec84fSDave May } 3210e2ec84fSDave May } 3229566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 3239566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 3249566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pos, &_coor)); 3250e2ec84fSDave May 3260e2ec84fSDave May if (redundant) { 32748a46eb9SPierre Jolivet if (rank > 0) PetscCall(PetscFree(my_coor)); 3280e2ec84fSDave May } 3299566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfcell)); 3309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pos)); 3310e2ec84fSDave May PetscFunctionReturn(0); 3320e2ec84fSDave May } 3330e2ec84fSDave May 3340e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM, DM, DMSwarmPICLayoutType, PetscInt); 3350e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM, DM, DMSwarmPICLayoutType, PetscInt); 3360e2ec84fSDave May 3370e2ec84fSDave May /*@C 3380e2ec84fSDave May DMSwarmInsertPointsUsingCellDM - Insert point coordinates within each cell 3390e2ec84fSDave May 3400e2ec84fSDave May Not collective 3410e2ec84fSDave May 3420e2ec84fSDave May Input parameters: 3430e2ec84fSDave May + dm - the DMSwarm 3440e2ec84fSDave May . layout_type - method used to fill each cell with the cell DM 3450e2ec84fSDave May - fill_param - parameter controlling how many points per cell are added (the meaning of this parameter is dependent on the layout type) 3460e2ec84fSDave May 3470e2ec84fSDave May Level: beginner 3480e2ec84fSDave May 3490e2ec84fSDave May Notes: 350e7af74c8SDave May 351e7af74c8SDave May The insert method will reset any previous defined points within the DMSwarm. 352e7af74c8SDave May 353e7af74c8SDave May When using a DMDA both 2D and 3D are supported for all layout types provided you are using DMDA_ELEMENT_Q1. 354e7af74c8SDave May 355e7af74c8SDave May When using a DMPLEX the following case are supported: 356ea3d7275SDave May (i) DMSWARMPIC_LAYOUT_REGULAR: 2D (triangle), 357ea3d7275SDave May (ii) DMSWARMPIC_LAYOUT_GAUSS: 2D and 3D provided the cell is a tri/tet or a quad/hex, 358ea3d7275SDave May (iii) DMSWARMPIC_LAYOUT_SUBDIVISION: 2D and 3D for quad/hex and 2D tri. 3590e2ec84fSDave May 360db781477SPatrick Sanan .seealso: `DMSwarmPICLayoutType`, `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType` 3610e2ec84fSDave May @*/ 362*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSwarmInsertPointsUsingCellDM(DM dm, DMSwarmPICLayoutType layout_type, PetscInt fill_param) 363*d71ae5a4SJacob Faibussowitsch { 3640e2ec84fSDave May DM celldm; 3650e2ec84fSDave May PetscBool isDA, isPLEX; 3660e2ec84fSDave May 3670e2ec84fSDave May PetscFunctionBegin; 3680e2ec84fSDave May DMSWARMPICVALID(dm); 3699566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 3709566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA)); 3719566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX)); 3720e2ec84fSDave May if (isDA) { 3739566063dSJacob Faibussowitsch PetscCall(private_DMSwarmInsertPointsUsingCellDM_DA(dm, celldm, layout_type, fill_param)); 3740e2ec84fSDave May } else if (isPLEX) { 3759566063dSJacob Faibussowitsch PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX(dm, celldm, layout_type, fill_param)); 3760e2ec84fSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX"); 3770e2ec84fSDave May PetscFunctionReturn(0); 3780e2ec84fSDave May } 3790e2ec84fSDave May 380431382f2SDave May extern PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM, DM, PetscInt, PetscReal *); 381431382f2SDave May 382431382f2SDave May /*@C 383431382f2SDave May DMSwarmSetPointCoordinatesCellwise - Insert point coordinates (defined over the reference cell) within each cell 384431382f2SDave May 385431382f2SDave May Not collective 386431382f2SDave May 387431382f2SDave May Input parameters: 388431382f2SDave May + dm - the DMSwarm 389431382f2SDave May . celldm - the cell DM 390431382f2SDave May . npoints - the number of points to insert in each cell 391431382f2SDave May - xi - the coordinates (defined in the local coordinate system for each cell) to insert 392431382f2SDave May 393431382f2SDave May Level: beginner 394431382f2SDave May 395431382f2SDave May Notes: 396431382f2SDave May The method will reset any previous defined points within the DMSwarm. 397ea3d7275SDave May Only supported for DMPLEX. If you are using a DMDA it is recommended to either use 398e7af74c8SDave May DMSwarmInsertPointsUsingCellDM(), or extract and set the coordinates yourself the following code 399e7af74c8SDave May 400e7af74c8SDave May $ PetscReal *coor; 401e7af74c8SDave May $ DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&coor); 402e7af74c8SDave May $ // user code to define the coordinates here 403e7af74c8SDave May $ DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&coor); 404431382f2SDave May 405db781477SPatrick Sanan .seealso: `DMSwarmSetCellDM()`, `DMSwarmInsertPointsUsingCellDM()` 406431382f2SDave May @*/ 407*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinatesCellwise(DM dm, PetscInt npoints, PetscReal xi[]) 408*d71ae5a4SJacob Faibussowitsch { 409431382f2SDave May DM celldm; 410431382f2SDave May PetscBool isDA, isPLEX; 411431382f2SDave May 4120e2ec84fSDave May PetscFunctionBegin; 413431382f2SDave May DMSWARMPICVALID(dm); 4149566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 4159566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA)); 4169566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX)); 41728b400f6SJacob Faibussowitsch PetscCheck(!isDA, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMPLEX. Recommended you use DMSwarmInsertPointsUsingCellDM()"); 418f7d195e4SLawrence Mitchell if (isPLEX) { 4199566063dSJacob Faibussowitsch PetscCall(private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm, celldm, npoints, xi)); 420431382f2SDave May } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX"); 4210e2ec84fSDave May PetscFunctionReturn(0); 4220e2ec84fSDave May } 423431382f2SDave May 4240e2ec84fSDave May /* Field projection API */ 42577048351SPatrick Sanan extern PetscErrorCode private_DMSwarmProjectFields_DA(DM swarm, DM celldm, PetscInt project_type, PetscInt nfields, DMSwarmDataField dfield[], Vec vecs[]); 42677048351SPatrick Sanan extern PetscErrorCode private_DMSwarmProjectFields_PLEX(DM swarm, DM celldm, PetscInt project_type, PetscInt nfields, DMSwarmDataField dfield[], Vec vecs[]); 42794f7d2dcSDave May 42894f7d2dcSDave May /*@C 42994f7d2dcSDave May DMSwarmProjectFields - Project a set of swarm fields onto the cell DM 43094f7d2dcSDave May 431d083f849SBarry Smith Collective on dm 43294f7d2dcSDave May 43394f7d2dcSDave May Input parameters: 43494f7d2dcSDave May + dm - the DMSwarm 43594f7d2dcSDave May . nfields - the number of swarm fields to project 43694f7d2dcSDave May . fieldnames - the textual names of the swarm fields to project 43794f7d2dcSDave May . fields - an array of Vec's of length nfields 43894f7d2dcSDave May - reuse - flag indicating whether the array and contents of fields should be re-used or internally allocated 43994f7d2dcSDave May 44094f7d2dcSDave May Currently, the only available projection method consists of 44194f7d2dcSDave May phi_i = \sum_{p=0}^{np} N_i(x_p) phi_p dJ / \sum_{p=0}^{np} N_i(x_p) dJ 44294f7d2dcSDave May where phi_p is the swarm field at point p, 44394f7d2dcSDave May N_i() is the cell DM basis function at vertex i, 44494f7d2dcSDave May dJ is the determinant of the cell Jacobian and 44594f7d2dcSDave May phi_i is the projected vertex value of the field phi. 44694f7d2dcSDave May 44794f7d2dcSDave May Level: beginner 44894f7d2dcSDave May 44994f7d2dcSDave May Notes: 450e7af74c8SDave May 451e7af74c8SDave May If reuse = PETSC_FALSE, this function will allocate the array of Vec's, and each individual Vec. 45294f7d2dcSDave May The user is responsible for destroying both the array and the individual Vec objects. 453e7af74c8SDave May 454e7af74c8SDave May Only swarm fields registered with data type = PETSC_REAL can be projected onto the cell DM. 455e7af74c8SDave May 456e7af74c8SDave May Only swarm fields of block size = 1 can currently be projected. 457e7af74c8SDave May 458e7af74c8SDave May The only projection methods currently only support the DA (2D) and PLEX (triangles 2D). 45994f7d2dcSDave May 460db781477SPatrick Sanan .seealso: `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType` 46194f7d2dcSDave May @*/ 462*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM dm, PetscInt nfields, const char *fieldnames[], Vec **fields, PetscBool reuse) 463*d71ae5a4SJacob Faibussowitsch { 46494f7d2dcSDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 46577048351SPatrick Sanan DMSwarmDataField *gfield; 46694f7d2dcSDave May DM celldm; 46794f7d2dcSDave May PetscBool isDA, isPLEX; 46894f7d2dcSDave May Vec *vecs; 46994f7d2dcSDave May PetscInt f, nvecs; 47094f7d2dcSDave May PetscInt project_type = 0; 47194f7d2dcSDave May 4720e2ec84fSDave May PetscFunctionBegin; 47394f7d2dcSDave May DMSWARMPICVALID(dm); 4749566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 4759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nfields, &gfield)); 47694f7d2dcSDave May nvecs = 0; 47794f7d2dcSDave May for (f = 0; f < nfields; f++) { 4789566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldnames[f], &gfield[f])); 4795f80ce2aSJacob Faibussowitsch PetscCheck(gfield[f]->petsc_type == PETSC_REAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Projection only valid for fields using a data type = PETSC_REAL"); 4805f80ce2aSJacob Faibussowitsch PetscCheck(gfield[f]->bs == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Projection only valid for fields with block size = 1"); 48194f7d2dcSDave May nvecs += gfield[f]->bs; 48294f7d2dcSDave May } 48394f7d2dcSDave May if (!reuse) { 4849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nvecs, &vecs)); 48594f7d2dcSDave May for (f = 0; f < nvecs; f++) { 4869566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(celldm, &vecs[f])); 4879566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)vecs[f], gfield[f]->name)); 48894f7d2dcSDave May } 48994f7d2dcSDave May } else { 49094f7d2dcSDave May vecs = *fields; 49194f7d2dcSDave May } 49294f7d2dcSDave May 4939566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA)); 4949566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX)); 49594f7d2dcSDave May if (isDA) { 4969566063dSJacob Faibussowitsch PetscCall(private_DMSwarmProjectFields_DA(dm, celldm, project_type, nfields, gfield, vecs)); 49794f7d2dcSDave May } else if (isPLEX) { 4989566063dSJacob Faibussowitsch PetscCall(private_DMSwarmProjectFields_PLEX(dm, celldm, project_type, nfields, gfield, vecs)); 49994f7d2dcSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX"); 50094f7d2dcSDave May 5019566063dSJacob Faibussowitsch PetscCall(PetscFree(gfield)); 5025f80ce2aSJacob Faibussowitsch if (!reuse) *fields = vecs; 5030e2ec84fSDave May PetscFunctionReturn(0); 5040e2ec84fSDave May } 5050e2ec84fSDave May 5060e2ec84fSDave May /*@C 507b799feefSDave May DMSwarmCreatePointPerCellCount - Count the number of points within all cells in the cell DM 5080e2ec84fSDave May 5090e2ec84fSDave May Not collective 5100e2ec84fSDave May 5110e2ec84fSDave May Input parameter: 5120e2ec84fSDave May . dm - the DMSwarm 5130e2ec84fSDave May 5140e2ec84fSDave May Output parameters: 5150e2ec84fSDave May + ncells - the number of cells in the cell DM (optional argument, pass NULL to ignore) 516b799feefSDave May - count - array of length ncells containing the number of points per cell 5170e2ec84fSDave May 5180e2ec84fSDave May Level: beginner 5190e2ec84fSDave May 5200e2ec84fSDave May Notes: 5210e2ec84fSDave May The array count is allocated internally and must be free'd by the user. 5220e2ec84fSDave May 523db781477SPatrick Sanan .seealso: `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType` 5240e2ec84fSDave May @*/ 525*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSwarmCreatePointPerCellCount(DM dm, PetscInt *ncells, PetscInt **count) 526*d71ae5a4SJacob Faibussowitsch { 527b799feefSDave May PetscBool isvalid; 528b799feefSDave May PetscInt nel; 529b799feefSDave May PetscInt *sum; 530b799feefSDave May 5310e2ec84fSDave May PetscFunctionBegin; 5329566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetIsValid(dm, &isvalid)); 533b799feefSDave May nel = 0; 534b799feefSDave May if (isvalid) { 535b799feefSDave May PetscInt e; 536b799feefSDave May 5379566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetSizes(dm, &nel, NULL)); 538b799feefSDave May 5399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nel, &sum)); 54048a46eb9SPierre Jolivet for (e = 0; e < nel; e++) PetscCall(DMSwarmSortGetNumberOfPointsPerCell(dm, e, &sum[e])); 541b799feefSDave May } else { 542b799feefSDave May DM celldm; 543b799feefSDave May PetscBool isda, isplex, isshell; 544b799feefSDave May PetscInt p, npoints; 545b799feefSDave May PetscInt *swarm_cellid; 546b799feefSDave May 547b799feefSDave May /* get the number of cells */ 5489566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 5499566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isda)); 5509566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isplex)); 5519566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMSHELL, &isshell)); 552b799feefSDave May if (isda) { 553b799feefSDave May PetscInt _nel, _npe; 554b799feefSDave May const PetscInt *_element; 555b799feefSDave May 5569566063dSJacob Faibussowitsch PetscCall(DMDAGetElements(celldm, &_nel, &_npe, &_element)); 557b799feefSDave May nel = _nel; 5589566063dSJacob Faibussowitsch PetscCall(DMDARestoreElements(celldm, &_nel, &_npe, &_element)); 559b799feefSDave May } else if (isplex) { 560b799feefSDave May PetscInt ps, pe; 561b799feefSDave May 5629566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(celldm, 0, &ps, &pe)); 563b799feefSDave May nel = pe - ps; 564b799feefSDave May } else if (isshell) { 565b799feefSDave May PetscErrorCode (*method_DMShellGetNumberOfCells)(DM, PetscInt *); 566b799feefSDave May 5679566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)celldm, "DMGetNumberOfCells_C", &method_DMShellGetNumberOfCells)); 568b799feefSDave May if (method_DMShellGetNumberOfCells) { 5699566063dSJacob Faibussowitsch PetscCall(method_DMShellGetNumberOfCells(celldm, &nel)); 5709371c9d4SSatish Balay } else 5719371c9d4SSatish 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);"); 572b799feefSDave 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"); 573b799feefSDave May 5749566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nel, &sum)); 5759566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(sum, nel)); 5769566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dm, &npoints)); 5779566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 578b799feefSDave May for (p = 0; p < npoints; p++) { 579ad540459SPierre Jolivet if (swarm_cellid[p] != DMLOCATEPOINT_POINT_NOT_FOUND) sum[swarm_cellid[p]]++; 580b799feefSDave May } 5819566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 582b799feefSDave May } 583ad540459SPierre Jolivet if (ncells) *ncells = nel; 584b799feefSDave May *count = sum; 5850e2ec84fSDave May PetscFunctionReturn(0); 5860e2ec84fSDave May } 58735b38c60SMatthew G. Knepley 58835b38c60SMatthew G. Knepley /*@ 58935b38c60SMatthew G. Knepley DMSwarmGetNumSpecies - Get the number of particle species 59035b38c60SMatthew G. Knepley 59135b38c60SMatthew G. Knepley Not collective 59235b38c60SMatthew G. Knepley 59335b38c60SMatthew G. Knepley Input parameter: 59435b38c60SMatthew G. Knepley . dm - the DMSwarm 59535b38c60SMatthew G. Knepley 59635b38c60SMatthew G. Knepley Output parameters: 59735b38c60SMatthew G. Knepley . Ns - the number of species 59835b38c60SMatthew G. Knepley 59935b38c60SMatthew G. Knepley Level: intermediate 60035b38c60SMatthew G. Knepley 601db781477SPatrick Sanan .seealso: `DMSwarmSetNumSpecies()`, `DMSwarmSetType()`, `DMSwarmType` 60235b38c60SMatthew G. Knepley @*/ 603*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetNumSpecies(DM sw, PetscInt *Ns) 604*d71ae5a4SJacob Faibussowitsch { 60535b38c60SMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 60635b38c60SMatthew G. Knepley 60735b38c60SMatthew G. Knepley PetscFunctionBegin; 60835b38c60SMatthew G. Knepley *Ns = swarm->Ns; 60935b38c60SMatthew G. Knepley PetscFunctionReturn(0); 61035b38c60SMatthew G. Knepley } 61135b38c60SMatthew G. Knepley 61235b38c60SMatthew G. Knepley /*@ 61335b38c60SMatthew G. Knepley DMSwarmSetNumSpecies - Set the number of particle species 61435b38c60SMatthew G. Knepley 61535b38c60SMatthew G. Knepley Not collective 61635b38c60SMatthew G. Knepley 61735b38c60SMatthew G. Knepley Input parameter: 61835b38c60SMatthew G. Knepley + dm - the DMSwarm 61935b38c60SMatthew G. Knepley - Ns - the number of species 62035b38c60SMatthew G. Knepley 62135b38c60SMatthew G. Knepley Level: intermediate 62235b38c60SMatthew G. Knepley 623db781477SPatrick Sanan .seealso: `DMSwarmGetNumSpecies()`, `DMSwarmSetType()`, `DMSwarmType` 62435b38c60SMatthew G. Knepley @*/ 625*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetNumSpecies(DM sw, PetscInt Ns) 626*d71ae5a4SJacob Faibussowitsch { 62735b38c60SMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 62835b38c60SMatthew G. Knepley 62935b38c60SMatthew G. Knepley PetscFunctionBegin; 63035b38c60SMatthew G. Knepley swarm->Ns = Ns; 63135b38c60SMatthew G. Knepley PetscFunctionReturn(0); 63235b38c60SMatthew G. Knepley } 63335b38c60SMatthew G. Knepley 63435b38c60SMatthew G. Knepley /*@C 6356c5a79ebSMatthew G. Knepley DMSwarmGetCoordinateFunction - Get the function setting initial particle positions, if it exists 6366c5a79ebSMatthew G. Knepley 6376c5a79ebSMatthew G. Knepley Not collective 6386c5a79ebSMatthew G. Knepley 6396c5a79ebSMatthew G. Knepley Input parameter: 6406c5a79ebSMatthew G. Knepley . dm - the DMSwarm 6416c5a79ebSMatthew G. Knepley 6426c5a79ebSMatthew G. Knepley Output Parameter: 6436c5a79ebSMatthew G. Knepley . coordFunc - the function setting initial particle positions, or NULL 6446c5a79ebSMatthew G. Knepley 6456c5a79ebSMatthew G. Knepley Level: intermediate 6466c5a79ebSMatthew G. Knepley 647db781477SPatrick Sanan .seealso: `DMSwarmSetCoordinateFunction()`, `DMSwarmGetVelocityFunction()`, `DMSwarmInitializeCoordinates()` 6486c5a79ebSMatthew G. Knepley @*/ 649*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetCoordinateFunction(DM sw, PetscSimplePointFunc *coordFunc) 650*d71ae5a4SJacob Faibussowitsch { 6516c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 6526c5a79ebSMatthew G. Knepley 6536c5a79ebSMatthew G. Knepley PetscFunctionBegin; 6546c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 6556c5a79ebSMatthew G. Knepley PetscValidPointer(coordFunc, 2); 6566c5a79ebSMatthew G. Knepley *coordFunc = swarm->coordFunc; 6576c5a79ebSMatthew G. Knepley PetscFunctionReturn(0); 6586c5a79ebSMatthew G. Knepley } 6596c5a79ebSMatthew G. Knepley 6606c5a79ebSMatthew G. Knepley /*@C 6616c5a79ebSMatthew G. Knepley DMSwarmSetCoordinateFunction - Set the function setting initial particle positions 6626c5a79ebSMatthew G. Knepley 6636c5a79ebSMatthew G. Knepley Not collective 6646c5a79ebSMatthew G. Knepley 6656c5a79ebSMatthew G. Knepley Input parameters: 6666c5a79ebSMatthew G. Knepley + dm - the DMSwarm 6676c5a79ebSMatthew G. Knepley - coordFunc - the function setting initial particle positions 6686c5a79ebSMatthew G. Knepley 6696c5a79ebSMatthew G. Knepley Level: intermediate 6706c5a79ebSMatthew G. Knepley 671db781477SPatrick Sanan .seealso: `DMSwarmGetCoordinateFunction()`, `DMSwarmSetVelocityFunction()`, `DMSwarmInitializeCoordinates()` 6726c5a79ebSMatthew G. Knepley @*/ 673*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetCoordinateFunction(DM sw, PetscSimplePointFunc coordFunc) 674*d71ae5a4SJacob Faibussowitsch { 6756c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 6766c5a79ebSMatthew G. Knepley 6776c5a79ebSMatthew G. Knepley PetscFunctionBegin; 6786c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 6796c5a79ebSMatthew G. Knepley PetscValidFunction(coordFunc, 2); 6806c5a79ebSMatthew G. Knepley swarm->coordFunc = coordFunc; 6816c5a79ebSMatthew G. Knepley PetscFunctionReturn(0); 6826c5a79ebSMatthew G. Knepley } 6836c5a79ebSMatthew G. Knepley 6846c5a79ebSMatthew G. Knepley /*@C 6856c5a79ebSMatthew G. Knepley DMSwarmGetCoordinateFunction - Get the function setting initial particle velocities, if it exists 6866c5a79ebSMatthew G. Knepley 6876c5a79ebSMatthew G. Knepley Not collective 6886c5a79ebSMatthew G. Knepley 6896c5a79ebSMatthew G. Knepley Input parameter: 6906c5a79ebSMatthew G. Knepley . dm - the DMSwarm 6916c5a79ebSMatthew G. Knepley 6926c5a79ebSMatthew G. Knepley Output Parameter: 6936c5a79ebSMatthew G. Knepley . velFunc - the function setting initial particle velocities, or NULL 6946c5a79ebSMatthew G. Knepley 6956c5a79ebSMatthew G. Knepley Level: intermediate 6966c5a79ebSMatthew G. Knepley 697db781477SPatrick Sanan .seealso: `DMSwarmSetVelocityFunction()`, `DMSwarmGetCoordinateFunction()`, `DMSwarmInitializeVelocities()` 6986c5a79ebSMatthew G. Knepley @*/ 699*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmGetVelocityFunction(DM sw, PetscSimplePointFunc *velFunc) 700*d71ae5a4SJacob Faibussowitsch { 7016c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 7026c5a79ebSMatthew G. Knepley 7036c5a79ebSMatthew G. Knepley PetscFunctionBegin; 7046c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 7056c5a79ebSMatthew G. Knepley PetscValidPointer(velFunc, 2); 7066c5a79ebSMatthew G. Knepley *velFunc = swarm->velFunc; 7076c5a79ebSMatthew G. Knepley PetscFunctionReturn(0); 7086c5a79ebSMatthew G. Knepley } 7096c5a79ebSMatthew G. Knepley 7106c5a79ebSMatthew G. Knepley /*@C 7116c5a79ebSMatthew G. Knepley DMSwarmSetVelocityFunction - Set the function setting initial particle velocities 7126c5a79ebSMatthew G. Knepley 7136c5a79ebSMatthew G. Knepley Not collective 7146c5a79ebSMatthew G. Knepley 7156c5a79ebSMatthew G. Knepley Input parameters: 7166c5a79ebSMatthew G. Knepley + dm - the DMSwarm 7176c5a79ebSMatthew G. Knepley - coordFunc - the function setting initial particle velocities 7186c5a79ebSMatthew G. Knepley 7196c5a79ebSMatthew G. Knepley Level: intermediate 7206c5a79ebSMatthew G. Knepley 721db781477SPatrick Sanan .seealso: `DMSwarmGetVelocityFunction()`, `DMSwarmSetCoordinateFunction()`, `DMSwarmInitializeVelocities()` 7226c5a79ebSMatthew G. Knepley @*/ 723*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmSetVelocityFunction(DM sw, PetscSimplePointFunc velFunc) 724*d71ae5a4SJacob Faibussowitsch { 7256c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 7266c5a79ebSMatthew G. Knepley 7276c5a79ebSMatthew G. Knepley PetscFunctionBegin; 7286c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 7296c5a79ebSMatthew G. Knepley PetscValidFunction(velFunc, 2); 7306c5a79ebSMatthew G. Knepley swarm->velFunc = velFunc; 7316c5a79ebSMatthew G. Knepley PetscFunctionReturn(0); 7326c5a79ebSMatthew G. Knepley } 7336c5a79ebSMatthew G. Knepley 7346c5a79ebSMatthew G. Knepley /*@C 73535b38c60SMatthew G. Knepley DMSwarmComputeLocalSize - Compute the local number and distribution of particles based upon a density function 73635b38c60SMatthew G. Knepley 73735b38c60SMatthew G. Knepley Not collective 73835b38c60SMatthew G. Knepley 73935b38c60SMatthew G. Knepley Input Parameters: 74035b38c60SMatthew G. Knepley + sw - The DMSwarm 74135b38c60SMatthew G. Knepley . N - The target number of particles 74235b38c60SMatthew G. Knepley - density - The density field for the particle layout, normalized to unity 74335b38c60SMatthew G. Knepley 74435b38c60SMatthew G. Knepley Note: One particle will be created for each species. 74535b38c60SMatthew G. Knepley 74635b38c60SMatthew G. Knepley Level: advanced 74735b38c60SMatthew G. Knepley 748db781477SPatrick Sanan .seealso: `DMSwarmComputeLocalSizeFromOptions()` 74935b38c60SMatthew G. Knepley @*/ 750*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmComputeLocalSize(DM sw, PetscInt N, PetscProbFunc density) 751*d71ae5a4SJacob Faibussowitsch { 75235b38c60SMatthew G. Knepley DM dm; 75335b38c60SMatthew G. Knepley PetscQuadrature quad; 75435b38c60SMatthew G. Knepley const PetscReal *xq, *wq; 75535b38c60SMatthew G. Knepley PetscInt *npc, *cellid; 7566c5a79ebSMatthew G. Knepley PetscReal xi0[3]; 75735b38c60SMatthew G. Knepley PetscInt Ns, cStart, cEnd, c, dim, d, Nq, q, Np = 0, p; 75835b38c60SMatthew G. Knepley PetscBool simplex; 75935b38c60SMatthew G. Knepley 76035b38c60SMatthew G. Knepley PetscFunctionBegin; 7619566063dSJacob Faibussowitsch PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 7629566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sw, &dm)); 7639566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 7649566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 7659566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 7666858538eSMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(dm)); 7679566063dSJacob Faibussowitsch if (simplex) PetscCall(PetscDTStroudConicalQuadrature(dim, 1, 5, -1.0, 1.0, &quad)); 7689566063dSJacob Faibussowitsch else PetscCall(PetscDTGaussTensorQuadrature(dim, 1, 5, -1.0, 1.0, &quad)); 7699566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &xq, &wq)); 7709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(cEnd - cStart, &npc)); 77135b38c60SMatthew G. Knepley /* Integrate the density function to get the number of particles in each cell */ 77235b38c60SMatthew G. Knepley for (d = 0; d < dim; ++d) xi0[d] = -1.0; 77335b38c60SMatthew G. Knepley for (c = 0; c < cEnd - cStart; ++c) { 77435b38c60SMatthew G. Knepley const PetscInt cell = c + cStart; 77535b38c60SMatthew G. Knepley PetscReal v0[3], J[9], invJ[9], detJ; 77635b38c60SMatthew G. Knepley PetscReal n_int = 0.; 77735b38c60SMatthew G. Knepley 7789566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, NULL, v0, J, invJ, &detJ)); 77935b38c60SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 78035b38c60SMatthew G. Knepley PetscReal xr[3], den[3]; 78135b38c60SMatthew G. Knepley 78235b38c60SMatthew G. Knepley CoordinatesRefToReal(dim, dim, xi0, v0, J, &xq[q * dim], xr); 7836c5a79ebSMatthew G. Knepley PetscCall(density(xr, NULL, den)); 7846c5a79ebSMatthew G. Knepley n_int += den[0] * wq[q]; 78535b38c60SMatthew G. Knepley } 7866c5a79ebSMatthew G. Knepley npc[c] = (PetscInt)(N * n_int); 78735b38c60SMatthew G. Knepley npc[c] *= Ns; 78835b38c60SMatthew G. Knepley Np += npc[c]; 78935b38c60SMatthew G. Knepley } 7909566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quad)); 7919566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(sw, Np, 0)); 79235b38c60SMatthew G. Knepley 7939566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 79435b38c60SMatthew G. Knepley for (c = 0, p = 0; c < cEnd - cStart; ++c) { 79535b38c60SMatthew G. Knepley for (q = 0; q < npc[c]; ++q, ++p) cellid[p] = c; 79635b38c60SMatthew G. Knepley } 7979566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 7989566063dSJacob Faibussowitsch PetscCall(PetscFree(npc)); 79935b38c60SMatthew G. Knepley PetscFunctionReturn(0); 80035b38c60SMatthew G. Knepley } 80135b38c60SMatthew G. Knepley 80235b38c60SMatthew G. Knepley /*@ 80335b38c60SMatthew G. Knepley DMSwarmComputeLocalSizeFromOptions - Compute the local number and distribution of particles based upon a density function determined by options 80435b38c60SMatthew G. Knepley 80535b38c60SMatthew G. Knepley Not collective 80635b38c60SMatthew G. Knepley 80735b38c60SMatthew G. Knepley Input Parameters: 80835b38c60SMatthew G. Knepley , sw - The DMSwarm 80935b38c60SMatthew G. Knepley 81035b38c60SMatthew G. Knepley Level: advanced 81135b38c60SMatthew G. Knepley 812db781477SPatrick Sanan .seealso: `DMSwarmComputeLocalSize()` 81335b38c60SMatthew G. Knepley @*/ 814*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmComputeLocalSizeFromOptions(DM sw) 815*d71ae5a4SJacob Faibussowitsch { 81635b38c60SMatthew G. Knepley PetscProbFunc pdf; 81735b38c60SMatthew G. Knepley const char *prefix; 8186c5a79ebSMatthew G. Knepley char funcname[PETSC_MAX_PATH_LEN]; 8196c5a79ebSMatthew G. Knepley PetscInt *N, Ns, dim, n; 8206c5a79ebSMatthew G. Knepley PetscBool flg; 8216c5a79ebSMatthew G. Knepley PetscMPIInt size, rank; 82235b38c60SMatthew G. Knepley 82335b38c60SMatthew G. Knepley PetscFunctionBegin; 8246c5a79ebSMatthew G. Knepley PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sw), &size)); 8256c5a79ebSMatthew G. Knepley PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sw), &rank)); 8266c5a79ebSMatthew G. Knepley PetscCall(PetscCalloc1(size, &N)); 827d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)sw), "", "DMSwarm Options", "DMSWARM"); 8286c5a79ebSMatthew G. Knepley n = size; 8296c5a79ebSMatthew G. Knepley PetscCall(PetscOptionsIntArray("-dm_swarm_num_particles", "The target number of particles", "", N, &n, NULL)); 8306c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 8319566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-dm_swarm_num_species", "The number of species", "DMSwarmSetNumSpecies", Ns, &Ns, &flg)); 8329566063dSJacob Faibussowitsch if (flg) PetscCall(DMSwarmSetNumSpecies(sw, Ns)); 8336c5a79ebSMatthew G. Knepley PetscCall(PetscOptionsString("-dm_swarm_coordinate_function", "Function to determine particle coordinates", "DMSwarmSetCoordinateFunction", funcname, funcname, sizeof(funcname), &flg)); 834d0609cedSBarry Smith PetscOptionsEnd(); 8356c5a79ebSMatthew G. Knepley if (flg) { 8366c5a79ebSMatthew G. Knepley PetscSimplePointFunc coordFunc; 83735b38c60SMatthew G. Knepley 8386c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 8396c5a79ebSMatthew G. Knepley PetscCall(PetscDLSym(NULL, funcname, (void **)&coordFunc)); 8406c5a79ebSMatthew G. Knepley PetscCheck(coordFunc, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname); 8416c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 8426c5a79ebSMatthew G. Knepley PetscCall(DMSwarmSetLocalSizes(sw, N[rank] * Ns, 0)); 8436c5a79ebSMatthew G. Knepley PetscCall(DMSwarmSetCoordinateFunction(sw, coordFunc)); 8446c5a79ebSMatthew G. Knepley } else { 8459566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 8469566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sw, &prefix)); 8479566063dSJacob Faibussowitsch PetscCall(PetscProbCreateFromOptions(dim, prefix, "-dm_swarm_coordinate_density", &pdf, NULL, NULL)); 8486c5a79ebSMatthew G. Knepley PetscCall(DMSwarmComputeLocalSize(sw, N[rank], pdf)); 8496c5a79ebSMatthew G. Knepley } 8506c5a79ebSMatthew G. Knepley PetscCall(PetscFree(N)); 85135b38c60SMatthew G. Knepley PetscFunctionReturn(0); 85235b38c60SMatthew G. Knepley } 85335b38c60SMatthew G. Knepley 85435b38c60SMatthew G. Knepley /*@ 85535b38c60SMatthew G. Knepley DMSwarmInitializeCoordinates - Determine the initial coordinates of particles for a PIC method 85635b38c60SMatthew G. Knepley 85735b38c60SMatthew G. Knepley Not collective 85835b38c60SMatthew G. Knepley 85935b38c60SMatthew G. Knepley Input Parameters: 86035b38c60SMatthew G. Knepley , sw - The DMSwarm 86135b38c60SMatthew G. Knepley 86235b38c60SMatthew G. Knepley Note: Currently, we randomly place particles in their assigned cell 86335b38c60SMatthew G. Knepley 86435b38c60SMatthew G. Knepley Level: advanced 86535b38c60SMatthew G. Knepley 866db781477SPatrick Sanan .seealso: `DMSwarmComputeLocalSize()`, `DMSwarmInitializeVelocities()` 86735b38c60SMatthew G. Knepley @*/ 868*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmInitializeCoordinates(DM sw) 869*d71ae5a4SJacob Faibussowitsch { 8706c5a79ebSMatthew G. Knepley PetscSimplePointFunc coordFunc; 87135b38c60SMatthew G. Knepley PetscScalar *weight; 8726c5a79ebSMatthew G. Knepley PetscReal *x; 87335b38c60SMatthew G. Knepley PetscInt *species; 8746c5a79ebSMatthew G. Knepley void *ctx; 87535b38c60SMatthew G. Knepley PetscBool removePoints = PETSC_TRUE; 87635b38c60SMatthew G. Knepley PetscDataType dtype; 8776c5a79ebSMatthew G. Knepley PetscInt Np, p, Ns, dim, d, bs; 87835b38c60SMatthew G. Knepley 87935b38c60SMatthew G. Knepley PetscFunctionBeginUser; 8806c5a79ebSMatthew G. Knepley PetscCall(DMGetDimension(sw, &dim)); 8816c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetLocalSize(sw, &Np)); 8829566063dSJacob Faibussowitsch PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 8836c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetCoordinateFunction(sw, &coordFunc)); 8846c5a79ebSMatthew G. Knepley 8856c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, &bs, &dtype, (void **)&x)); 8866c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&weight)); 8876c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 8886c5a79ebSMatthew G. Knepley if (coordFunc) { 8896c5a79ebSMatthew G. Knepley PetscCall(DMGetApplicationContext(sw, &ctx)); 8906c5a79ebSMatthew G. Knepley for (p = 0; p < Np; ++p) { 8916c5a79ebSMatthew G. Knepley PetscScalar X[3]; 8926c5a79ebSMatthew G. Knepley 8936c5a79ebSMatthew G. Knepley PetscCall((*coordFunc)(dim, 0., NULL, p, X, ctx)); 8946c5a79ebSMatthew G. Knepley for (d = 0; d < dim; ++d) x[p * dim + d] = PetscRealPart(X[d]); 8956c5a79ebSMatthew G. Knepley weight[p] = 1.0; 8966c5a79ebSMatthew G. Knepley species[p] = p % Ns; 8976c5a79ebSMatthew G. Knepley } 8986c5a79ebSMatthew G. Knepley } else { 8996c5a79ebSMatthew G. Knepley DM dm; 9006c5a79ebSMatthew G. Knepley PetscRandom rnd; 9016c5a79ebSMatthew G. Knepley PetscReal xi0[3]; 9026c5a79ebSMatthew G. Knepley PetscInt cStart, cEnd, c; 9036c5a79ebSMatthew G. Knepley 9049566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sw, &dm)); 9059566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 90635b38c60SMatthew G. Knepley 90735b38c60SMatthew G. Knepley /* Set particle position randomly in cell, set weights to 1 */ 9089566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd)); 9099566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0)); 9109566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rnd)); 9119566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetAccess(sw)); 91235b38c60SMatthew G. Knepley for (d = 0; d < dim; ++d) xi0[d] = -1.0; 91335b38c60SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 91435b38c60SMatthew G. Knepley PetscReal v0[3], J[9], invJ[9], detJ; 91535b38c60SMatthew G. Knepley PetscInt *pidx, Npc, q; 91635b38c60SMatthew G. Knepley 9179566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx)); 9189566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); 91935b38c60SMatthew G. Knepley for (q = 0; q < Npc; ++q) { 92035b38c60SMatthew G. Knepley const PetscInt p = pidx[q]; 92135b38c60SMatthew G. Knepley PetscReal xref[3]; 92235b38c60SMatthew G. Knepley 9239566063dSJacob Faibussowitsch for (d = 0; d < dim; ++d) PetscCall(PetscRandomGetValueReal(rnd, &xref[d])); 92435b38c60SMatthew G. Knepley CoordinatesRefToReal(dim, dim, xi0, v0, J, xref, &x[p * dim]); 92535b38c60SMatthew G. Knepley 92635b38c60SMatthew G. Knepley weight[p] = 1.0; 9276c5a79ebSMatthew G. Knepley species[p] = p % Ns; 92835b38c60SMatthew G. Knepley } 9299566063dSJacob Faibussowitsch PetscCall(PetscFree(pidx)); 93035b38c60SMatthew G. Knepley } 9319566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rnd)); 9329566063dSJacob Faibussowitsch PetscCall(DMSwarmSortRestoreAccess(sw)); 9336c5a79ebSMatthew G. Knepley } 9349566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x)); 9356c5a79ebSMatthew G. Knepley PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight)); 9369566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 9376c5a79ebSMatthew G. Knepley 9389566063dSJacob Faibussowitsch PetscCall(DMSwarmMigrate(sw, removePoints)); 9399566063dSJacob Faibussowitsch PetscCall(DMLocalizeCoordinates(sw)); 94035b38c60SMatthew G. Knepley PetscFunctionReturn(0); 94135b38c60SMatthew G. Knepley } 94235b38c60SMatthew G. Knepley 94335b38c60SMatthew G. Knepley /*@C 94435b38c60SMatthew G. Knepley DMSwarmInitializeVelocities - Set the initial velocities of particles using a distribution. 94535b38c60SMatthew G. Knepley 94635b38c60SMatthew G. Knepley Collective on dm 94735b38c60SMatthew G. Knepley 94835b38c60SMatthew G. Knepley Input Parameters: 94935b38c60SMatthew G. Knepley + sw - The DMSwarm object 95035b38c60SMatthew G. Knepley . sampler - A function which uniformly samples the velocity PDF 95135b38c60SMatthew G. Knepley - v0 - The velocity scale for nondimensionalization for each species 95235b38c60SMatthew G. Knepley 9536c5a79ebSMatthew G. Knepley Note: 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. 9546c5a79ebSMatthew G. Knepley 95535b38c60SMatthew G. Knepley Level: advanced 95635b38c60SMatthew G. Knepley 957db781477SPatrick Sanan .seealso: `DMSwarmComputeLocalSize()`, `DMSwarmInitializeCoordinates()`, `DMSwarmInitializeVelocitiesFromOptions()` 95835b38c60SMatthew G. Knepley @*/ 959*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmInitializeVelocities(DM sw, PetscProbFunc sampler, const PetscReal v0[]) 960*d71ae5a4SJacob Faibussowitsch { 9616c5a79ebSMatthew G. Knepley PetscSimplePointFunc velFunc; 96235b38c60SMatthew G. Knepley PetscReal *v; 96335b38c60SMatthew G. Knepley PetscInt *species; 9646c5a79ebSMatthew G. Knepley void *ctx; 96535b38c60SMatthew G. Knepley PetscInt dim, Np, p; 96635b38c60SMatthew G. Knepley 96735b38c60SMatthew G. Knepley PetscFunctionBegin; 9686c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetVelocityFunction(sw, &velFunc)); 96935b38c60SMatthew G. Knepley 9709566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 9719566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(sw, &Np)); 9729566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&v)); 9739566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 9746c5a79ebSMatthew G. Knepley if (v0[0] == 0.) { 9756c5a79ebSMatthew G. Knepley PetscCall(PetscArrayzero(v, Np * dim)); 9766c5a79ebSMatthew G. Knepley } else if (velFunc) { 9776c5a79ebSMatthew G. Knepley PetscCall(DMGetApplicationContext(sw, &ctx)); 9786c5a79ebSMatthew G. Knepley for (p = 0; p < Np; ++p) { 9796c5a79ebSMatthew G. Knepley PetscInt s = species[p], d; 9806c5a79ebSMatthew G. Knepley PetscScalar vel[3]; 9816c5a79ebSMatthew G. Knepley 9826c5a79ebSMatthew G. Knepley PetscCall((*velFunc)(dim, 0., NULL, p, vel, ctx)); 9836c5a79ebSMatthew G. Knepley for (d = 0; d < dim; ++d) v[p * dim + d] = (v0[s] / v0[0]) * PetscRealPart(vel[d]); 9846c5a79ebSMatthew G. Knepley } 9856c5a79ebSMatthew G. Knepley } else { 9866c5a79ebSMatthew G. Knepley PetscRandom rnd; 9876c5a79ebSMatthew G. Knepley 9886c5a79ebSMatthew G. Knepley PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rnd)); 9896c5a79ebSMatthew G. Knepley PetscCall(PetscRandomSetInterval(rnd, 0, 1.)); 9906c5a79ebSMatthew G. Knepley PetscCall(PetscRandomSetFromOptions(rnd)); 9916c5a79ebSMatthew G. Knepley 99235b38c60SMatthew G. Knepley for (p = 0; p < Np; ++p) { 99335b38c60SMatthew G. Knepley PetscInt s = species[p], d; 99435b38c60SMatthew G. Knepley PetscReal a[3], vel[3]; 99535b38c60SMatthew G. Knepley 9969566063dSJacob Faibussowitsch for (d = 0; d < dim; ++d) PetscCall(PetscRandomGetValueReal(rnd, &a[d])); 9979566063dSJacob Faibussowitsch PetscCall(sampler(a, NULL, vel)); 998ad540459SPierre Jolivet for (d = 0; d < dim; ++d) v[p * dim + d] = (v0[s] / v0[0]) * vel[d]; 99935b38c60SMatthew G. Knepley } 10006c5a79ebSMatthew G. Knepley PetscCall(PetscRandomDestroy(&rnd)); 10016c5a79ebSMatthew G. Knepley } 10029566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&v)); 10039566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 100435b38c60SMatthew G. Knepley PetscFunctionReturn(0); 100535b38c60SMatthew G. Knepley } 100635b38c60SMatthew G. Knepley 100735b38c60SMatthew G. Knepley /*@ 100835b38c60SMatthew G. Knepley DMSwarmInitializeVelocitiesFromOptions - Set the initial velocities of particles using a distribution determined from options. 100935b38c60SMatthew G. Knepley 101035b38c60SMatthew G. Knepley Collective on dm 101135b38c60SMatthew G. Knepley 101235b38c60SMatthew G. Knepley Input Parameters: 101335b38c60SMatthew G. Knepley + sw - The DMSwarm object 101435b38c60SMatthew G. Knepley - v0 - The velocity scale for nondimensionalization for each species 101535b38c60SMatthew G. Knepley 101635b38c60SMatthew G. Knepley Level: advanced 101735b38c60SMatthew G. Knepley 1018db781477SPatrick Sanan .seealso: `DMSwarmComputeLocalSize()`, `DMSwarmInitializeCoordinates()`, `DMSwarmInitializeVelocities()` 101935b38c60SMatthew G. Knepley @*/ 1020*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSwarmInitializeVelocitiesFromOptions(DM sw, const PetscReal v0[]) 1021*d71ae5a4SJacob Faibussowitsch { 102235b38c60SMatthew G. Knepley PetscProbFunc sampler; 102335b38c60SMatthew G. Knepley PetscInt dim; 102435b38c60SMatthew G. Knepley const char *prefix; 10256c5a79ebSMatthew G. Knepley char funcname[PETSC_MAX_PATH_LEN]; 10266c5a79ebSMatthew G. Knepley PetscBool flg; 102735b38c60SMatthew G. Knepley 102835b38c60SMatthew G. Knepley PetscFunctionBegin; 1029d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)sw), "", "DMSwarm Options", "DMSWARM"); 10306c5a79ebSMatthew G. Knepley PetscCall(PetscOptionsString("-dm_swarm_velocity_function", "Function to determine particle velocities", "DMSwarmSetVelocityFunction", funcname, funcname, sizeof(funcname), &flg)); 1031d0609cedSBarry Smith PetscOptionsEnd(); 10326c5a79ebSMatthew G. Knepley if (flg) { 10336c5a79ebSMatthew G. Knepley PetscSimplePointFunc velFunc; 10346c5a79ebSMatthew G. Knepley 10356c5a79ebSMatthew G. Knepley PetscCall(PetscDLSym(NULL, funcname, (void **)&velFunc)); 10366c5a79ebSMatthew G. Knepley PetscCheck(velFunc, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname); 10376c5a79ebSMatthew G. Knepley PetscCall(DMSwarmSetVelocityFunction(sw, velFunc)); 10386c5a79ebSMatthew G. Knepley } 10399566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 10409566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sw, &prefix)); 10419566063dSJacob Faibussowitsch PetscCall(PetscProbCreateFromOptions(dim, prefix, "-dm_swarm_velocity_density", NULL, NULL, &sampler)); 10429566063dSJacob Faibussowitsch PetscCall(DMSwarmInitializeVelocities(sw, sampler, v0)); 104335b38c60SMatthew G. Knepley PetscFunctionReturn(0); 104435b38c60SMatthew G. Knepley } 1045