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