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