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