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