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