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 @*/ 439371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMSwarmSetPointsUniformCoordinates(DM dm, PetscReal min[], PetscReal max[], PetscInt npoints[], InsertMode mode) { 440e2ec84fSDave May PetscReal gmin[] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL}; 450e2ec84fSDave May PetscReal gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 460e2ec84fSDave May PetscInt i, j, k, N, bs, b, n_estimate, n_curr, n_new_est, p, n_found; 470e2ec84fSDave May Vec coorlocal; 480e2ec84fSDave May const PetscScalar *_coor; 490e2ec84fSDave May DM celldm; 500e2ec84fSDave May PetscReal dx[3]; 512e3d3761SDave May PetscInt _npoints[] = {0, 0, 1}; 520e2ec84fSDave May Vec pos; 530e2ec84fSDave May PetscScalar *_pos; 540e2ec84fSDave May PetscReal *swarm_coor; 550e2ec84fSDave May PetscInt *swarm_cellid; 560e2ec84fSDave May PetscSF sfcell = NULL; 570e2ec84fSDave May const PetscSFNode *LA_sfcell; 580e2ec84fSDave May 590e2ec84fSDave May PetscFunctionBegin; 600e2ec84fSDave May DMSWARMPICVALID(dm); 619566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 629566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(celldm, &coorlocal)); 639566063dSJacob Faibussowitsch PetscCall(VecGetSize(coorlocal, &N)); 649566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(coorlocal, &bs)); 650e2ec84fSDave May N = N / bs; 669566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coorlocal, &_coor)); 670e2ec84fSDave May for (i = 0; i < N; i++) { 680e2ec84fSDave May for (b = 0; b < bs; b++) { 69a5f152d1SDave May gmin[b] = PetscMin(gmin[b], PetscRealPart(_coor[bs * i + b])); 70a5f152d1SDave May gmax[b] = PetscMax(gmax[b], PetscRealPart(_coor[bs * i + b])); 710e2ec84fSDave May } 720e2ec84fSDave May } 739566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coorlocal, &_coor)); 740e2ec84fSDave May 750e2ec84fSDave May for (b = 0; b < bs; b++) { 7671844bb8SDave May if (npoints[b] > 1) { 770e2ec84fSDave May dx[b] = (max[b] - min[b]) / ((PetscReal)(npoints[b] - 1)); 7871844bb8SDave May } else { 7971844bb8SDave May dx[b] = 0.0; 8071844bb8SDave May } 812e3d3761SDave May _npoints[b] = npoints[b]; 820e2ec84fSDave May } 830e2ec84fSDave May 840e2ec84fSDave May /* determine number of points living in the bounding box */ 850e2ec84fSDave May n_estimate = 0; 862e3d3761SDave May for (k = 0; k < _npoints[2]; k++) { 872e3d3761SDave May for (j = 0; j < _npoints[1]; j++) { 882e3d3761SDave May for (i = 0; i < _npoints[0]; i++) { 890e2ec84fSDave May PetscReal xp[] = {0.0, 0.0, 0.0}; 900e2ec84fSDave May PetscInt ijk[3]; 910e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 920e2ec84fSDave May 930e2ec84fSDave May ijk[0] = i; 940e2ec84fSDave May ijk[1] = j; 950e2ec84fSDave May ijk[2] = k; 969371c9d4SSatish Balay for (b = 0; b < bs; b++) { xp[b] = min[b] + ijk[b] * dx[b]; } 970e2ec84fSDave May for (b = 0; b < bs; b++) { 980e2ec84fSDave May if (xp[b] < gmin[b]) { point_inside = PETSC_FALSE; } 990e2ec84fSDave May if (xp[b] > gmax[b]) { point_inside = PETSC_FALSE; } 1000e2ec84fSDave May } 1010e2ec84fSDave May if (point_inside) { n_estimate++; } 1020e2ec84fSDave May } 1030e2ec84fSDave May } 1040e2ec84fSDave May } 1050e2ec84fSDave May 1060e2ec84fSDave May /* create candidate list */ 1079566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &pos)); 1089566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pos, bs * n_estimate, PETSC_DECIDE)); 1099566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(pos, bs)); 1109566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pos)); 1119566063dSJacob Faibussowitsch PetscCall(VecGetArray(pos, &_pos)); 1120e2ec84fSDave May 1130e2ec84fSDave May n_estimate = 0; 1142e3d3761SDave May for (k = 0; k < _npoints[2]; k++) { 1152e3d3761SDave May for (j = 0; j < _npoints[1]; j++) { 1162e3d3761SDave May for (i = 0; i < _npoints[0]; i++) { 1170e2ec84fSDave May PetscReal xp[] = {0.0, 0.0, 0.0}; 1180e2ec84fSDave May PetscInt ijk[3]; 1190e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 1200e2ec84fSDave May 1210e2ec84fSDave May ijk[0] = i; 1220e2ec84fSDave May ijk[1] = j; 1230e2ec84fSDave May ijk[2] = k; 1249371c9d4SSatish Balay for (b = 0; b < bs; b++) { xp[b] = min[b] + ijk[b] * dx[b]; } 1250e2ec84fSDave May for (b = 0; b < bs; b++) { 1260e2ec84fSDave May if (xp[b] < gmin[b]) { point_inside = PETSC_FALSE; } 1270e2ec84fSDave May if (xp[b] > gmax[b]) { point_inside = PETSC_FALSE; } 1280e2ec84fSDave May } 1290e2ec84fSDave May if (point_inside) { 1309371c9d4SSatish Balay for (b = 0; b < bs; b++) { _pos[bs * n_estimate + b] = xp[b]; } 1310e2ec84fSDave May n_estimate++; 1320e2ec84fSDave May } 1330e2ec84fSDave May } 1340e2ec84fSDave May } 1350e2ec84fSDave May } 1369566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pos, &_pos)); 1370e2ec84fSDave May 1380e2ec84fSDave May /* locate points */ 1399566063dSJacob Faibussowitsch PetscCall(DMLocatePoints(celldm, pos, DM_POINTLOCATION_NONE, &sfcell)); 1409566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 1410e2ec84fSDave May n_found = 0; 1420e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 1439371c9d4SSatish Balay if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { n_found++; } 1440e2ec84fSDave May } 1450e2ec84fSDave May 1460e2ec84fSDave May /* adjust size */ 1470e2ec84fSDave May if (mode == ADD_VALUES) { 1489566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dm, &n_curr)); 1490e2ec84fSDave May n_new_est = n_curr + n_found; 1509566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 1510e2ec84fSDave May } 1520e2ec84fSDave May if (mode == INSERT_VALUES) { 1530e2ec84fSDave May n_curr = 0; 1540e2ec84fSDave May n_new_est = n_found; 1559566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 1560e2ec84fSDave May } 1570e2ec84fSDave May 1580e2ec84fSDave May /* initialize new coords, cell owners, pid */ 1599566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pos, &_coor)); 1609566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 1619566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 1620e2ec84fSDave May n_found = 0; 1630e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 1640e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 1659371c9d4SSatish Balay for (b = 0; b < bs; b++) { swarm_coor[bs * (n_curr + n_found) + b] = PetscRealPart(_coor[bs * p + b]); } 1660e2ec84fSDave May swarm_cellid[n_curr + n_found] = LA_sfcell[p].index; 1670e2ec84fSDave May n_found++; 1680e2ec84fSDave May } 1690e2ec84fSDave May } 1709566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 1719566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 1729566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pos, &_coor)); 1730e2ec84fSDave May 1749566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfcell)); 1759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pos)); 1760e2ec84fSDave May PetscFunctionReturn(0); 1770e2ec84fSDave May } 1780e2ec84fSDave May 1790e2ec84fSDave May /*@C 1800e2ec84fSDave May DMSwarmSetPointCoordinates - Set point coordinates in a DMSwarm from a user defined list 1810e2ec84fSDave May 182d083f849SBarry Smith Collective on dm 1830e2ec84fSDave May 1840e2ec84fSDave May Input parameters: 1850e2ec84fSDave May + dm - the DMSwarm 1860e2ec84fSDave May . npoints - the number of points to insert 1870e2ec84fSDave May . coor - the coordinate values 1880e2ec84fSDave 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 1890e2ec84fSDave May - mode - indicates whether to append points to the swarm (ADD_VALUES), or over-ride existing points (INSERT_VALUES) 1900e2ec84fSDave May 1910e2ec84fSDave May Level: beginner 1920e2ec84fSDave May 1930e2ec84fSDave May Notes: 1940e2ec84fSDave May If the user has specified redundant = PETSC_FALSE, the cell DM will attempt to locate the coordinates provided by coor[] within 1950e2ec84fSDave 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 1960e2ec84fSDave May added to the DMSwarm. 1970e2ec84fSDave May 198db781477SPatrick Sanan .seealso: `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType`, `DMSwarmSetPointsUniformCoordinates()` 1990e2ec84fSDave May @*/ 2009371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinates(DM dm, PetscInt npoints, PetscReal coor[], PetscBool redundant, InsertMode mode) { 2010e2ec84fSDave May PetscReal gmin[] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL}; 2020e2ec84fSDave May PetscReal gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 2030e2ec84fSDave May PetscInt i, N, bs, b, n_estimate, n_curr, n_new_est, p, n_found; 2040e2ec84fSDave May Vec coorlocal; 2050e2ec84fSDave May const PetscScalar *_coor; 2060e2ec84fSDave May DM celldm; 2070e2ec84fSDave May Vec pos; 2080e2ec84fSDave May PetscScalar *_pos; 2090e2ec84fSDave May PetscReal *swarm_coor; 2100e2ec84fSDave May PetscInt *swarm_cellid; 2110e2ec84fSDave May PetscSF sfcell = NULL; 2120e2ec84fSDave May const PetscSFNode *LA_sfcell; 2130e2ec84fSDave May PetscReal *my_coor; 2140e2ec84fSDave May PetscInt my_npoints; 2150e2ec84fSDave May PetscMPIInt rank; 2160e2ec84fSDave May MPI_Comm comm; 2170e2ec84fSDave May 2180e2ec84fSDave May PetscFunctionBegin; 2190e2ec84fSDave May DMSWARMPICVALID(dm); 2209566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 2219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2220e2ec84fSDave May 2239566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 2249566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(celldm, &coorlocal)); 2259566063dSJacob Faibussowitsch PetscCall(VecGetSize(coorlocal, &N)); 2269566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(coorlocal, &bs)); 2270e2ec84fSDave May N = N / bs; 2289566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coorlocal, &_coor)); 2290e2ec84fSDave May for (i = 0; i < N; i++) { 2300e2ec84fSDave May for (b = 0; b < bs; b++) { 231a5f152d1SDave May gmin[b] = PetscMin(gmin[b], PetscRealPart(_coor[bs * i + b])); 232a5f152d1SDave May gmax[b] = PetscMax(gmax[b], PetscRealPart(_coor[bs * i + b])); 2330e2ec84fSDave May } 2340e2ec84fSDave May } 2359566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coorlocal, &_coor)); 2360e2ec84fSDave May 2370e2ec84fSDave May /* broadcast points from rank 0 if requested */ 2380e2ec84fSDave May if (redundant) { 2390e2ec84fSDave May my_npoints = npoints; 2409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&my_npoints, 1, MPIU_INT, 0, comm)); 2410e2ec84fSDave May 2420e2ec84fSDave May if (rank > 0) { /* allocate space */ 2439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * my_npoints, &my_coor)); 2440e2ec84fSDave May } else { 2450e2ec84fSDave May my_coor = coor; 2460e2ec84fSDave May } 2479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(my_coor, bs * my_npoints, MPIU_REAL, 0, comm)); 2480e2ec84fSDave May } else { 2490e2ec84fSDave May my_npoints = npoints; 2500e2ec84fSDave May my_coor = coor; 2510e2ec84fSDave May } 2520e2ec84fSDave May 2530e2ec84fSDave May /* determine the number of points living in the bounding box */ 2540e2ec84fSDave May n_estimate = 0; 2550e2ec84fSDave May for (i = 0; i < my_npoints; i++) { 2560e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 2570e2ec84fSDave May 2580e2ec84fSDave May for (b = 0; b < bs; b++) { 2590e2ec84fSDave May if (my_coor[bs * i + b] < gmin[b]) { point_inside = PETSC_FALSE; } 2600e2ec84fSDave May if (my_coor[bs * i + b] > gmax[b]) { point_inside = PETSC_FALSE; } 2610e2ec84fSDave May } 2620e2ec84fSDave May if (point_inside) { n_estimate++; } 2630e2ec84fSDave May } 2640e2ec84fSDave May 2650e2ec84fSDave May /* create candidate list */ 2669566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &pos)); 2679566063dSJacob Faibussowitsch PetscCall(VecSetSizes(pos, bs * n_estimate, PETSC_DECIDE)); 2689566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(pos, bs)); 2699566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(pos)); 2709566063dSJacob Faibussowitsch PetscCall(VecGetArray(pos, &_pos)); 2710e2ec84fSDave May 2720e2ec84fSDave May n_estimate = 0; 2730e2ec84fSDave May for (i = 0; i < my_npoints; i++) { 2740e2ec84fSDave May PetscBool point_inside = PETSC_TRUE; 2750e2ec84fSDave May 2760e2ec84fSDave May for (b = 0; b < bs; b++) { 2770e2ec84fSDave May if (my_coor[bs * i + b] < gmin[b]) { point_inside = PETSC_FALSE; } 2780e2ec84fSDave May if (my_coor[bs * i + b] > gmax[b]) { point_inside = PETSC_FALSE; } 2790e2ec84fSDave May } 2800e2ec84fSDave May if (point_inside) { 2819371c9d4SSatish Balay for (b = 0; b < bs; b++) { _pos[bs * n_estimate + b] = my_coor[bs * i + b]; } 2820e2ec84fSDave May n_estimate++; 2830e2ec84fSDave May } 2840e2ec84fSDave May } 2859566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pos, &_pos)); 2860e2ec84fSDave May 2870e2ec84fSDave May /* locate points */ 2889566063dSJacob Faibussowitsch PetscCall(DMLocatePoints(celldm, pos, DM_POINTLOCATION_NONE, &sfcell)); 2890e2ec84fSDave May 2909566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 2910e2ec84fSDave May n_found = 0; 2920e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 2939371c9d4SSatish Balay if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { n_found++; } 2940e2ec84fSDave May } 2950e2ec84fSDave May 2960e2ec84fSDave May /* adjust size */ 2970e2ec84fSDave May if (mode == ADD_VALUES) { 2989566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dm, &n_curr)); 2990e2ec84fSDave May n_new_est = n_curr + n_found; 3009566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 3010e2ec84fSDave May } 3020e2ec84fSDave May if (mode == INSERT_VALUES) { 3030e2ec84fSDave May n_curr = 0; 3040e2ec84fSDave May n_new_est = n_found; 3059566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 3060e2ec84fSDave May } 3070e2ec84fSDave May 3080e2ec84fSDave May /* initialize new coords, cell owners, pid */ 3099566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pos, &_coor)); 3109566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 3119566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 3120e2ec84fSDave May n_found = 0; 3130e2ec84fSDave May for (p = 0; p < n_estimate; p++) { 3140e2ec84fSDave May if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 3159371c9d4SSatish Balay for (b = 0; b < bs; b++) { swarm_coor[bs * (n_curr + n_found) + b] = PetscRealPart(_coor[bs * p + b]); } 3160e2ec84fSDave May swarm_cellid[n_curr + n_found] = LA_sfcell[p].index; 3170e2ec84fSDave May n_found++; 3180e2ec84fSDave May } 3190e2ec84fSDave May } 3209566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 3219566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 3229566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pos, &_coor)); 3230e2ec84fSDave May 3240e2ec84fSDave May if (redundant) { 325*48a46eb9SPierre Jolivet if (rank > 0) PetscCall(PetscFree(my_coor)); 3260e2ec84fSDave May } 3279566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfcell)); 3289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pos)); 3290e2ec84fSDave May PetscFunctionReturn(0); 3300e2ec84fSDave May } 3310e2ec84fSDave May 3320e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM, DM, DMSwarmPICLayoutType, PetscInt); 3330e2ec84fSDave May extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM, DM, DMSwarmPICLayoutType, PetscInt); 3340e2ec84fSDave May 3350e2ec84fSDave May /*@C 3360e2ec84fSDave May DMSwarmInsertPointsUsingCellDM - Insert point coordinates within each cell 3370e2ec84fSDave May 3380e2ec84fSDave May Not collective 3390e2ec84fSDave May 3400e2ec84fSDave May Input parameters: 3410e2ec84fSDave May + dm - the DMSwarm 3420e2ec84fSDave May . layout_type - method used to fill each cell with the cell DM 3430e2ec84fSDave May - fill_param - parameter controlling how many points per cell are added (the meaning of this parameter is dependent on the layout type) 3440e2ec84fSDave May 3450e2ec84fSDave May Level: beginner 3460e2ec84fSDave May 3470e2ec84fSDave May Notes: 348e7af74c8SDave May 349e7af74c8SDave May The insert method will reset any previous defined points within the DMSwarm. 350e7af74c8SDave May 351e7af74c8SDave May When using a DMDA both 2D and 3D are supported for all layout types provided you are using DMDA_ELEMENT_Q1. 352e7af74c8SDave May 353e7af74c8SDave May When using a DMPLEX the following case are supported: 354ea3d7275SDave May (i) DMSWARMPIC_LAYOUT_REGULAR: 2D (triangle), 355ea3d7275SDave May (ii) DMSWARMPIC_LAYOUT_GAUSS: 2D and 3D provided the cell is a tri/tet or a quad/hex, 356ea3d7275SDave May (iii) DMSWARMPIC_LAYOUT_SUBDIVISION: 2D and 3D for quad/hex and 2D tri. 3570e2ec84fSDave May 358db781477SPatrick Sanan .seealso: `DMSwarmPICLayoutType`, `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType` 3590e2ec84fSDave May @*/ 3609371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMSwarmInsertPointsUsingCellDM(DM dm, DMSwarmPICLayoutType layout_type, PetscInt fill_param) { 3610e2ec84fSDave May DM celldm; 3620e2ec84fSDave May PetscBool isDA, isPLEX; 3630e2ec84fSDave May 3640e2ec84fSDave May PetscFunctionBegin; 3650e2ec84fSDave May DMSWARMPICVALID(dm); 3669566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 3679566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA)); 3689566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX)); 3690e2ec84fSDave May if (isDA) { 3709566063dSJacob Faibussowitsch PetscCall(private_DMSwarmInsertPointsUsingCellDM_DA(dm, celldm, layout_type, fill_param)); 3710e2ec84fSDave May } else if (isPLEX) { 3729566063dSJacob Faibussowitsch PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX(dm, celldm, layout_type, fill_param)); 3730e2ec84fSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX"); 3740e2ec84fSDave May PetscFunctionReturn(0); 3750e2ec84fSDave May } 3760e2ec84fSDave May 377431382f2SDave May extern PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM, DM, PetscInt, PetscReal *); 378431382f2SDave May 379431382f2SDave May /*@C 380431382f2SDave May DMSwarmSetPointCoordinatesCellwise - Insert point coordinates (defined over the reference cell) within each cell 381431382f2SDave May 382431382f2SDave May Not collective 383431382f2SDave May 384431382f2SDave May Input parameters: 385431382f2SDave May + dm - the DMSwarm 386431382f2SDave May . celldm - the cell DM 387431382f2SDave May . npoints - the number of points to insert in each cell 388431382f2SDave May - xi - the coordinates (defined in the local coordinate system for each cell) to insert 389431382f2SDave May 390431382f2SDave May Level: beginner 391431382f2SDave May 392431382f2SDave May Notes: 393431382f2SDave May The method will reset any previous defined points within the DMSwarm. 394ea3d7275SDave May Only supported for DMPLEX. If you are using a DMDA it is recommended to either use 395e7af74c8SDave May DMSwarmInsertPointsUsingCellDM(), or extract and set the coordinates yourself the following code 396e7af74c8SDave May 397e7af74c8SDave May $ PetscReal *coor; 398e7af74c8SDave May $ DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&coor); 399e7af74c8SDave May $ // user code to define the coordinates here 400e7af74c8SDave May $ DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&coor); 401431382f2SDave May 402db781477SPatrick Sanan .seealso: `DMSwarmSetCellDM()`, `DMSwarmInsertPointsUsingCellDM()` 403431382f2SDave May @*/ 4049371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinatesCellwise(DM dm, PetscInt npoints, PetscReal xi[]) { 405431382f2SDave May DM celldm; 406431382f2SDave May PetscBool isDA, isPLEX; 407431382f2SDave May 4080e2ec84fSDave May PetscFunctionBegin; 409431382f2SDave May DMSWARMPICVALID(dm); 4109566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 4119566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA)); 4129566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX)); 41328b400f6SJacob Faibussowitsch PetscCheck(!isDA, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMPLEX. Recommended you use DMSwarmInsertPointsUsingCellDM()"); 414f7d195e4SLawrence Mitchell if (isPLEX) { 4159566063dSJacob Faibussowitsch PetscCall(private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm, celldm, npoints, xi)); 416431382f2SDave May } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX"); 4170e2ec84fSDave May PetscFunctionReturn(0); 4180e2ec84fSDave May } 419431382f2SDave May 4200e2ec84fSDave May /* Field projection API */ 42177048351SPatrick Sanan extern PetscErrorCode private_DMSwarmProjectFields_DA(DM swarm, DM celldm, PetscInt project_type, PetscInt nfields, DMSwarmDataField dfield[], Vec vecs[]); 42277048351SPatrick Sanan extern PetscErrorCode private_DMSwarmProjectFields_PLEX(DM swarm, DM celldm, PetscInt project_type, PetscInt nfields, DMSwarmDataField dfield[], Vec vecs[]); 42394f7d2dcSDave May 42494f7d2dcSDave May /*@C 42594f7d2dcSDave May DMSwarmProjectFields - Project a set of swarm fields onto the cell DM 42694f7d2dcSDave May 427d083f849SBarry Smith Collective on dm 42894f7d2dcSDave May 42994f7d2dcSDave May Input parameters: 43094f7d2dcSDave May + dm - the DMSwarm 43194f7d2dcSDave May . nfields - the number of swarm fields to project 43294f7d2dcSDave May . fieldnames - the textual names of the swarm fields to project 43394f7d2dcSDave May . fields - an array of Vec's of length nfields 43494f7d2dcSDave May - reuse - flag indicating whether the array and contents of fields should be re-used or internally allocated 43594f7d2dcSDave May 43694f7d2dcSDave May Currently, the only available projection method consists of 43794f7d2dcSDave May phi_i = \sum_{p=0}^{np} N_i(x_p) phi_p dJ / \sum_{p=0}^{np} N_i(x_p) dJ 43894f7d2dcSDave May where phi_p is the swarm field at point p, 43994f7d2dcSDave May N_i() is the cell DM basis function at vertex i, 44094f7d2dcSDave May dJ is the determinant of the cell Jacobian and 44194f7d2dcSDave May phi_i is the projected vertex value of the field phi. 44294f7d2dcSDave May 44394f7d2dcSDave May Level: beginner 44494f7d2dcSDave May 44594f7d2dcSDave May Notes: 446e7af74c8SDave May 447e7af74c8SDave May If reuse = PETSC_FALSE, this function will allocate the array of Vec's, and each individual Vec. 44894f7d2dcSDave May The user is responsible for destroying both the array and the individual Vec objects. 449e7af74c8SDave May 450e7af74c8SDave May Only swarm fields registered with data type = PETSC_REAL can be projected onto the cell DM. 451e7af74c8SDave May 452e7af74c8SDave May Only swarm fields of block size = 1 can currently be projected. 453e7af74c8SDave May 454e7af74c8SDave May The only projection methods currently only support the DA (2D) and PLEX (triangles 2D). 45594f7d2dcSDave May 456db781477SPatrick Sanan .seealso: `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType` 45794f7d2dcSDave May @*/ 4589371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM dm, PetscInt nfields, const char *fieldnames[], Vec **fields, PetscBool reuse) { 45994f7d2dcSDave May DM_Swarm *swarm = (DM_Swarm *)dm->data; 46077048351SPatrick Sanan DMSwarmDataField *gfield; 46194f7d2dcSDave May DM celldm; 46294f7d2dcSDave May PetscBool isDA, isPLEX; 46394f7d2dcSDave May Vec *vecs; 46494f7d2dcSDave May PetscInt f, nvecs; 46594f7d2dcSDave May PetscInt project_type = 0; 46694f7d2dcSDave May 4670e2ec84fSDave May PetscFunctionBegin; 46894f7d2dcSDave May DMSWARMPICVALID(dm); 4699566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 4709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nfields, &gfield)); 47194f7d2dcSDave May nvecs = 0; 47294f7d2dcSDave May for (f = 0; f < nfields; f++) { 4739566063dSJacob Faibussowitsch PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldnames[f], &gfield[f])); 4745f80ce2aSJacob 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"); 4755f80ce2aSJacob Faibussowitsch PetscCheck(gfield[f]->bs == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Projection only valid for fields with block size = 1"); 47694f7d2dcSDave May nvecs += gfield[f]->bs; 47794f7d2dcSDave May } 47894f7d2dcSDave May if (!reuse) { 4799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nvecs, &vecs)); 48094f7d2dcSDave May for (f = 0; f < nvecs; f++) { 4819566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(celldm, &vecs[f])); 4829566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)vecs[f], gfield[f]->name)); 48394f7d2dcSDave May } 48494f7d2dcSDave May } else { 48594f7d2dcSDave May vecs = *fields; 48694f7d2dcSDave May } 48794f7d2dcSDave May 4889566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA)); 4899566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX)); 49094f7d2dcSDave May if (isDA) { 4919566063dSJacob Faibussowitsch PetscCall(private_DMSwarmProjectFields_DA(dm, celldm, project_type, nfields, gfield, vecs)); 49294f7d2dcSDave May } else if (isPLEX) { 4939566063dSJacob Faibussowitsch PetscCall(private_DMSwarmProjectFields_PLEX(dm, celldm, project_type, nfields, gfield, vecs)); 49494f7d2dcSDave May } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX"); 49594f7d2dcSDave May 4969566063dSJacob Faibussowitsch PetscCall(PetscFree(gfield)); 4975f80ce2aSJacob Faibussowitsch if (!reuse) *fields = vecs; 4980e2ec84fSDave May PetscFunctionReturn(0); 4990e2ec84fSDave May } 5000e2ec84fSDave May 5010e2ec84fSDave May /*@C 502b799feefSDave May DMSwarmCreatePointPerCellCount - Count the number of points within all cells in the cell DM 5030e2ec84fSDave May 5040e2ec84fSDave May Not collective 5050e2ec84fSDave May 5060e2ec84fSDave May Input parameter: 5070e2ec84fSDave May . dm - the DMSwarm 5080e2ec84fSDave May 5090e2ec84fSDave May Output parameters: 5100e2ec84fSDave May + ncells - the number of cells in the cell DM (optional argument, pass NULL to ignore) 511b799feefSDave May - count - array of length ncells containing the number of points per cell 5120e2ec84fSDave May 5130e2ec84fSDave May Level: beginner 5140e2ec84fSDave May 5150e2ec84fSDave May Notes: 5160e2ec84fSDave May The array count is allocated internally and must be free'd by the user. 5170e2ec84fSDave May 518db781477SPatrick Sanan .seealso: `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType` 5190e2ec84fSDave May @*/ 5209371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMSwarmCreatePointPerCellCount(DM dm, PetscInt *ncells, PetscInt **count) { 521b799feefSDave May PetscBool isvalid; 522b799feefSDave May PetscInt nel; 523b799feefSDave May PetscInt *sum; 524b799feefSDave May 5250e2ec84fSDave May PetscFunctionBegin; 5269566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetIsValid(dm, &isvalid)); 527b799feefSDave May nel = 0; 528b799feefSDave May if (isvalid) { 529b799feefSDave May PetscInt e; 530b799feefSDave May 5319566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetSizes(dm, &nel, NULL)); 532b799feefSDave May 5339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nel, &sum)); 534*48a46eb9SPierre Jolivet for (e = 0; e < nel; e++) PetscCall(DMSwarmSortGetNumberOfPointsPerCell(dm, e, &sum[e])); 535b799feefSDave May } else { 536b799feefSDave May DM celldm; 537b799feefSDave May PetscBool isda, isplex, isshell; 538b799feefSDave May PetscInt p, npoints; 539b799feefSDave May PetscInt *swarm_cellid; 540b799feefSDave May 541b799feefSDave May /* get the number of cells */ 5429566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(dm, &celldm)); 5439566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isda)); 5449566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isplex)); 5459566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMSHELL, &isshell)); 546b799feefSDave May if (isda) { 547b799feefSDave May PetscInt _nel, _npe; 548b799feefSDave May const PetscInt *_element; 549b799feefSDave May 5509566063dSJacob Faibussowitsch PetscCall(DMDAGetElements(celldm, &_nel, &_npe, &_element)); 551b799feefSDave May nel = _nel; 5529566063dSJacob Faibussowitsch PetscCall(DMDARestoreElements(celldm, &_nel, &_npe, &_element)); 553b799feefSDave May } else if (isplex) { 554b799feefSDave May PetscInt ps, pe; 555b799feefSDave May 5569566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(celldm, 0, &ps, &pe)); 557b799feefSDave May nel = pe - ps; 558b799feefSDave May } else if (isshell) { 559b799feefSDave May PetscErrorCode (*method_DMShellGetNumberOfCells)(DM, PetscInt *); 560b799feefSDave May 5619566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)celldm, "DMGetNumberOfCells_C", &method_DMShellGetNumberOfCells)); 562b799feefSDave May if (method_DMShellGetNumberOfCells) { 5639566063dSJacob Faibussowitsch PetscCall(method_DMShellGetNumberOfCells(celldm, &nel)); 5649371c9d4SSatish Balay } else 5659371c9d4SSatish 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);"); 566b799feefSDave 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"); 567b799feefSDave May 5689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nel, &sum)); 5699566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(sum, nel)); 5709566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dm, &npoints)); 5719566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 572b799feefSDave May for (p = 0; p < npoints; p++) { 5739371c9d4SSatish Balay if (swarm_cellid[p] != DMLOCATEPOINT_POINT_NOT_FOUND) { sum[swarm_cellid[p]]++; } 574b799feefSDave May } 5759566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 576b799feefSDave May } 577b799feefSDave May if (ncells) { *ncells = nel; } 578b799feefSDave May *count = sum; 5790e2ec84fSDave May PetscFunctionReturn(0); 5800e2ec84fSDave May } 58135b38c60SMatthew G. Knepley 58235b38c60SMatthew G. Knepley /*@ 58335b38c60SMatthew G. Knepley DMSwarmGetNumSpecies - Get the number of particle species 58435b38c60SMatthew G. Knepley 58535b38c60SMatthew G. Knepley Not collective 58635b38c60SMatthew G. Knepley 58735b38c60SMatthew G. Knepley Input parameter: 58835b38c60SMatthew G. Knepley . dm - the DMSwarm 58935b38c60SMatthew G. Knepley 59035b38c60SMatthew G. Knepley Output parameters: 59135b38c60SMatthew G. Knepley . Ns - the number of species 59235b38c60SMatthew G. Knepley 59335b38c60SMatthew G. Knepley Level: intermediate 59435b38c60SMatthew G. Knepley 595db781477SPatrick Sanan .seealso: `DMSwarmSetNumSpecies()`, `DMSwarmSetType()`, `DMSwarmType` 59635b38c60SMatthew G. Knepley @*/ 5979371c9d4SSatish Balay PetscErrorCode DMSwarmGetNumSpecies(DM sw, PetscInt *Ns) { 59835b38c60SMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 59935b38c60SMatthew G. Knepley 60035b38c60SMatthew G. Knepley PetscFunctionBegin; 60135b38c60SMatthew G. Knepley *Ns = swarm->Ns; 60235b38c60SMatthew G. Knepley PetscFunctionReturn(0); 60335b38c60SMatthew G. Knepley } 60435b38c60SMatthew G. Knepley 60535b38c60SMatthew G. Knepley /*@ 60635b38c60SMatthew G. Knepley DMSwarmSetNumSpecies - Set the number of particle species 60735b38c60SMatthew G. Knepley 60835b38c60SMatthew G. Knepley Not collective 60935b38c60SMatthew G. Knepley 61035b38c60SMatthew G. Knepley Input parameter: 61135b38c60SMatthew G. Knepley + dm - the DMSwarm 61235b38c60SMatthew G. Knepley - Ns - the number of species 61335b38c60SMatthew G. Knepley 61435b38c60SMatthew G. Knepley Level: intermediate 61535b38c60SMatthew G. Knepley 616db781477SPatrick Sanan .seealso: `DMSwarmGetNumSpecies()`, `DMSwarmSetType()`, `DMSwarmType` 61735b38c60SMatthew G. Knepley @*/ 6189371c9d4SSatish Balay PetscErrorCode DMSwarmSetNumSpecies(DM sw, PetscInt Ns) { 61935b38c60SMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 62035b38c60SMatthew G. Knepley 62135b38c60SMatthew G. Knepley PetscFunctionBegin; 62235b38c60SMatthew G. Knepley swarm->Ns = Ns; 62335b38c60SMatthew G. Knepley PetscFunctionReturn(0); 62435b38c60SMatthew G. Knepley } 62535b38c60SMatthew G. Knepley 62635b38c60SMatthew G. Knepley /*@C 6276c5a79ebSMatthew G. Knepley DMSwarmGetCoordinateFunction - Get the function setting initial particle positions, if it exists 6286c5a79ebSMatthew G. Knepley 6296c5a79ebSMatthew G. Knepley Not collective 6306c5a79ebSMatthew G. Knepley 6316c5a79ebSMatthew G. Knepley Input parameter: 6326c5a79ebSMatthew G. Knepley . dm - the DMSwarm 6336c5a79ebSMatthew G. Knepley 6346c5a79ebSMatthew G. Knepley Output Parameter: 6356c5a79ebSMatthew G. Knepley . coordFunc - the function setting initial particle positions, or NULL 6366c5a79ebSMatthew G. Knepley 6376c5a79ebSMatthew G. Knepley Level: intermediate 6386c5a79ebSMatthew G. Knepley 639db781477SPatrick Sanan .seealso: `DMSwarmSetCoordinateFunction()`, `DMSwarmGetVelocityFunction()`, `DMSwarmInitializeCoordinates()` 6406c5a79ebSMatthew G. Knepley @*/ 6419371c9d4SSatish Balay PetscErrorCode DMSwarmGetCoordinateFunction(DM sw, PetscSimplePointFunc *coordFunc) { 6426c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 6436c5a79ebSMatthew G. Knepley 6446c5a79ebSMatthew G. Knepley PetscFunctionBegin; 6456c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 6466c5a79ebSMatthew G. Knepley PetscValidPointer(coordFunc, 2); 6476c5a79ebSMatthew G. Knepley *coordFunc = swarm->coordFunc; 6486c5a79ebSMatthew G. Knepley PetscFunctionReturn(0); 6496c5a79ebSMatthew G. Knepley } 6506c5a79ebSMatthew G. Knepley 6516c5a79ebSMatthew G. Knepley /*@C 6526c5a79ebSMatthew G. Knepley DMSwarmSetCoordinateFunction - Set the function setting initial particle positions 6536c5a79ebSMatthew G. Knepley 6546c5a79ebSMatthew G. Knepley Not collective 6556c5a79ebSMatthew G. Knepley 6566c5a79ebSMatthew G. Knepley Input parameters: 6576c5a79ebSMatthew G. Knepley + dm - the DMSwarm 6586c5a79ebSMatthew G. Knepley - coordFunc - the function setting initial particle positions 6596c5a79ebSMatthew G. Knepley 6606c5a79ebSMatthew G. Knepley Level: intermediate 6616c5a79ebSMatthew G. Knepley 662db781477SPatrick Sanan .seealso: `DMSwarmGetCoordinateFunction()`, `DMSwarmSetVelocityFunction()`, `DMSwarmInitializeCoordinates()` 6636c5a79ebSMatthew G. Knepley @*/ 6649371c9d4SSatish Balay PetscErrorCode DMSwarmSetCoordinateFunction(DM sw, PetscSimplePointFunc coordFunc) { 6656c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 6666c5a79ebSMatthew G. Knepley 6676c5a79ebSMatthew G. Knepley PetscFunctionBegin; 6686c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 6696c5a79ebSMatthew G. Knepley PetscValidFunction(coordFunc, 2); 6706c5a79ebSMatthew G. Knepley swarm->coordFunc = coordFunc; 6716c5a79ebSMatthew G. Knepley PetscFunctionReturn(0); 6726c5a79ebSMatthew G. Knepley } 6736c5a79ebSMatthew G. Knepley 6746c5a79ebSMatthew G. Knepley /*@C 6756c5a79ebSMatthew G. Knepley DMSwarmGetCoordinateFunction - Get the function setting initial particle velocities, if it exists 6766c5a79ebSMatthew G. Knepley 6776c5a79ebSMatthew G. Knepley Not collective 6786c5a79ebSMatthew G. Knepley 6796c5a79ebSMatthew G. Knepley Input parameter: 6806c5a79ebSMatthew G. Knepley . dm - the DMSwarm 6816c5a79ebSMatthew G. Knepley 6826c5a79ebSMatthew G. Knepley Output Parameter: 6836c5a79ebSMatthew G. Knepley . velFunc - the function setting initial particle velocities, or NULL 6846c5a79ebSMatthew G. Knepley 6856c5a79ebSMatthew G. Knepley Level: intermediate 6866c5a79ebSMatthew G. Knepley 687db781477SPatrick Sanan .seealso: `DMSwarmSetVelocityFunction()`, `DMSwarmGetCoordinateFunction()`, `DMSwarmInitializeVelocities()` 6886c5a79ebSMatthew G. Knepley @*/ 6899371c9d4SSatish Balay PetscErrorCode DMSwarmGetVelocityFunction(DM sw, PetscSimplePointFunc *velFunc) { 6906c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 6916c5a79ebSMatthew G. Knepley 6926c5a79ebSMatthew G. Knepley PetscFunctionBegin; 6936c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 6946c5a79ebSMatthew G. Knepley PetscValidPointer(velFunc, 2); 6956c5a79ebSMatthew G. Knepley *velFunc = swarm->velFunc; 6966c5a79ebSMatthew G. Knepley PetscFunctionReturn(0); 6976c5a79ebSMatthew G. Knepley } 6986c5a79ebSMatthew G. Knepley 6996c5a79ebSMatthew G. Knepley /*@C 7006c5a79ebSMatthew G. Knepley DMSwarmSetVelocityFunction - Set the function setting initial particle velocities 7016c5a79ebSMatthew G. Knepley 7026c5a79ebSMatthew G. Knepley Not collective 7036c5a79ebSMatthew G. Knepley 7046c5a79ebSMatthew G. Knepley Input parameters: 7056c5a79ebSMatthew G. Knepley + dm - the DMSwarm 7066c5a79ebSMatthew G. Knepley - coordFunc - the function setting initial particle velocities 7076c5a79ebSMatthew G. Knepley 7086c5a79ebSMatthew G. Knepley Level: intermediate 7096c5a79ebSMatthew G. Knepley 710db781477SPatrick Sanan .seealso: `DMSwarmGetVelocityFunction()`, `DMSwarmSetCoordinateFunction()`, `DMSwarmInitializeVelocities()` 7116c5a79ebSMatthew G. Knepley @*/ 7129371c9d4SSatish Balay PetscErrorCode DMSwarmSetVelocityFunction(DM sw, PetscSimplePointFunc velFunc) { 7136c5a79ebSMatthew G. Knepley DM_Swarm *swarm = (DM_Swarm *)sw->data; 7146c5a79ebSMatthew G. Knepley 7156c5a79ebSMatthew G. Knepley PetscFunctionBegin; 7166c5a79ebSMatthew G. Knepley PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 7176c5a79ebSMatthew G. Knepley PetscValidFunction(velFunc, 2); 7186c5a79ebSMatthew G. Knepley swarm->velFunc = velFunc; 7196c5a79ebSMatthew G. Knepley PetscFunctionReturn(0); 7206c5a79ebSMatthew G. Knepley } 7216c5a79ebSMatthew G. Knepley 7226c5a79ebSMatthew G. Knepley /*@C 72335b38c60SMatthew G. Knepley DMSwarmComputeLocalSize - Compute the local number and distribution of particles based upon a density function 72435b38c60SMatthew G. Knepley 72535b38c60SMatthew G. Knepley Not collective 72635b38c60SMatthew G. Knepley 72735b38c60SMatthew G. Knepley Input Parameters: 72835b38c60SMatthew G. Knepley + sw - The DMSwarm 72935b38c60SMatthew G. Knepley . N - The target number of particles 73035b38c60SMatthew G. Knepley - density - The density field for the particle layout, normalized to unity 73135b38c60SMatthew G. Knepley 73235b38c60SMatthew G. Knepley Note: One particle will be created for each species. 73335b38c60SMatthew G. Knepley 73435b38c60SMatthew G. Knepley Level: advanced 73535b38c60SMatthew G. Knepley 736db781477SPatrick Sanan .seealso: `DMSwarmComputeLocalSizeFromOptions()` 73735b38c60SMatthew G. Knepley @*/ 7389371c9d4SSatish Balay PetscErrorCode DMSwarmComputeLocalSize(DM sw, PetscInt N, PetscProbFunc density) { 73935b38c60SMatthew G. Knepley DM dm; 74035b38c60SMatthew G. Knepley PetscQuadrature quad; 74135b38c60SMatthew G. Knepley const PetscReal *xq, *wq; 74235b38c60SMatthew G. Knepley PetscInt *npc, *cellid; 7436c5a79ebSMatthew G. Knepley PetscReal xi0[3]; 74435b38c60SMatthew G. Knepley PetscInt Ns, cStart, cEnd, c, dim, d, Nq, q, Np = 0, p; 74535b38c60SMatthew G. Knepley PetscBool simplex; 74635b38c60SMatthew G. Knepley 74735b38c60SMatthew G. Knepley PetscFunctionBegin; 7489566063dSJacob Faibussowitsch PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 7499566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sw, &dm)); 7509566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 7519566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 7529566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 7536858538eSMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(dm)); 7549566063dSJacob Faibussowitsch if (simplex) PetscCall(PetscDTStroudConicalQuadrature(dim, 1, 5, -1.0, 1.0, &quad)); 7559566063dSJacob Faibussowitsch else PetscCall(PetscDTGaussTensorQuadrature(dim, 1, 5, -1.0, 1.0, &quad)); 7569566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &xq, &wq)); 7579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(cEnd - cStart, &npc)); 75835b38c60SMatthew G. Knepley /* Integrate the density function to get the number of particles in each cell */ 75935b38c60SMatthew G. Knepley for (d = 0; d < dim; ++d) xi0[d] = -1.0; 76035b38c60SMatthew G. Knepley for (c = 0; c < cEnd - cStart; ++c) { 76135b38c60SMatthew G. Knepley const PetscInt cell = c + cStart; 76235b38c60SMatthew G. Knepley PetscReal v0[3], J[9], invJ[9], detJ; 76335b38c60SMatthew G. Knepley PetscReal n_int = 0.; 76435b38c60SMatthew G. Knepley 7659566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, NULL, v0, J, invJ, &detJ)); 76635b38c60SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 76735b38c60SMatthew G. Knepley PetscReal xr[3], den[3]; 76835b38c60SMatthew G. Knepley 76935b38c60SMatthew G. Knepley CoordinatesRefToReal(dim, dim, xi0, v0, J, &xq[q * dim], xr); 7706c5a79ebSMatthew G. Knepley PetscCall(density(xr, NULL, den)); 7716c5a79ebSMatthew G. Knepley n_int += den[0] * wq[q]; 77235b38c60SMatthew G. Knepley } 7736c5a79ebSMatthew G. Knepley npc[c] = (PetscInt)(N * n_int); 77435b38c60SMatthew G. Knepley npc[c] *= Ns; 77535b38c60SMatthew G. Knepley Np += npc[c]; 77635b38c60SMatthew G. Knepley } 7779566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quad)); 7789566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(sw, Np, 0)); 77935b38c60SMatthew G. Knepley 7809566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 78135b38c60SMatthew G. Knepley for (c = 0, p = 0; c < cEnd - cStart; ++c) { 78235b38c60SMatthew G. Knepley for (q = 0; q < npc[c]; ++q, ++p) cellid[p] = c; 78335b38c60SMatthew G. Knepley } 7849566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 7859566063dSJacob Faibussowitsch PetscCall(PetscFree(npc)); 78635b38c60SMatthew G. Knepley PetscFunctionReturn(0); 78735b38c60SMatthew G. Knepley } 78835b38c60SMatthew G. Knepley 78935b38c60SMatthew G. Knepley /*@ 79035b38c60SMatthew G. Knepley DMSwarmComputeLocalSizeFromOptions - Compute the local number and distribution of particles based upon a density function determined by options 79135b38c60SMatthew G. Knepley 79235b38c60SMatthew G. Knepley Not collective 79335b38c60SMatthew G. Knepley 79435b38c60SMatthew G. Knepley Input Parameters: 79535b38c60SMatthew G. Knepley , sw - The DMSwarm 79635b38c60SMatthew G. Knepley 79735b38c60SMatthew G. Knepley Level: advanced 79835b38c60SMatthew G. Knepley 799db781477SPatrick Sanan .seealso: `DMSwarmComputeLocalSize()` 80035b38c60SMatthew G. Knepley @*/ 8019371c9d4SSatish Balay PetscErrorCode DMSwarmComputeLocalSizeFromOptions(DM sw) { 80235b38c60SMatthew G. Knepley PetscProbFunc pdf; 80335b38c60SMatthew G. Knepley const char *prefix; 8046c5a79ebSMatthew G. Knepley char funcname[PETSC_MAX_PATH_LEN]; 8056c5a79ebSMatthew G. Knepley PetscInt *N, Ns, dim, n; 8066c5a79ebSMatthew G. Knepley PetscBool flg; 8076c5a79ebSMatthew G. Knepley PetscMPIInt size, rank; 80835b38c60SMatthew G. Knepley 80935b38c60SMatthew G. Knepley PetscFunctionBegin; 8106c5a79ebSMatthew G. Knepley PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sw), &size)); 8116c5a79ebSMatthew G. Knepley PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sw), &rank)); 8126c5a79ebSMatthew G. Knepley PetscCall(PetscCalloc1(size, &N)); 813d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)sw), "", "DMSwarm Options", "DMSWARM"); 8146c5a79ebSMatthew G. Knepley n = size; 8156c5a79ebSMatthew G. Knepley PetscCall(PetscOptionsIntArray("-dm_swarm_num_particles", "The target number of particles", "", N, &n, NULL)); 8166c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 8179566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-dm_swarm_num_species", "The number of species", "DMSwarmSetNumSpecies", Ns, &Ns, &flg)); 8189566063dSJacob Faibussowitsch if (flg) PetscCall(DMSwarmSetNumSpecies(sw, Ns)); 8196c5a79ebSMatthew G. Knepley PetscCall(PetscOptionsString("-dm_swarm_coordinate_function", "Function to determine particle coordinates", "DMSwarmSetCoordinateFunction", funcname, funcname, sizeof(funcname), &flg)); 820d0609cedSBarry Smith PetscOptionsEnd(); 8216c5a79ebSMatthew G. Knepley if (flg) { 8226c5a79ebSMatthew G. Knepley PetscSimplePointFunc coordFunc; 82335b38c60SMatthew G. Knepley 8246c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 8256c5a79ebSMatthew G. Knepley PetscCall(PetscDLSym(NULL, funcname, (void **)&coordFunc)); 8266c5a79ebSMatthew G. Knepley PetscCheck(coordFunc, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname); 8276c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 8286c5a79ebSMatthew G. Knepley PetscCall(DMSwarmSetLocalSizes(sw, N[rank] * Ns, 0)); 8296c5a79ebSMatthew G. Knepley PetscCall(DMSwarmSetCoordinateFunction(sw, coordFunc)); 8306c5a79ebSMatthew G. Knepley } else { 8319566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 8329566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sw, &prefix)); 8339566063dSJacob Faibussowitsch PetscCall(PetscProbCreateFromOptions(dim, prefix, "-dm_swarm_coordinate_density", &pdf, NULL, NULL)); 8346c5a79ebSMatthew G. Knepley PetscCall(DMSwarmComputeLocalSize(sw, N[rank], pdf)); 8356c5a79ebSMatthew G. Knepley } 8366c5a79ebSMatthew G. Knepley PetscCall(PetscFree(N)); 83735b38c60SMatthew G. Knepley PetscFunctionReturn(0); 83835b38c60SMatthew G. Knepley } 83935b38c60SMatthew G. Knepley 84035b38c60SMatthew G. Knepley /*@ 84135b38c60SMatthew G. Knepley DMSwarmInitializeCoordinates - Determine the initial coordinates of particles for a PIC method 84235b38c60SMatthew G. Knepley 84335b38c60SMatthew G. Knepley Not collective 84435b38c60SMatthew G. Knepley 84535b38c60SMatthew G. Knepley Input Parameters: 84635b38c60SMatthew G. Knepley , sw - The DMSwarm 84735b38c60SMatthew G. Knepley 84835b38c60SMatthew G. Knepley Note: Currently, we randomly place particles in their assigned cell 84935b38c60SMatthew G. Knepley 85035b38c60SMatthew G. Knepley Level: advanced 85135b38c60SMatthew G. Knepley 852db781477SPatrick Sanan .seealso: `DMSwarmComputeLocalSize()`, `DMSwarmInitializeVelocities()` 85335b38c60SMatthew G. Knepley @*/ 8549371c9d4SSatish Balay PetscErrorCode DMSwarmInitializeCoordinates(DM sw) { 8556c5a79ebSMatthew G. Knepley PetscSimplePointFunc coordFunc; 85635b38c60SMatthew G. Knepley PetscScalar *weight; 8576c5a79ebSMatthew G. Knepley PetscReal *x; 85835b38c60SMatthew G. Knepley PetscInt *species; 8596c5a79ebSMatthew G. Knepley void *ctx; 86035b38c60SMatthew G. Knepley PetscBool removePoints = PETSC_TRUE; 86135b38c60SMatthew G. Knepley PetscDataType dtype; 8626c5a79ebSMatthew G. Knepley PetscInt Np, p, Ns, dim, d, bs; 86335b38c60SMatthew G. Knepley 86435b38c60SMatthew G. Knepley PetscFunctionBeginUser; 8656c5a79ebSMatthew G. Knepley PetscCall(DMGetDimension(sw, &dim)); 8666c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetLocalSize(sw, &Np)); 8679566063dSJacob Faibussowitsch PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 8686c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetCoordinateFunction(sw, &coordFunc)); 8696c5a79ebSMatthew G. Knepley 8706c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, &bs, &dtype, (void **)&x)); 8716c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&weight)); 8726c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 8736c5a79ebSMatthew G. Knepley if (coordFunc) { 8746c5a79ebSMatthew G. Knepley PetscCall(DMGetApplicationContext(sw, &ctx)); 8756c5a79ebSMatthew G. Knepley for (p = 0; p < Np; ++p) { 8766c5a79ebSMatthew G. Knepley PetscScalar X[3]; 8776c5a79ebSMatthew G. Knepley 8786c5a79ebSMatthew G. Knepley PetscCall((*coordFunc)(dim, 0., NULL, p, X, ctx)); 8796c5a79ebSMatthew G. Knepley for (d = 0; d < dim; ++d) x[p * dim + d] = PetscRealPart(X[d]); 8806c5a79ebSMatthew G. Knepley weight[p] = 1.0; 8816c5a79ebSMatthew G. Knepley species[p] = p % Ns; 8826c5a79ebSMatthew G. Knepley } 8836c5a79ebSMatthew G. Knepley } else { 8846c5a79ebSMatthew G. Knepley DM dm; 8856c5a79ebSMatthew G. Knepley PetscRandom rnd; 8866c5a79ebSMatthew G. Knepley PetscReal xi0[3]; 8876c5a79ebSMatthew G. Knepley PetscInt cStart, cEnd, c; 8886c5a79ebSMatthew G. Knepley 8899566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sw, &dm)); 8909566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 89135b38c60SMatthew G. Knepley 89235b38c60SMatthew G. Knepley /* Set particle position randomly in cell, set weights to 1 */ 8939566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd)); 8949566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0)); 8959566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rnd)); 8969566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetAccess(sw)); 89735b38c60SMatthew G. Knepley for (d = 0; d < dim; ++d) xi0[d] = -1.0; 89835b38c60SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 89935b38c60SMatthew G. Knepley PetscReal v0[3], J[9], invJ[9], detJ; 90035b38c60SMatthew G. Knepley PetscInt *pidx, Npc, q; 90135b38c60SMatthew G. Knepley 9029566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx)); 9039566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); 90435b38c60SMatthew G. Knepley for (q = 0; q < Npc; ++q) { 90535b38c60SMatthew G. Knepley const PetscInt p = pidx[q]; 90635b38c60SMatthew G. Knepley PetscReal xref[3]; 90735b38c60SMatthew G. Knepley 9089566063dSJacob Faibussowitsch for (d = 0; d < dim; ++d) PetscCall(PetscRandomGetValueReal(rnd, &xref[d])); 90935b38c60SMatthew G. Knepley CoordinatesRefToReal(dim, dim, xi0, v0, J, xref, &x[p * dim]); 91035b38c60SMatthew G. Knepley 91135b38c60SMatthew G. Knepley weight[p] = 1.0; 9126c5a79ebSMatthew G. Knepley species[p] = p % Ns; 91335b38c60SMatthew G. Knepley } 9149566063dSJacob Faibussowitsch PetscCall(PetscFree(pidx)); 91535b38c60SMatthew G. Knepley } 9169566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rnd)); 9179566063dSJacob Faibussowitsch PetscCall(DMSwarmSortRestoreAccess(sw)); 9186c5a79ebSMatthew G. Knepley } 9199566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x)); 9206c5a79ebSMatthew G. Knepley PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight)); 9219566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 9226c5a79ebSMatthew G. Knepley 9239566063dSJacob Faibussowitsch PetscCall(DMSwarmMigrate(sw, removePoints)); 9249566063dSJacob Faibussowitsch PetscCall(DMLocalizeCoordinates(sw)); 92535b38c60SMatthew G. Knepley PetscFunctionReturn(0); 92635b38c60SMatthew G. Knepley } 92735b38c60SMatthew G. Knepley 92835b38c60SMatthew G. Knepley /*@C 92935b38c60SMatthew G. Knepley DMSwarmInitializeVelocities - Set the initial velocities of particles using a distribution. 93035b38c60SMatthew G. Knepley 93135b38c60SMatthew G. Knepley Collective on dm 93235b38c60SMatthew G. Knepley 93335b38c60SMatthew G. Knepley Input Parameters: 93435b38c60SMatthew G. Knepley + sw - The DMSwarm object 93535b38c60SMatthew G. Knepley . sampler - A function which uniformly samples the velocity PDF 93635b38c60SMatthew G. Knepley - v0 - The velocity scale for nondimensionalization for each species 93735b38c60SMatthew G. Knepley 9386c5a79ebSMatthew 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. 9396c5a79ebSMatthew G. Knepley 94035b38c60SMatthew G. Knepley Level: advanced 94135b38c60SMatthew G. Knepley 942db781477SPatrick Sanan .seealso: `DMSwarmComputeLocalSize()`, `DMSwarmInitializeCoordinates()`, `DMSwarmInitializeVelocitiesFromOptions()` 94335b38c60SMatthew G. Knepley @*/ 9449371c9d4SSatish Balay PetscErrorCode DMSwarmInitializeVelocities(DM sw, PetscProbFunc sampler, const PetscReal v0[]) { 9456c5a79ebSMatthew G. Knepley PetscSimplePointFunc velFunc; 94635b38c60SMatthew G. Knepley PetscReal *v; 94735b38c60SMatthew G. Knepley PetscInt *species; 9486c5a79ebSMatthew G. Knepley void *ctx; 94935b38c60SMatthew G. Knepley PetscInt dim, Np, p; 95035b38c60SMatthew G. Knepley 95135b38c60SMatthew G. Knepley PetscFunctionBegin; 9526c5a79ebSMatthew G. Knepley PetscCall(DMSwarmGetVelocityFunction(sw, &velFunc)); 95335b38c60SMatthew G. Knepley 9549566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 9559566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(sw, &Np)); 9569566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&v)); 9579566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 9586c5a79ebSMatthew G. Knepley if (v0[0] == 0.) { 9596c5a79ebSMatthew G. Knepley PetscCall(PetscArrayzero(v, Np * dim)); 9606c5a79ebSMatthew G. Knepley } else if (velFunc) { 9616c5a79ebSMatthew G. Knepley PetscCall(DMGetApplicationContext(sw, &ctx)); 9626c5a79ebSMatthew G. Knepley for (p = 0; p < Np; ++p) { 9636c5a79ebSMatthew G. Knepley PetscInt s = species[p], d; 9646c5a79ebSMatthew G. Knepley PetscScalar vel[3]; 9656c5a79ebSMatthew G. Knepley 9666c5a79ebSMatthew G. Knepley PetscCall((*velFunc)(dim, 0., NULL, p, vel, ctx)); 9676c5a79ebSMatthew G. Knepley for (d = 0; d < dim; ++d) v[p * dim + d] = (v0[s] / v0[0]) * PetscRealPart(vel[d]); 9686c5a79ebSMatthew G. Knepley } 9696c5a79ebSMatthew G. Knepley } else { 9706c5a79ebSMatthew G. Knepley PetscRandom rnd; 9716c5a79ebSMatthew G. Knepley 9726c5a79ebSMatthew G. Knepley PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rnd)); 9736c5a79ebSMatthew G. Knepley PetscCall(PetscRandomSetInterval(rnd, 0, 1.)); 9746c5a79ebSMatthew G. Knepley PetscCall(PetscRandomSetFromOptions(rnd)); 9756c5a79ebSMatthew G. Knepley 97635b38c60SMatthew G. Knepley for (p = 0; p < Np; ++p) { 97735b38c60SMatthew G. Knepley PetscInt s = species[p], d; 97835b38c60SMatthew G. Knepley PetscReal a[3], vel[3]; 97935b38c60SMatthew G. Knepley 9809566063dSJacob Faibussowitsch for (d = 0; d < dim; ++d) PetscCall(PetscRandomGetValueReal(rnd, &a[d])); 9819566063dSJacob Faibussowitsch PetscCall(sampler(a, NULL, vel)); 98235b38c60SMatthew G. Knepley for (d = 0; d < dim; ++d) { v[p * dim + d] = (v0[s] / v0[0]) * vel[d]; } 98335b38c60SMatthew G. Knepley } 9846c5a79ebSMatthew G. Knepley PetscCall(PetscRandomDestroy(&rnd)); 9856c5a79ebSMatthew G. Knepley } 9869566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&v)); 9879566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 98835b38c60SMatthew G. Knepley PetscFunctionReturn(0); 98935b38c60SMatthew G. Knepley } 99035b38c60SMatthew G. Knepley 99135b38c60SMatthew G. Knepley /*@ 99235b38c60SMatthew G. Knepley DMSwarmInitializeVelocitiesFromOptions - Set the initial velocities of particles using a distribution determined from options. 99335b38c60SMatthew G. Knepley 99435b38c60SMatthew G. Knepley Collective on dm 99535b38c60SMatthew G. Knepley 99635b38c60SMatthew G. Knepley Input Parameters: 99735b38c60SMatthew G. Knepley + sw - The DMSwarm object 99835b38c60SMatthew G. Knepley - v0 - The velocity scale for nondimensionalization for each species 99935b38c60SMatthew G. Knepley 100035b38c60SMatthew G. Knepley Level: advanced 100135b38c60SMatthew G. Knepley 1002db781477SPatrick Sanan .seealso: `DMSwarmComputeLocalSize()`, `DMSwarmInitializeCoordinates()`, `DMSwarmInitializeVelocities()` 100335b38c60SMatthew G. Knepley @*/ 10049371c9d4SSatish Balay PetscErrorCode DMSwarmInitializeVelocitiesFromOptions(DM sw, const PetscReal v0[]) { 100535b38c60SMatthew G. Knepley PetscProbFunc sampler; 100635b38c60SMatthew G. Knepley PetscInt dim; 100735b38c60SMatthew G. Knepley const char *prefix; 10086c5a79ebSMatthew G. Knepley char funcname[PETSC_MAX_PATH_LEN]; 10096c5a79ebSMatthew G. Knepley PetscBool flg; 101035b38c60SMatthew G. Knepley 101135b38c60SMatthew G. Knepley PetscFunctionBegin; 1012d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)sw), "", "DMSwarm Options", "DMSWARM"); 10136c5a79ebSMatthew G. Knepley PetscCall(PetscOptionsString("-dm_swarm_velocity_function", "Function to determine particle velocities", "DMSwarmSetVelocityFunction", funcname, funcname, sizeof(funcname), &flg)); 1014d0609cedSBarry Smith PetscOptionsEnd(); 10156c5a79ebSMatthew G. Knepley if (flg) { 10166c5a79ebSMatthew G. Knepley PetscSimplePointFunc velFunc; 10176c5a79ebSMatthew G. Knepley 10186c5a79ebSMatthew G. Knepley PetscCall(PetscDLSym(NULL, funcname, (void **)&velFunc)); 10196c5a79ebSMatthew G. Knepley PetscCheck(velFunc, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname); 10206c5a79ebSMatthew G. Knepley PetscCall(DMSwarmSetVelocityFunction(sw, velFunc)); 10216c5a79ebSMatthew G. Knepley } 10229566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 10239566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sw, &prefix)); 10249566063dSJacob Faibussowitsch PetscCall(PetscProbCreateFromOptions(dim, prefix, "-dm_swarm_velocity_density", NULL, NULL, &sampler)); 10259566063dSJacob Faibussowitsch PetscCall(DMSwarmInitializeVelocities(sw, sampler, v0)); 102635b38c60SMatthew G. Knepley PetscFunctionReturn(0); 102735b38c60SMatthew G. Knepley } 1028