1 #define PETSCDM_DLL 2 #include <petsc/private/dmswarmimpl.h> /*I "petscdmswarm.h" I*/ 3 #include <petsc/private/hashsetij.h> 4 #include <petsc/private/petscfeimpl.h> 5 #include <petscviewer.h> 6 #include <petscdraw.h> 7 #include <petscdmplex.h> 8 #include <petscblaslapack.h> 9 #include "../src/dm/impls/swarm/data_bucket.h" 10 #include <petscdmlabel.h> 11 #include <petscsection.h> 12 13 PetscLogEvent DMSWARM_Migrate, DMSWARM_SetSizes, DMSWARM_AddPoints, DMSWARM_RemovePoints, DMSWARM_Sort; 14 PetscLogEvent DMSWARM_DataExchangerTopologySetup, DMSWARM_DataExchangerBegin, DMSWARM_DataExchangerEnd; 15 PetscLogEvent DMSWARM_DataExchangerSendCount, DMSWARM_DataExchangerPack; 16 17 const char *DMSwarmTypeNames[] = {"basic", "pic", NULL}; 18 const char *DMSwarmMigrateTypeNames[] = {"basic", "dmcellnscatter", "dmcellexact", "user", NULL}; 19 const char *DMSwarmCollectTypeNames[] = {"basic", "boundingbox", "general", "user", NULL}; 20 const char *DMSwarmPICLayoutTypeNames[] = {"regular", "gauss", "subdivision", NULL}; 21 22 const char DMSwarmField_pid[] = "DMSwarm_pid"; 23 const char DMSwarmField_rank[] = "DMSwarm_rank"; 24 const char DMSwarmPICField_coor[] = "DMSwarmPIC_coor"; 25 const char DMSwarmPICField_cellid[] = "DMSwarm_cellid"; 26 27 PetscInt SwarmDataFieldId = -1; 28 29 #if defined(PETSC_HAVE_HDF5) 30 #include <petscviewerhdf5.h> 31 32 PetscErrorCode VecView_Swarm_HDF5_Internal(Vec v, PetscViewer viewer) 33 { 34 DM dm; 35 PetscReal seqval; 36 PetscInt seqnum, bs; 37 PetscBool isseq; 38 39 PetscFunctionBegin; 40 PetscCall(VecGetDM(v, &dm)); 41 PetscCall(VecGetBlockSize(v, &bs)); 42 PetscCall(PetscViewerHDF5PushGroup(viewer, "/particle_fields")); 43 PetscCall(PetscObjectTypeCompare((PetscObject)v, VECSEQ, &isseq)); 44 PetscCall(DMGetOutputSequenceNumber(dm, &seqnum, &seqval)); 45 PetscCall(PetscViewerHDF5SetTimestep(viewer, seqnum)); 46 /* PetscCall(DMSequenceView_HDF5(dm, "time", seqnum, (PetscScalar) seqval, viewer)); */ 47 PetscCall(VecViewNative(v, viewer)); 48 PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)v, "Nc", PETSC_INT, (void *)&bs)); 49 PetscCall(PetscViewerHDF5PopGroup(viewer)); 50 PetscFunctionReturn(PETSC_SUCCESS); 51 } 52 53 PetscErrorCode DMSwarmView_HDF5(DM dm, PetscViewer viewer) 54 { 55 Vec coordinates; 56 PetscInt Np; 57 PetscBool isseq; 58 59 PetscFunctionBegin; 60 PetscCall(DMSwarmGetSize(dm, &Np)); 61 PetscCall(DMSwarmCreateGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates)); 62 PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 63 PetscCall(PetscViewerHDF5PushGroup(viewer, "/particles")); 64 PetscCall(PetscObjectTypeCompare((PetscObject)coordinates, VECSEQ, &isseq)); 65 PetscCall(VecViewNative(coordinates, viewer)); 66 PetscCall(PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject)coordinates, "Np", PETSC_INT, (void *)&Np)); 67 PetscCall(PetscViewerHDF5PopGroup(viewer)); 68 PetscCall(DMSwarmDestroyGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates)); 69 PetscFunctionReturn(PETSC_SUCCESS); 70 } 71 #endif 72 73 PetscErrorCode VecView_Swarm(Vec v, PetscViewer viewer) 74 { 75 DM dm; 76 #if defined(PETSC_HAVE_HDF5) 77 PetscBool ishdf5; 78 #endif 79 80 PetscFunctionBegin; 81 PetscCall(VecGetDM(v, &dm)); 82 PetscCheck(dm, PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM"); 83 #if defined(PETSC_HAVE_HDF5) 84 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5)); 85 if (ishdf5) { 86 PetscCall(VecView_Swarm_HDF5_Internal(v, viewer)); 87 PetscFunctionReturn(PETSC_SUCCESS); 88 } 89 #endif 90 PetscCall(VecViewNative(v, viewer)); 91 PetscFunctionReturn(PETSC_SUCCESS); 92 } 93 94 /*@C 95 DMSwarmVectorDefineField - Sets the field from which to define a Vec object 96 when DMCreateLocalVector(), or DMCreateGlobalVector() is called 97 98 Collective on dm 99 100 Input parameters: 101 + dm - a DMSwarm 102 - fieldname - the textual name given to a registered field 103 104 Level: beginner 105 106 Notes: 107 108 The field with name fieldname must be defined as having a data type of PetscScalar. 109 110 This function must be called prior to calling DMCreateLocalVector(), DMCreateGlobalVector(). 111 Multiple calls to DMSwarmVectorDefineField() are permitted. 112 113 .seealso: `DMSwarmRegisterPetscDatatypeField()`, `DMCreateGlobalVector()`, `DMCreateLocalVector()` 114 @*/ 115 PetscErrorCode DMSwarmVectorDefineField(DM dm, const char fieldname[]) 116 { 117 DM_Swarm *swarm = (DM_Swarm *)dm->data; 118 PetscInt bs, n; 119 PetscScalar *array; 120 PetscDataType type; 121 122 PetscFunctionBegin; 123 if (!swarm->issetup) PetscCall(DMSetUp(dm)); 124 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL)); 125 PetscCall(DMSwarmGetField(dm, fieldname, &bs, &type, (void **)&array)); 126 127 /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */ 128 PetscCheck(type == PETSC_REAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL"); 129 PetscCall(PetscSNPrintf(swarm->vec_field_name, PETSC_MAX_PATH_LEN - 1, "%s", fieldname)); 130 swarm->vec_field_set = PETSC_TRUE; 131 swarm->vec_field_bs = bs; 132 swarm->vec_field_nlocal = n; 133 PetscCall(DMSwarmRestoreField(dm, fieldname, &bs, &type, (void **)&array)); 134 PetscFunctionReturn(PETSC_SUCCESS); 135 } 136 137 /* requires DMSwarmDefineFieldVector has been called */ 138 PetscErrorCode DMCreateGlobalVector_Swarm(DM dm, Vec *vec) 139 { 140 DM_Swarm *swarm = (DM_Swarm *)dm->data; 141 Vec x; 142 char name[PETSC_MAX_PATH_LEN]; 143 144 PetscFunctionBegin; 145 if (!swarm->issetup) PetscCall(DMSetUp(dm)); 146 PetscCheck(swarm->vec_field_set, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmVectorDefineField first"); 147 PetscCheck(swarm->vec_field_nlocal == swarm->db->L, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSwarm sizes have changed since last call to VectorDefineField first"); /* Stale data */ 148 149 PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmField_%s", swarm->vec_field_name)); 150 PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), &x)); 151 PetscCall(PetscObjectSetName((PetscObject)x, name)); 152 PetscCall(VecSetSizes(x, swarm->db->L * swarm->vec_field_bs, PETSC_DETERMINE)); 153 PetscCall(VecSetBlockSize(x, swarm->vec_field_bs)); 154 PetscCall(VecSetDM(x, dm)); 155 PetscCall(VecSetFromOptions(x)); 156 PetscCall(VecSetDM(x, dm)); 157 PetscCall(VecSetOperation(x, VECOP_VIEW, (void (*)(void))VecView_Swarm)); 158 *vec = x; 159 PetscFunctionReturn(PETSC_SUCCESS); 160 } 161 162 /* requires DMSwarmDefineFieldVector has been called */ 163 PetscErrorCode DMCreateLocalVector_Swarm(DM dm, Vec *vec) 164 { 165 DM_Swarm *swarm = (DM_Swarm *)dm->data; 166 Vec x; 167 char name[PETSC_MAX_PATH_LEN]; 168 169 PetscFunctionBegin; 170 if (!swarm->issetup) PetscCall(DMSetUp(dm)); 171 PetscCheck(swarm->vec_field_set, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmVectorDefineField first"); 172 PetscCheck(swarm->vec_field_nlocal == swarm->db->L, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSwarm sizes have changed since last call to VectorDefineField first"); 173 174 PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmField_%s", swarm->vec_field_name)); 175 PetscCall(VecCreate(PETSC_COMM_SELF, &x)); 176 PetscCall(PetscObjectSetName((PetscObject)x, name)); 177 PetscCall(VecSetSizes(x, swarm->db->L * swarm->vec_field_bs, PETSC_DETERMINE)); 178 PetscCall(VecSetBlockSize(x, swarm->vec_field_bs)); 179 PetscCall(VecSetDM(x, dm)); 180 PetscCall(VecSetFromOptions(x)); 181 *vec = x; 182 PetscFunctionReturn(PETSC_SUCCESS); 183 } 184 185 static PetscErrorCode DMSwarmDestroyVectorFromField_Private(DM dm, const char fieldname[], Vec *vec) 186 { 187 DM_Swarm *swarm = (DM_Swarm *)dm->data; 188 DMSwarmDataField gfield; 189 PetscInt bs, nlocal, fid = -1, cfid = -2; 190 PetscBool flg; 191 192 PetscFunctionBegin; 193 /* check vector is an inplace array */ 194 PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldname, &fid)); 195 PetscCall(PetscObjectComposedDataGetInt((PetscObject)*vec, SwarmDataFieldId, cfid, flg)); 196 PetscCheck(cfid == fid, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Vector being destroyed was not created from DMSwarm field(%s)! %" PetscInt_FMT " != %" PetscInt_FMT "", fieldname, cfid, fid); 197 PetscCall(VecGetLocalSize(*vec, &nlocal)); 198 PetscCall(VecGetBlockSize(*vec, &bs)); 199 PetscCheck(nlocal / bs == swarm->db->L, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSwarm sizes have changed since vector was created - cannot ensure pointers are valid"); 200 PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield)); 201 PetscCall(DMSwarmDataFieldRestoreAccess(gfield)); 202 PetscCall(VecResetArray(*vec)); 203 PetscCall(VecDestroy(vec)); 204 PetscFunctionReturn(PETSC_SUCCESS); 205 } 206 207 static PetscErrorCode DMSwarmCreateVectorFromField_Private(DM dm, const char fieldname[], MPI_Comm comm, Vec *vec) 208 { 209 DM_Swarm *swarm = (DM_Swarm *)dm->data; 210 PetscDataType type; 211 PetscScalar *array; 212 PetscInt bs, n, fid; 213 char name[PETSC_MAX_PATH_LEN]; 214 PetscMPIInt size; 215 PetscBool iscuda, iskokkos; 216 217 PetscFunctionBegin; 218 if (!swarm->issetup) PetscCall(DMSetUp(dm)); 219 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL)); 220 PetscCall(DMSwarmGetField(dm, fieldname, &bs, &type, (void **)&array)); 221 PetscCheck(type == PETSC_REAL, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL"); 222 223 PetscCallMPI(MPI_Comm_size(comm, &size)); 224 PetscCall(PetscStrcmp(dm->vectype, VECKOKKOS, &iskokkos)); 225 PetscCall(PetscStrcmp(dm->vectype, VECCUDA, &iscuda)); 226 PetscCall(VecCreate(comm, vec)); 227 PetscCall(VecSetSizes(*vec, n * bs, PETSC_DETERMINE)); 228 PetscCall(VecSetBlockSize(*vec, bs)); 229 if (iskokkos) PetscCall(VecSetType(*vec, VECKOKKOS)); 230 else if (iscuda) PetscCall(VecSetType(*vec, VECCUDA)); 231 else PetscCall(VecSetType(*vec, VECSTANDARD)); 232 PetscCall(VecPlaceArray(*vec, array)); 233 234 PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "DMSwarmSharedField_%s", fieldname)); 235 PetscCall(PetscObjectSetName((PetscObject)*vec, name)); 236 237 /* Set guard */ 238 PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldIdByName(swarm->db, fieldname, &fid)); 239 PetscCall(PetscObjectComposedDataSetInt((PetscObject)*vec, SwarmDataFieldId, fid)); 240 241 PetscCall(VecSetDM(*vec, dm)); 242 PetscCall(VecSetOperation(*vec, VECOP_VIEW, (void (*)(void))VecView_Swarm)); 243 PetscFunctionReturn(PETSC_SUCCESS); 244 } 245 246 /* This creates a "mass matrix" between a finite element and particle space. If a finite element interpolant is given by 247 248 \hat f = \sum_i f_i \phi_i 249 250 and a particle function is given by 251 252 f = \sum_i w_i \delta(x - x_i) 253 254 then we want to require that 255 256 M \hat f = M_p f 257 258 where the particle mass matrix is given by 259 260 (M_p)_{ij} = \int \phi_i \delta(x - x_j) 261 262 The way Dave May does particles, they amount to quadratue weights rather than delta functions, so he has |J| is in 263 his integral. We allow this with the boolean flag. 264 */ 265 static PetscErrorCode DMSwarmComputeMassMatrix_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx) 266 { 267 const char *name = "Mass Matrix"; 268 MPI_Comm comm; 269 PetscDS prob; 270 PetscSection fsection, globalFSection; 271 PetscHSetIJ ht; 272 PetscLayout rLayout, colLayout; 273 PetscInt *dnz, *onz; 274 PetscInt locRows, locCols, rStart, colStart, colEnd, *rowIDXs; 275 PetscReal *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0}; 276 PetscScalar *elemMat; 277 PetscInt dim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0; 278 279 PetscFunctionBegin; 280 PetscCall(PetscObjectGetComm((PetscObject)mass, &comm)); 281 PetscCall(DMGetCoordinateDim(dmf, &dim)); 282 PetscCall(DMGetDS(dmf, &prob)); 283 PetscCall(PetscDSGetNumFields(prob, &Nf)); 284 PetscCall(PetscDSGetTotalDimension(prob, &totDim)); 285 PetscCall(PetscMalloc3(dim, &v0, dim * dim, &J, dim * dim, &invJ)); 286 PetscCall(DMGetLocalSection(dmf, &fsection)); 287 PetscCall(DMGetGlobalSection(dmf, &globalFSection)); 288 PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd)); 289 PetscCall(MatGetLocalSize(mass, &locRows, &locCols)); 290 291 PetscCall(PetscLayoutCreate(comm, &colLayout)); 292 PetscCall(PetscLayoutSetLocalSize(colLayout, locCols)); 293 PetscCall(PetscLayoutSetBlockSize(colLayout, 1)); 294 PetscCall(PetscLayoutSetUp(colLayout)); 295 PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd)); 296 PetscCall(PetscLayoutDestroy(&colLayout)); 297 298 PetscCall(PetscLayoutCreate(comm, &rLayout)); 299 PetscCall(PetscLayoutSetLocalSize(rLayout, locRows)); 300 PetscCall(PetscLayoutSetBlockSize(rLayout, 1)); 301 PetscCall(PetscLayoutSetUp(rLayout)); 302 PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL)); 303 PetscCall(PetscLayoutDestroy(&rLayout)); 304 305 PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz)); 306 PetscCall(PetscHSetIJCreate(&ht)); 307 308 PetscCall(PetscSynchronizedFlush(comm, NULL)); 309 /* count non-zeros */ 310 PetscCall(DMSwarmSortGetAccess(dmc)); 311 for (field = 0; field < Nf; ++field) { 312 for (cell = cStart; cell < cEnd; ++cell) { 313 PetscInt c, i; 314 PetscInt *findices, *cindices; /* fine is vertices, coarse is particles */ 315 PetscInt numFIndices, numCIndices; 316 317 PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL)); 318 PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices)); 319 maxC = PetscMax(maxC, numCIndices); 320 { 321 PetscHashIJKey key; 322 PetscBool missing; 323 for (i = 0; i < numFIndices; ++i) { 324 key.j = findices[i]; /* global column (from Plex) */ 325 if (key.j >= 0) { 326 /* Get indices for coarse elements */ 327 for (c = 0; c < numCIndices; ++c) { 328 key.i = cindices[c] + rStart; /* global cols (from Swarm) */ 329 if (key.i < 0) continue; 330 PetscCall(PetscHSetIJQueryAdd(ht, key, &missing)); 331 if (missing) { 332 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart]; 333 else ++onz[key.i - rStart]; 334 } else SETERRQ(PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Set new value at %" PetscInt_FMT ",%" PetscInt_FMT, key.i, key.j); 335 } 336 } 337 } 338 PetscCall(PetscFree(cindices)); 339 } 340 PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL)); 341 } 342 } 343 PetscCall(PetscHSetIJDestroy(&ht)); 344 PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL)); 345 PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE)); 346 PetscCall(PetscFree2(dnz, onz)); 347 PetscCall(PetscMalloc3(maxC * totDim, &elemMat, maxC, &rowIDXs, maxC * dim, &xi)); 348 for (field = 0; field < Nf; ++field) { 349 PetscTabulation Tcoarse; 350 PetscObject obj; 351 PetscReal *fieldVals; 352 PetscInt Nc, i; 353 354 PetscCall(PetscDSGetDiscretization(prob, field, &obj)); 355 PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc)); 356 PetscCheck(Nc == 1, PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %" PetscInt_FMT, Nc); 357 358 PetscCall(DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&fieldVals)); 359 for (cell = cStart; cell < cEnd; ++cell) { 360 PetscInt *findices, *cindices; 361 PetscInt numFIndices, numCIndices; 362 PetscInt p, c; 363 364 /* TODO: Use DMField instead of assuming affine */ 365 PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ)); 366 PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL)); 367 PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices)); 368 for (p = 0; p < numCIndices; ++p) CoordinatesRealToRef(dim, dim, v0ref, v0, invJ, &fieldVals[cindices[p] * dim], &xi[p * dim]); 369 PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, numCIndices, xi, 0, &Tcoarse)); 370 /* Get elemMat entries by multiplying by weight */ 371 PetscCall(PetscArrayzero(elemMat, numCIndices * totDim)); 372 for (i = 0; i < numFIndices; ++i) { 373 for (p = 0; p < numCIndices; ++p) { 374 for (c = 0; c < Nc; ++c) { 375 /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */ 376 elemMat[p * numFIndices + i] += Tcoarse->T[0][(p * numFIndices + i) * Nc + c] * (useDeltaFunction ? 1.0 : detJ); 377 } 378 } 379 } 380 for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart; 381 if (0) PetscCall(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat)); 382 PetscCall(MatSetValues(mass, numCIndices, rowIDXs, numFIndices, findices, elemMat, ADD_VALUES)); 383 PetscCall(PetscFree(cindices)); 384 PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL)); 385 PetscCall(PetscTabulationDestroy(&Tcoarse)); 386 } 387 PetscCall(DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&fieldVals)); 388 } 389 PetscCall(PetscFree3(elemMat, rowIDXs, xi)); 390 PetscCall(DMSwarmSortRestoreAccess(dmc)); 391 PetscCall(PetscFree3(v0, J, invJ)); 392 PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY)); 393 PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY)); 394 PetscFunctionReturn(PETSC_SUCCESS); 395 } 396 397 /* Returns empty matrix for use with SNES FD */ 398 static PetscErrorCode DMCreateMatrix_Swarm(DM sw, Mat *m) 399 { 400 Vec field; 401 PetscInt size; 402 403 PetscFunctionBegin; 404 PetscCall(DMGetGlobalVector(sw, &field)); 405 PetscCall(VecGetLocalSize(field, &size)); 406 PetscCall(DMRestoreGlobalVector(sw, &field)); 407 PetscCall(MatCreate(PETSC_COMM_WORLD, m)); 408 PetscCall(MatSetFromOptions(*m)); 409 PetscCall(MatSetSizes(*m, PETSC_DECIDE, PETSC_DECIDE, size, size)); 410 PetscCall(MatSeqAIJSetPreallocation(*m, 1, NULL)); 411 PetscCall(MatZeroEntries(*m)); 412 PetscCall(MatAssemblyBegin(*m, MAT_FINAL_ASSEMBLY)); 413 PetscCall(MatAssemblyEnd(*m, MAT_FINAL_ASSEMBLY)); 414 PetscCall(MatShift(*m, 1.0)); 415 PetscCall(MatSetDM(*m, sw)); 416 PetscFunctionReturn(PETSC_SUCCESS); 417 } 418 419 /* FEM cols, Particle rows */ 420 static PetscErrorCode DMCreateMassMatrix_Swarm(DM dmCoarse, DM dmFine, Mat *mass) 421 { 422 PetscSection gsf; 423 PetscInt m, n; 424 void *ctx; 425 426 PetscFunctionBegin; 427 PetscCall(DMGetGlobalSection(dmFine, &gsf)); 428 PetscCall(PetscSectionGetConstrainedStorageSize(gsf, &m)); 429 PetscCall(DMSwarmGetLocalSize(dmCoarse, &n)); 430 PetscCall(MatCreate(PetscObjectComm((PetscObject)dmCoarse), mass)); 431 PetscCall(MatSetSizes(*mass, n, m, PETSC_DETERMINE, PETSC_DETERMINE)); 432 PetscCall(MatSetType(*mass, dmCoarse->mattype)); 433 PetscCall(DMGetApplicationContext(dmFine, &ctx)); 434 435 PetscCall(DMSwarmComputeMassMatrix_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx)); 436 PetscCall(MatViewFromOptions(*mass, NULL, "-mass_mat_view")); 437 PetscFunctionReturn(PETSC_SUCCESS); 438 } 439 440 static PetscErrorCode DMSwarmComputeMassMatrixSquare_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx) 441 { 442 const char *name = "Mass Matrix Square"; 443 MPI_Comm comm; 444 PetscDS prob; 445 PetscSection fsection, globalFSection; 446 PetscHSetIJ ht; 447 PetscLayout rLayout, colLayout; 448 PetscInt *dnz, *onz, *adj, depth, maxConeSize, maxSupportSize, maxAdjSize; 449 PetscInt locRows, locCols, rStart, colStart, colEnd, *rowIDXs; 450 PetscReal *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0}; 451 PetscScalar *elemMat, *elemMatSq; 452 PetscInt cdim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0; 453 454 PetscFunctionBegin; 455 PetscCall(PetscObjectGetComm((PetscObject)mass, &comm)); 456 PetscCall(DMGetCoordinateDim(dmf, &cdim)); 457 PetscCall(DMGetDS(dmf, &prob)); 458 PetscCall(PetscDSGetNumFields(prob, &Nf)); 459 PetscCall(PetscDSGetTotalDimension(prob, &totDim)); 460 PetscCall(PetscMalloc3(cdim, &v0, cdim * cdim, &J, cdim * cdim, &invJ)); 461 PetscCall(DMGetLocalSection(dmf, &fsection)); 462 PetscCall(DMGetGlobalSection(dmf, &globalFSection)); 463 PetscCall(DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd)); 464 PetscCall(MatGetLocalSize(mass, &locRows, &locCols)); 465 466 PetscCall(PetscLayoutCreate(comm, &colLayout)); 467 PetscCall(PetscLayoutSetLocalSize(colLayout, locCols)); 468 PetscCall(PetscLayoutSetBlockSize(colLayout, 1)); 469 PetscCall(PetscLayoutSetUp(colLayout)); 470 PetscCall(PetscLayoutGetRange(colLayout, &colStart, &colEnd)); 471 PetscCall(PetscLayoutDestroy(&colLayout)); 472 473 PetscCall(PetscLayoutCreate(comm, &rLayout)); 474 PetscCall(PetscLayoutSetLocalSize(rLayout, locRows)); 475 PetscCall(PetscLayoutSetBlockSize(rLayout, 1)); 476 PetscCall(PetscLayoutSetUp(rLayout)); 477 PetscCall(PetscLayoutGetRange(rLayout, &rStart, NULL)); 478 PetscCall(PetscLayoutDestroy(&rLayout)); 479 480 PetscCall(DMPlexGetDepth(dmf, &depth)); 481 PetscCall(DMPlexGetMaxSizes(dmf, &maxConeSize, &maxSupportSize)); 482 maxAdjSize = PetscPowInt(maxConeSize * maxSupportSize, depth); 483 PetscCall(PetscMalloc1(maxAdjSize, &adj)); 484 485 PetscCall(PetscCalloc2(locRows, &dnz, locRows, &onz)); 486 PetscCall(PetscHSetIJCreate(&ht)); 487 /* Count nonzeros 488 This is just FVM++, but we cannot use the Plex P0 allocation since unknowns in a cell will not be contiguous 489 */ 490 PetscCall(DMSwarmSortGetAccess(dmc)); 491 for (cell = cStart; cell < cEnd; ++cell) { 492 PetscInt i; 493 PetscInt *cindices; 494 PetscInt numCIndices; 495 #if 0 496 PetscInt adjSize = maxAdjSize, a, j; 497 #endif 498 499 PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices)); 500 maxC = PetscMax(maxC, numCIndices); 501 /* Diagonal block */ 502 for (i = 0; i < numCIndices; ++i) dnz[cindices[i]] += numCIndices; 503 #if 0 504 /* Off-diagonal blocks */ 505 PetscCall(DMPlexGetAdjacency(dmf, cell, &adjSize, &adj)); 506 for (a = 0; a < adjSize; ++a) { 507 if (adj[a] >= cStart && adj[a] < cEnd && adj[a] != cell) { 508 const PetscInt ncell = adj[a]; 509 PetscInt *ncindices; 510 PetscInt numNCIndices; 511 512 PetscCall(DMSwarmSortGetPointsPerCell(dmc, ncell, &numNCIndices, &ncindices)); 513 { 514 PetscHashIJKey key; 515 PetscBool missing; 516 517 for (i = 0; i < numCIndices; ++i) { 518 key.i = cindices[i] + rStart; /* global rows (from Swarm) */ 519 if (key.i < 0) continue; 520 for (j = 0; j < numNCIndices; ++j) { 521 key.j = ncindices[j] + rStart; /* global column (from Swarm) */ 522 if (key.j < 0) continue; 523 PetscCall(PetscHSetIJQueryAdd(ht, key, &missing)); 524 if (missing) { 525 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart]; 526 else ++onz[key.i - rStart]; 527 } 528 } 529 } 530 } 531 PetscCall(PetscFree(ncindices)); 532 } 533 } 534 #endif 535 PetscCall(PetscFree(cindices)); 536 } 537 PetscCall(PetscHSetIJDestroy(&ht)); 538 PetscCall(MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL)); 539 PetscCall(MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE)); 540 PetscCall(PetscFree2(dnz, onz)); 541 PetscCall(PetscMalloc4(maxC * totDim, &elemMat, maxC * maxC, &elemMatSq, maxC, &rowIDXs, maxC * cdim, &xi)); 542 /* Fill in values 543 Each entry is a sum of terms \phi_i(x_p) \phi_i(x_q) 544 Start just by producing block diagonal 545 Could loop over adjacent cells 546 Produce neighboring element matrix 547 TODO Determine which columns and rows correspond to shared dual vector 548 Do MatMatMult with rectangular matrices 549 Insert block 550 */ 551 for (field = 0; field < Nf; ++field) { 552 PetscTabulation Tcoarse; 553 PetscObject obj; 554 PetscReal *coords; 555 PetscInt Nc, i; 556 557 PetscCall(PetscDSGetDiscretization(prob, field, &obj)); 558 PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc)); 559 PetscCheck(Nc == 1, PetscObjectComm((PetscObject)dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %" PetscInt_FMT, Nc); 560 PetscCall(DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 561 for (cell = cStart; cell < cEnd; ++cell) { 562 PetscInt *findices, *cindices; 563 PetscInt numFIndices, numCIndices; 564 PetscInt p, c; 565 566 /* TODO: Use DMField instead of assuming affine */ 567 PetscCall(DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ)); 568 PetscCall(DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL)); 569 PetscCall(DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices)); 570 for (p = 0; p < numCIndices; ++p) CoordinatesRealToRef(cdim, cdim, v0ref, v0, invJ, &coords[cindices[p] * cdim], &xi[p * cdim]); 571 PetscCall(PetscFECreateTabulation((PetscFE)obj, 1, numCIndices, xi, 0, &Tcoarse)); 572 /* Get elemMat entries by multiplying by weight */ 573 PetscCall(PetscArrayzero(elemMat, numCIndices * totDim)); 574 for (i = 0; i < numFIndices; ++i) { 575 for (p = 0; p < numCIndices; ++p) { 576 for (c = 0; c < Nc; ++c) { 577 /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */ 578 elemMat[p * numFIndices + i] += Tcoarse->T[0][(p * numFIndices + i) * Nc + c] * (useDeltaFunction ? 1.0 : detJ); 579 } 580 } 581 } 582 PetscCall(PetscTabulationDestroy(&Tcoarse)); 583 for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart; 584 if (0) PetscCall(DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat)); 585 /* Block diagonal */ 586 if (numCIndices) { 587 PetscBLASInt blasn, blask; 588 PetscScalar one = 1.0, zero = 0.0; 589 590 PetscCall(PetscBLASIntCast(numCIndices, &blasn)); 591 PetscCall(PetscBLASIntCast(numFIndices, &blask)); 592 PetscCallBLAS("BLASgemm", BLASgemm_("T", "N", &blasn, &blasn, &blask, &one, elemMat, &blask, elemMat, &blask, &zero, elemMatSq, &blasn)); 593 } 594 PetscCall(MatSetValues(mass, numCIndices, rowIDXs, numCIndices, rowIDXs, elemMatSq, ADD_VALUES)); 595 /* TODO Off-diagonal */ 596 PetscCall(PetscFree(cindices)); 597 PetscCall(DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL)); 598 } 599 PetscCall(DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 600 } 601 PetscCall(PetscFree4(elemMat, elemMatSq, rowIDXs, xi)); 602 PetscCall(PetscFree(adj)); 603 PetscCall(DMSwarmSortRestoreAccess(dmc)); 604 PetscCall(PetscFree3(v0, J, invJ)); 605 PetscCall(MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY)); 606 PetscCall(MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY)); 607 PetscFunctionReturn(PETSC_SUCCESS); 608 } 609 610 /*@C 611 DMSwarmCreateMassMatrixSquare - Creates the block-diagonal of the square, M^T_p M_p, of the particle mass matrix M_p 612 613 Collective on dmCoarse 614 615 Input parameters: 616 + dmCoarse - a DMSwarm 617 - dmFine - a DMPlex 618 619 Output parameter: 620 . mass - the square of the particle mass matrix 621 622 Level: advanced 623 624 Notes: 625 We only compute the block diagonal since this provides a good preconditioner and is completely local. It would be possible in the 626 future to compute the full normal equations. 627 628 .seealso: `DMCreateMassMatrix()` 629 @*/ 630 PetscErrorCode DMSwarmCreateMassMatrixSquare(DM dmCoarse, DM dmFine, Mat *mass) 631 { 632 PetscInt n; 633 void *ctx; 634 635 PetscFunctionBegin; 636 PetscCall(DMSwarmGetLocalSize(dmCoarse, &n)); 637 PetscCall(MatCreate(PetscObjectComm((PetscObject)dmCoarse), mass)); 638 PetscCall(MatSetSizes(*mass, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 639 PetscCall(MatSetType(*mass, dmCoarse->mattype)); 640 PetscCall(DMGetApplicationContext(dmFine, &ctx)); 641 642 PetscCall(DMSwarmComputeMassMatrixSquare_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx)); 643 PetscCall(MatViewFromOptions(*mass, NULL, "-mass_sq_mat_view")); 644 PetscFunctionReturn(PETSC_SUCCESS); 645 } 646 647 /*@C 648 DMSwarmCreateGlobalVectorFromField - Creates a Vec object sharing the array associated with a given field 649 650 Collective on dm 651 652 Input parameters: 653 + dm - a DMSwarm 654 - fieldname - the textual name given to a registered field 655 656 Output parameter: 657 . vec - the vector 658 659 Level: beginner 660 661 Notes: 662 The vector must be returned using a matching call to DMSwarmDestroyGlobalVectorFromField(). 663 664 .seealso: `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyGlobalVectorFromField()` 665 @*/ 666 PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm, const char fieldname[], Vec *vec) 667 { 668 MPI_Comm comm = PetscObjectComm((PetscObject)dm); 669 670 PetscFunctionBegin; 671 PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM); 672 PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec)); 673 PetscFunctionReturn(PETSC_SUCCESS); 674 } 675 676 /*@C 677 DMSwarmDestroyGlobalVectorFromField - Destroys the Vec object which share the array associated with a given field 678 679 Collective on dm 680 681 Input parameters: 682 + dm - a DMSwarm 683 - fieldname - the textual name given to a registered field 684 685 Output parameter: 686 . vec - the vector 687 688 Level: beginner 689 690 .seealso: `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateGlobalVectorFromField()` 691 @*/ 692 PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm, const char fieldname[], Vec *vec) 693 { 694 PetscFunctionBegin; 695 PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM); 696 PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec)); 697 PetscFunctionReturn(PETSC_SUCCESS); 698 } 699 700 /*@C 701 DMSwarmCreateLocalVectorFromField - Creates a Vec object sharing the array associated with a given field 702 703 Collective on dm 704 705 Input parameters: 706 + dm - a DMSwarm 707 - fieldname - the textual name given to a registered field 708 709 Output parameter: 710 . vec - the vector 711 712 Level: beginner 713 714 Notes: 715 The vector must be returned using a matching call to DMSwarmDestroyLocalVectorFromField(). 716 717 .seealso: `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmDestroyLocalVectorFromField()` 718 @*/ 719 PetscErrorCode DMSwarmCreateLocalVectorFromField(DM dm, const char fieldname[], Vec *vec) 720 { 721 MPI_Comm comm = PETSC_COMM_SELF; 722 723 PetscFunctionBegin; 724 PetscCall(DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec)); 725 PetscFunctionReturn(PETSC_SUCCESS); 726 } 727 728 /*@C 729 DMSwarmDestroyLocalVectorFromField - Destroys the Vec object which share the array associated with a given field 730 731 Collective on dm 732 733 Input parameters: 734 + dm - a DMSwarm 735 - fieldname - the textual name given to a registered field 736 737 Output parameter: 738 . vec - the vector 739 740 Level: beginner 741 742 .seealso: `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmCreateLocalVectorFromField()` 743 @*/ 744 PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM dm, const char fieldname[], Vec *vec) 745 { 746 PetscFunctionBegin; 747 PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM); 748 PetscCall(DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec)); 749 PetscFunctionReturn(PETSC_SUCCESS); 750 } 751 752 /*@C 753 DMSwarmInitializeFieldRegister - Initiates the registration of fields to a DMSwarm 754 755 Collective on dm 756 757 Input parameter: 758 . dm - a DMSwarm 759 760 Level: beginner 761 762 Notes: 763 After all fields have been registered, you must call DMSwarmFinalizeFieldRegister(). 764 765 .seealso: `DMSwarmFinalizeFieldRegister()`, `DMSwarmRegisterPetscDatatypeField()`, 766 `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()` 767 @*/ 768 PetscErrorCode DMSwarmInitializeFieldRegister(DM dm) 769 { 770 DM_Swarm *swarm = (DM_Swarm *)dm->data; 771 772 PetscFunctionBegin; 773 if (!swarm->field_registration_initialized) { 774 swarm->field_registration_initialized = PETSC_TRUE; 775 PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmField_pid, 1, PETSC_INT64)); /* unique identifier */ 776 PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmField_rank, 1, PETSC_INT)); /* used for communication */ 777 } 778 PetscFunctionReturn(PETSC_SUCCESS); 779 } 780 781 /*@ 782 DMSwarmFinalizeFieldRegister - Finalizes the registration of fields to a DMSwarm 783 784 Collective on dm 785 786 Input parameter: 787 . dm - a DMSwarm 788 789 Level: beginner 790 791 Notes: 792 After DMSwarmFinalizeFieldRegister() has been called, no new fields can be defined on the DMSwarm. 793 794 .seealso: `DMSwarmInitializeFieldRegister()`, `DMSwarmRegisterPetscDatatypeField()`, 795 `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()` 796 @*/ 797 PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm) 798 { 799 DM_Swarm *swarm = (DM_Swarm *)dm->data; 800 801 PetscFunctionBegin; 802 if (!swarm->field_registration_finalized) PetscCall(DMSwarmDataBucketFinalize(swarm->db)); 803 swarm->field_registration_finalized = PETSC_TRUE; 804 PetscFunctionReturn(PETSC_SUCCESS); 805 } 806 807 /*@ 808 DMSwarmSetLocalSizes - Sets the length of all registered fields on the DMSwarm 809 810 Not collective 811 812 Input parameters: 813 + dm - a DMSwarm 814 . nlocal - the length of each registered field 815 - buffer - the length of the buffer used to efficient dynamic re-sizing 816 817 Level: beginner 818 819 .seealso: `DMSwarmGetLocalSize()` 820 @*/ 821 PetscErrorCode DMSwarmSetLocalSizes(DM dm, PetscInt nlocal, PetscInt buffer) 822 { 823 DM_Swarm *swarm = (DM_Swarm *)dm->data; 824 825 PetscFunctionBegin; 826 PetscCall(PetscLogEventBegin(DMSWARM_SetSizes, 0, 0, 0, 0)); 827 PetscCall(DMSwarmDataBucketSetSizes(swarm->db, nlocal, buffer)); 828 PetscCall(PetscLogEventEnd(DMSWARM_SetSizes, 0, 0, 0, 0)); 829 PetscFunctionReturn(PETSC_SUCCESS); 830 } 831 832 /*@ 833 DMSwarmSetCellDM - Attaches a DM to a DMSwarm 834 835 Collective on dm 836 837 Input parameters: 838 + dm - a DMSwarm 839 - dmcell - the DM to attach to the DMSwarm 840 841 Level: beginner 842 843 Notes: 844 The attached DM (dmcell) will be queried for point location and 845 neighbor MPI-rank information if DMSwarmMigrate() is called. 846 847 .seealso: `DMSwarmSetType()`, `DMSwarmGetCellDM()`, `DMSwarmMigrate()` 848 @*/ 849 PetscErrorCode DMSwarmSetCellDM(DM dm, DM dmcell) 850 { 851 DM_Swarm *swarm = (DM_Swarm *)dm->data; 852 853 PetscFunctionBegin; 854 swarm->dmcell = dmcell; 855 PetscFunctionReturn(PETSC_SUCCESS); 856 } 857 858 /*@ 859 DMSwarmGetCellDM - Fetches the attached cell DM 860 861 Collective on dm 862 863 Input parameter: 864 . dm - a DMSwarm 865 866 Output parameter: 867 . dmcell - the DM which was attached to the DMSwarm 868 869 Level: beginner 870 871 .seealso: `DMSwarmSetCellDM()` 872 @*/ 873 PetscErrorCode DMSwarmGetCellDM(DM dm, DM *dmcell) 874 { 875 DM_Swarm *swarm = (DM_Swarm *)dm->data; 876 877 PetscFunctionBegin; 878 *dmcell = swarm->dmcell; 879 PetscFunctionReturn(PETSC_SUCCESS); 880 } 881 882 /*@ 883 DMSwarmGetLocalSize - Retrieves the local length of fields registered 884 885 Not collective 886 887 Input parameter: 888 . dm - a DMSwarm 889 890 Output parameter: 891 . nlocal - the length of each registered field 892 893 Level: beginner 894 895 .seealso: `DMSwarmGetSize()`, `DMSwarmSetLocalSizes()` 896 @*/ 897 PetscErrorCode DMSwarmGetLocalSize(DM dm, PetscInt *nlocal) 898 { 899 DM_Swarm *swarm = (DM_Swarm *)dm->data; 900 901 PetscFunctionBegin; 902 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, nlocal, NULL, NULL)); 903 PetscFunctionReturn(PETSC_SUCCESS); 904 } 905 906 /*@ 907 DMSwarmGetSize - Retrieves the total length of fields registered 908 909 Collective on dm 910 911 Input parameter: 912 . dm - a DMSwarm 913 914 Output parameter: 915 . n - the total length of each registered field 916 917 Level: beginner 918 919 Note: 920 This calls MPI_Allreduce upon each call (inefficient but safe) 921 922 .seealso: `DMSwarmGetLocalSize()`, `DMSwarmSetLocalSizes()` 923 @*/ 924 PetscErrorCode DMSwarmGetSize(DM dm, PetscInt *n) 925 { 926 DM_Swarm *swarm = (DM_Swarm *)dm->data; 927 PetscInt nlocal; 928 929 PetscFunctionBegin; 930 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &nlocal, NULL, NULL)); 931 PetscCallMPI(MPI_Allreduce(&nlocal, n, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm))); 932 PetscFunctionReturn(PETSC_SUCCESS); 933 } 934 935 /*@C 936 DMSwarmRegisterPetscDatatypeField - Register a field to a DMSwarm with a native PETSc data type 937 938 Collective on dm 939 940 Input parameters: 941 + dm - a DMSwarm 942 . fieldname - the textual name to identify this field 943 . blocksize - the number of each data type 944 - type - a valid PETSc data type (PETSC_CHAR, PETSC_SHORT, PETSC_INT, PETSC_FLOAT, PETSC_REAL, PETSC_LONG) 945 946 Level: beginner 947 948 Notes: 949 The textual name for each registered field must be unique. 950 951 .seealso: `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()` 952 @*/ 953 PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm, const char fieldname[], PetscInt blocksize, PetscDataType type) 954 { 955 DM_Swarm *swarm = (DM_Swarm *)dm->data; 956 size_t size; 957 958 PetscFunctionBegin; 959 PetscCheck(swarm->field_registration_initialized, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Must call DMSwarmInitializeFieldRegister() first"); 960 PetscCheck(!swarm->field_registration_finalized, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first"); 961 962 PetscCheck(type != PETSC_OBJECT, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}"); 963 PetscCheck(type != PETSC_FUNCTION, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}"); 964 PetscCheck(type != PETSC_STRING, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}"); 965 PetscCheck(type != PETSC_STRUCT, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}"); 966 PetscCheck(type != PETSC_DATATYPE_UNKNOWN, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Valid for {char,short,int,long,float,double}"); 967 968 PetscCall(PetscDataTypeGetSize(type, &size)); 969 /* Load a specific data type into data bucket, specifying textual name and its size in bytes */ 970 PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterPetscDatatypeField", fieldname, blocksize * size, NULL)); 971 { 972 DMSwarmDataField gfield; 973 974 PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield)); 975 PetscCall(DMSwarmDataFieldSetBlockSize(gfield, blocksize)); 976 } 977 swarm->db->field[swarm->db->nfields - 1]->petsc_type = type; 978 PetscFunctionReturn(PETSC_SUCCESS); 979 } 980 981 /*@C 982 DMSwarmRegisterUserStructField - Register a user defined struct to a DMSwarm 983 984 Collective on dm 985 986 Input parameters: 987 + dm - a DMSwarm 988 . fieldname - the textual name to identify this field 989 - size - the size in bytes of the user struct of each data type 990 991 Level: beginner 992 993 Notes: 994 The textual name for each registered field must be unique. 995 996 .seealso: `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmRegisterUserDatatypeField()` 997 @*/ 998 PetscErrorCode DMSwarmRegisterUserStructField(DM dm, const char fieldname[], size_t size) 999 { 1000 DM_Swarm *swarm = (DM_Swarm *)dm->data; 1001 1002 PetscFunctionBegin; 1003 PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterUserStructField", fieldname, size, NULL)); 1004 swarm->db->field[swarm->db->nfields - 1]->petsc_type = PETSC_STRUCT; 1005 PetscFunctionReturn(PETSC_SUCCESS); 1006 } 1007 1008 /*@C 1009 DMSwarmRegisterUserDatatypeField - Register a user defined data type to a DMSwarm 1010 1011 Collective on dm 1012 1013 Input parameters: 1014 + dm - a DMSwarm 1015 . fieldname - the textual name to identify this field 1016 . size - the size in bytes of the user data type 1017 - blocksize - the number of each data type 1018 1019 Level: beginner 1020 1021 Notes: 1022 The textual name for each registered field must be unique. 1023 1024 .seealso: `DMSwarmRegisterPetscDatatypeField()`, `DMSwarmRegisterUserStructField()`, `DMSwarmRegisterUserDatatypeField()` 1025 @*/ 1026 PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm, const char fieldname[], size_t size, PetscInt blocksize) 1027 { 1028 DM_Swarm *swarm = (DM_Swarm *)dm->data; 1029 1030 PetscFunctionBegin; 1031 PetscCall(DMSwarmDataBucketRegisterField(swarm->db, "DMSwarmRegisterUserDatatypeField", fieldname, blocksize * size, NULL)); 1032 { 1033 DMSwarmDataField gfield; 1034 1035 PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield)); 1036 PetscCall(DMSwarmDataFieldSetBlockSize(gfield, blocksize)); 1037 } 1038 swarm->db->field[swarm->db->nfields - 1]->petsc_type = PETSC_DATATYPE_UNKNOWN; 1039 PetscFunctionReturn(PETSC_SUCCESS); 1040 } 1041 1042 /*@C 1043 DMSwarmGetField - Get access to the underlying array storing all entries associated with a registered field 1044 1045 Not collective 1046 1047 Input parameters: 1048 + dm - a DMSwarm 1049 - fieldname - the textual name to identify this field 1050 1051 Output parameters: 1052 + blocksize - the number of each data type 1053 . type - the data type 1054 - data - pointer to raw array 1055 1056 Level: beginner 1057 1058 Notes: 1059 The array must be returned using a matching call to DMSwarmRestoreField(). 1060 1061 .seealso: `DMSwarmRestoreField()` 1062 @*/ 1063 PetscErrorCode DMSwarmGetField(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type, void **data) 1064 { 1065 DM_Swarm *swarm = (DM_Swarm *)dm->data; 1066 DMSwarmDataField gfield; 1067 1068 PetscFunctionBegin; 1069 PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM); 1070 if (!swarm->issetup) PetscCall(DMSetUp(dm)); 1071 PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield)); 1072 PetscCall(DMSwarmDataFieldGetAccess(gfield)); 1073 PetscCall(DMSwarmDataFieldGetEntries(gfield, data)); 1074 if (blocksize) *blocksize = gfield->bs; 1075 if (type) *type = gfield->petsc_type; 1076 PetscFunctionReturn(PETSC_SUCCESS); 1077 } 1078 1079 /*@C 1080 DMSwarmRestoreField - Restore access to the underlying array storing all entries associated with a registered field 1081 1082 Not collective 1083 1084 Input parameters: 1085 + dm - a DMSwarm 1086 - fieldname - the textual name to identify this field 1087 1088 Output parameters: 1089 + blocksize - the number of each data type 1090 . type - the data type 1091 - data - pointer to raw array 1092 1093 Level: beginner 1094 1095 Notes: 1096 The user must call DMSwarmGetField() prior to calling DMSwarmRestoreField(). 1097 1098 .seealso: `DMSwarmGetField()` 1099 @*/ 1100 PetscErrorCode DMSwarmRestoreField(DM dm, const char fieldname[], PetscInt *blocksize, PetscDataType *type, void **data) 1101 { 1102 DM_Swarm *swarm = (DM_Swarm *)dm->data; 1103 DMSwarmDataField gfield; 1104 1105 PetscFunctionBegin; 1106 PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMSWARM); 1107 PetscCall(DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield)); 1108 PetscCall(DMSwarmDataFieldRestoreAccess(gfield)); 1109 if (data) *data = NULL; 1110 PetscFunctionReturn(PETSC_SUCCESS); 1111 } 1112 1113 /*@ 1114 DMSwarmAddPoint - Add space for one new point in the DMSwarm 1115 1116 Not collective 1117 1118 Input parameter: 1119 . dm - a DMSwarm 1120 1121 Level: beginner 1122 1123 Notes: 1124 The new point will have all fields initialized to zero. 1125 1126 .seealso: `DMSwarmAddNPoints()` 1127 @*/ 1128 PetscErrorCode DMSwarmAddPoint(DM dm) 1129 { 1130 DM_Swarm *swarm = (DM_Swarm *)dm->data; 1131 1132 PetscFunctionBegin; 1133 if (!swarm->issetup) PetscCall(DMSetUp(dm)); 1134 PetscCall(PetscLogEventBegin(DMSWARM_AddPoints, 0, 0, 0, 0)); 1135 PetscCall(DMSwarmDataBucketAddPoint(swarm->db)); 1136 PetscCall(PetscLogEventEnd(DMSWARM_AddPoints, 0, 0, 0, 0)); 1137 PetscFunctionReturn(PETSC_SUCCESS); 1138 } 1139 1140 /*@ 1141 DMSwarmAddNPoints - Add space for a number of new points in the DMSwarm 1142 1143 Not collective 1144 1145 Input parameters: 1146 + dm - a DMSwarm 1147 - npoints - the number of new points to add 1148 1149 Level: beginner 1150 1151 Notes: 1152 The new point will have all fields initialized to zero. 1153 1154 .seealso: `DMSwarmAddPoint()` 1155 @*/ 1156 PetscErrorCode DMSwarmAddNPoints(DM dm, PetscInt npoints) 1157 { 1158 DM_Swarm *swarm = (DM_Swarm *)dm->data; 1159 PetscInt nlocal; 1160 1161 PetscFunctionBegin; 1162 PetscCall(PetscLogEventBegin(DMSWARM_AddPoints, 0, 0, 0, 0)); 1163 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &nlocal, NULL, NULL)); 1164 nlocal = nlocal + npoints; 1165 PetscCall(DMSwarmDataBucketSetSizes(swarm->db, nlocal, DMSWARM_DATA_BUCKET_BUFFER_DEFAULT)); 1166 PetscCall(PetscLogEventEnd(DMSWARM_AddPoints, 0, 0, 0, 0)); 1167 PetscFunctionReturn(PETSC_SUCCESS); 1168 } 1169 1170 /*@ 1171 DMSwarmRemovePoint - Remove the last point from the DMSwarm 1172 1173 Not collective 1174 1175 Input parameter: 1176 . dm - a DMSwarm 1177 1178 Level: beginner 1179 1180 .seealso: `DMSwarmRemovePointAtIndex()` 1181 @*/ 1182 PetscErrorCode DMSwarmRemovePoint(DM dm) 1183 { 1184 DM_Swarm *swarm = (DM_Swarm *)dm->data; 1185 1186 PetscFunctionBegin; 1187 PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints, 0, 0, 0, 0)); 1188 PetscCall(DMSwarmDataBucketRemovePoint(swarm->db)); 1189 PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints, 0, 0, 0, 0)); 1190 PetscFunctionReturn(PETSC_SUCCESS); 1191 } 1192 1193 /*@ 1194 DMSwarmRemovePointAtIndex - Removes a specific point from the DMSwarm 1195 1196 Not collective 1197 1198 Input parameters: 1199 + dm - a DMSwarm 1200 - idx - index of point to remove 1201 1202 Level: beginner 1203 1204 .seealso: `DMSwarmRemovePoint()` 1205 @*/ 1206 PetscErrorCode DMSwarmRemovePointAtIndex(DM dm, PetscInt idx) 1207 { 1208 DM_Swarm *swarm = (DM_Swarm *)dm->data; 1209 1210 PetscFunctionBegin; 1211 PetscCall(PetscLogEventBegin(DMSWARM_RemovePoints, 0, 0, 0, 0)); 1212 PetscCall(DMSwarmDataBucketRemovePointAtIndex(swarm->db, idx)); 1213 PetscCall(PetscLogEventEnd(DMSWARM_RemovePoints, 0, 0, 0, 0)); 1214 PetscFunctionReturn(PETSC_SUCCESS); 1215 } 1216 1217 /*@ 1218 DMSwarmCopyPoint - Copy point pj to point pi in the DMSwarm 1219 1220 Not collective 1221 1222 Input parameters: 1223 + dm - a DMSwarm 1224 . pi - the index of the point to copy 1225 - pj - the point index where the copy should be located 1226 1227 Level: beginner 1228 1229 .seealso: `DMSwarmRemovePoint()` 1230 @*/ 1231 PetscErrorCode DMSwarmCopyPoint(DM dm, PetscInt pi, PetscInt pj) 1232 { 1233 DM_Swarm *swarm = (DM_Swarm *)dm->data; 1234 1235 PetscFunctionBegin; 1236 if (!swarm->issetup) PetscCall(DMSetUp(dm)); 1237 PetscCall(DMSwarmDataBucketCopyPoint(swarm->db, pi, swarm->db, pj)); 1238 PetscFunctionReturn(PETSC_SUCCESS); 1239 } 1240 1241 PetscErrorCode DMSwarmMigrate_Basic(DM dm, PetscBool remove_sent_points) 1242 { 1243 PetscFunctionBegin; 1244 PetscCall(DMSwarmMigrate_Push_Basic(dm, remove_sent_points)); 1245 PetscFunctionReturn(PETSC_SUCCESS); 1246 } 1247 1248 /*@ 1249 DMSwarmMigrate - Relocates points defined in the DMSwarm to other MPI-ranks 1250 1251 Collective on dm 1252 1253 Input parameters: 1254 + dm - the DMSwarm 1255 - remove_sent_points - flag indicating if sent points should be removed from the current MPI-rank 1256 1257 Notes: 1258 The DM will be modified to accommodate received points. 1259 If remove_sent_points = PETSC_TRUE, any points that were sent will be removed from the DM. 1260 Different styles of migration are supported. See DMSwarmSetMigrateType(). 1261 1262 Level: advanced 1263 1264 .seealso: `DMSwarmSetMigrateType()` 1265 @*/ 1266 PetscErrorCode DMSwarmMigrate(DM dm, PetscBool remove_sent_points) 1267 { 1268 DM_Swarm *swarm = (DM_Swarm *)dm->data; 1269 1270 PetscFunctionBegin; 1271 PetscCall(PetscLogEventBegin(DMSWARM_Migrate, 0, 0, 0, 0)); 1272 switch (swarm->migrate_type) { 1273 case DMSWARM_MIGRATE_BASIC: 1274 PetscCall(DMSwarmMigrate_Basic(dm, remove_sent_points)); 1275 break; 1276 case DMSWARM_MIGRATE_DMCELLNSCATTER: 1277 PetscCall(DMSwarmMigrate_CellDMScatter(dm, remove_sent_points)); 1278 break; 1279 case DMSWARM_MIGRATE_DMCELLEXACT: 1280 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE_DMCELLEXACT not implemented"); 1281 case DMSWARM_MIGRATE_USER: 1282 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE_USER not implemented"); 1283 default: 1284 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_MIGRATE type unknown"); 1285 } 1286 PetscCall(PetscLogEventEnd(DMSWARM_Migrate, 0, 0, 0, 0)); 1287 PetscCall(DMClearGlobalVectors(dm)); 1288 PetscFunctionReturn(PETSC_SUCCESS); 1289 } 1290 1291 PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm, PetscInt *globalsize); 1292 1293 /* 1294 DMSwarmCollectViewCreate 1295 1296 * Applies a collection method and gathers point neighbour points into dm 1297 1298 Notes: 1299 Users should call DMSwarmCollectViewDestroy() after 1300 they have finished computations associated with the collected points 1301 */ 1302 1303 /*@ 1304 DMSwarmCollectViewCreate - Applies a collection method and gathers points 1305 in neighbour ranks into the DMSwarm 1306 1307 Collective on dm 1308 1309 Input parameter: 1310 . dm - the DMSwarm 1311 1312 Notes: 1313 Users should call DMSwarmCollectViewDestroy() after 1314 they have finished computations associated with the collected points 1315 Different collect methods are supported. See DMSwarmSetCollectType(). 1316 1317 Level: advanced 1318 1319 .seealso: `DMSwarmCollectViewDestroy()`, `DMSwarmSetCollectType()` 1320 @*/ 1321 PetscErrorCode DMSwarmCollectViewCreate(DM dm) 1322 { 1323 DM_Swarm *swarm = (DM_Swarm *)dm->data; 1324 PetscInt ng; 1325 1326 PetscFunctionBegin; 1327 PetscCheck(!swarm->collect_view_active, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "CollectView currently active"); 1328 PetscCall(DMSwarmGetLocalSize(dm, &ng)); 1329 switch (swarm->collect_type) { 1330 case DMSWARM_COLLECT_BASIC: 1331 PetscCall(DMSwarmMigrate_GlobalToLocal_Basic(dm, &ng)); 1332 break; 1333 case DMSWARM_COLLECT_DMDABOUNDINGBOX: 1334 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT_DMDABOUNDINGBOX not implemented"); 1335 case DMSWARM_COLLECT_GENERAL: 1336 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT_GENERAL not implemented"); 1337 default: 1338 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "DMSWARM_COLLECT type unknown"); 1339 } 1340 swarm->collect_view_active = PETSC_TRUE; 1341 swarm->collect_view_reset_nlocal = ng; 1342 PetscFunctionReturn(PETSC_SUCCESS); 1343 } 1344 1345 /*@ 1346 DMSwarmCollectViewDestroy - Resets the DMSwarm to the size prior to calling DMSwarmCollectViewCreate() 1347 1348 Collective on dm 1349 1350 Input parameters: 1351 . dm - the DMSwarm 1352 1353 Notes: 1354 Users should call DMSwarmCollectViewCreate() before this function is called. 1355 1356 Level: advanced 1357 1358 .seealso: `DMSwarmCollectViewCreate()`, `DMSwarmSetCollectType()` 1359 @*/ 1360 PetscErrorCode DMSwarmCollectViewDestroy(DM dm) 1361 { 1362 DM_Swarm *swarm = (DM_Swarm *)dm->data; 1363 1364 PetscFunctionBegin; 1365 PetscCheck(swarm->collect_view_active, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "CollectView is currently not active"); 1366 PetscCall(DMSwarmSetLocalSizes(dm, swarm->collect_view_reset_nlocal, -1)); 1367 swarm->collect_view_active = PETSC_FALSE; 1368 PetscFunctionReturn(PETSC_SUCCESS); 1369 } 1370 1371 PetscErrorCode DMSwarmSetUpPIC(DM dm) 1372 { 1373 PetscInt dim; 1374 1375 PetscFunctionBegin; 1376 PetscCall(DMSwarmSetNumSpecies(dm, 1)); 1377 PetscCall(DMGetDimension(dm, &dim)); 1378 PetscCheck(dim >= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Dimension must be 1,2,3 - found %" PetscInt_FMT, dim); 1379 PetscCheck(dim <= 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Dimension must be 1,2,3 - found %" PetscInt_FMT, dim); 1380 PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmPICField_coor, dim, PETSC_DOUBLE)); 1381 PetscCall(DMSwarmRegisterPetscDatatypeField(dm, DMSwarmPICField_cellid, 1, PETSC_INT)); 1382 PetscFunctionReturn(PETSC_SUCCESS); 1383 } 1384 1385 /*@ 1386 DMSwarmSetPointCoordinatesRandom - Sets initial coordinates for particles in each cell 1387 1388 Collective on dm 1389 1390 Input parameters: 1391 + dm - the DMSwarm 1392 - Npc - The number of particles per cell in the cell DM 1393 1394 Notes: 1395 The user must use DMSwarmSetCellDM() to set the cell DM first. The particles are placed randomly inside each cell. If only 1396 one particle is in each cell, it is placed at the centroid. 1397 1398 Level: intermediate 1399 1400 .seealso: `DMSwarmSetCellDM()` 1401 @*/ 1402 PetscErrorCode DMSwarmSetPointCoordinatesRandom(DM dm, PetscInt Npc) 1403 { 1404 DM cdm; 1405 PetscRandom rnd; 1406 DMPolytopeType ct; 1407 PetscBool simplex; 1408 PetscReal *centroid, *coords, *xi0, *v0, *J, *invJ, detJ; 1409 PetscInt dim, d, cStart, cEnd, c, p; 1410 1411 PetscFunctionBeginUser; 1412 PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd)); 1413 PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0)); 1414 PetscCall(PetscRandomSetType(rnd, PETSCRAND48)); 1415 1416 PetscCall(DMSwarmGetCellDM(dm, &cdm)); 1417 PetscCall(DMGetDimension(cdm, &dim)); 1418 PetscCall(DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd)); 1419 PetscCall(DMPlexGetCellType(cdm, cStart, &ct)); 1420 simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 1421 1422 PetscCall(PetscMalloc5(dim, ¢roid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ)); 1423 for (d = 0; d < dim; ++d) xi0[d] = -1.0; 1424 PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 1425 for (c = cStart; c < cEnd; ++c) { 1426 if (Npc == 1) { 1427 PetscCall(DMPlexComputeCellGeometryFVM(cdm, c, NULL, centroid, NULL)); 1428 for (d = 0; d < dim; ++d) coords[c * dim + d] = centroid[d]; 1429 } else { 1430 PetscCall(DMPlexComputeCellGeometryFEM(cdm, c, NULL, v0, J, invJ, &detJ)); /* affine */ 1431 for (p = 0; p < Npc; ++p) { 1432 const PetscInt n = c * Npc + p; 1433 PetscReal sum = 0.0, refcoords[3]; 1434 1435 for (d = 0; d < dim; ++d) { 1436 PetscCall(PetscRandomGetValueReal(rnd, &refcoords[d])); 1437 sum += refcoords[d]; 1438 } 1439 if (simplex && sum > 0.0) 1440 for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum; 1441 CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n * dim]); 1442 } 1443 } 1444 } 1445 PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 1446 PetscCall(PetscFree5(centroid, xi0, v0, J, invJ)); 1447 PetscCall(PetscRandomDestroy(&rnd)); 1448 PetscFunctionReturn(PETSC_SUCCESS); 1449 } 1450 1451 /*@ 1452 DMSwarmSetType - Set particular flavor of DMSwarm 1453 1454 Collective on dm 1455 1456 Input parameters: 1457 + dm - the DMSwarm 1458 - stype - the DMSwarm type (e.g. DMSWARM_PIC) 1459 1460 Level: advanced 1461 1462 .seealso: `DMSwarmSetMigrateType()`, `DMSwarmSetCollectType()`, `DMSwarmType`, `DMSWARM_PIC`, `DMSWARM_BASIC` 1463 @*/ 1464 PetscErrorCode DMSwarmSetType(DM dm, DMSwarmType stype) 1465 { 1466 DM_Swarm *swarm = (DM_Swarm *)dm->data; 1467 1468 PetscFunctionBegin; 1469 swarm->swarm_type = stype; 1470 if (swarm->swarm_type == DMSWARM_PIC) PetscCall(DMSwarmSetUpPIC(dm)); 1471 PetscFunctionReturn(PETSC_SUCCESS); 1472 } 1473 1474 PetscErrorCode DMSetup_Swarm(DM dm) 1475 { 1476 DM_Swarm *swarm = (DM_Swarm *)dm->data; 1477 PetscMPIInt rank; 1478 PetscInt p, npoints, *rankval; 1479 1480 PetscFunctionBegin; 1481 if (swarm->issetup) PetscFunctionReturn(PETSC_SUCCESS); 1482 swarm->issetup = PETSC_TRUE; 1483 1484 if (swarm->swarm_type == DMSWARM_PIC) { 1485 /* check dmcell exists */ 1486 PetscCheck(swarm->dmcell, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires you call DMSwarmSetCellDM"); 1487 1488 if (swarm->dmcell->ops->locatepointssubdomain) { 1489 /* check methods exists for exact ownership identificiation */ 1490 PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->ops->LocatePointsSubdomain\n")); 1491 swarm->migrate_type = DMSWARM_MIGRATE_DMCELLEXACT; 1492 } else { 1493 /* check methods exist for point location AND rank neighbor identification */ 1494 if (swarm->dmcell->ops->locatepoints) { 1495 PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->LocatePoints\n")); 1496 } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires the method CellDM->ops->locatepoints be defined"); 1497 1498 if (swarm->dmcell->ops->getneighbors) { 1499 PetscCall(PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->GetNeigbors\n")); 1500 } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "DMSWARM_PIC requires the method CellDM->ops->getneighbors be defined"); 1501 1502 swarm->migrate_type = DMSWARM_MIGRATE_DMCELLNSCATTER; 1503 } 1504 } 1505 1506 PetscCall(DMSwarmFinalizeFieldRegister(dm)); 1507 1508 /* check some fields were registered */ 1509 PetscCheck(swarm->db->nfields > 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "At least one field user must be registered via DMSwarmRegisterXXX()"); 1510 1511 /* check local sizes were set */ 1512 PetscCheck(swarm->db->L != -1, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Local sizes must be set via DMSwarmSetLocalSizes()"); 1513 1514 /* initialize values in pid and rank placeholders */ 1515 /* TODO: [pid - use MPI_Scan] */ 1516 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 1517 PetscCall(DMSwarmDataBucketGetSizes(swarm->db, &npoints, NULL, NULL)); 1518 PetscCall(DMSwarmGetField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 1519 for (p = 0; p < npoints; p++) rankval[p] = (PetscInt)rank; 1520 PetscCall(DMSwarmRestoreField(dm, DMSwarmField_rank, NULL, NULL, (void **)&rankval)); 1521 PetscFunctionReturn(PETSC_SUCCESS); 1522 } 1523 1524 extern PetscErrorCode DMSwarmSortDestroy(DMSwarmSort *_ctx); 1525 1526 PetscErrorCode DMDestroy_Swarm(DM dm) 1527 { 1528 DM_Swarm *swarm = (DM_Swarm *)dm->data; 1529 1530 PetscFunctionBegin; 1531 if (--swarm->refct > 0) PetscFunctionReturn(PETSC_SUCCESS); 1532 PetscCall(DMSwarmDataBucketDestroy(&swarm->db)); 1533 if (swarm->sort_context) PetscCall(DMSwarmSortDestroy(&swarm->sort_context)); 1534 PetscCall(PetscFree(swarm)); 1535 PetscFunctionReturn(PETSC_SUCCESS); 1536 } 1537 1538 PetscErrorCode DMSwarmView_Draw(DM dm, PetscViewer viewer) 1539 { 1540 DM cdm; 1541 PetscDraw draw; 1542 PetscReal *coords, oldPause, radius = 0.01; 1543 PetscInt Np, p, bs; 1544 1545 PetscFunctionBegin; 1546 PetscCall(PetscOptionsGetReal(NULL, ((PetscObject)dm)->prefix, "-dm_view_swarm_radius", &radius, NULL)); 1547 PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 1548 PetscCall(DMSwarmGetCellDM(dm, &cdm)); 1549 PetscCall(PetscDrawGetPause(draw, &oldPause)); 1550 PetscCall(PetscDrawSetPause(draw, 0.0)); 1551 PetscCall(DMView(cdm, viewer)); 1552 PetscCall(PetscDrawSetPause(draw, oldPause)); 1553 1554 PetscCall(DMSwarmGetLocalSize(dm, &Np)); 1555 PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, &bs, NULL, (void **)&coords)); 1556 for (p = 0; p < Np; ++p) { 1557 const PetscInt i = p * bs; 1558 1559 PetscCall(PetscDrawEllipse(draw, coords[i], coords[i + 1], radius, radius, PETSC_DRAW_BLUE)); 1560 } 1561 PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, &bs, NULL, (void **)&coords)); 1562 PetscCall(PetscDrawFlush(draw)); 1563 PetscCall(PetscDrawPause(draw)); 1564 PetscCall(PetscDrawSave(draw)); 1565 PetscFunctionReturn(PETSC_SUCCESS); 1566 } 1567 1568 static PetscErrorCode DMView_Swarm_Ascii(DM dm, PetscViewer viewer) 1569 { 1570 PetscViewerFormat format; 1571 PetscInt *sizes; 1572 PetscInt dim, Np, maxSize = 17; 1573 MPI_Comm comm; 1574 PetscMPIInt rank, size, p; 1575 const char *name; 1576 1577 PetscFunctionBegin; 1578 PetscCall(PetscViewerGetFormat(viewer, &format)); 1579 PetscCall(DMGetDimension(dm, &dim)); 1580 PetscCall(DMSwarmGetLocalSize(dm, &Np)); 1581 PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 1582 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1583 PetscCallMPI(MPI_Comm_size(comm, &size)); 1584 PetscCall(PetscObjectGetName((PetscObject)dm, &name)); 1585 if (name) PetscCall(PetscViewerASCIIPrintf(viewer, "%s in %" PetscInt_FMT " dimension%s:\n", name, dim, dim == 1 ? "" : "s")); 1586 else PetscCall(PetscViewerASCIIPrintf(viewer, "Swarm in %" PetscInt_FMT " dimension%s:\n", dim, dim == 1 ? "" : "s")); 1587 if (size < maxSize) PetscCall(PetscCalloc1(size, &sizes)); 1588 else PetscCall(PetscCalloc1(3, &sizes)); 1589 if (size < maxSize) { 1590 PetscCallMPI(MPI_Gather(&Np, 1, MPIU_INT, sizes, 1, MPIU_INT, 0, comm)); 1591 PetscCall(PetscViewerASCIIPrintf(viewer, " Number of particles per rank:")); 1592 for (p = 0; p < size; ++p) { 1593 if (rank == 0) PetscCall(PetscViewerASCIIPrintf(viewer, " %" PetscInt_FMT, sizes[p])); 1594 } 1595 } else { 1596 PetscInt locMinMax[2] = {Np, Np}; 1597 1598 PetscCall(PetscGlobalMinMaxInt(comm, locMinMax, sizes)); 1599 PetscCall(PetscViewerASCIIPrintf(viewer, " Min/Max of particles per rank: %" PetscInt_FMT "/%" PetscInt_FMT, sizes[0], sizes[1])); 1600 } 1601 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 1602 PetscCall(PetscFree(sizes)); 1603 if (format == PETSC_VIEWER_ASCII_INFO) { 1604 PetscInt *cell; 1605 1606 PetscCall(PetscViewerASCIIPrintf(viewer, " Cells containing each particle:\n")); 1607 PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 1608 PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&cell)); 1609 for (p = 0; p < Np; ++p) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " p%d: %" PetscInt_FMT "\n", p, cell[p])); 1610 PetscCall(PetscViewerFlush(viewer)); 1611 PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 1612 PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&cell)); 1613 } 1614 PetscFunctionReturn(PETSC_SUCCESS); 1615 } 1616 1617 PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer) 1618 { 1619 DM_Swarm *swarm = (DM_Swarm *)dm->data; 1620 PetscBool iascii, ibinary, isvtk, isdraw; 1621 #if defined(PETSC_HAVE_HDF5) 1622 PetscBool ishdf5; 1623 #endif 1624 1625 PetscFunctionBegin; 1626 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1627 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1628 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1629 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary)); 1630 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk)); 1631 #if defined(PETSC_HAVE_HDF5) 1632 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5)); 1633 #endif 1634 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 1635 if (iascii) { 1636 PetscViewerFormat format; 1637 1638 PetscCall(PetscViewerGetFormat(viewer, &format)); 1639 switch (format) { 1640 case PETSC_VIEWER_ASCII_INFO_DETAIL: 1641 PetscCall(DMSwarmDataBucketView(PetscObjectComm((PetscObject)dm), swarm->db, NULL, DATABUCKET_VIEW_STDOUT)); 1642 break; 1643 default: 1644 PetscCall(DMView_Swarm_Ascii(dm, viewer)); 1645 } 1646 } else { 1647 #if defined(PETSC_HAVE_HDF5) 1648 if (ishdf5) PetscCall(DMSwarmView_HDF5(dm, viewer)); 1649 #endif 1650 if (isdraw) PetscCall(DMSwarmView_Draw(dm, viewer)); 1651 } 1652 PetscFunctionReturn(PETSC_SUCCESS); 1653 } 1654 1655 /*@C 1656 DMSwarmGetCellSwarm - Extracts a single cell from the DMSwarm object, returns it as a single cell DMSwarm. 1657 The cell DM is filtered for fields of that cell, and the filtered DM is used as the cell DM of the new swarm object. 1658 1659 Important: Changes to this cell of the swarm will be lost if they are made prior to restoring this cell. 1660 1661 Noncollective 1662 1663 Input parameters: 1664 + sw - the DMSwarm 1665 . cellID - the integer id of the cell to be extracted and filtered 1666 - cellswarm - The DMSwarm to receive the cell 1667 1668 Level: beginner 1669 1670 Notes: 1671 This presently only supports DMSWARM_PIC type 1672 1673 Should be restored with DMSwarmRestoreCellSwarm() 1674 1675 .seealso: `DMSwarmRestoreCellSwarm()` 1676 @*/ 1677 PETSC_EXTERN PetscErrorCode DMSwarmGetCellSwarm(DM sw, PetscInt cellID, DM cellswarm) 1678 { 1679 DM_Swarm *original = (DM_Swarm *)sw->data; 1680 DMLabel label; 1681 DM dmc, subdmc; 1682 PetscInt *pids, particles, dim; 1683 1684 PetscFunctionBegin; 1685 /* Configure new swarm */ 1686 PetscCall(DMSetType(cellswarm, DMSWARM)); 1687 PetscCall(DMGetDimension(sw, &dim)); 1688 PetscCall(DMSetDimension(cellswarm, dim)); 1689 PetscCall(DMSwarmSetType(cellswarm, DMSWARM_PIC)); 1690 /* Destroy the unused, unconfigured data bucket to prevent stragglers in memory */ 1691 PetscCall(DMSwarmDataBucketDestroy(&((DM_Swarm *)cellswarm->data)->db)); 1692 PetscCall(DMSwarmSortGetAccess(sw)); 1693 PetscCall(DMSwarmSortGetNumberOfPointsPerCell(sw, cellID, &particles)); 1694 PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids)); 1695 PetscCall(DMSwarmDataBucketCreateFromSubset(original->db, particles, pids, &((DM_Swarm *)cellswarm->data)->db)); 1696 PetscCall(DMSwarmSortRestoreAccess(sw)); 1697 PetscCall(PetscFree(pids)); 1698 PetscCall(DMSwarmGetCellDM(sw, &dmc)); 1699 PetscCall(DMLabelCreate(PetscObjectComm((PetscObject)sw), "singlecell", &label)); 1700 PetscCall(DMAddLabel(dmc, label)); 1701 PetscCall(DMLabelSetValue(label, cellID, 1)); 1702 PetscCall(DMPlexFilter(dmc, label, 1, &subdmc)); 1703 PetscCall(DMSwarmSetCellDM(cellswarm, subdmc)); 1704 PetscCall(DMLabelDestroy(&label)); 1705 PetscFunctionReturn(PETSC_SUCCESS); 1706 } 1707 1708 /*@C 1709 DMSwarmRestoreCellSwarm - Restores a DMSwarm object obtained with DMSwarmGetCellSwarm(). All fields are copied back into the parent swarm. 1710 1711 Noncollective 1712 1713 Input parameters: 1714 + sw - the parent DMSwarm 1715 . cellID - the integer id of the cell to be copied back into the parent swarm 1716 - cellswarm - the cell swarm object 1717 1718 Level: beginner 1719 1720 Note: 1721 This only supports DMSWARM_PIC types of DMSwarms 1722 1723 .seealso: `DMSwarmGetCellSwarm()` 1724 @*/ 1725 PETSC_EXTERN PetscErrorCode DMSwarmRestoreCellSwarm(DM sw, PetscInt cellID, DM cellswarm) 1726 { 1727 DM dmc; 1728 PetscInt *pids, particles, p; 1729 1730 PetscFunctionBegin; 1731 PetscCall(DMSwarmSortGetAccess(sw)); 1732 PetscCall(DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids)); 1733 PetscCall(DMSwarmSortRestoreAccess(sw)); 1734 /* Pointwise copy of each particle based on pid. The parent swarm may not be altered during this process. */ 1735 for (p = 0; p < particles; ++p) PetscCall(DMSwarmDataBucketCopyPoint(((DM_Swarm *)cellswarm->data)->db, pids[p], ((DM_Swarm *)sw->data)->db, pids[p])); 1736 /* Free memory, destroy cell dm */ 1737 PetscCall(DMSwarmGetCellDM(cellswarm, &dmc)); 1738 PetscCall(DMDestroy(&dmc)); 1739 PetscCall(PetscFree(pids)); 1740 PetscFunctionReturn(PETSC_SUCCESS); 1741 } 1742 1743 PETSC_INTERN PetscErrorCode DMClone_Swarm(DM, DM *); 1744 1745 static PetscErrorCode DMInitialize_Swarm(DM sw) 1746 { 1747 PetscFunctionBegin; 1748 sw->dim = 0; 1749 sw->ops->view = DMView_Swarm; 1750 sw->ops->load = NULL; 1751 sw->ops->setfromoptions = NULL; 1752 sw->ops->clone = DMClone_Swarm; 1753 sw->ops->setup = DMSetup_Swarm; 1754 sw->ops->createlocalsection = NULL; 1755 sw->ops->createdefaultconstraints = NULL; 1756 sw->ops->createglobalvector = DMCreateGlobalVector_Swarm; 1757 sw->ops->createlocalvector = DMCreateLocalVector_Swarm; 1758 sw->ops->getlocaltoglobalmapping = NULL; 1759 sw->ops->createfieldis = NULL; 1760 sw->ops->createcoordinatedm = NULL; 1761 sw->ops->getcoloring = NULL; 1762 sw->ops->creatematrix = DMCreateMatrix_Swarm; 1763 sw->ops->createinterpolation = NULL; 1764 sw->ops->createinjection = NULL; 1765 sw->ops->createmassmatrix = DMCreateMassMatrix_Swarm; 1766 sw->ops->refine = NULL; 1767 sw->ops->coarsen = NULL; 1768 sw->ops->refinehierarchy = NULL; 1769 sw->ops->coarsenhierarchy = NULL; 1770 sw->ops->globaltolocalbegin = NULL; 1771 sw->ops->globaltolocalend = NULL; 1772 sw->ops->localtoglobalbegin = NULL; 1773 sw->ops->localtoglobalend = NULL; 1774 sw->ops->destroy = DMDestroy_Swarm; 1775 sw->ops->createsubdm = NULL; 1776 sw->ops->getdimpoints = NULL; 1777 sw->ops->locatepoints = NULL; 1778 PetscFunctionReturn(PETSC_SUCCESS); 1779 } 1780 1781 PETSC_INTERN PetscErrorCode DMClone_Swarm(DM dm, DM *newdm) 1782 { 1783 DM_Swarm *swarm = (DM_Swarm *)dm->data; 1784 1785 PetscFunctionBegin; 1786 swarm->refct++; 1787 (*newdm)->data = swarm; 1788 PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, DMSWARM)); 1789 PetscCall(DMInitialize_Swarm(*newdm)); 1790 (*newdm)->dim = dm->dim; 1791 PetscFunctionReturn(PETSC_SUCCESS); 1792 } 1793 1794 /*MC 1795 1796 DMSWARM = "swarm" - A DM object used to represent arrays of data (fields) of arbitrary data type. 1797 This implementation was designed for particle methods in which the underlying 1798 data required to be represented is both (i) dynamic in length, (ii) and of arbitrary data type. 1799 1800 User data can be represented by DMSwarm through a registering "fields". 1801 To register a field, the user must provide; 1802 (a) a unique name; 1803 (b) the data type (or size in bytes); 1804 (c) the block size of the data. 1805 1806 For example, suppose the application requires a unique id, energy, momentum and density to be stored 1807 on a set of particles. Then the following code could be used 1808 1809 $ DMSwarmInitializeFieldRegister(dm) 1810 $ DMSwarmRegisterPetscDatatypeField(dm,"uid",1,PETSC_LONG); 1811 $ DMSwarmRegisterPetscDatatypeField(dm,"energy",1,PETSC_REAL); 1812 $ DMSwarmRegisterPetscDatatypeField(dm,"momentum",3,PETSC_REAL); 1813 $ DMSwarmRegisterPetscDatatypeField(dm,"density",1,PETSC_FLOAT); 1814 $ DMSwarmFinalizeFieldRegister(dm) 1815 1816 The fields represented by DMSwarm are dynamic and can be re-sized at any time. 1817 The only restriction imposed by DMSwarm is that all fields contain the same number of points. 1818 1819 To support particle methods, "migration" techniques are provided. These methods migrate data 1820 between ranks. 1821 1822 DMSwarm supports the methods DMCreateGlobalVector() and DMCreateLocalVector(). 1823 As a DMSwarm may internally define and store values of different data types, 1824 before calling DMCreateGlobalVector() or DMCreateLocalVector(), the user must inform DMSwarm which 1825 fields should be used to define a Vec object via 1826 DMSwarmVectorDefineField() 1827 The specified field can be changed at any time - thereby permitting vectors 1828 compatible with different fields to be created. 1829 1830 A dual representation of fields in the DMSwarm and a Vec object is permitted via 1831 DMSwarmCreateGlobalVectorFromField() 1832 Here the data defining the field in the DMSwarm is shared with a Vec. 1833 This is inherently unsafe if you alter the size of the field at any time between 1834 calls to DMSwarmCreateGlobalVectorFromField() and DMSwarmDestroyGlobalVectorFromField(). 1835 If the local size of the DMSwarm does not match the local size of the global vector 1836 when DMSwarmDestroyGlobalVectorFromField() is called, an error is thrown. 1837 1838 Additional high-level support is provided for Particle-In-Cell methods. 1839 Please refer to the man page for DMSwarmSetType(). 1840 1841 Level: beginner 1842 1843 .seealso: `DMType`, `DMCreate()`, `DMSetType()` 1844 M*/ 1845 PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm) 1846 { 1847 DM_Swarm *swarm; 1848 1849 PetscFunctionBegin; 1850 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1851 PetscCall(PetscNew(&swarm)); 1852 dm->data = swarm; 1853 PetscCall(DMSwarmDataBucketCreate(&swarm->db)); 1854 PetscCall(DMSwarmInitializeFieldRegister(dm)); 1855 swarm->refct = 1; 1856 swarm->vec_field_set = PETSC_FALSE; 1857 swarm->issetup = PETSC_FALSE; 1858 swarm->swarm_type = DMSWARM_BASIC; 1859 swarm->migrate_type = DMSWARM_MIGRATE_BASIC; 1860 swarm->collect_type = DMSWARM_COLLECT_BASIC; 1861 swarm->migrate_error_on_missing_point = PETSC_FALSE; 1862 swarm->dmcell = NULL; 1863 swarm->collect_view_active = PETSC_FALSE; 1864 swarm->collect_view_reset_nlocal = -1; 1865 PetscCall(DMInitialize_Swarm(dm)); 1866 if (SwarmDataFieldId == -1) PetscCall(PetscObjectComposedDataRegister(&SwarmDataFieldId)); 1867 PetscFunctionReturn(PETSC_SUCCESS); 1868 } 1869