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