1 #define PETSCDM_DLL 2 #include <petsc/private/dmswarmimpl.h> /*I "petscdmswarm.h" I*/ 3 #include <petscsf.h> 4 #include <petscdmda.h> 5 #include <petscdmplex.h> 6 #include <petscdt.h> 7 #include "../src/dm/impls/swarm/data_bucket.h" 8 9 #include <petsc/private/petscfeimpl.h> /* For CoordinatesRefToReal() */ 10 11 /* 12 Error checking to ensure the swarm type is correct and that a cell DM has been set 13 */ 14 #define DMSWARMPICVALID(dm) \ 15 do { \ 16 DM_Swarm *_swarm = (DM_Swarm *)(dm)->data; \ 17 PetscCheck(_swarm->swarm_type == DMSWARM_PIC, PetscObjectComm((PetscObject)(dm)), PETSC_ERR_SUP, "Valid only for DMSwarm-PIC. You must call DMSwarmSetType(dm,DMSWARM_PIC)"); \ 18 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)"); \ 19 } while (0) 20 21 /* Coordinate insertition/addition API */ 22 /*@C 23 DMSwarmSetPointsUniformCoordinates - Set point coordinates in a DMSwarm on a regular (ijk) grid 24 25 Collective on dm 26 27 Input parameters: 28 + dm - the DMSwarm 29 . min - minimum coordinate values in the x, y, z directions (array of length dim) 30 . max - maximum coordinate values in the x, y, z directions (array of length dim) 31 . npoints - number of points in each spatial direction (array of length dim) 32 - mode - indicates whether to append points to the swarm (ADD_VALUES), or over-ride existing points (INSERT_VALUES) 33 34 Level: beginner 35 36 Notes: 37 When using mode = INSERT_VALUES, this method will reset the number of particles in the DMSwarm 38 to be npoints[0]*npoints[1] (2D) or npoints[0]*npoints[1]*npoints[2] (3D). When using mode = ADD_VALUES, 39 new points will be appended to any already existing in the DMSwarm 40 41 .seealso: `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType` 42 @*/ 43 PETSC_EXTERN PetscErrorCode DMSwarmSetPointsUniformCoordinates(DM dm, PetscReal min[], PetscReal max[], PetscInt npoints[], InsertMode mode) { 44 PetscReal gmin[] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL}; 45 PetscReal gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 46 PetscInt i, j, k, N, bs, b, n_estimate, n_curr, n_new_est, p, n_found; 47 Vec coorlocal; 48 const PetscScalar *_coor; 49 DM celldm; 50 PetscReal dx[3]; 51 PetscInt _npoints[] = {0, 0, 1}; 52 Vec pos; 53 PetscScalar *_pos; 54 PetscReal *swarm_coor; 55 PetscInt *swarm_cellid; 56 PetscSF sfcell = NULL; 57 const PetscSFNode *LA_sfcell; 58 59 PetscFunctionBegin; 60 DMSWARMPICVALID(dm); 61 PetscCall(DMSwarmGetCellDM(dm, &celldm)); 62 PetscCall(DMGetCoordinatesLocal(celldm, &coorlocal)); 63 PetscCall(VecGetSize(coorlocal, &N)); 64 PetscCall(VecGetBlockSize(coorlocal, &bs)); 65 N = N / bs; 66 PetscCall(VecGetArrayRead(coorlocal, &_coor)); 67 for (i = 0; i < N; i++) { 68 for (b = 0; b < bs; b++) { 69 gmin[b] = PetscMin(gmin[b], PetscRealPart(_coor[bs * i + b])); 70 gmax[b] = PetscMax(gmax[b], PetscRealPart(_coor[bs * i + b])); 71 } 72 } 73 PetscCall(VecRestoreArrayRead(coorlocal, &_coor)); 74 75 for (b = 0; b < bs; b++) { 76 if (npoints[b] > 1) { 77 dx[b] = (max[b] - min[b]) / ((PetscReal)(npoints[b] - 1)); 78 } else { 79 dx[b] = 0.0; 80 } 81 _npoints[b] = npoints[b]; 82 } 83 84 /* determine number of points living in the bounding box */ 85 n_estimate = 0; 86 for (k = 0; k < _npoints[2]; k++) { 87 for (j = 0; j < _npoints[1]; j++) { 88 for (i = 0; i < _npoints[0]; i++) { 89 PetscReal xp[] = {0.0, 0.0, 0.0}; 90 PetscInt ijk[3]; 91 PetscBool point_inside = PETSC_TRUE; 92 93 ijk[0] = i; 94 ijk[1] = j; 95 ijk[2] = k; 96 for (b = 0; b < bs; b++) { xp[b] = min[b] + ijk[b] * dx[b]; } 97 for (b = 0; b < bs; b++) { 98 if (xp[b] < gmin[b]) { point_inside = PETSC_FALSE; } 99 if (xp[b] > gmax[b]) { point_inside = PETSC_FALSE; } 100 } 101 if (point_inside) { n_estimate++; } 102 } 103 } 104 } 105 106 /* create candidate list */ 107 PetscCall(VecCreate(PETSC_COMM_SELF, &pos)); 108 PetscCall(VecSetSizes(pos, bs * n_estimate, PETSC_DECIDE)); 109 PetscCall(VecSetBlockSize(pos, bs)); 110 PetscCall(VecSetFromOptions(pos)); 111 PetscCall(VecGetArray(pos, &_pos)); 112 113 n_estimate = 0; 114 for (k = 0; k < _npoints[2]; k++) { 115 for (j = 0; j < _npoints[1]; j++) { 116 for (i = 0; i < _npoints[0]; i++) { 117 PetscReal xp[] = {0.0, 0.0, 0.0}; 118 PetscInt ijk[3]; 119 PetscBool point_inside = PETSC_TRUE; 120 121 ijk[0] = i; 122 ijk[1] = j; 123 ijk[2] = k; 124 for (b = 0; b < bs; b++) { xp[b] = min[b] + ijk[b] * dx[b]; } 125 for (b = 0; b < bs; b++) { 126 if (xp[b] < gmin[b]) { point_inside = PETSC_FALSE; } 127 if (xp[b] > gmax[b]) { point_inside = PETSC_FALSE; } 128 } 129 if (point_inside) { 130 for (b = 0; b < bs; b++) { _pos[bs * n_estimate + b] = xp[b]; } 131 n_estimate++; 132 } 133 } 134 } 135 } 136 PetscCall(VecRestoreArray(pos, &_pos)); 137 138 /* locate points */ 139 PetscCall(DMLocatePoints(celldm, pos, DM_POINTLOCATION_NONE, &sfcell)); 140 PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 141 n_found = 0; 142 for (p = 0; p < n_estimate; p++) { 143 if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { n_found++; } 144 } 145 146 /* adjust size */ 147 if (mode == ADD_VALUES) { 148 PetscCall(DMSwarmGetLocalSize(dm, &n_curr)); 149 n_new_est = n_curr + n_found; 150 PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 151 } 152 if (mode == INSERT_VALUES) { 153 n_curr = 0; 154 n_new_est = n_found; 155 PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 156 } 157 158 /* initialize new coords, cell owners, pid */ 159 PetscCall(VecGetArrayRead(pos, &_coor)); 160 PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 161 PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 162 n_found = 0; 163 for (p = 0; p < n_estimate; p++) { 164 if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 165 for (b = 0; b < bs; b++) { swarm_coor[bs * (n_curr + n_found) + b] = PetscRealPart(_coor[bs * p + b]); } 166 swarm_cellid[n_curr + n_found] = LA_sfcell[p].index; 167 n_found++; 168 } 169 } 170 PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 171 PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 172 PetscCall(VecRestoreArrayRead(pos, &_coor)); 173 174 PetscCall(PetscSFDestroy(&sfcell)); 175 PetscCall(VecDestroy(&pos)); 176 PetscFunctionReturn(0); 177 } 178 179 /*@C 180 DMSwarmSetPointCoordinates - Set point coordinates in a DMSwarm from a user defined list 181 182 Collective on dm 183 184 Input parameters: 185 + dm - the DMSwarm 186 . npoints - the number of points to insert 187 . coor - the coordinate values 188 . 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 189 - mode - indicates whether to append points to the swarm (ADD_VALUES), or over-ride existing points (INSERT_VALUES) 190 191 Level: beginner 192 193 Notes: 194 If the user has specified redundant = PETSC_FALSE, the cell DM will attempt to locate the coordinates provided by coor[] within 195 its sub-domain. If they any values within coor[] are not located in the sub-domain, they will be ignored and will not get 196 added to the DMSwarm. 197 198 .seealso: `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType`, `DMSwarmSetPointsUniformCoordinates()` 199 @*/ 200 PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinates(DM dm, PetscInt npoints, PetscReal coor[], PetscBool redundant, InsertMode mode) { 201 PetscReal gmin[] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL}; 202 PetscReal gmax[] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL}; 203 PetscInt i, N, bs, b, n_estimate, n_curr, n_new_est, p, n_found; 204 Vec coorlocal; 205 const PetscScalar *_coor; 206 DM celldm; 207 Vec pos; 208 PetscScalar *_pos; 209 PetscReal *swarm_coor; 210 PetscInt *swarm_cellid; 211 PetscSF sfcell = NULL; 212 const PetscSFNode *LA_sfcell; 213 PetscReal *my_coor; 214 PetscInt my_npoints; 215 PetscMPIInt rank; 216 MPI_Comm comm; 217 218 PetscFunctionBegin; 219 DMSWARMPICVALID(dm); 220 PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 221 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 222 223 PetscCall(DMSwarmGetCellDM(dm, &celldm)); 224 PetscCall(DMGetCoordinatesLocal(celldm, &coorlocal)); 225 PetscCall(VecGetSize(coorlocal, &N)); 226 PetscCall(VecGetBlockSize(coorlocal, &bs)); 227 N = N / bs; 228 PetscCall(VecGetArrayRead(coorlocal, &_coor)); 229 for (i = 0; i < N; i++) { 230 for (b = 0; b < bs; b++) { 231 gmin[b] = PetscMin(gmin[b], PetscRealPart(_coor[bs * i + b])); 232 gmax[b] = PetscMax(gmax[b], PetscRealPart(_coor[bs * i + b])); 233 } 234 } 235 PetscCall(VecRestoreArrayRead(coorlocal, &_coor)); 236 237 /* broadcast points from rank 0 if requested */ 238 if (redundant) { 239 my_npoints = npoints; 240 PetscCallMPI(MPI_Bcast(&my_npoints, 1, MPIU_INT, 0, comm)); 241 242 if (rank > 0) { /* allocate space */ 243 PetscCall(PetscMalloc1(bs * my_npoints, &my_coor)); 244 } else { 245 my_coor = coor; 246 } 247 PetscCallMPI(MPI_Bcast(my_coor, bs * my_npoints, MPIU_REAL, 0, comm)); 248 } else { 249 my_npoints = npoints; 250 my_coor = coor; 251 } 252 253 /* determine the number of points living in the bounding box */ 254 n_estimate = 0; 255 for (i = 0; i < my_npoints; i++) { 256 PetscBool point_inside = PETSC_TRUE; 257 258 for (b = 0; b < bs; b++) { 259 if (my_coor[bs * i + b] < gmin[b]) { point_inside = PETSC_FALSE; } 260 if (my_coor[bs * i + b] > gmax[b]) { point_inside = PETSC_FALSE; } 261 } 262 if (point_inside) { n_estimate++; } 263 } 264 265 /* create candidate list */ 266 PetscCall(VecCreate(PETSC_COMM_SELF, &pos)); 267 PetscCall(VecSetSizes(pos, bs * n_estimate, PETSC_DECIDE)); 268 PetscCall(VecSetBlockSize(pos, bs)); 269 PetscCall(VecSetFromOptions(pos)); 270 PetscCall(VecGetArray(pos, &_pos)); 271 272 n_estimate = 0; 273 for (i = 0; i < my_npoints; i++) { 274 PetscBool point_inside = PETSC_TRUE; 275 276 for (b = 0; b < bs; b++) { 277 if (my_coor[bs * i + b] < gmin[b]) { point_inside = PETSC_FALSE; } 278 if (my_coor[bs * i + b] > gmax[b]) { point_inside = PETSC_FALSE; } 279 } 280 if (point_inside) { 281 for (b = 0; b < bs; b++) { _pos[bs * n_estimate + b] = my_coor[bs * i + b]; } 282 n_estimate++; 283 } 284 } 285 PetscCall(VecRestoreArray(pos, &_pos)); 286 287 /* locate points */ 288 PetscCall(DMLocatePoints(celldm, pos, DM_POINTLOCATION_NONE, &sfcell)); 289 290 PetscCall(PetscSFGetGraph(sfcell, NULL, NULL, NULL, &LA_sfcell)); 291 n_found = 0; 292 for (p = 0; p < n_estimate; p++) { 293 if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { n_found++; } 294 } 295 296 /* adjust size */ 297 if (mode == ADD_VALUES) { 298 PetscCall(DMSwarmGetLocalSize(dm, &n_curr)); 299 n_new_est = n_curr + n_found; 300 PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 301 } 302 if (mode == INSERT_VALUES) { 303 n_curr = 0; 304 n_new_est = n_found; 305 PetscCall(DMSwarmSetLocalSizes(dm, n_new_est, -1)); 306 } 307 308 /* initialize new coords, cell owners, pid */ 309 PetscCall(VecGetArrayRead(pos, &_coor)); 310 PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 311 PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 312 n_found = 0; 313 for (p = 0; p < n_estimate; p++) { 314 if (LA_sfcell[p].index != DMLOCATEPOINT_POINT_NOT_FOUND) { 315 for (b = 0; b < bs; b++) { swarm_coor[bs * (n_curr + n_found) + b] = PetscRealPart(_coor[bs * p + b]); } 316 swarm_cellid[n_curr + n_found] = LA_sfcell[p].index; 317 n_found++; 318 } 319 } 320 PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 321 PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 322 PetscCall(VecRestoreArrayRead(pos, &_coor)); 323 324 if (redundant) { 325 if (rank > 0) { PetscCall(PetscFree(my_coor)); } 326 } 327 PetscCall(PetscSFDestroy(&sfcell)); 328 PetscCall(VecDestroy(&pos)); 329 PetscFunctionReturn(0); 330 } 331 332 extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM, DM, DMSwarmPICLayoutType, PetscInt); 333 extern PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM, DM, DMSwarmPICLayoutType, PetscInt); 334 335 /*@C 336 DMSwarmInsertPointsUsingCellDM - Insert point coordinates within each cell 337 338 Not collective 339 340 Input parameters: 341 + dm - the DMSwarm 342 . layout_type - method used to fill each cell with the cell DM 343 - fill_param - parameter controlling how many points per cell are added (the meaning of this parameter is dependent on the layout type) 344 345 Level: beginner 346 347 Notes: 348 349 The insert method will reset any previous defined points within the DMSwarm. 350 351 When using a DMDA both 2D and 3D are supported for all layout types provided you are using DMDA_ELEMENT_Q1. 352 353 When using a DMPLEX the following case are supported: 354 (i) DMSWARMPIC_LAYOUT_REGULAR: 2D (triangle), 355 (ii) DMSWARMPIC_LAYOUT_GAUSS: 2D and 3D provided the cell is a tri/tet or a quad/hex, 356 (iii) DMSWARMPIC_LAYOUT_SUBDIVISION: 2D and 3D for quad/hex and 2D tri. 357 358 .seealso: `DMSwarmPICLayoutType`, `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType` 359 @*/ 360 PETSC_EXTERN PetscErrorCode DMSwarmInsertPointsUsingCellDM(DM dm, DMSwarmPICLayoutType layout_type, PetscInt fill_param) { 361 DM celldm; 362 PetscBool isDA, isPLEX; 363 364 PetscFunctionBegin; 365 DMSWARMPICVALID(dm); 366 PetscCall(DMSwarmGetCellDM(dm, &celldm)); 367 PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA)); 368 PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX)); 369 if (isDA) { 370 PetscCall(private_DMSwarmInsertPointsUsingCellDM_DA(dm, celldm, layout_type, fill_param)); 371 } else if (isPLEX) { 372 PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX(dm, celldm, layout_type, fill_param)); 373 } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX"); 374 PetscFunctionReturn(0); 375 } 376 377 extern PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM, DM, PetscInt, PetscReal *); 378 379 /*@C 380 DMSwarmSetPointCoordinatesCellwise - Insert point coordinates (defined over the reference cell) within each cell 381 382 Not collective 383 384 Input parameters: 385 + dm - the DMSwarm 386 . celldm - the cell DM 387 . npoints - the number of points to insert in each cell 388 - xi - the coordinates (defined in the local coordinate system for each cell) to insert 389 390 Level: beginner 391 392 Notes: 393 The method will reset any previous defined points within the DMSwarm. 394 Only supported for DMPLEX. If you are using a DMDA it is recommended to either use 395 DMSwarmInsertPointsUsingCellDM(), or extract and set the coordinates yourself the following code 396 397 $ PetscReal *coor; 398 $ DMSwarmGetField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&coor); 399 $ // user code to define the coordinates here 400 $ DMSwarmRestoreField(dm,DMSwarmPICField_coor,NULL,NULL,(void**)&coor); 401 402 .seealso: `DMSwarmSetCellDM()`, `DMSwarmInsertPointsUsingCellDM()` 403 @*/ 404 PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinatesCellwise(DM dm, PetscInt npoints, PetscReal xi[]) { 405 DM celldm; 406 PetscBool isDA, isPLEX; 407 408 PetscFunctionBegin; 409 DMSWARMPICVALID(dm); 410 PetscCall(DMSwarmGetCellDM(dm, &celldm)); 411 PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA)); 412 PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX)); 413 PetscCheck(!isDA, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMPLEX. Recommended you use DMSwarmInsertPointsUsingCellDM()"); 414 if (isPLEX) { 415 PetscCall(private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm, celldm, npoints, xi)); 416 } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX"); 417 PetscFunctionReturn(0); 418 } 419 420 /* Field projection API */ 421 extern PetscErrorCode private_DMSwarmProjectFields_DA(DM swarm, DM celldm, PetscInt project_type, PetscInt nfields, DMSwarmDataField dfield[], Vec vecs[]); 422 extern PetscErrorCode private_DMSwarmProjectFields_PLEX(DM swarm, DM celldm, PetscInt project_type, PetscInt nfields, DMSwarmDataField dfield[], Vec vecs[]); 423 424 /*@C 425 DMSwarmProjectFields - Project a set of swarm fields onto the cell DM 426 427 Collective on dm 428 429 Input parameters: 430 + dm - the DMSwarm 431 . nfields - the number of swarm fields to project 432 . fieldnames - the textual names of the swarm fields to project 433 . fields - an array of Vec's of length nfields 434 - reuse - flag indicating whether the array and contents of fields should be re-used or internally allocated 435 436 Currently, the only available projection method consists of 437 phi_i = \sum_{p=0}^{np} N_i(x_p) phi_p dJ / \sum_{p=0}^{np} N_i(x_p) dJ 438 where phi_p is the swarm field at point p, 439 N_i() is the cell DM basis function at vertex i, 440 dJ is the determinant of the cell Jacobian and 441 phi_i is the projected vertex value of the field phi. 442 443 Level: beginner 444 445 Notes: 446 447 If reuse = PETSC_FALSE, this function will allocate the array of Vec's, and each individual Vec. 448 The user is responsible for destroying both the array and the individual Vec objects. 449 450 Only swarm fields registered with data type = PETSC_REAL can be projected onto the cell DM. 451 452 Only swarm fields of block size = 1 can currently be projected. 453 454 The only projection methods currently only support the DA (2D) and PLEX (triangles 2D). 455 456 .seealso: `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType` 457 @*/ 458 PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM dm, PetscInt nfields, const char *fieldnames[], Vec **fields, PetscBool reuse) { 459 DM_Swarm *swarm = (DM_Swarm *)dm->data; 460 DMSwarmDataField *gfield; 461 DM celldm; 462 PetscBool isDA, isPLEX; 463 Vec *vecs; 464 PetscInt f, nvecs; 465 PetscInt project_type = 0; 466 467 PetscFunctionBegin; 468 DMSWARMPICVALID(dm); 469 PetscCall(DMSwarmGetCellDM(dm, &celldm)); 470 PetscCall(PetscMalloc1(nfields, &gfield)); 471 nvecs = 0; 472 for (f = 0; f < nfields; f++) { 473 PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldnames[f], &gfield[f])); 474 PetscCheck(gfield[f]->petsc_type == PETSC_REAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Projection only valid for fields using a data type = PETSC_REAL"); 475 PetscCheck(gfield[f]->bs == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Projection only valid for fields with block size = 1"); 476 nvecs += gfield[f]->bs; 477 } 478 if (!reuse) { 479 PetscCall(PetscMalloc1(nvecs, &vecs)); 480 for (f = 0; f < nvecs; f++) { 481 PetscCall(DMCreateGlobalVector(celldm, &vecs[f])); 482 PetscCall(PetscObjectSetName((PetscObject)vecs[f], gfield[f]->name)); 483 } 484 } else { 485 vecs = *fields; 486 } 487 488 PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isDA)); 489 PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isPLEX)); 490 if (isDA) { 491 PetscCall(private_DMSwarmProjectFields_DA(dm, celldm, project_type, nfields, gfield, vecs)); 492 } else if (isPLEX) { 493 PetscCall(private_DMSwarmProjectFields_PLEX(dm, celldm, project_type, nfields, gfield, vecs)); 494 } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only supported for cell DMs of type DMDA and DMPLEX"); 495 496 PetscCall(PetscFree(gfield)); 497 if (!reuse) *fields = vecs; 498 PetscFunctionReturn(0); 499 } 500 501 /*@C 502 DMSwarmCreatePointPerCellCount - Count the number of points within all cells in the cell DM 503 504 Not collective 505 506 Input parameter: 507 . dm - the DMSwarm 508 509 Output parameters: 510 + ncells - the number of cells in the cell DM (optional argument, pass NULL to ignore) 511 - count - array of length ncells containing the number of points per cell 512 513 Level: beginner 514 515 Notes: 516 The array count is allocated internally and must be free'd by the user. 517 518 .seealso: `DMSwarmSetType()`, `DMSwarmSetCellDM()`, `DMSwarmType` 519 @*/ 520 PETSC_EXTERN PetscErrorCode DMSwarmCreatePointPerCellCount(DM dm, PetscInt *ncells, PetscInt **count) { 521 PetscBool isvalid; 522 PetscInt nel; 523 PetscInt *sum; 524 525 PetscFunctionBegin; 526 PetscCall(DMSwarmSortGetIsValid(dm, &isvalid)); 527 nel = 0; 528 if (isvalid) { 529 PetscInt e; 530 531 PetscCall(DMSwarmSortGetSizes(dm, &nel, NULL)); 532 533 PetscCall(PetscMalloc1(nel, &sum)); 534 for (e = 0; e < nel; e++) { PetscCall(DMSwarmSortGetNumberOfPointsPerCell(dm, e, &sum[e])); } 535 } else { 536 DM celldm; 537 PetscBool isda, isplex, isshell; 538 PetscInt p, npoints; 539 PetscInt *swarm_cellid; 540 541 /* get the number of cells */ 542 PetscCall(DMSwarmGetCellDM(dm, &celldm)); 543 PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMDA, &isda)); 544 PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMPLEX, &isplex)); 545 PetscCall(PetscObjectTypeCompare((PetscObject)celldm, DMSHELL, &isshell)); 546 if (isda) { 547 PetscInt _nel, _npe; 548 const PetscInt *_element; 549 550 PetscCall(DMDAGetElements(celldm, &_nel, &_npe, &_element)); 551 nel = _nel; 552 PetscCall(DMDARestoreElements(celldm, &_nel, &_npe, &_element)); 553 } else if (isplex) { 554 PetscInt ps, pe; 555 556 PetscCall(DMPlexGetHeightStratum(celldm, 0, &ps, &pe)); 557 nel = pe - ps; 558 } else if (isshell) { 559 PetscErrorCode (*method_DMShellGetNumberOfCells)(DM, PetscInt *); 560 561 PetscCall(PetscObjectQueryFunction((PetscObject)celldm, "DMGetNumberOfCells_C", &method_DMShellGetNumberOfCells)); 562 if (method_DMShellGetNumberOfCells) { 563 PetscCall(method_DMShellGetNumberOfCells(celldm, &nel)); 564 } else 565 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);"); 566 } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot determine the number of cells for a DM not of type DA, PLEX or SHELL"); 567 568 PetscCall(PetscMalloc1(nel, &sum)); 569 PetscCall(PetscArrayzero(sum, nel)); 570 PetscCall(DMSwarmGetLocalSize(dm, &npoints)); 571 PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 572 for (p = 0; p < npoints; p++) { 573 if (swarm_cellid[p] != DMLOCATEPOINT_POINT_NOT_FOUND) { sum[swarm_cellid[p]]++; } 574 } 575 PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 576 } 577 if (ncells) { *ncells = nel; } 578 *count = sum; 579 PetscFunctionReturn(0); 580 } 581 582 /*@ 583 DMSwarmGetNumSpecies - Get the number of particle species 584 585 Not collective 586 587 Input parameter: 588 . dm - the DMSwarm 589 590 Output parameters: 591 . Ns - the number of species 592 593 Level: intermediate 594 595 .seealso: `DMSwarmSetNumSpecies()`, `DMSwarmSetType()`, `DMSwarmType` 596 @*/ 597 PetscErrorCode DMSwarmGetNumSpecies(DM sw, PetscInt *Ns) { 598 DM_Swarm *swarm = (DM_Swarm *)sw->data; 599 600 PetscFunctionBegin; 601 *Ns = swarm->Ns; 602 PetscFunctionReturn(0); 603 } 604 605 /*@ 606 DMSwarmSetNumSpecies - Set the number of particle species 607 608 Not collective 609 610 Input parameter: 611 + dm - the DMSwarm 612 - Ns - the number of species 613 614 Level: intermediate 615 616 .seealso: `DMSwarmGetNumSpecies()`, `DMSwarmSetType()`, `DMSwarmType` 617 @*/ 618 PetscErrorCode DMSwarmSetNumSpecies(DM sw, PetscInt Ns) { 619 DM_Swarm *swarm = (DM_Swarm *)sw->data; 620 621 PetscFunctionBegin; 622 swarm->Ns = Ns; 623 PetscFunctionReturn(0); 624 } 625 626 /*@C 627 DMSwarmGetCoordinateFunction - Get the function setting initial particle positions, if it exists 628 629 Not collective 630 631 Input parameter: 632 . dm - the DMSwarm 633 634 Output Parameter: 635 . coordFunc - the function setting initial particle positions, or NULL 636 637 Level: intermediate 638 639 .seealso: `DMSwarmSetCoordinateFunction()`, `DMSwarmGetVelocityFunction()`, `DMSwarmInitializeCoordinates()` 640 @*/ 641 PetscErrorCode DMSwarmGetCoordinateFunction(DM sw, PetscSimplePointFunc *coordFunc) { 642 DM_Swarm *swarm = (DM_Swarm *)sw->data; 643 644 PetscFunctionBegin; 645 PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 646 PetscValidPointer(coordFunc, 2); 647 *coordFunc = swarm->coordFunc; 648 PetscFunctionReturn(0); 649 } 650 651 /*@C 652 DMSwarmSetCoordinateFunction - Set the function setting initial particle positions 653 654 Not collective 655 656 Input parameters: 657 + dm - the DMSwarm 658 - coordFunc - the function setting initial particle positions 659 660 Level: intermediate 661 662 .seealso: `DMSwarmGetCoordinateFunction()`, `DMSwarmSetVelocityFunction()`, `DMSwarmInitializeCoordinates()` 663 @*/ 664 PetscErrorCode DMSwarmSetCoordinateFunction(DM sw, PetscSimplePointFunc coordFunc) { 665 DM_Swarm *swarm = (DM_Swarm *)sw->data; 666 667 PetscFunctionBegin; 668 PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 669 PetscValidFunction(coordFunc, 2); 670 swarm->coordFunc = coordFunc; 671 PetscFunctionReturn(0); 672 } 673 674 /*@C 675 DMSwarmGetCoordinateFunction - Get the function setting initial particle velocities, if it exists 676 677 Not collective 678 679 Input parameter: 680 . dm - the DMSwarm 681 682 Output Parameter: 683 . velFunc - the function setting initial particle velocities, or NULL 684 685 Level: intermediate 686 687 .seealso: `DMSwarmSetVelocityFunction()`, `DMSwarmGetCoordinateFunction()`, `DMSwarmInitializeVelocities()` 688 @*/ 689 PetscErrorCode DMSwarmGetVelocityFunction(DM sw, PetscSimplePointFunc *velFunc) { 690 DM_Swarm *swarm = (DM_Swarm *)sw->data; 691 692 PetscFunctionBegin; 693 PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 694 PetscValidPointer(velFunc, 2); 695 *velFunc = swarm->velFunc; 696 PetscFunctionReturn(0); 697 } 698 699 /*@C 700 DMSwarmSetVelocityFunction - Set the function setting initial particle velocities 701 702 Not collective 703 704 Input parameters: 705 + dm - the DMSwarm 706 - coordFunc - the function setting initial particle velocities 707 708 Level: intermediate 709 710 .seealso: `DMSwarmGetVelocityFunction()`, `DMSwarmSetCoordinateFunction()`, `DMSwarmInitializeVelocities()` 711 @*/ 712 PetscErrorCode DMSwarmSetVelocityFunction(DM sw, PetscSimplePointFunc velFunc) { 713 DM_Swarm *swarm = (DM_Swarm *)sw->data; 714 715 PetscFunctionBegin; 716 PetscValidHeaderSpecific(sw, DM_CLASSID, 1); 717 PetscValidFunction(velFunc, 2); 718 swarm->velFunc = velFunc; 719 PetscFunctionReturn(0); 720 } 721 722 /*@C 723 DMSwarmComputeLocalSize - Compute the local number and distribution of particles based upon a density function 724 725 Not collective 726 727 Input Parameters: 728 + sw - The DMSwarm 729 . N - The target number of particles 730 - density - The density field for the particle layout, normalized to unity 731 732 Note: One particle will be created for each species. 733 734 Level: advanced 735 736 .seealso: `DMSwarmComputeLocalSizeFromOptions()` 737 @*/ 738 PetscErrorCode DMSwarmComputeLocalSize(DM sw, PetscInt N, PetscProbFunc density) { 739 DM dm; 740 PetscQuadrature quad; 741 const PetscReal *xq, *wq; 742 PetscInt *npc, *cellid; 743 PetscReal xi0[3]; 744 PetscInt Ns, cStart, cEnd, c, dim, d, Nq, q, Np = 0, p; 745 PetscBool simplex; 746 747 PetscFunctionBegin; 748 PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 749 PetscCall(DMSwarmGetCellDM(sw, &dm)); 750 PetscCall(DMGetDimension(dm, &dim)); 751 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 752 PetscCall(DMPlexIsSimplex(dm, &simplex)); 753 PetscCall(DMGetCoordinatesLocalSetUp(dm)); 754 if (simplex) PetscCall(PetscDTStroudConicalQuadrature(dim, 1, 5, -1.0, 1.0, &quad)); 755 else PetscCall(PetscDTGaussTensorQuadrature(dim, 1, 5, -1.0, 1.0, &quad)); 756 PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &xq, &wq)); 757 PetscCall(PetscMalloc1(cEnd - cStart, &npc)); 758 /* Integrate the density function to get the number of particles in each cell */ 759 for (d = 0; d < dim; ++d) xi0[d] = -1.0; 760 for (c = 0; c < cEnd - cStart; ++c) { 761 const PetscInt cell = c + cStart; 762 PetscReal v0[3], J[9], invJ[9], detJ; 763 PetscReal n_int = 0.; 764 765 PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, NULL, v0, J, invJ, &detJ)); 766 for (q = 0; q < Nq; ++q) { 767 PetscReal xr[3], den[3]; 768 769 CoordinatesRefToReal(dim, dim, xi0, v0, J, &xq[q * dim], xr); 770 PetscCall(density(xr, NULL, den)); 771 n_int += den[0] * wq[q]; 772 } 773 npc[c] = (PetscInt)(N * n_int); 774 npc[c] *= Ns; 775 Np += npc[c]; 776 } 777 PetscCall(PetscQuadratureDestroy(&quad)); 778 PetscCall(DMSwarmSetLocalSizes(sw, Np, 0)); 779 780 PetscCall(DMSwarmGetField(sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 781 for (c = 0, p = 0; c < cEnd - cStart; ++c) { 782 for (q = 0; q < npc[c]; ++q, ++p) cellid[p] = c; 783 } 784 PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 785 PetscCall(PetscFree(npc)); 786 PetscFunctionReturn(0); 787 } 788 789 /*@ 790 DMSwarmComputeLocalSizeFromOptions - Compute the local number and distribution of particles based upon a density function determined by options 791 792 Not collective 793 794 Input Parameters: 795 , sw - The DMSwarm 796 797 Level: advanced 798 799 .seealso: `DMSwarmComputeLocalSize()` 800 @*/ 801 PetscErrorCode DMSwarmComputeLocalSizeFromOptions(DM sw) { 802 PetscProbFunc pdf; 803 const char *prefix; 804 char funcname[PETSC_MAX_PATH_LEN]; 805 PetscInt *N, Ns, dim, n; 806 PetscBool flg; 807 PetscMPIInt size, rank; 808 809 PetscFunctionBegin; 810 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sw), &size)); 811 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sw), &rank)); 812 PetscCall(PetscCalloc1(size, &N)); 813 PetscOptionsBegin(PetscObjectComm((PetscObject)sw), "", "DMSwarm Options", "DMSWARM"); 814 n = size; 815 PetscCall(PetscOptionsIntArray("-dm_swarm_num_particles", "The target number of particles", "", N, &n, NULL)); 816 PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 817 PetscCall(PetscOptionsInt("-dm_swarm_num_species", "The number of species", "DMSwarmSetNumSpecies", Ns, &Ns, &flg)); 818 if (flg) PetscCall(DMSwarmSetNumSpecies(sw, Ns)); 819 PetscCall(PetscOptionsString("-dm_swarm_coordinate_function", "Function to determine particle coordinates", "DMSwarmSetCoordinateFunction", funcname, funcname, sizeof(funcname), &flg)); 820 PetscOptionsEnd(); 821 if (flg) { 822 PetscSimplePointFunc coordFunc; 823 824 PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 825 PetscCall(PetscDLSym(NULL, funcname, (void **)&coordFunc)); 826 PetscCheck(coordFunc, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname); 827 PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 828 PetscCall(DMSwarmSetLocalSizes(sw, N[rank] * Ns, 0)); 829 PetscCall(DMSwarmSetCoordinateFunction(sw, coordFunc)); 830 } else { 831 PetscCall(DMGetDimension(sw, &dim)); 832 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sw, &prefix)); 833 PetscCall(PetscProbCreateFromOptions(dim, prefix, "-dm_swarm_coordinate_density", &pdf, NULL, NULL)); 834 PetscCall(DMSwarmComputeLocalSize(sw, N[rank], pdf)); 835 } 836 PetscCall(PetscFree(N)); 837 PetscFunctionReturn(0); 838 } 839 840 /*@ 841 DMSwarmInitializeCoordinates - Determine the initial coordinates of particles for a PIC method 842 843 Not collective 844 845 Input Parameters: 846 , sw - The DMSwarm 847 848 Note: Currently, we randomly place particles in their assigned cell 849 850 Level: advanced 851 852 .seealso: `DMSwarmComputeLocalSize()`, `DMSwarmInitializeVelocities()` 853 @*/ 854 PetscErrorCode DMSwarmInitializeCoordinates(DM sw) { 855 PetscSimplePointFunc coordFunc; 856 PetscScalar *weight; 857 PetscReal *x; 858 PetscInt *species; 859 void *ctx; 860 PetscBool removePoints = PETSC_TRUE; 861 PetscDataType dtype; 862 PetscInt Np, p, Ns, dim, d, bs; 863 864 PetscFunctionBeginUser; 865 PetscCall(DMGetDimension(sw, &dim)); 866 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 867 PetscCall(DMSwarmGetNumSpecies(sw, &Ns)); 868 PetscCall(DMSwarmGetCoordinateFunction(sw, &coordFunc)); 869 870 PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, &bs, &dtype, (void **)&x)); 871 PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&weight)); 872 PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 873 if (coordFunc) { 874 PetscCall(DMGetApplicationContext(sw, &ctx)); 875 for (p = 0; p < Np; ++p) { 876 PetscScalar X[3]; 877 878 PetscCall((*coordFunc)(dim, 0., NULL, p, X, ctx)); 879 for (d = 0; d < dim; ++d) x[p * dim + d] = PetscRealPart(X[d]); 880 weight[p] = 1.0; 881 species[p] = p % Ns; 882 } 883 } else { 884 DM dm; 885 PetscRandom rnd; 886 PetscReal xi0[3]; 887 PetscInt cStart, cEnd, c; 888 889 PetscCall(DMSwarmGetCellDM(sw, &dm)); 890 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 891 892 /* Set particle position randomly in cell, set weights to 1 */ 893 PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd)); 894 PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0)); 895 PetscCall(PetscRandomSetFromOptions(rnd)); 896 PetscCall(DMSwarmSortGetAccess(sw)); 897 for (d = 0; d < dim; ++d) xi0[d] = -1.0; 898 for (c = cStart; c < cEnd; ++c) { 899 PetscReal v0[3], J[9], invJ[9], detJ; 900 PetscInt *pidx, Npc, q; 901 902 PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Npc, &pidx)); 903 PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); 904 for (q = 0; q < Npc; ++q) { 905 const PetscInt p = pidx[q]; 906 PetscReal xref[3]; 907 908 for (d = 0; d < dim; ++d) PetscCall(PetscRandomGetValueReal(rnd, &xref[d])); 909 CoordinatesRefToReal(dim, dim, xi0, v0, J, xref, &x[p * dim]); 910 911 weight[p] = 1.0; 912 species[p] = p % Ns; 913 } 914 PetscCall(PetscFree(pidx)); 915 } 916 PetscCall(PetscRandomDestroy(&rnd)); 917 PetscCall(DMSwarmSortRestoreAccess(sw)); 918 } 919 PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&x)); 920 PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight)); 921 PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 922 923 PetscCall(DMSwarmMigrate(sw, removePoints)); 924 PetscCall(DMLocalizeCoordinates(sw)); 925 PetscFunctionReturn(0); 926 } 927 928 /*@C 929 DMSwarmInitializeVelocities - Set the initial velocities of particles using a distribution. 930 931 Collective on dm 932 933 Input Parameters: 934 + sw - The DMSwarm object 935 . sampler - A function which uniformly samples the velocity PDF 936 - v0 - The velocity scale for nondimensionalization for each species 937 938 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. 939 940 Level: advanced 941 942 .seealso: `DMSwarmComputeLocalSize()`, `DMSwarmInitializeCoordinates()`, `DMSwarmInitializeVelocitiesFromOptions()` 943 @*/ 944 PetscErrorCode DMSwarmInitializeVelocities(DM sw, PetscProbFunc sampler, const PetscReal v0[]) { 945 PetscSimplePointFunc velFunc; 946 PetscReal *v; 947 PetscInt *species; 948 void *ctx; 949 PetscInt dim, Np, p; 950 951 PetscFunctionBegin; 952 PetscCall(DMSwarmGetVelocityFunction(sw, &velFunc)); 953 954 PetscCall(DMGetDimension(sw, &dim)); 955 PetscCall(DMSwarmGetLocalSize(sw, &Np)); 956 PetscCall(DMSwarmGetField(sw, "velocity", NULL, NULL, (void **)&v)); 957 PetscCall(DMSwarmGetField(sw, "species", NULL, NULL, (void **)&species)); 958 if (v0[0] == 0.) { 959 PetscCall(PetscArrayzero(v, Np * dim)); 960 } else if (velFunc) { 961 PetscCall(DMGetApplicationContext(sw, &ctx)); 962 for (p = 0; p < Np; ++p) { 963 PetscInt s = species[p], d; 964 PetscScalar vel[3]; 965 966 PetscCall((*velFunc)(dim, 0., NULL, p, vel, ctx)); 967 for (d = 0; d < dim; ++d) v[p * dim + d] = (v0[s] / v0[0]) * PetscRealPart(vel[d]); 968 } 969 } else { 970 PetscRandom rnd; 971 972 PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)sw), &rnd)); 973 PetscCall(PetscRandomSetInterval(rnd, 0, 1.)); 974 PetscCall(PetscRandomSetFromOptions(rnd)); 975 976 for (p = 0; p < Np; ++p) { 977 PetscInt s = species[p], d; 978 PetscReal a[3], vel[3]; 979 980 for (d = 0; d < dim; ++d) PetscCall(PetscRandomGetValueReal(rnd, &a[d])); 981 PetscCall(sampler(a, NULL, vel)); 982 for (d = 0; d < dim; ++d) { v[p * dim + d] = (v0[s] / v0[0]) * vel[d]; } 983 } 984 PetscCall(PetscRandomDestroy(&rnd)); 985 } 986 PetscCall(DMSwarmRestoreField(sw, "velocity", NULL, NULL, (void **)&v)); 987 PetscCall(DMSwarmRestoreField(sw, "species", NULL, NULL, (void **)&species)); 988 PetscFunctionReturn(0); 989 } 990 991 /*@ 992 DMSwarmInitializeVelocitiesFromOptions - Set the initial velocities of particles using a distribution determined from options. 993 994 Collective on dm 995 996 Input Parameters: 997 + sw - The DMSwarm object 998 - v0 - The velocity scale for nondimensionalization for each species 999 1000 Level: advanced 1001 1002 .seealso: `DMSwarmComputeLocalSize()`, `DMSwarmInitializeCoordinates()`, `DMSwarmInitializeVelocities()` 1003 @*/ 1004 PetscErrorCode DMSwarmInitializeVelocitiesFromOptions(DM sw, const PetscReal v0[]) { 1005 PetscProbFunc sampler; 1006 PetscInt dim; 1007 const char *prefix; 1008 char funcname[PETSC_MAX_PATH_LEN]; 1009 PetscBool flg; 1010 1011 PetscFunctionBegin; 1012 PetscOptionsBegin(PetscObjectComm((PetscObject)sw), "", "DMSwarm Options", "DMSWARM"); 1013 PetscCall(PetscOptionsString("-dm_swarm_velocity_function", "Function to determine particle velocities", "DMSwarmSetVelocityFunction", funcname, funcname, sizeof(funcname), &flg)); 1014 PetscOptionsEnd(); 1015 if (flg) { 1016 PetscSimplePointFunc velFunc; 1017 1018 PetscCall(PetscDLSym(NULL, funcname, (void **)&velFunc)); 1019 PetscCheck(velFunc, PetscObjectComm((PetscObject)sw), PETSC_ERR_ARG_WRONG, "Could not locate function %s", funcname); 1020 PetscCall(DMSwarmSetVelocityFunction(sw, velFunc)); 1021 } 1022 PetscCall(DMGetDimension(sw, &dim)); 1023 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sw, &prefix)); 1024 PetscCall(PetscProbCreateFromOptions(dim, prefix, "-dm_swarm_velocity_density", NULL, NULL, &sampler)); 1025 PetscCall(DMSwarmInitializeVelocities(sw, sampler, v0)); 1026 PetscFunctionReturn(0); 1027 } 1028