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