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