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 "../src/dm/impls/swarm/data_bucket.h" 9 10 PetscLogEvent DMSWARM_Migrate, DMSWARM_SetSizes, DMSWARM_AddPoints, DMSWARM_RemovePoints, DMSWARM_Sort; 11 PetscLogEvent DMSWARM_DataExchangerTopologySetup, DMSWARM_DataExchangerBegin, DMSWARM_DataExchangerEnd; 12 PetscLogEvent DMSWARM_DataExchangerSendCount, DMSWARM_DataExchangerPack; 13 14 const char* DMSwarmTypeNames[] = { "basic", "pic", 0 }; 15 const char* DMSwarmMigrateTypeNames[] = { "basic", "dmcellnscatter", "dmcellexact", "user", 0 }; 16 const char* DMSwarmCollectTypeNames[] = { "basic", "boundingbox", "general", "user", 0 }; 17 const char* DMSwarmPICLayoutTypeNames[] = { "regular", "gauss", "subdivision", 0 }; 18 19 const char DMSwarmField_pid[] = "DMSwarm_pid"; 20 const char DMSwarmField_rank[] = "DMSwarm_rank"; 21 const char DMSwarmPICField_coor[] = "DMSwarmPIC_coor"; 22 const char DMSwarmPICField_cellid[] = "DMSwarm_cellid"; 23 24 /*@C 25 DMSwarmVectorDefineField - Sets the field from which to define a Vec object 26 when DMCreateLocalVector(), or DMCreateGlobalVector() is called 27 28 Collective on dm 29 30 Input parameters: 31 + dm - a DMSwarm 32 - fieldname - the textual name given to a registered field 33 34 Level: beginner 35 36 Notes: 37 38 The field with name fieldname must be defined as having a data type of PetscScalar. 39 40 This function must be called prior to calling DMCreateLocalVector(), DMCreateGlobalVector(). 41 Mutiple calls to DMSwarmVectorDefineField() are permitted. 42 43 .seealso: DMSwarmRegisterPetscDatatypeField(), DMCreateGlobalVector(), DMCreateLocalVector() 44 @*/ 45 PETSC_EXTERN PetscErrorCode DMSwarmVectorDefineField(DM dm,const char fieldname[]) 46 { 47 DM_Swarm *swarm = (DM_Swarm*)dm->data; 48 PetscErrorCode ierr; 49 PetscInt bs,n; 50 PetscScalar *array; 51 PetscDataType type; 52 53 PetscFunctionBegin; 54 if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); } 55 ierr = DMSwarmDataBucketGetSizes(swarm->db,&n,NULL,NULL);CHKERRQ(ierr); 56 ierr = DMSwarmGetField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr); 57 58 /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */ 59 if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid for PETSC_REAL"); 60 ierr = PetscSNPrintf(swarm->vec_field_name,PETSC_MAX_PATH_LEN-1,"%s",fieldname);CHKERRQ(ierr); 61 swarm->vec_field_set = PETSC_TRUE; 62 swarm->vec_field_bs = bs; 63 swarm->vec_field_nlocal = n; 64 ierr = DMSwarmRestoreField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr); 65 PetscFunctionReturn(0); 66 } 67 68 /* requires DMSwarmDefineFieldVector has been called */ 69 PetscErrorCode DMCreateGlobalVector_Swarm(DM dm,Vec *vec) 70 { 71 DM_Swarm *swarm = (DM_Swarm*)dm->data; 72 PetscErrorCode ierr; 73 Vec x; 74 char name[PETSC_MAX_PATH_LEN]; 75 76 PetscFunctionBegin; 77 if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); } 78 if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first"); 79 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 */ 80 81 ierr = PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name);CHKERRQ(ierr); 82 ierr = VecCreate(PetscObjectComm((PetscObject)dm),&x);CHKERRQ(ierr); 83 ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr); 84 ierr = VecSetSizes(x,swarm->db->L*swarm->vec_field_bs,PETSC_DETERMINE);CHKERRQ(ierr); 85 ierr = VecSetBlockSize(x,swarm->vec_field_bs);CHKERRQ(ierr); 86 ierr = VecSetFromOptions(x);CHKERRQ(ierr); 87 *vec = x; 88 PetscFunctionReturn(0); 89 } 90 91 /* requires DMSwarmDefineFieldVector has been called */ 92 PetscErrorCode DMCreateLocalVector_Swarm(DM dm,Vec *vec) 93 { 94 DM_Swarm *swarm = (DM_Swarm*)dm->data; 95 PetscErrorCode ierr; 96 Vec x; 97 char name[PETSC_MAX_PATH_LEN]; 98 99 PetscFunctionBegin; 100 if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); } 101 if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first"); 102 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 */ 103 104 ierr = PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name);CHKERRQ(ierr); 105 ierr = VecCreate(PETSC_COMM_SELF,&x);CHKERRQ(ierr); 106 ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr); 107 ierr = VecSetSizes(x,swarm->db->L*swarm->vec_field_bs,PETSC_DETERMINE);CHKERRQ(ierr); 108 ierr = VecSetBlockSize(x,swarm->vec_field_bs);CHKERRQ(ierr); 109 ierr = VecSetFromOptions(x);CHKERRQ(ierr); 110 *vec = x; 111 PetscFunctionReturn(0); 112 } 113 114 static PetscErrorCode DMSwarmDestroyVectorFromField_Private(DM dm, const char fieldname[], Vec *vec) 115 { 116 DM_Swarm *swarm = (DM_Swarm *) dm->data; 117 DMSwarmDataField gfield; 118 void (*fptr)(void); 119 PetscInt bs, nlocal; 120 char name[PETSC_MAX_PATH_LEN]; 121 PetscErrorCode ierr; 122 123 PetscFunctionBegin; 124 ierr = VecGetLocalSize(*vec, &nlocal);CHKERRQ(ierr); 125 ierr = VecGetBlockSize(*vec, &bs);CHKERRQ(ierr); 126 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 */ 127 ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield);CHKERRQ(ierr); 128 /* check vector is an inplace array */ 129 ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarm_VecFieldInPlace_%s", fieldname);CHKERRQ(ierr); 130 ierr = PetscObjectQueryFunction((PetscObject) *vec, name, &fptr);CHKERRQ(ierr); 131 if (!fptr) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_USER, "Vector being destroyed was not created from DMSwarm field(%s)", fieldname); 132 ierr = DMSwarmDataFieldRestoreAccess(gfield);CHKERRQ(ierr); 133 ierr = VecDestroy(vec);CHKERRQ(ierr); 134 PetscFunctionReturn(0); 135 } 136 137 static PetscErrorCode DMSwarmCreateVectorFromField_Private(DM dm, const char fieldname[], MPI_Comm comm, Vec *vec) 138 { 139 DM_Swarm *swarm = (DM_Swarm *) dm->data; 140 PetscDataType type; 141 PetscScalar *array; 142 PetscInt bs, n; 143 char name[PETSC_MAX_PATH_LEN]; 144 PetscMPIInt size; 145 PetscErrorCode ierr; 146 147 PetscFunctionBegin; 148 if (!swarm->issetup) {ierr = DMSetUp(dm);CHKERRQ(ierr);} 149 ierr = DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL);CHKERRQ(ierr); 150 ierr = DMSwarmGetField(dm, fieldname, &bs, &type, (void **) &array);CHKERRQ(ierr); 151 if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL"); 152 153 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 154 if (size == 1) { 155 ierr = VecCreateSeqWithArray(comm, bs, n*bs, array, vec);CHKERRQ(ierr); 156 } else { 157 ierr = VecCreateMPIWithArray(comm, bs, n*bs, PETSC_DETERMINE, array, vec);CHKERRQ(ierr); 158 } 159 ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarmSharedField_%s", fieldname);CHKERRQ(ierr); 160 ierr = PetscObjectSetName((PetscObject) *vec, name);CHKERRQ(ierr); 161 162 /* Set guard */ 163 ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarm_VecFieldInPlace_%s", fieldname);CHKERRQ(ierr); 164 ierr = PetscObjectComposeFunction((PetscObject) *vec, name, DMSwarmDestroyVectorFromField_Private);CHKERRQ(ierr); 165 PetscFunctionReturn(0); 166 } 167 168 /* This creates a "mass matrix" between a finite element and particle space. If a finite element interpolant is given by 169 170 \hat f = \sum_i f_i \phi_i 171 172 and a particle function is given by 173 174 f = \sum_i w_i \delta(x - x_i) 175 176 then we want to require that 177 178 M \hat f = M_p f 179 180 where the particle mass matrix is given by 181 182 (M_p)_{ij} = \int \phi_i \delta(x - x_j) 183 184 The way Dave May does particles, they amount to quadratue weights rather than delta functions, so he has |J| is in 185 his integral. We allow this with the boolean flag. 186 */ 187 static PetscErrorCode DMSwarmComputeMassMatrix_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx) 188 { 189 const char *name = "Mass Matrix"; 190 MPI_Comm comm; 191 PetscDS prob; 192 PetscSection fsection, globalFSection; 193 PetscHSetIJ ht; 194 PetscLayout rLayout, colLayout; 195 PetscInt *dnz, *onz; 196 PetscInt locRows, locCols, rStart, colStart, colEnd, *rowIDXs; 197 PetscReal *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0}; 198 PetscScalar *elemMat; 199 PetscInt dim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0; 200 PetscErrorCode ierr; 201 202 PetscFunctionBegin; 203 ierr = PetscObjectGetComm((PetscObject) mass, &comm);CHKERRQ(ierr); 204 ierr = DMGetCoordinateDim(dmf, &dim);CHKERRQ(ierr); 205 ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr); 206 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 207 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 208 ierr = PetscMalloc3(dim, &v0, dim*dim, &J, dim*dim,&invJ);CHKERRQ(ierr); 209 ierr = DMGetSection(dmf, &fsection);CHKERRQ(ierr); 210 ierr = DMGetGlobalSection(dmf, &globalFSection);CHKERRQ(ierr); 211 ierr = DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);CHKERRQ(ierr); 212 ierr = MatGetLocalSize(mass, &locRows, &locCols);CHKERRQ(ierr); 213 214 ierr = PetscLayoutCreate(comm, &colLayout);CHKERRQ(ierr); 215 ierr = PetscLayoutSetLocalSize(colLayout, locCols);CHKERRQ(ierr); 216 ierr = PetscLayoutSetBlockSize(colLayout, 1);CHKERRQ(ierr); 217 ierr = PetscLayoutSetUp(colLayout);CHKERRQ(ierr); 218 ierr = PetscLayoutGetRange(colLayout, &colStart, &colEnd);CHKERRQ(ierr); 219 ierr = PetscLayoutDestroy(&colLayout);CHKERRQ(ierr); 220 221 ierr = PetscLayoutCreate(comm, &rLayout);CHKERRQ(ierr); 222 ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr); 223 ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr); 224 ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr); 225 ierr = PetscLayoutGetRange(rLayout, &rStart, NULL);CHKERRQ(ierr); 226 ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr); 227 228 ierr = PetscCalloc2(locRows, &dnz, locRows, &onz);CHKERRQ(ierr); 229 ierr = PetscHSetIJCreate(&ht);CHKERRQ(ierr); 230 231 ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr); 232 /* count non-zeros */ 233 ierr = DMSwarmSortGetAccess(dmc);CHKERRQ(ierr); 234 for (field = 0; field < Nf; ++field) { 235 for (cell = cStart; cell < cEnd; ++cell) { 236 PetscInt c, i; 237 PetscInt *findices, *cindices; /* fine is vertices, coarse is particles */ 238 PetscInt numFIndices, numCIndices; 239 240 ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr); 241 ierr = DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices);CHKERRQ(ierr); 242 maxC = PetscMax(maxC, numCIndices); 243 { 244 PetscHashIJKey key; 245 PetscBool missing; 246 for (i = 0; i < numFIndices; ++i) { 247 key.j = findices[i]; /* global column (from Plex) */ 248 if (key.j >= 0) { 249 /* Get indices for coarse elements */ 250 for (c = 0; c < numCIndices; ++c) { 251 key.i = cindices[c] + rStart; /* global cols (from Swarm) */ 252 if (key.i < 0) continue; 253 ierr = PetscHSetIJQueryAdd(ht, key, &missing);CHKERRQ(ierr); 254 if (missing) { 255 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart]; 256 else ++onz[key.i - rStart]; 257 } else SETERRQ2(PetscObjectComm((PetscObject) dmf), PETSC_ERR_SUP, "Set new value at %D,%D", key.i, key.j); 258 } 259 } 260 } 261 ierr = PetscFree(cindices);CHKERRQ(ierr); 262 } 263 ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr); 264 } 265 } 266 ierr = PetscHSetIJDestroy(&ht);CHKERRQ(ierr); 267 ierr = MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr); 268 ierr = MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE);CHKERRQ(ierr); 269 ierr = PetscFree2(dnz, onz);CHKERRQ(ierr); 270 ierr = PetscMalloc3(maxC*totDim, &elemMat, maxC, &rowIDXs, maxC*dim, &xi);CHKERRQ(ierr); 271 for (field = 0; field < Nf; ++field) { 272 PetscObject obj; 273 PetscReal *Bcoarse, *coords; 274 PetscInt Nc, i; 275 276 ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr); 277 ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr); 278 if (Nc != 1) SETERRQ1(PetscObjectComm((PetscObject) dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %D", Nc); 279 ierr = DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 280 for (cell = cStart; cell < cEnd; ++cell) { 281 PetscInt *findices , *cindices; 282 PetscInt numFIndices, numCIndices; 283 PetscInt p, c; 284 285 /* TODO: Use DMField instead of assuming affine */ 286 ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 287 ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr); 288 ierr = DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices);CHKERRQ(ierr); 289 for (p = 0; p < numCIndices; ++p) { 290 CoordinatesRealToRef(dim, dim, v0ref, v0, invJ, &coords[cindices[p]*dim], &xi[p*dim]); 291 } 292 ierr = PetscFEGetTabulation((PetscFE) obj, numCIndices, xi, &Bcoarse, NULL, NULL);CHKERRQ(ierr); 293 /* Get elemMat entries by multiplying by weight */ 294 ierr = PetscArrayzero(elemMat, numCIndices*totDim);CHKERRQ(ierr); 295 for (i = 0; i < numFIndices; ++i) { 296 for (p = 0; p < numCIndices; ++p) { 297 for (c = 0; c < Nc; ++c) { 298 /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */ 299 elemMat[p*numFIndices+i] += Bcoarse[(p*numFIndices + i)*Nc + c]*(useDeltaFunction ? 1.0 : detJ); 300 } 301 } 302 } 303 for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart; 304 if (0) {ierr = DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);CHKERRQ(ierr);} 305 ierr = MatSetValues(mass, numCIndices, rowIDXs, numFIndices, findices, elemMat, ADD_VALUES);CHKERRQ(ierr); 306 ierr = PetscFree(cindices);CHKERRQ(ierr); 307 ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr); 308 ierr = PetscFERestoreTabulation((PetscFE) obj, numCIndices, xi, &Bcoarse, NULL, NULL);CHKERRQ(ierr); 309 } 310 ierr = DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 311 } 312 ierr = PetscFree3(elemMat, rowIDXs, xi);CHKERRQ(ierr); 313 ierr = DMSwarmSortRestoreAccess(dmc);CHKERRQ(ierr); 314 ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr); 315 ierr = MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 316 ierr = MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 317 PetscFunctionReturn(0); 318 } 319 320 /* FEM cols, Particle rows */ 321 static PetscErrorCode DMCreateMassMatrix_Swarm(DM dmCoarse, DM dmFine, Mat *mass) 322 { 323 PetscSection gsf; 324 PetscInt m, n; 325 void *ctx; 326 PetscErrorCode ierr; 327 328 PetscFunctionBegin; 329 ierr = DMGetGlobalSection(dmFine, &gsf);CHKERRQ(ierr); 330 ierr = PetscSectionGetConstrainedStorageSize(gsf, &m);CHKERRQ(ierr); 331 ierr = DMSwarmGetLocalSize(dmCoarse, &n);CHKERRQ(ierr); 332 ierr = MatCreate(PetscObjectComm((PetscObject) dmCoarse), mass);CHKERRQ(ierr); 333 ierr = MatSetSizes(*mass, n, m, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 334 ierr = MatSetType(*mass, dmCoarse->mattype);CHKERRQ(ierr); 335 ierr = DMGetApplicationContext(dmFine, &ctx);CHKERRQ(ierr); 336 337 ierr = DMSwarmComputeMassMatrix_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx);CHKERRQ(ierr); 338 ierr = MatViewFromOptions(*mass, NULL, "-mass_mat_view");CHKERRQ(ierr); 339 PetscFunctionReturn(0); 340 } 341 342 /*@C 343 DMSwarmCreateGlobalVectorFromField - Creates a Vec object sharing the array associated with a given field 344 345 Collective on dm 346 347 Input parameters: 348 + dm - a DMSwarm 349 - fieldname - the textual name given to a registered field 350 351 Output parameter: 352 . vec - the vector 353 354 Level: beginner 355 356 Notes: 357 The vector must be returned using a matching call to DMSwarmDestroyGlobalVectorFromField(). 358 359 .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmDestroyGlobalVectorFromField() 360 @*/ 361 PETSC_EXTERN PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec) 362 { 363 MPI_Comm comm = PetscObjectComm((PetscObject) dm); 364 PetscErrorCode ierr; 365 366 PetscFunctionBegin; 367 ierr = DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec);CHKERRQ(ierr); 368 PetscFunctionReturn(0); 369 } 370 371 /*@C 372 DMSwarmDestroyGlobalVectorFromField - Destroys the Vec object which share the array associated with a given field 373 374 Collective on dm 375 376 Input parameters: 377 + dm - a DMSwarm 378 - fieldname - the textual name given to a registered field 379 380 Output parameter: 381 . vec - the vector 382 383 Level: beginner 384 385 .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmCreateGlobalVectorFromField() 386 @*/ 387 PETSC_EXTERN PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec) 388 { 389 PetscErrorCode ierr; 390 391 PetscFunctionBegin; 392 ierr = DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec);CHKERRQ(ierr); 393 PetscFunctionReturn(0); 394 } 395 396 /*@C 397 DMSwarmCreateLocalVectorFromField - Creates a Vec object sharing the array associated with a given field 398 399 Collective on dm 400 401 Input parameters: 402 + dm - a DMSwarm 403 - fieldname - the textual name given to a registered field 404 405 Output parameter: 406 . vec - the vector 407 408 Level: beginner 409 410 Notes: 411 The vector must be returned using a matching call to DMSwarmDestroyLocalVectorFromField(). 412 413 .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmDestroyLocalVectorFromField() 414 @*/ 415 PETSC_EXTERN PetscErrorCode DMSwarmCreateLocalVectorFromField(DM dm,const char fieldname[],Vec *vec) 416 { 417 MPI_Comm comm = PETSC_COMM_SELF; 418 PetscErrorCode ierr; 419 420 PetscFunctionBegin; 421 ierr = DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec);CHKERRQ(ierr); 422 PetscFunctionReturn(0); 423 } 424 425 /*@C 426 DMSwarmDestroyLocalVectorFromField - Destroys the Vec object which share the array associated with a given field 427 428 Collective on dm 429 430 Input parameters: 431 + dm - a DMSwarm 432 - fieldname - the textual name given to a registered field 433 434 Output parameter: 435 . vec - the vector 436 437 Level: beginner 438 439 .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmCreateLocalVectorFromField() 440 @*/ 441 PETSC_EXTERN PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM dm,const char fieldname[],Vec *vec) 442 { 443 PetscErrorCode ierr; 444 445 PetscFunctionBegin; 446 ierr = DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec);CHKERRQ(ierr); 447 PetscFunctionReturn(0); 448 } 449 450 /* 451 PETSC_EXTERN PetscErrorCode DMSwarmCreateGlobalVectorFromFields(DM dm,const PetscInt nf,const char *fieldnames[],Vec *vec) 452 { 453 PetscFunctionReturn(0); 454 } 455 456 PETSC_EXTERN PetscErrorCode DMSwarmRestoreGlobalVectorFromFields(DM dm,Vec *vec) 457 { 458 PetscFunctionReturn(0); 459 } 460 */ 461 462 /*@C 463 DMSwarmInitializeFieldRegister - Initiates the registration of fields to a DMSwarm 464 465 Collective on dm 466 467 Input parameter: 468 . dm - a DMSwarm 469 470 Level: beginner 471 472 Notes: 473 After all fields have been registered, you must call DMSwarmFinalizeFieldRegister(). 474 475 .seealso: DMSwarmFinalizeFieldRegister(), DMSwarmRegisterPetscDatatypeField(), 476 DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField() 477 @*/ 478 PETSC_EXTERN PetscErrorCode DMSwarmInitializeFieldRegister(DM dm) 479 { 480 DM_Swarm *swarm = (DM_Swarm *) dm->data; 481 PetscErrorCode ierr; 482 483 PetscFunctionBegin; 484 if (!swarm->field_registration_initialized) { 485 swarm->field_registration_initialized = PETSC_TRUE; 486 ierr = DMSwarmRegisterPetscDatatypeField(dm,DMSwarmField_pid,1,PETSC_LONG);CHKERRQ(ierr); /* unique identifer */ 487 ierr = DMSwarmRegisterPetscDatatypeField(dm,DMSwarmField_rank,1,PETSC_INT);CHKERRQ(ierr); /* used for communication */ 488 } 489 PetscFunctionReturn(0); 490 } 491 492 /*@C 493 DMSwarmFinalizeFieldRegister - Finalizes the registration of fields to a DMSwarm 494 495 Collective on dm 496 497 Input parameter: 498 . dm - a DMSwarm 499 500 Level: beginner 501 502 Notes: 503 After DMSwarmFinalizeFieldRegister() has been called, no new fields can be defined on the DMSwarm. 504 505 .seealso: DMSwarmInitializeFieldRegister(), DMSwarmRegisterPetscDatatypeField(), 506 DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField() 507 @*/ 508 PETSC_EXTERN PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm) 509 { 510 DM_Swarm *swarm = (DM_Swarm*)dm->data; 511 PetscErrorCode ierr; 512 513 PetscFunctionBegin; 514 if (!swarm->field_registration_finalized) { 515 ierr = DMSwarmDataBucketFinalize(swarm->db);CHKERRQ(ierr); 516 } 517 swarm->field_registration_finalized = PETSC_TRUE; 518 PetscFunctionReturn(0); 519 } 520 521 /*@C 522 DMSwarmSetLocalSizes - Sets the length of all registered fields on the DMSwarm 523 524 Not collective 525 526 Input parameters: 527 + dm - a DMSwarm 528 . nlocal - the length of each registered field 529 - buffer - the length of the buffer used to efficient dynamic re-sizing 530 531 Level: beginner 532 533 .seealso: DMSwarmGetLocalSize() 534 @*/ 535 PETSC_EXTERN PetscErrorCode DMSwarmSetLocalSizes(DM dm,PetscInt nlocal,PetscInt buffer) 536 { 537 DM_Swarm *swarm = (DM_Swarm*)dm->data; 538 PetscErrorCode ierr; 539 540 PetscFunctionBegin; 541 ierr = PetscLogEventBegin(DMSWARM_SetSizes,0,0,0,0);CHKERRQ(ierr); 542 ierr = DMSwarmDataBucketSetSizes(swarm->db,nlocal,buffer);CHKERRQ(ierr); 543 ierr = PetscLogEventEnd(DMSWARM_SetSizes,0,0,0,0);CHKERRQ(ierr); 544 PetscFunctionReturn(0); 545 } 546 547 /*@C 548 DMSwarmSetCellDM - Attachs a DM to a DMSwarm 549 550 Collective on dm 551 552 Input parameters: 553 + dm - a DMSwarm 554 - dmcell - the DM to attach to the DMSwarm 555 556 Level: beginner 557 558 Notes: 559 The attached DM (dmcell) will be queried for point location and 560 neighbor MPI-rank information if DMSwarmMigrate() is called. 561 562 .seealso: DMSwarmSetType(), DMSwarmGetCellDM(), DMSwarmMigrate() 563 @*/ 564 PETSC_EXTERN PetscErrorCode DMSwarmSetCellDM(DM dm,DM dmcell) 565 { 566 DM_Swarm *swarm = (DM_Swarm*)dm->data; 567 568 PetscFunctionBegin; 569 swarm->dmcell = dmcell; 570 PetscFunctionReturn(0); 571 } 572 573 /*@C 574 DMSwarmGetCellDM - Fetches the attached cell DM 575 576 Collective on dm 577 578 Input parameter: 579 . dm - a DMSwarm 580 581 Output parameter: 582 . dmcell - the DM which was attached to the DMSwarm 583 584 Level: beginner 585 586 .seealso: DMSwarmSetCellDM() 587 @*/ 588 PETSC_EXTERN PetscErrorCode DMSwarmGetCellDM(DM dm,DM *dmcell) 589 { 590 DM_Swarm *swarm = (DM_Swarm*)dm->data; 591 592 PetscFunctionBegin; 593 *dmcell = swarm->dmcell; 594 PetscFunctionReturn(0); 595 } 596 597 /*@C 598 DMSwarmGetLocalSize - Retrives the local length of fields registered 599 600 Not collective 601 602 Input parameter: 603 . dm - a DMSwarm 604 605 Output parameter: 606 . nlocal - the length of each registered field 607 608 Level: beginner 609 610 .seealso: DMSwarmGetSize(), DMSwarmSetLocalSizes() 611 @*/ 612 PETSC_EXTERN PetscErrorCode DMSwarmGetLocalSize(DM dm,PetscInt *nlocal) 613 { 614 DM_Swarm *swarm = (DM_Swarm*)dm->data; 615 PetscErrorCode ierr; 616 617 PetscFunctionBegin; 618 if (nlocal) {ierr = DMSwarmDataBucketGetSizes(swarm->db,nlocal,NULL,NULL);CHKERRQ(ierr);} 619 PetscFunctionReturn(0); 620 } 621 622 /*@C 623 DMSwarmGetSize - Retrives the total length of fields registered 624 625 Collective on dm 626 627 Input parameter: 628 . dm - a DMSwarm 629 630 Output parameter: 631 . n - the total length of each registered field 632 633 Level: beginner 634 635 Note: 636 This calls MPI_Allreduce upon each call (inefficient but safe) 637 638 .seealso: DMSwarmGetLocalSize(), DMSwarmSetLocalSizes() 639 @*/ 640 PETSC_EXTERN PetscErrorCode DMSwarmGetSize(DM dm,PetscInt *n) 641 { 642 DM_Swarm *swarm = (DM_Swarm*)dm->data; 643 PetscErrorCode ierr; 644 PetscInt nlocal,ng; 645 646 PetscFunctionBegin; 647 ierr = DMSwarmDataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);CHKERRQ(ierr); 648 ierr = MPI_Allreduce(&nlocal,&ng,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 649 if (n) { *n = ng; } 650 PetscFunctionReturn(0); 651 } 652 653 /*@C 654 DMSwarmRegisterPetscDatatypeField - Register a field to a DMSwarm with a native PETSc data type 655 656 Collective on dm 657 658 Input parameters: 659 + dm - a DMSwarm 660 . fieldname - the textual name to identify this field 661 . blocksize - the number of each data type 662 - type - a valid PETSc data type (PETSC_CHAR, PETSC_SHORT, PETSC_INT, PETSC_FLOAT, PETSC_REAL, PETSC_LONG) 663 664 Level: beginner 665 666 Notes: 667 The textual name for each registered field must be unique. 668 669 .seealso: DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField() 670 @*/ 671 PETSC_EXTERN PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm,const char fieldname[],PetscInt blocksize,PetscDataType type) 672 { 673 PetscErrorCode ierr; 674 DM_Swarm *swarm = (DM_Swarm*)dm->data; 675 size_t size; 676 677 PetscFunctionBegin; 678 if (!swarm->field_registration_initialized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmInitializeFieldRegister() first"); 679 if (swarm->field_registration_finalized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first"); 680 681 if (type == PETSC_OBJECT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 682 if (type == PETSC_FUNCTION) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 683 if (type == PETSC_STRING) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 684 if (type == PETSC_STRUCT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 685 if (type == PETSC_DATATYPE_UNKNOWN) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 686 687 ierr = PetscDataTypeGetSize(type, &size);CHKERRQ(ierr); 688 /* Load a specific data type into data bucket, specifying textual name and its size in bytes */ 689 ierr = DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterPetscDatatypeField",fieldname,blocksize*size,NULL);CHKERRQ(ierr); 690 { 691 DMSwarmDataField gfield; 692 693 ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr); 694 ierr = DMSwarmDataFieldSetBlockSize(gfield,blocksize);CHKERRQ(ierr); 695 } 696 swarm->db->field[swarm->db->nfields-1]->petsc_type = type; 697 PetscFunctionReturn(0); 698 } 699 700 /*@C 701 DMSwarmRegisterUserStructField - Register a user defined struct to a DMSwarm 702 703 Collective on dm 704 705 Input parameters: 706 + dm - a DMSwarm 707 . fieldname - the textual name to identify this field 708 - size - the size in bytes of the user struct of each data type 709 710 Level: beginner 711 712 Notes: 713 The textual name for each registered field must be unique. 714 715 .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmRegisterUserDatatypeField() 716 @*/ 717 PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserStructField(DM dm,const char fieldname[],size_t size) 718 { 719 PetscErrorCode ierr; 720 DM_Swarm *swarm = (DM_Swarm*)dm->data; 721 722 PetscFunctionBegin; 723 ierr = DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterUserStructField",fieldname,size,NULL);CHKERRQ(ierr); 724 swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_STRUCT ; 725 PetscFunctionReturn(0); 726 } 727 728 /*@C 729 DMSwarmRegisterUserDatatypeField - Register a user defined data type to a DMSwarm 730 731 Collective on dm 732 733 Input parameters: 734 + dm - a DMSwarm 735 . fieldname - the textual name to identify this field 736 . size - the size in bytes of the user data type 737 - blocksize - the number of each data type 738 739 Level: beginner 740 741 Notes: 742 The textual name for each registered field must be unique. 743 744 .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField() 745 @*/ 746 PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm,const char fieldname[],size_t size,PetscInt blocksize) 747 { 748 DM_Swarm *swarm = (DM_Swarm*)dm->data; 749 PetscErrorCode ierr; 750 751 PetscFunctionBegin; 752 ierr = DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterUserDatatypeField",fieldname,blocksize*size,NULL);CHKERRQ(ierr); 753 { 754 DMSwarmDataField gfield; 755 756 ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr); 757 ierr = DMSwarmDataFieldSetBlockSize(gfield,blocksize);CHKERRQ(ierr); 758 } 759 swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_DATATYPE_UNKNOWN; 760 PetscFunctionReturn(0); 761 } 762 763 /*@C 764 DMSwarmGetField - Get access to the underlying array storing all entries associated with a registered field 765 766 Not collective 767 768 Input parameters: 769 + dm - a DMSwarm 770 - fieldname - the textual name to identify this field 771 772 Output parameters: 773 + blocksize - the number of each data type 774 . type - the data type 775 - data - pointer to raw array 776 777 Level: beginner 778 779 Notes: 780 The array must be returned using a matching call to DMSwarmRestoreField(). 781 782 .seealso: DMSwarmRestoreField() 783 @*/ 784 PETSC_EXTERN PetscErrorCode DMSwarmGetField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data) 785 { 786 DM_Swarm *swarm = (DM_Swarm*)dm->data; 787 DMSwarmDataField gfield; 788 PetscErrorCode ierr; 789 790 PetscFunctionBegin; 791 if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); } 792 ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr); 793 ierr = DMSwarmDataFieldGetAccess(gfield);CHKERRQ(ierr); 794 ierr = DMSwarmDataFieldGetEntries(gfield,data);CHKERRQ(ierr); 795 if (blocksize) {*blocksize = gfield->bs; } 796 if (type) { *type = gfield->petsc_type; } 797 PetscFunctionReturn(0); 798 } 799 800 /*@C 801 DMSwarmRestoreField - Restore access to the underlying array storing all entries associated with a registered field 802 803 Not collective 804 805 Input parameters: 806 + dm - a DMSwarm 807 - fieldname - the textual name to identify this field 808 809 Output parameters: 810 + blocksize - the number of each data type 811 . type - the data type 812 - data - pointer to raw array 813 814 Level: beginner 815 816 Notes: 817 The user must call DMSwarmGetField() prior to calling DMSwarmRestoreField(). 818 819 .seealso: DMSwarmGetField() 820 @*/ 821 PETSC_EXTERN PetscErrorCode DMSwarmRestoreField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data) 822 { 823 DM_Swarm *swarm = (DM_Swarm*)dm->data; 824 DMSwarmDataField gfield; 825 PetscErrorCode ierr; 826 827 PetscFunctionBegin; 828 ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr); 829 ierr = DMSwarmDataFieldRestoreAccess(gfield);CHKERRQ(ierr); 830 if (data) *data = NULL; 831 PetscFunctionReturn(0); 832 } 833 834 /*@C 835 DMSwarmAddPoint - Add space for one new point in the DMSwarm 836 837 Not collective 838 839 Input parameter: 840 . dm - a DMSwarm 841 842 Level: beginner 843 844 Notes: 845 The new point will have all fields initialized to zero. 846 847 .seealso: DMSwarmAddNPoints() 848 @*/ 849 PETSC_EXTERN PetscErrorCode DMSwarmAddPoint(DM dm) 850 { 851 DM_Swarm *swarm = (DM_Swarm*)dm->data; 852 PetscErrorCode ierr; 853 854 PetscFunctionBegin; 855 if (!swarm->issetup) {ierr = DMSetUp(dm);CHKERRQ(ierr);} 856 ierr = PetscLogEventBegin(DMSWARM_AddPoints,0,0,0,0);CHKERRQ(ierr); 857 ierr = DMSwarmDataBucketAddPoint(swarm->db);CHKERRQ(ierr); 858 ierr = PetscLogEventEnd(DMSWARM_AddPoints,0,0,0,0);CHKERRQ(ierr); 859 PetscFunctionReturn(0); 860 } 861 862 /*@C 863 DMSwarmAddNPoints - Add space for a number of new points in the DMSwarm 864 865 Not collective 866 867 Input parameters: 868 + dm - a DMSwarm 869 - npoints - the number of new points to add 870 871 Level: beginner 872 873 Notes: 874 The new point will have all fields initialized to zero. 875 876 .seealso: DMSwarmAddPoint() 877 @*/ 878 PETSC_EXTERN PetscErrorCode DMSwarmAddNPoints(DM dm,PetscInt npoints) 879 { 880 DM_Swarm *swarm = (DM_Swarm*)dm->data; 881 PetscErrorCode ierr; 882 PetscInt nlocal; 883 884 PetscFunctionBegin; 885 ierr = PetscLogEventBegin(DMSWARM_AddPoints,0,0,0,0);CHKERRQ(ierr); 886 ierr = DMSwarmDataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);CHKERRQ(ierr); 887 nlocal = nlocal + npoints; 888 ierr = DMSwarmDataBucketSetSizes(swarm->db,nlocal,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr); 889 ierr = PetscLogEventEnd(DMSWARM_AddPoints,0,0,0,0);CHKERRQ(ierr); 890 PetscFunctionReturn(0); 891 } 892 893 /*@C 894 DMSwarmRemovePoint - Remove the last point from the DMSwarm 895 896 Not collective 897 898 Input parameter: 899 . dm - a DMSwarm 900 901 Level: beginner 902 903 .seealso: DMSwarmRemovePointAtIndex() 904 @*/ 905 PETSC_EXTERN PetscErrorCode DMSwarmRemovePoint(DM dm) 906 { 907 DM_Swarm *swarm = (DM_Swarm*)dm->data; 908 PetscErrorCode ierr; 909 910 PetscFunctionBegin; 911 ierr = PetscLogEventBegin(DMSWARM_RemovePoints,0,0,0,0);CHKERRQ(ierr); 912 ierr = DMSwarmDataBucketRemovePoint(swarm->db);CHKERRQ(ierr); 913 ierr = PetscLogEventEnd(DMSWARM_RemovePoints,0,0,0,0);CHKERRQ(ierr); 914 PetscFunctionReturn(0); 915 } 916 917 /*@C 918 DMSwarmRemovePointAtIndex - Removes a specific point from the DMSwarm 919 920 Not collective 921 922 Input parameters: 923 + dm - a DMSwarm 924 - idx - index of point to remove 925 926 Level: beginner 927 928 .seealso: DMSwarmRemovePoint() 929 @*/ 930 PETSC_EXTERN PetscErrorCode DMSwarmRemovePointAtIndex(DM dm,PetscInt idx) 931 { 932 DM_Swarm *swarm = (DM_Swarm*)dm->data; 933 PetscErrorCode ierr; 934 935 PetscFunctionBegin; 936 ierr = PetscLogEventBegin(DMSWARM_RemovePoints,0,0,0,0);CHKERRQ(ierr); 937 ierr = DMSwarmDataBucketRemovePointAtIndex(swarm->db,idx);CHKERRQ(ierr); 938 ierr = PetscLogEventEnd(DMSWARM_RemovePoints,0,0,0,0);CHKERRQ(ierr); 939 PetscFunctionReturn(0); 940 } 941 942 /*@C 943 DMSwarmCopyPoint - Copy point pj to point pi in the DMSwarm 944 945 Not collective 946 947 Input parameters: 948 + dm - a DMSwarm 949 . pi - the index of the point to copy 950 - pj - the point index where the copy should be located 951 952 Level: beginner 953 954 .seealso: DMSwarmRemovePoint() 955 @*/ 956 PETSC_EXTERN PetscErrorCode DMSwarmCopyPoint(DM dm,PetscInt pi,PetscInt pj) 957 { 958 DM_Swarm *swarm = (DM_Swarm*)dm->data; 959 PetscErrorCode ierr; 960 961 PetscFunctionBegin; 962 if (!swarm->issetup) {ierr = DMSetUp(dm);CHKERRQ(ierr);} 963 ierr = DMSwarmDataBucketCopyPoint(swarm->db,pi,swarm->db,pj);CHKERRQ(ierr); 964 PetscFunctionReturn(0); 965 } 966 967 PetscErrorCode DMSwarmMigrate_Basic(DM dm,PetscBool remove_sent_points) 968 { 969 PetscErrorCode ierr; 970 971 PetscFunctionBegin; 972 ierr = DMSwarmMigrate_Push_Basic(dm,remove_sent_points);CHKERRQ(ierr); 973 PetscFunctionReturn(0); 974 } 975 976 /*@C 977 DMSwarmMigrate - Relocates points defined in the DMSwarm to other MPI-ranks 978 979 Collective on dm 980 981 Input parameters: 982 + dm - the DMSwarm 983 - remove_sent_points - flag indicating if sent points should be removed from the current MPI-rank 984 985 Notes: 986 The DM will be modified to accomodate received points. 987 If remove_sent_points = PETSC_TRUE, any points that were sent will be removed from the DM. 988 Different styles of migration are supported. See DMSwarmSetMigrateType(). 989 990 Level: advanced 991 992 .seealso: DMSwarmSetMigrateType() 993 @*/ 994 PETSC_EXTERN PetscErrorCode DMSwarmMigrate(DM dm,PetscBool remove_sent_points) 995 { 996 DM_Swarm *swarm = (DM_Swarm*)dm->data; 997 PetscErrorCode ierr; 998 999 PetscFunctionBegin; 1000 ierr = PetscLogEventBegin(DMSWARM_Migrate,0,0,0,0);CHKERRQ(ierr); 1001 switch (swarm->migrate_type) { 1002 case DMSWARM_MIGRATE_BASIC: 1003 ierr = DMSwarmMigrate_Basic(dm,remove_sent_points);CHKERRQ(ierr); 1004 break; 1005 case DMSWARM_MIGRATE_DMCELLNSCATTER: 1006 ierr = DMSwarmMigrate_CellDMScatter(dm,remove_sent_points);CHKERRQ(ierr); 1007 break; 1008 case DMSWARM_MIGRATE_DMCELLEXACT: 1009 SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE_DMCELLEXACT not implemented"); 1010 break; 1011 case DMSWARM_MIGRATE_USER: 1012 SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE_USER not implemented"); 1013 break; 1014 default: 1015 SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE type unknown"); 1016 break; 1017 } 1018 ierr = PetscLogEventEnd(DMSWARM_Migrate,0,0,0,0);CHKERRQ(ierr); 1019 PetscFunctionReturn(0); 1020 } 1021 1022 PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize); 1023 1024 /* 1025 DMSwarmCollectViewCreate 1026 1027 * Applies a collection method and gathers point neighbour points into dm 1028 1029 Notes: 1030 Users should call DMSwarmCollectViewDestroy() after 1031 they have finished computations associated with the collected points 1032 */ 1033 1034 /*@C 1035 DMSwarmCollectViewCreate - Applies a collection method and gathers points 1036 in neighbour MPI-ranks into the DMSwarm 1037 1038 Collective on dm 1039 1040 Input parameter: 1041 . dm - the DMSwarm 1042 1043 Notes: 1044 Users should call DMSwarmCollectViewDestroy() after 1045 they have finished computations associated with the collected points 1046 Different collect methods are supported. See DMSwarmSetCollectType(). 1047 1048 Level: advanced 1049 1050 .seealso: DMSwarmCollectViewDestroy(), DMSwarmSetCollectType() 1051 @*/ 1052 PETSC_EXTERN PetscErrorCode DMSwarmCollectViewCreate(DM dm) 1053 { 1054 PetscErrorCode ierr; 1055 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1056 PetscInt ng; 1057 1058 PetscFunctionBegin; 1059 if (swarm->collect_view_active) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView currently active"); 1060 ierr = DMSwarmGetLocalSize(dm,&ng);CHKERRQ(ierr); 1061 switch (swarm->collect_type) { 1062 1063 case DMSWARM_COLLECT_BASIC: 1064 ierr = DMSwarmMigrate_GlobalToLocal_Basic(dm,&ng);CHKERRQ(ierr); 1065 break; 1066 case DMSWARM_COLLECT_DMDABOUNDINGBOX: 1067 SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT_DMDABOUNDINGBOX not implemented"); 1068 break; 1069 case DMSWARM_COLLECT_GENERAL: 1070 SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT_GENERAL not implemented"); 1071 break; 1072 default: 1073 SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT type unknown"); 1074 break; 1075 } 1076 swarm->collect_view_active = PETSC_TRUE; 1077 swarm->collect_view_reset_nlocal = ng; 1078 PetscFunctionReturn(0); 1079 } 1080 1081 /*@C 1082 DMSwarmCollectViewDestroy - Resets the DMSwarm to the size prior to calling DMSwarmCollectViewCreate() 1083 1084 Collective on dm 1085 1086 Input parameters: 1087 . dm - the DMSwarm 1088 1089 Notes: 1090 Users should call DMSwarmCollectViewCreate() before this function is called. 1091 1092 Level: advanced 1093 1094 .seealso: DMSwarmCollectViewCreate(), DMSwarmSetCollectType() 1095 @*/ 1096 PETSC_EXTERN PetscErrorCode DMSwarmCollectViewDestroy(DM dm) 1097 { 1098 PetscErrorCode ierr; 1099 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1100 1101 PetscFunctionBegin; 1102 if (!swarm->collect_view_active) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView is currently not active"); 1103 ierr = DMSwarmSetLocalSizes(dm,swarm->collect_view_reset_nlocal,-1);CHKERRQ(ierr); 1104 swarm->collect_view_active = PETSC_FALSE; 1105 PetscFunctionReturn(0); 1106 } 1107 1108 PetscErrorCode DMSwarmSetUpPIC(DM dm) 1109 { 1110 PetscInt dim; 1111 PetscErrorCode ierr; 1112 1113 PetscFunctionBegin; 1114 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 1115 if (dim < 1) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Dimension must be 1,2,3 - found %D",dim); 1116 if (dim > 3) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Dimension must be 1,2,3 - found %D",dim); 1117 ierr = DMSwarmRegisterPetscDatatypeField(dm,DMSwarmPICField_coor,dim,PETSC_DOUBLE);CHKERRQ(ierr); 1118 ierr = DMSwarmRegisterPetscDatatypeField(dm,DMSwarmPICField_cellid,1,PETSC_INT);CHKERRQ(ierr); 1119 PetscFunctionReturn(0); 1120 } 1121 1122 /*@C 1123 DMSwarmSetType - Set particular flavor of DMSwarm 1124 1125 Collective on dm 1126 1127 Input parameters: 1128 + dm - the DMSwarm 1129 - stype - the DMSwarm type (e.g. DMSWARM_PIC) 1130 1131 Level: advanced 1132 1133 .seealso: DMSwarmSetMigrateType(), DMSwarmSetCollectType() 1134 @*/ 1135 PETSC_EXTERN PetscErrorCode DMSwarmSetType(DM dm,DMSwarmType stype) 1136 { 1137 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1138 PetscErrorCode ierr; 1139 1140 PetscFunctionBegin; 1141 swarm->swarm_type = stype; 1142 if (swarm->swarm_type == DMSWARM_PIC) { 1143 ierr = DMSwarmSetUpPIC(dm);CHKERRQ(ierr); 1144 } 1145 PetscFunctionReturn(0); 1146 } 1147 1148 PetscErrorCode DMSetup_Swarm(DM dm) 1149 { 1150 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1151 PetscErrorCode ierr; 1152 PetscMPIInt rank; 1153 PetscInt p,npoints,*rankval; 1154 1155 PetscFunctionBegin; 1156 if (swarm->issetup) PetscFunctionReturn(0); 1157 1158 swarm->issetup = PETSC_TRUE; 1159 1160 if (swarm->swarm_type == DMSWARM_PIC) { 1161 /* check dmcell exists */ 1162 if (!swarm->dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires you call DMSwarmSetCellDM"); 1163 1164 if (swarm->dmcell->ops->locatepointssubdomain) { 1165 /* check methods exists for exact ownership identificiation */ 1166 ierr = PetscPrintf(PetscObjectComm((PetscObject)dm)," DMSWARM_PIC: Using method CellDM->ops->LocatePointsSubdomain\n");CHKERRQ(ierr); 1167 swarm->migrate_type = DMSWARM_MIGRATE_DMCELLEXACT; 1168 } else { 1169 /* check methods exist for point location AND rank neighbor identification */ 1170 if (swarm->dmcell->ops->locatepoints) { 1171 ierr = PetscPrintf(PetscObjectComm((PetscObject)dm)," DMSWARM_PIC: Using method CellDM->LocatePoints\n");CHKERRQ(ierr); 1172 } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires the method CellDM->ops->locatepoints be defined"); 1173 1174 if (swarm->dmcell->ops->getneighbors) { 1175 ierr = PetscPrintf(PetscObjectComm((PetscObject)dm)," DMSWARM_PIC: Using method CellDM->GetNeigbors\n");CHKERRQ(ierr); 1176 } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires the method CellDM->ops->getneighbors be defined"); 1177 1178 swarm->migrate_type = DMSWARM_MIGRATE_DMCELLNSCATTER; 1179 } 1180 } 1181 1182 ierr = DMSwarmFinalizeFieldRegister(dm);CHKERRQ(ierr); 1183 1184 /* check some fields were registered */ 1185 if (swarm->db->nfields <= 2) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"At least one field user must be registered via DMSwarmRegisterXXX()"); 1186 1187 /* check local sizes were set */ 1188 if (swarm->db->L == -1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Local sizes must be set via DMSwarmSetLocalSizes()"); 1189 1190 /* initialize values in pid and rank placeholders */ 1191 /* TODO: [pid - use MPI_Scan] */ 1192 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 1193 ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 1194 ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 1195 for (p=0; p<npoints; p++) { 1196 rankval[p] = (PetscInt)rank; 1197 } 1198 ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 1199 PetscFunctionReturn(0); 1200 } 1201 1202 extern PetscErrorCode DMSwarmSortDestroy(DMSwarmSort *_ctx); 1203 1204 PetscErrorCode DMDestroy_Swarm(DM dm) 1205 { 1206 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1207 PetscErrorCode ierr; 1208 1209 PetscFunctionBegin; 1210 ierr = DMSwarmDataBucketDestroy(&swarm->db);CHKERRQ(ierr); 1211 if (swarm->sort_context) { 1212 ierr = DMSwarmSortDestroy(&swarm->sort_context);CHKERRQ(ierr); 1213 } 1214 ierr = PetscFree(swarm);CHKERRQ(ierr); 1215 PetscFunctionReturn(0); 1216 } 1217 1218 PetscErrorCode DMSwarmView_Draw(DM dm, PetscViewer viewer) 1219 { 1220 DM cdm; 1221 PetscDraw draw; 1222 PetscReal *coords, oldPause; 1223 PetscInt Np, p, bs; 1224 PetscErrorCode ierr; 1225 1226 PetscFunctionBegin; 1227 ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr); 1228 ierr = DMSwarmGetCellDM(dm, &cdm);CHKERRQ(ierr); 1229 ierr = PetscDrawGetPause(draw, &oldPause);CHKERRQ(ierr); 1230 ierr = PetscDrawSetPause(draw, 0.0);CHKERRQ(ierr); 1231 ierr = DMView(cdm, viewer);CHKERRQ(ierr); 1232 ierr = PetscDrawSetPause(draw, oldPause);CHKERRQ(ierr); 1233 1234 ierr = DMSwarmGetLocalSize(dm, &Np);CHKERRQ(ierr); 1235 ierr = DMSwarmGetField(dm, DMSwarmPICField_coor, &bs, NULL, (void **) &coords);CHKERRQ(ierr); 1236 for (p = 0; p < Np; ++p) { 1237 const PetscInt i = p*bs; 1238 1239 ierr = PetscDrawEllipse(draw, coords[i], coords[i+1], 0.01, 0.01, PETSC_DRAW_BLUE);CHKERRQ(ierr); 1240 } 1241 ierr = DMSwarmRestoreField(dm, DMSwarmPICField_coor, &bs, NULL, (void **) &coords);CHKERRQ(ierr); 1242 ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 1243 ierr = PetscDrawPause(draw);CHKERRQ(ierr); 1244 ierr = PetscDrawSave(draw);CHKERRQ(ierr); 1245 PetscFunctionReturn(0); 1246 } 1247 1248 PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer) 1249 { 1250 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1251 PetscBool iascii,ibinary,ishdf5,isvtk,isdraw; 1252 PetscErrorCode ierr; 1253 1254 PetscFunctionBegin; 1255 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1256 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 1257 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 1258 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY,&ibinary);CHKERRQ(ierr); 1259 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk);CHKERRQ(ierr); 1260 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5);CHKERRQ(ierr); 1261 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw);CHKERRQ(ierr); 1262 if (iascii) { 1263 ierr = DMSwarmDataBucketView(PetscObjectComm((PetscObject)dm),swarm->db,NULL,DATABUCKET_VIEW_STDOUT);CHKERRQ(ierr); 1264 } else if (ibinary) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO Binary support"); 1265 #if defined(PETSC_HAVE_HDF5) 1266 else if (ishdf5) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO HDF5 support"); 1267 #else 1268 else if (ishdf5) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"HDF5 not supported. Please reconfigure using --download-hdf5"); 1269 #endif 1270 else if (isvtk) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support"); 1271 else if (isdraw) { 1272 ierr = DMSwarmView_Draw(dm, viewer);CHKERRQ(ierr); 1273 } 1274 PetscFunctionReturn(0); 1275 } 1276 1277 /*MC 1278 1279 DMSWARM = "swarm" - A DM object used to represent arrays of data (fields) of arbitrary data type. 1280 This implementation was designed for particle methods in which the underlying 1281 data required to be represented is both (i) dynamic in length, (ii) and of arbitrary data type. 1282 1283 User data can be represented by DMSwarm through a registering "fields". 1284 To register a field, the user must provide; 1285 (a) a unique name; 1286 (b) the data type (or size in bytes); 1287 (c) the block size of the data. 1288 1289 For example, suppose the application requires a unique id, energy, momentum and density to be stored 1290 on a set of of particles. Then the following code could be used 1291 1292 $ DMSwarmInitializeFieldRegister(dm) 1293 $ DMSwarmRegisterPetscDatatypeField(dm,"uid",1,PETSC_LONG); 1294 $ DMSwarmRegisterPetscDatatypeField(dm,"energy",1,PETSC_REAL); 1295 $ DMSwarmRegisterPetscDatatypeField(dm,"momentum",3,PETSC_REAL); 1296 $ DMSwarmRegisterPetscDatatypeField(dm,"density",1,PETSC_FLOAT); 1297 $ DMSwarmFinalizeFieldRegister(dm) 1298 1299 The fields represented by DMSwarm are dynamic and can be re-sized at any time. 1300 The only restriction imposed by DMSwarm is that all fields contain the same number of points. 1301 1302 To support particle methods, "migration" techniques are provided. These methods migrate data 1303 between MPI-ranks. 1304 1305 DMSwarm supports the methods DMCreateGlobalVector() and DMCreateLocalVector(). 1306 As a DMSwarm may internally define and store values of different data types, 1307 before calling DMCreateGlobalVector() or DMCreateLocalVector(), the user must inform DMSwarm which 1308 fields should be used to define a Vec object via 1309 DMSwarmVectorDefineField() 1310 The specified field can can changed be changed at any time - thereby permitting vectors 1311 compatable with different fields to be created. 1312 1313 A dual representation of fields in the DMSwarm and a Vec object is permitted via 1314 DMSwarmCreateGlobalVectorFromField() 1315 Here the data defining the field in the DMSwarm is shared with a Vec. 1316 This is inherently unsafe if you alter the size of the field at any time between 1317 calls to DMSwarmCreateGlobalVectorFromField() and DMSwarmDestroyGlobalVectorFromField(). 1318 If the local size of the DMSwarm does not match the local size of the global vector 1319 when DMSwarmDestroyGlobalVectorFromField() is called, an error is thrown. 1320 1321 Additional high-level support is provided for Particle-In-Cell methods. 1322 Please refer to the man page for DMSwarmSetType(). 1323 1324 Level: beginner 1325 1326 .seealso: DMType, DMCreate(), DMSetType() 1327 M*/ 1328 PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm) 1329 { 1330 DM_Swarm *swarm; 1331 PetscErrorCode ierr; 1332 1333 PetscFunctionBegin; 1334 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1335 ierr = PetscNewLog(dm,&swarm);CHKERRQ(ierr); 1336 dm->data = swarm; 1337 1338 ierr = DMSwarmDataBucketCreate(&swarm->db);CHKERRQ(ierr); 1339 ierr = DMSwarmInitializeFieldRegister(dm);CHKERRQ(ierr); 1340 1341 swarm->vec_field_set = PETSC_FALSE; 1342 swarm->issetup = PETSC_FALSE; 1343 swarm->swarm_type = DMSWARM_BASIC; 1344 swarm->migrate_type = DMSWARM_MIGRATE_BASIC; 1345 swarm->collect_type = DMSWARM_COLLECT_BASIC; 1346 swarm->migrate_error_on_missing_point = PETSC_FALSE; 1347 1348 swarm->dmcell = NULL; 1349 swarm->collect_view_active = PETSC_FALSE; 1350 swarm->collect_view_reset_nlocal = -1; 1351 1352 dm->dim = 0; 1353 dm->ops->view = DMView_Swarm; 1354 dm->ops->load = NULL; 1355 dm->ops->setfromoptions = NULL; 1356 dm->ops->clone = NULL; 1357 dm->ops->setup = DMSetup_Swarm; 1358 dm->ops->createdefaultsection = NULL; 1359 dm->ops->createdefaultconstraints = NULL; 1360 dm->ops->createglobalvector = DMCreateGlobalVector_Swarm; 1361 dm->ops->createlocalvector = DMCreateLocalVector_Swarm; 1362 dm->ops->getlocaltoglobalmapping = NULL; 1363 dm->ops->createfieldis = NULL; 1364 dm->ops->createcoordinatedm = NULL; 1365 dm->ops->getcoloring = NULL; 1366 dm->ops->creatematrix = NULL; 1367 dm->ops->createinterpolation = NULL; 1368 dm->ops->createinjection = NULL; 1369 dm->ops->createmassmatrix = DMCreateMassMatrix_Swarm; 1370 dm->ops->refine = NULL; 1371 dm->ops->coarsen = NULL; 1372 dm->ops->refinehierarchy = NULL; 1373 dm->ops->coarsenhierarchy = NULL; 1374 dm->ops->globaltolocalbegin = NULL; 1375 dm->ops->globaltolocalend = NULL; 1376 dm->ops->localtoglobalbegin = NULL; 1377 dm->ops->localtoglobalend = NULL; 1378 dm->ops->destroy = DMDestroy_Swarm; 1379 dm->ops->createsubdm = NULL; 1380 dm->ops->getdimpoints = NULL; 1381 dm->ops->locatepoints = NULL; 1382 PetscFunctionReturn(0); 1383 } 1384