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 <petscblaslapack.h> 9 #include "../src/dm/impls/swarm/data_bucket.h" 10 #include <petscdmlabel.h> 11 #include <petscsection.h> 12 13 PetscLogEvent DMSWARM_Migrate, DMSWARM_SetSizes, DMSWARM_AddPoints, DMSWARM_RemovePoints, DMSWARM_Sort; 14 PetscLogEvent DMSWARM_DataExchangerTopologySetup, DMSWARM_DataExchangerBegin, DMSWARM_DataExchangerEnd; 15 PetscLogEvent DMSWARM_DataExchangerSendCount, DMSWARM_DataExchangerPack; 16 17 const char* DMSwarmTypeNames[] = { "basic", "pic", NULL }; 18 const char* DMSwarmMigrateTypeNames[] = { "basic", "dmcellnscatter", "dmcellexact", "user", NULL }; 19 const char* DMSwarmCollectTypeNames[] = { "basic", "boundingbox", "general", "user", NULL }; 20 const char* DMSwarmPICLayoutTypeNames[] = { "regular", "gauss", "subdivision", NULL }; 21 22 const char DMSwarmField_pid[] = "DMSwarm_pid"; 23 const char DMSwarmField_rank[] = "DMSwarm_rank"; 24 const char DMSwarmPICField_coor[] = "DMSwarmPIC_coor"; 25 const char DMSwarmPICField_cellid[] = "DMSwarm_cellid"; 26 27 #if defined(PETSC_HAVE_HDF5) 28 #include <petscviewerhdf5.h> 29 30 PetscErrorCode VecView_Swarm_HDF5_Internal(Vec v, PetscViewer viewer) 31 { 32 DM dm; 33 PetscReal seqval; 34 PetscInt seqnum, bs; 35 PetscBool isseq; 36 PetscErrorCode ierr; 37 38 PetscFunctionBegin; 39 ierr = VecGetDM(v, &dm);CHKERRQ(ierr); 40 ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr); 41 ierr = PetscViewerHDF5PushGroup(viewer, "/particle_fields");CHKERRQ(ierr); 42 ierr = PetscObjectTypeCompare((PetscObject) v, VECSEQ, &isseq);CHKERRQ(ierr); 43 ierr = DMGetOutputSequenceNumber(dm, &seqnum, &seqval);CHKERRQ(ierr); 44 ierr = PetscViewerHDF5SetTimestep(viewer, seqnum);CHKERRQ(ierr); 45 /* ierr = DMSequenceView_HDF5(dm, "time", seqnum, (PetscScalar) seqval, viewer);CHKERRQ(ierr); */ 46 ierr = VecViewNative(v, viewer);CHKERRQ(ierr); 47 ierr = PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) v, "Nc", PETSC_INT, (void *) &bs);CHKERRQ(ierr); 48 ierr = PetscViewerHDF5PopGroup(viewer);CHKERRQ(ierr); 49 PetscFunctionReturn(0); 50 } 51 52 PetscErrorCode DMSwarmView_HDF5(DM dm, PetscViewer viewer) 53 { 54 Vec coordinates; 55 PetscInt Np; 56 PetscBool isseq; 57 PetscErrorCode ierr; 58 59 PetscFunctionBegin; 60 ierr = DMSwarmGetSize(dm, &Np);CHKERRQ(ierr); 61 ierr = DMSwarmCreateGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates);CHKERRQ(ierr); 62 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 63 ierr = PetscViewerHDF5PushGroup(viewer, "/particles");CHKERRQ(ierr); 64 ierr = PetscObjectTypeCompare((PetscObject) coordinates, VECSEQ, &isseq);CHKERRQ(ierr); 65 ierr = VecViewNative(coordinates, viewer);CHKERRQ(ierr); 66 ierr = PetscViewerHDF5WriteObjectAttribute(viewer, (PetscObject) coordinates, "Np", PETSC_INT, (void *) &Np);CHKERRQ(ierr); 67 ierr = PetscViewerHDF5PopGroup(viewer);CHKERRQ(ierr); 68 ierr = DMSwarmDestroyGlobalVectorFromField(dm, DMSwarmPICField_coor, &coordinates);CHKERRQ(ierr); 69 PetscFunctionReturn(0); 70 } 71 #endif 72 73 PetscErrorCode VecView_Swarm(Vec v, PetscViewer viewer) 74 { 75 DM dm; 76 #if defined(PETSC_HAVE_HDF5) 77 PetscBool ishdf5; 78 #endif 79 PetscErrorCode ierr; 80 81 PetscFunctionBegin; 82 ierr = VecGetDM(v, &dm);CHKERRQ(ierr); 83 if (!dm) SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM"); 84 #if defined(PETSC_HAVE_HDF5) 85 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);CHKERRQ(ierr); 86 if (ishdf5) { 87 ierr = VecView_Swarm_HDF5_Internal(v, viewer);CHKERRQ(ierr); 88 PetscFunctionReturn(0); 89 } 90 #endif 91 ierr = VecViewNative(v, viewer);CHKERRQ(ierr); 92 PetscFunctionReturn(0); 93 } 94 95 /*@C 96 DMSwarmVectorDefineField - Sets the field from which to define a Vec object 97 when DMCreateLocalVector(), or DMCreateGlobalVector() is called 98 99 Collective on dm 100 101 Input parameters: 102 + dm - a DMSwarm 103 - fieldname - the textual name given to a registered field 104 105 Level: beginner 106 107 Notes: 108 109 The field with name fieldname must be defined as having a data type of PetscScalar. 110 111 This function must be called prior to calling DMCreateLocalVector(), DMCreateGlobalVector(). 112 Mutiple calls to DMSwarmVectorDefineField() are permitted. 113 114 .seealso: DMSwarmRegisterPetscDatatypeField(), DMCreateGlobalVector(), DMCreateLocalVector() 115 @*/ 116 PetscErrorCode DMSwarmVectorDefineField(DM dm,const char fieldname[]) 117 { 118 DM_Swarm *swarm = (DM_Swarm*)dm->data; 119 PetscErrorCode ierr; 120 PetscInt bs,n; 121 PetscScalar *array; 122 PetscDataType type; 123 124 PetscFunctionBegin; 125 if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); } 126 ierr = DMSwarmDataBucketGetSizes(swarm->db,&n,NULL,NULL);CHKERRQ(ierr); 127 ierr = DMSwarmGetField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr); 128 129 /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */ 130 if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid for PETSC_REAL"); 131 ierr = PetscSNPrintf(swarm->vec_field_name,PETSC_MAX_PATH_LEN-1,"%s",fieldname);CHKERRQ(ierr); 132 swarm->vec_field_set = PETSC_TRUE; 133 swarm->vec_field_bs = bs; 134 swarm->vec_field_nlocal = n; 135 ierr = DMSwarmRestoreField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr); 136 PetscFunctionReturn(0); 137 } 138 139 /* requires DMSwarmDefineFieldVector has been called */ 140 PetscErrorCode DMCreateGlobalVector_Swarm(DM dm,Vec *vec) 141 { 142 DM_Swarm *swarm = (DM_Swarm*)dm->data; 143 PetscErrorCode ierr; 144 Vec x; 145 char name[PETSC_MAX_PATH_LEN]; 146 147 PetscFunctionBegin; 148 if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); } 149 if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first"); 150 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 */ 151 152 ierr = PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name);CHKERRQ(ierr); 153 ierr = VecCreate(PetscObjectComm((PetscObject)dm),&x);CHKERRQ(ierr); 154 ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr); 155 ierr = VecSetSizes(x,swarm->db->L*swarm->vec_field_bs,PETSC_DETERMINE);CHKERRQ(ierr); 156 ierr = VecSetBlockSize(x,swarm->vec_field_bs);CHKERRQ(ierr); 157 ierr = VecSetDM(x,dm);CHKERRQ(ierr); 158 ierr = VecSetFromOptions(x);CHKERRQ(ierr); 159 ierr = VecSetDM(x, dm);CHKERRQ(ierr); 160 ierr = VecSetOperation(x, VECOP_VIEW, (void (*)(void)) VecView_Swarm);CHKERRQ(ierr); 161 *vec = x; 162 PetscFunctionReturn(0); 163 } 164 165 /* requires DMSwarmDefineFieldVector has been called */ 166 PetscErrorCode DMCreateLocalVector_Swarm(DM dm,Vec *vec) 167 { 168 DM_Swarm *swarm = (DM_Swarm*)dm->data; 169 PetscErrorCode ierr; 170 Vec x; 171 char name[PETSC_MAX_PATH_LEN]; 172 173 PetscFunctionBegin; 174 if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); } 175 if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first"); 176 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 */ 177 178 ierr = PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name);CHKERRQ(ierr); 179 ierr = VecCreate(PETSC_COMM_SELF,&x);CHKERRQ(ierr); 180 ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr); 181 ierr = VecSetSizes(x,swarm->db->L*swarm->vec_field_bs,PETSC_DETERMINE);CHKERRQ(ierr); 182 ierr = VecSetBlockSize(x,swarm->vec_field_bs);CHKERRQ(ierr); 183 ierr = VecSetDM(x,dm);CHKERRQ(ierr); 184 ierr = VecSetFromOptions(x);CHKERRQ(ierr); 185 *vec = x; 186 PetscFunctionReturn(0); 187 } 188 189 static PetscErrorCode DMSwarmDestroyVectorFromField_Private(DM dm, const char fieldname[], Vec *vec) 190 { 191 DM_Swarm *swarm = (DM_Swarm *) dm->data; 192 DMSwarmDataField gfield; 193 void (*fptr)(void); 194 PetscInt bs, nlocal; 195 char name[PETSC_MAX_PATH_LEN]; 196 PetscErrorCode ierr; 197 198 PetscFunctionBegin; 199 ierr = VecGetLocalSize(*vec, &nlocal);CHKERRQ(ierr); 200 ierr = VecGetBlockSize(*vec, &bs);CHKERRQ(ierr); 201 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 */ 202 ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db, fieldname, &gfield);CHKERRQ(ierr); 203 /* check vector is an inplace array */ 204 ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarm_VecFieldInPlace_%s", fieldname);CHKERRQ(ierr); 205 ierr = PetscObjectQueryFunction((PetscObject) *vec, name, &fptr);CHKERRQ(ierr); 206 if (!fptr) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_USER, "Vector being destroyed was not created from DMSwarm field(%s)", fieldname); 207 ierr = DMSwarmDataFieldRestoreAccess(gfield);CHKERRQ(ierr); 208 ierr = VecDestroy(vec);CHKERRQ(ierr); 209 PetscFunctionReturn(0); 210 } 211 212 static PetscErrorCode DMSwarmCreateVectorFromField_Private(DM dm, const char fieldname[], MPI_Comm comm, Vec *vec) 213 { 214 DM_Swarm *swarm = (DM_Swarm *) dm->data; 215 PetscDataType type; 216 PetscScalar *array; 217 PetscInt bs, n; 218 char name[PETSC_MAX_PATH_LEN]; 219 PetscMPIInt size; 220 PetscErrorCode ierr; 221 222 PetscFunctionBegin; 223 if (!swarm->issetup) {ierr = DMSetUp(dm);CHKERRQ(ierr);} 224 ierr = DMSwarmDataBucketGetSizes(swarm->db, &n, NULL, NULL);CHKERRQ(ierr); 225 ierr = DMSwarmGetField(dm, fieldname, &bs, &type, (void **) &array);CHKERRQ(ierr); 226 if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Only valid for PETSC_REAL"); 227 228 ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr); 229 if (size == 1) { 230 ierr = VecCreateSeqWithArray(comm, bs, n*bs, array, vec);CHKERRQ(ierr); 231 } else { 232 ierr = VecCreateMPIWithArray(comm, bs, n*bs, PETSC_DETERMINE, array, vec);CHKERRQ(ierr); 233 } 234 ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarmSharedField_%s", fieldname);CHKERRQ(ierr); 235 ierr = PetscObjectSetName((PetscObject) *vec, name);CHKERRQ(ierr); 236 237 /* Set guard */ 238 ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "DMSwarm_VecFieldInPlace_%s", fieldname);CHKERRQ(ierr); 239 ierr = PetscObjectComposeFunction((PetscObject) *vec, name, DMSwarmDestroyVectorFromField_Private);CHKERRQ(ierr); 240 241 ierr = VecSetDM(*vec, dm);CHKERRQ(ierr); 242 ierr = VecSetOperation(*vec, VECOP_VIEW, (void (*)(void)) VecView_Swarm);CHKERRQ(ierr); 243 PetscFunctionReturn(0); 244 } 245 246 /* This creates a "mass matrix" between a finite element and particle space. If a finite element interpolant is given by 247 248 \hat f = \sum_i f_i \phi_i 249 250 and a particle function is given by 251 252 f = \sum_i w_i \delta(x - x_i) 253 254 then we want to require that 255 256 M \hat f = M_p f 257 258 where the particle mass matrix is given by 259 260 (M_p)_{ij} = \int \phi_i \delta(x - x_j) 261 262 The way Dave May does particles, they amount to quadratue weights rather than delta functions, so he has |J| is in 263 his integral. We allow this with the boolean flag. 264 */ 265 static PetscErrorCode DMSwarmComputeMassMatrix_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx) 266 { 267 const char *name = "Mass Matrix"; 268 MPI_Comm comm; 269 PetscDS prob; 270 PetscSection fsection, globalFSection; 271 PetscHSetIJ ht; 272 PetscLayout rLayout, colLayout; 273 PetscInt *dnz, *onz; 274 PetscInt locRows, locCols, rStart, colStart, colEnd, *rowIDXs; 275 PetscReal *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0}; 276 PetscScalar *elemMat; 277 PetscInt dim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0; 278 PetscErrorCode ierr; 279 280 PetscFunctionBegin; 281 ierr = PetscObjectGetComm((PetscObject) mass, &comm);CHKERRQ(ierr); 282 ierr = DMGetCoordinateDim(dmf, &dim);CHKERRQ(ierr); 283 ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr); 284 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 285 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 286 ierr = PetscMalloc3(dim, &v0, dim*dim, &J, dim*dim,&invJ);CHKERRQ(ierr); 287 ierr = DMGetLocalSection(dmf, &fsection);CHKERRQ(ierr); 288 ierr = DMGetGlobalSection(dmf, &globalFSection);CHKERRQ(ierr); 289 ierr = DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);CHKERRQ(ierr); 290 ierr = MatGetLocalSize(mass, &locRows, &locCols);CHKERRQ(ierr); 291 292 ierr = PetscLayoutCreate(comm, &colLayout);CHKERRQ(ierr); 293 ierr = PetscLayoutSetLocalSize(colLayout, locCols);CHKERRQ(ierr); 294 ierr = PetscLayoutSetBlockSize(colLayout, 1);CHKERRQ(ierr); 295 ierr = PetscLayoutSetUp(colLayout);CHKERRQ(ierr); 296 ierr = PetscLayoutGetRange(colLayout, &colStart, &colEnd);CHKERRQ(ierr); 297 ierr = PetscLayoutDestroy(&colLayout);CHKERRQ(ierr); 298 299 ierr = PetscLayoutCreate(comm, &rLayout);CHKERRQ(ierr); 300 ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr); 301 ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr); 302 ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr); 303 ierr = PetscLayoutGetRange(rLayout, &rStart, NULL);CHKERRQ(ierr); 304 ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr); 305 306 ierr = PetscCalloc2(locRows, &dnz, locRows, &onz);CHKERRQ(ierr); 307 ierr = PetscHSetIJCreate(&ht);CHKERRQ(ierr); 308 309 ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr); 310 /* count non-zeros */ 311 ierr = DMSwarmSortGetAccess(dmc);CHKERRQ(ierr); 312 for (field = 0; field < Nf; ++field) { 313 for (cell = cStart; cell < cEnd; ++cell) { 314 PetscInt c, i; 315 PetscInt *findices, *cindices; /* fine is vertices, coarse is particles */ 316 PetscInt numFIndices, numCIndices; 317 318 ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);CHKERRQ(ierr); 319 ierr = DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices);CHKERRQ(ierr); 320 maxC = PetscMax(maxC, numCIndices); 321 { 322 PetscHashIJKey key; 323 PetscBool missing; 324 for (i = 0; i < numFIndices; ++i) { 325 key.j = findices[i]; /* global column (from Plex) */ 326 if (key.j >= 0) { 327 /* Get indices for coarse elements */ 328 for (c = 0; c < numCIndices; ++c) { 329 key.i = cindices[c] + rStart; /* global cols (from Swarm) */ 330 if (key.i < 0) continue; 331 ierr = PetscHSetIJQueryAdd(ht, key, &missing);CHKERRQ(ierr); 332 if (missing) { 333 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart]; 334 else ++onz[key.i - rStart]; 335 } else SETERRQ2(PetscObjectComm((PetscObject) dmf), PETSC_ERR_SUP, "Set new value at %D,%D", key.i, key.j); 336 } 337 } 338 } 339 ierr = PetscFree(cindices);CHKERRQ(ierr); 340 } 341 ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);CHKERRQ(ierr); 342 } 343 } 344 ierr = PetscHSetIJDestroy(&ht);CHKERRQ(ierr); 345 ierr = MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr); 346 ierr = MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE);CHKERRQ(ierr); 347 ierr = PetscFree2(dnz, onz);CHKERRQ(ierr); 348 ierr = PetscMalloc3(maxC*totDim, &elemMat, maxC, &rowIDXs, maxC*dim, &xi);CHKERRQ(ierr); 349 for (field = 0; field < Nf; ++field) { 350 PetscTabulation Tcoarse; 351 PetscObject obj; 352 PetscReal *coords; 353 PetscInt Nc, i; 354 355 ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr); 356 ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr); 357 if (Nc != 1) SETERRQ1(PetscObjectComm((PetscObject) dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %D", Nc); 358 ierr = DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 359 for (cell = cStart; cell < cEnd; ++cell) { 360 PetscInt *findices , *cindices; 361 PetscInt numFIndices, numCIndices; 362 PetscInt p, c; 363 364 /* TODO: Use DMField instead of assuming affine */ 365 ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 366 ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);CHKERRQ(ierr); 367 ierr = DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices);CHKERRQ(ierr); 368 for (p = 0; p < numCIndices; ++p) { 369 CoordinatesRealToRef(dim, dim, v0ref, v0, invJ, &coords[cindices[p]*dim], &xi[p*dim]); 370 } 371 ierr = PetscFECreateTabulation((PetscFE) obj, 1, numCIndices, xi, 0, &Tcoarse);CHKERRQ(ierr); 372 /* Get elemMat entries by multiplying by weight */ 373 ierr = PetscArrayzero(elemMat, numCIndices*totDim);CHKERRQ(ierr); 374 for (i = 0; i < numFIndices; ++i) { 375 for (p = 0; p < numCIndices; ++p) { 376 for (c = 0; c < Nc; ++c) { 377 /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */ 378 elemMat[p*numFIndices+i] += Tcoarse->T[0][(p*numFIndices + i)*Nc + c]*(useDeltaFunction ? 1.0 : detJ); 379 } 380 } 381 } 382 for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart; 383 if (0) {ierr = DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);CHKERRQ(ierr);} 384 ierr = MatSetValues(mass, numCIndices, rowIDXs, numFIndices, findices, elemMat, ADD_VALUES);CHKERRQ(ierr); 385 ierr = PetscFree(cindices);CHKERRQ(ierr); 386 ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);CHKERRQ(ierr); 387 ierr = PetscTabulationDestroy(&Tcoarse);CHKERRQ(ierr); 388 } 389 ierr = DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 390 } 391 ierr = PetscFree3(elemMat, rowIDXs, xi);CHKERRQ(ierr); 392 ierr = DMSwarmSortRestoreAccess(dmc);CHKERRQ(ierr); 393 ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr); 394 ierr = MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 395 ierr = MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 396 PetscFunctionReturn(0); 397 } 398 399 /* Returns empty matrix for use with SNES FD */ 400 static PetscErrorCode DMCreateMatrix_Swarm(DM sw, Mat* m) 401 { 402 Vec field; 403 PetscInt size; 404 PetscErrorCode ierr; 405 406 PetscFunctionBegin; 407 ierr = DMGetGlobalVector(sw, &field);CHKERRQ(ierr); 408 ierr = VecGetLocalSize(field, &size);CHKERRQ(ierr); 409 ierr = DMRestoreGlobalVector(sw, &field);CHKERRQ(ierr); 410 ierr = MatCreate(PETSC_COMM_WORLD, m);CHKERRQ(ierr); 411 ierr = MatSetFromOptions(*m);CHKERRQ(ierr); 412 ierr = MatSetSizes(*m, PETSC_DECIDE, PETSC_DECIDE, size, size);CHKERRQ(ierr); 413 ierr = MatSeqAIJSetPreallocation(*m, 1, NULL);CHKERRQ(ierr); 414 ierr = MatZeroEntries(*m);CHKERRQ(ierr); 415 ierr = MatAssemblyBegin(*m, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 416 ierr = MatAssemblyEnd(*m, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 417 ierr = MatShift(*m, 1.0);CHKERRQ(ierr); 418 ierr = MatSetDM(*m, sw);CHKERRQ(ierr); 419 PetscFunctionReturn(0); 420 } 421 422 /* FEM cols, Particle rows */ 423 static PetscErrorCode DMCreateMassMatrix_Swarm(DM dmCoarse, DM dmFine, Mat *mass) 424 { 425 PetscSection gsf; 426 PetscInt m, n; 427 void *ctx; 428 PetscErrorCode ierr; 429 430 PetscFunctionBegin; 431 ierr = DMGetGlobalSection(dmFine, &gsf);CHKERRQ(ierr); 432 ierr = PetscSectionGetConstrainedStorageSize(gsf, &m);CHKERRQ(ierr); 433 ierr = DMSwarmGetLocalSize(dmCoarse, &n);CHKERRQ(ierr); 434 ierr = MatCreate(PetscObjectComm((PetscObject) dmCoarse), mass);CHKERRQ(ierr); 435 ierr = MatSetSizes(*mass, n, m, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 436 ierr = MatSetType(*mass, dmCoarse->mattype);CHKERRQ(ierr); 437 ierr = DMGetApplicationContext(dmFine, &ctx);CHKERRQ(ierr); 438 439 ierr = DMSwarmComputeMassMatrix_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx);CHKERRQ(ierr); 440 ierr = MatViewFromOptions(*mass, NULL, "-mass_mat_view");CHKERRQ(ierr); 441 PetscFunctionReturn(0); 442 } 443 444 static PetscErrorCode DMSwarmComputeMassMatrixSquare_Private(DM dmc, DM dmf, Mat mass, PetscBool useDeltaFunction, void *ctx) 445 { 446 const char *name = "Mass Matrix Square"; 447 MPI_Comm comm; 448 PetscDS prob; 449 PetscSection fsection, globalFSection; 450 PetscHSetIJ ht; 451 PetscLayout rLayout, colLayout; 452 PetscInt *dnz, *onz, *adj, depth, maxConeSize, maxSupportSize, maxAdjSize; 453 PetscInt locRows, locCols, rStart, colStart, colEnd, *rowIDXs; 454 PetscReal *xi, *v0, *J, *invJ, detJ = 1.0, v0ref[3] = {-1.0, -1.0, -1.0}; 455 PetscScalar *elemMat, *elemMatSq; 456 PetscInt cdim, Nf, field, cStart, cEnd, cell, totDim, maxC = 0; 457 PetscErrorCode ierr; 458 459 PetscFunctionBegin; 460 ierr = PetscObjectGetComm((PetscObject) mass, &comm);CHKERRQ(ierr); 461 ierr = DMGetCoordinateDim(dmf, &cdim);CHKERRQ(ierr); 462 ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr); 463 ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 464 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 465 ierr = PetscMalloc3(cdim, &v0, cdim*cdim, &J, cdim*cdim,&invJ);CHKERRQ(ierr); 466 ierr = DMGetLocalSection(dmf, &fsection);CHKERRQ(ierr); 467 ierr = DMGetGlobalSection(dmf, &globalFSection);CHKERRQ(ierr); 468 ierr = DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);CHKERRQ(ierr); 469 ierr = MatGetLocalSize(mass, &locRows, &locCols);CHKERRQ(ierr); 470 471 ierr = PetscLayoutCreate(comm, &colLayout);CHKERRQ(ierr); 472 ierr = PetscLayoutSetLocalSize(colLayout, locCols);CHKERRQ(ierr); 473 ierr = PetscLayoutSetBlockSize(colLayout, 1);CHKERRQ(ierr); 474 ierr = PetscLayoutSetUp(colLayout);CHKERRQ(ierr); 475 ierr = PetscLayoutGetRange(colLayout, &colStart, &colEnd);CHKERRQ(ierr); 476 ierr = PetscLayoutDestroy(&colLayout);CHKERRQ(ierr); 477 478 ierr = PetscLayoutCreate(comm, &rLayout);CHKERRQ(ierr); 479 ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr); 480 ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr); 481 ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr); 482 ierr = PetscLayoutGetRange(rLayout, &rStart, NULL);CHKERRQ(ierr); 483 ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr); 484 485 ierr = DMPlexGetDepth(dmf, &depth);CHKERRQ(ierr); 486 ierr = DMPlexGetMaxSizes(dmf, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 487 maxAdjSize = PetscPowInt(maxConeSize*maxSupportSize, depth); 488 ierr = PetscMalloc1(maxAdjSize, &adj);CHKERRQ(ierr); 489 490 ierr = PetscCalloc2(locRows, &dnz, locRows, &onz);CHKERRQ(ierr); 491 ierr = PetscHSetIJCreate(&ht);CHKERRQ(ierr); 492 /* Count nonzeros 493 This is just FVM++, but we cannot use the Plex P0 allocation since unknowns in a cell will not be contiguous 494 */ 495 ierr = DMSwarmSortGetAccess(dmc);CHKERRQ(ierr); 496 for (cell = cStart; cell < cEnd; ++cell) { 497 PetscInt i; 498 PetscInt *cindices; 499 PetscInt numCIndices; 500 #if 0 501 PetscInt adjSize = maxAdjSize, a, j; 502 #endif 503 504 ierr = DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices);CHKERRQ(ierr); 505 maxC = PetscMax(maxC, numCIndices); 506 /* Diagonal block */ 507 for (i = 0; i < numCIndices; ++i) {dnz[cindices[i]] += numCIndices;} 508 #if 0 509 /* Off-diagonal blocks */ 510 ierr = DMPlexGetAdjacency(dmf, cell, &adjSize, &adj);CHKERRQ(ierr); 511 for (a = 0; a < adjSize; ++a) { 512 if (adj[a] >= cStart && adj[a] < cEnd && adj[a] != cell) { 513 const PetscInt ncell = adj[a]; 514 PetscInt *ncindices; 515 PetscInt numNCIndices; 516 517 ierr = DMSwarmSortGetPointsPerCell(dmc, ncell, &numNCIndices, &ncindices);CHKERRQ(ierr); 518 { 519 PetscHashIJKey key; 520 PetscBool missing; 521 522 for (i = 0; i < numCIndices; ++i) { 523 key.i = cindices[i] + rStart; /* global rows (from Swarm) */ 524 if (key.i < 0) continue; 525 for (j = 0; j < numNCIndices; ++j) { 526 key.j = ncindices[j] + rStart; /* global column (from Swarm) */ 527 if (key.j < 0) continue; 528 ierr = PetscHSetIJQueryAdd(ht, key, &missing);CHKERRQ(ierr); 529 if (missing) { 530 if ((key.j >= colStart) && (key.j < colEnd)) ++dnz[key.i - rStart]; 531 else ++onz[key.i - rStart]; 532 } 533 } 534 } 535 } 536 ierr = PetscFree(ncindices);CHKERRQ(ierr); 537 } 538 } 539 #endif 540 ierr = PetscFree(cindices);CHKERRQ(ierr); 541 } 542 ierr = PetscHSetIJDestroy(&ht);CHKERRQ(ierr); 543 ierr = MatXAIJSetPreallocation(mass, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr); 544 ierr = MatSetOption(mass, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE);CHKERRQ(ierr); 545 ierr = PetscFree2(dnz, onz);CHKERRQ(ierr); 546 ierr = PetscMalloc4(maxC*totDim, &elemMat, maxC*maxC, &elemMatSq, maxC, &rowIDXs, maxC*cdim, &xi);CHKERRQ(ierr); 547 /* Fill in values 548 Each entry is a sum of terms \phi_i(x_p) \phi_i(x_q) 549 Start just by producing block diagonal 550 Could loop over adjacent cells 551 Produce neighboring element matrix 552 TODO Determine which columns and rows correspond to shared dual vector 553 Do MatMatMult with rectangular matrices 554 Insert block 555 */ 556 for (field = 0; field < Nf; ++field) { 557 PetscTabulation Tcoarse; 558 PetscObject obj; 559 PetscReal *coords; 560 PetscInt Nc, i; 561 562 ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr); 563 ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr); 564 if (Nc != 1) SETERRQ1(PetscObjectComm((PetscObject) dmf), PETSC_ERR_SUP, "Can only interpolate a scalar field from particles, Nc = %D", Nc); 565 ierr = DMSwarmGetField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 566 for (cell = cStart; cell < cEnd; ++cell) { 567 PetscInt *findices , *cindices; 568 PetscInt numFIndices, numCIndices; 569 PetscInt p, c; 570 571 /* TODO: Use DMField instead of assuming affine */ 572 ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 573 ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);CHKERRQ(ierr); 574 ierr = DMSwarmSortGetPointsPerCell(dmc, cell, &numCIndices, &cindices);CHKERRQ(ierr); 575 for (p = 0; p < numCIndices; ++p) { 576 CoordinatesRealToRef(cdim, cdim, v0ref, v0, invJ, &coords[cindices[p]*cdim], &xi[p*cdim]); 577 } 578 ierr = PetscFECreateTabulation((PetscFE) obj, 1, numCIndices, xi, 0, &Tcoarse);CHKERRQ(ierr); 579 /* Get elemMat entries by multiplying by weight */ 580 ierr = PetscArrayzero(elemMat, numCIndices*totDim);CHKERRQ(ierr); 581 for (i = 0; i < numFIndices; ++i) { 582 for (p = 0; p < numCIndices; ++p) { 583 for (c = 0; c < Nc; ++c) { 584 /* B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c */ 585 elemMat[p*numFIndices+i] += Tcoarse->T[0][(p*numFIndices + i)*Nc + c]*(useDeltaFunction ? 1.0 : detJ); 586 } 587 } 588 } 589 ierr = PetscTabulationDestroy(&Tcoarse);CHKERRQ(ierr); 590 for (p = 0; p < numCIndices; ++p) rowIDXs[p] = cindices[p] + rStart; 591 if (0) {ierr = DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);CHKERRQ(ierr);} 592 /* Block diagonal */ 593 if (numCIndices) { 594 PetscBLASInt blasn, blask; 595 PetscScalar one = 1.0, zero = 0.0; 596 597 ierr = PetscBLASIntCast(numCIndices, &blasn);CHKERRQ(ierr); 598 ierr = PetscBLASIntCast(numFIndices, &blask);CHKERRQ(ierr); 599 PetscStackCallBLAS("BLASgemm",BLASgemm_("T","N",&blasn,&blasn,&blask,&one,elemMat,&blask,elemMat,&blask,&zero,elemMatSq,&blasn)); 600 } 601 ierr = MatSetValues(mass, numCIndices, rowIDXs, numCIndices, rowIDXs, elemMatSq, ADD_VALUES);CHKERRQ(ierr); 602 /* TODO Off-diagonal */ 603 ierr = PetscFree(cindices);CHKERRQ(ierr); 604 ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, PETSC_FALSE, &numFIndices, &findices, NULL, NULL);CHKERRQ(ierr); 605 } 606 ierr = DMSwarmRestoreField(dmc, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 607 } 608 ierr = PetscFree4(elemMat, elemMatSq, rowIDXs, xi);CHKERRQ(ierr); 609 ierr = PetscFree(adj);CHKERRQ(ierr); 610 ierr = DMSwarmSortRestoreAccess(dmc);CHKERRQ(ierr); 611 ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr); 612 ierr = MatAssemblyBegin(mass, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 613 ierr = MatAssemblyEnd(mass, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 614 PetscFunctionReturn(0); 615 } 616 617 /*@C 618 DMSwarmCreateMassMatrixSquare - Creates the block-diagonal of the square, M^T_p M_p, of the particle mass matrix M_p 619 620 Collective on dmCoarse 621 622 Input parameters: 623 + dmCoarse - a DMSwarm 624 - dmFine - a DMPlex 625 626 Output parameter: 627 . mass - the square of the particle mass matrix 628 629 Level: advanced 630 631 Notes: 632 We only compute the block diagonal since this provides a good preconditioner and is completely local. It would be possible in the 633 future to compute the full normal equations. 634 635 .seealso: DMCreateMassMatrix() 636 @*/ 637 PetscErrorCode DMSwarmCreateMassMatrixSquare(DM dmCoarse, DM dmFine, Mat *mass) 638 { 639 PetscInt n; 640 void *ctx; 641 PetscErrorCode ierr; 642 643 PetscFunctionBegin; 644 ierr = DMSwarmGetLocalSize(dmCoarse, &n);CHKERRQ(ierr); 645 ierr = MatCreate(PetscObjectComm((PetscObject) dmCoarse), mass);CHKERRQ(ierr); 646 ierr = MatSetSizes(*mass, n, n, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 647 ierr = MatSetType(*mass, dmCoarse->mattype);CHKERRQ(ierr); 648 ierr = DMGetApplicationContext(dmFine, &ctx);CHKERRQ(ierr); 649 650 ierr = DMSwarmComputeMassMatrixSquare_Private(dmCoarse, dmFine, *mass, PETSC_TRUE, ctx);CHKERRQ(ierr); 651 ierr = MatViewFromOptions(*mass, NULL, "-mass_sq_mat_view");CHKERRQ(ierr); 652 PetscFunctionReturn(0); 653 } 654 655 /*@C 656 DMSwarmCreateGlobalVectorFromField - Creates a Vec object sharing the array associated with a given field 657 658 Collective on dm 659 660 Input parameters: 661 + dm - a DMSwarm 662 - fieldname - the textual name given to a registered field 663 664 Output parameter: 665 . vec - the vector 666 667 Level: beginner 668 669 Notes: 670 The vector must be returned using a matching call to DMSwarmDestroyGlobalVectorFromField(). 671 672 .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmDestroyGlobalVectorFromField() 673 @*/ 674 PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec) 675 { 676 MPI_Comm comm = PetscObjectComm((PetscObject) dm); 677 PetscErrorCode ierr; 678 679 PetscFunctionBegin; 680 ierr = DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec);CHKERRQ(ierr); 681 PetscFunctionReturn(0); 682 } 683 684 /*@C 685 DMSwarmDestroyGlobalVectorFromField - Destroys the Vec object which share the array associated with a given field 686 687 Collective on dm 688 689 Input parameters: 690 + dm - a DMSwarm 691 - fieldname - the textual name given to a registered field 692 693 Output parameter: 694 . vec - the vector 695 696 Level: beginner 697 698 .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmCreateGlobalVectorFromField() 699 @*/ 700 PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec) 701 { 702 PetscErrorCode ierr; 703 704 PetscFunctionBegin; 705 ierr = DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec);CHKERRQ(ierr); 706 PetscFunctionReturn(0); 707 } 708 709 /*@C 710 DMSwarmCreateLocalVectorFromField - Creates a Vec object sharing the array associated with a given field 711 712 Collective on dm 713 714 Input parameters: 715 + dm - a DMSwarm 716 - fieldname - the textual name given to a registered field 717 718 Output parameter: 719 . vec - the vector 720 721 Level: beginner 722 723 Notes: 724 The vector must be returned using a matching call to DMSwarmDestroyLocalVectorFromField(). 725 726 .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmDestroyLocalVectorFromField() 727 @*/ 728 PetscErrorCode DMSwarmCreateLocalVectorFromField(DM dm,const char fieldname[],Vec *vec) 729 { 730 MPI_Comm comm = PETSC_COMM_SELF; 731 PetscErrorCode ierr; 732 733 PetscFunctionBegin; 734 ierr = DMSwarmCreateVectorFromField_Private(dm, fieldname, comm, vec);CHKERRQ(ierr); 735 PetscFunctionReturn(0); 736 } 737 738 /*@C 739 DMSwarmDestroyLocalVectorFromField - Destroys the Vec object which share the array associated with a given field 740 741 Collective on dm 742 743 Input parameters: 744 + dm - a DMSwarm 745 - fieldname - the textual name given to a registered field 746 747 Output parameter: 748 . vec - the vector 749 750 Level: beginner 751 752 .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmCreateLocalVectorFromField() 753 @*/ 754 PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM dm,const char fieldname[],Vec *vec) 755 { 756 PetscErrorCode ierr; 757 758 PetscFunctionBegin; 759 ierr = DMSwarmDestroyVectorFromField_Private(dm, fieldname, vec);CHKERRQ(ierr); 760 PetscFunctionReturn(0); 761 } 762 763 /* 764 PetscErrorCode DMSwarmCreateGlobalVectorFromFields(DM dm,const PetscInt nf,const char *fieldnames[],Vec *vec) 765 { 766 PetscFunctionReturn(0); 767 } 768 769 PetscErrorCode DMSwarmRestoreGlobalVectorFromFields(DM dm,Vec *vec) 770 { 771 PetscFunctionReturn(0); 772 } 773 */ 774 775 /*@C 776 DMSwarmInitializeFieldRegister - Initiates the registration of fields to a DMSwarm 777 778 Collective on dm 779 780 Input parameter: 781 . dm - a DMSwarm 782 783 Level: beginner 784 785 Notes: 786 After all fields have been registered, you must call DMSwarmFinalizeFieldRegister(). 787 788 .seealso: DMSwarmFinalizeFieldRegister(), DMSwarmRegisterPetscDatatypeField(), 789 DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField() 790 @*/ 791 PetscErrorCode DMSwarmInitializeFieldRegister(DM dm) 792 { 793 DM_Swarm *swarm = (DM_Swarm *) dm->data; 794 PetscErrorCode ierr; 795 796 PetscFunctionBegin; 797 if (!swarm->field_registration_initialized) { 798 swarm->field_registration_initialized = PETSC_TRUE; 799 ierr = DMSwarmRegisterPetscDatatypeField(dm,DMSwarmField_pid,1,PETSC_INT64);CHKERRQ(ierr); /* unique identifer */ 800 ierr = DMSwarmRegisterPetscDatatypeField(dm,DMSwarmField_rank,1,PETSC_INT);CHKERRQ(ierr); /* used for communication */ 801 } 802 PetscFunctionReturn(0); 803 } 804 805 /*@ 806 DMSwarmFinalizeFieldRegister - Finalizes the registration of fields to a DMSwarm 807 808 Collective on dm 809 810 Input parameter: 811 . dm - a DMSwarm 812 813 Level: beginner 814 815 Notes: 816 After DMSwarmFinalizeFieldRegister() has been called, no new fields can be defined on the DMSwarm. 817 818 .seealso: DMSwarmInitializeFieldRegister(), DMSwarmRegisterPetscDatatypeField(), 819 DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField() 820 @*/ 821 PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm) 822 { 823 DM_Swarm *swarm = (DM_Swarm*)dm->data; 824 PetscErrorCode ierr; 825 826 PetscFunctionBegin; 827 if (!swarm->field_registration_finalized) { 828 ierr = DMSwarmDataBucketFinalize(swarm->db);CHKERRQ(ierr); 829 } 830 swarm->field_registration_finalized = PETSC_TRUE; 831 PetscFunctionReturn(0); 832 } 833 834 /*@ 835 DMSwarmSetLocalSizes - Sets the length of all registered fields on the DMSwarm 836 837 Not collective 838 839 Input parameters: 840 + dm - a DMSwarm 841 . nlocal - the length of each registered field 842 - buffer - the length of the buffer used to efficient dynamic re-sizing 843 844 Level: beginner 845 846 .seealso: DMSwarmGetLocalSize() 847 @*/ 848 PetscErrorCode DMSwarmSetLocalSizes(DM dm,PetscInt nlocal,PetscInt buffer) 849 { 850 DM_Swarm *swarm = (DM_Swarm*)dm->data; 851 PetscErrorCode ierr; 852 853 PetscFunctionBegin; 854 ierr = PetscLogEventBegin(DMSWARM_SetSizes,0,0,0,0);CHKERRQ(ierr); 855 ierr = DMSwarmDataBucketSetSizes(swarm->db,nlocal,buffer);CHKERRQ(ierr); 856 ierr = PetscLogEventEnd(DMSWARM_SetSizes,0,0,0,0);CHKERRQ(ierr); 857 PetscFunctionReturn(0); 858 } 859 860 /*@ 861 DMSwarmSetCellDM - Attachs a DM to a DMSwarm 862 863 Collective on dm 864 865 Input parameters: 866 + dm - a DMSwarm 867 - dmcell - the DM to attach to the DMSwarm 868 869 Level: beginner 870 871 Notes: 872 The attached DM (dmcell) will be queried for point location and 873 neighbor MPI-rank information if DMSwarmMigrate() is called. 874 875 .seealso: DMSwarmSetType(), DMSwarmGetCellDM(), DMSwarmMigrate() 876 @*/ 877 PetscErrorCode DMSwarmSetCellDM(DM dm,DM dmcell) 878 { 879 DM_Swarm *swarm = (DM_Swarm*)dm->data; 880 881 PetscFunctionBegin; 882 swarm->dmcell = dmcell; 883 PetscFunctionReturn(0); 884 } 885 886 /*@ 887 DMSwarmGetCellDM - Fetches the attached cell DM 888 889 Collective on dm 890 891 Input parameter: 892 . dm - a DMSwarm 893 894 Output parameter: 895 . dmcell - the DM which was attached to the DMSwarm 896 897 Level: beginner 898 899 .seealso: DMSwarmSetCellDM() 900 @*/ 901 PetscErrorCode DMSwarmGetCellDM(DM dm,DM *dmcell) 902 { 903 DM_Swarm *swarm = (DM_Swarm*)dm->data; 904 905 PetscFunctionBegin; 906 *dmcell = swarm->dmcell; 907 PetscFunctionReturn(0); 908 } 909 910 /*@ 911 DMSwarmGetLocalSize - Retrives the local length of fields registered 912 913 Not collective 914 915 Input parameter: 916 . dm - a DMSwarm 917 918 Output parameter: 919 . nlocal - the length of each registered field 920 921 Level: beginner 922 923 .seealso: DMSwarmGetSize(), DMSwarmSetLocalSizes() 924 @*/ 925 PetscErrorCode DMSwarmGetLocalSize(DM dm,PetscInt *nlocal) 926 { 927 DM_Swarm *swarm = (DM_Swarm*)dm->data; 928 PetscErrorCode ierr; 929 930 PetscFunctionBegin; 931 if (nlocal) {ierr = DMSwarmDataBucketGetSizes(swarm->db,nlocal,NULL,NULL);CHKERRQ(ierr);} 932 PetscFunctionReturn(0); 933 } 934 935 /*@ 936 DMSwarmGetSize - Retrives the total length of fields registered 937 938 Collective on dm 939 940 Input parameter: 941 . dm - a DMSwarm 942 943 Output parameter: 944 . n - the total length of each registered field 945 946 Level: beginner 947 948 Note: 949 This calls MPI_Allreduce upon each call (inefficient but safe) 950 951 .seealso: DMSwarmGetLocalSize(), DMSwarmSetLocalSizes() 952 @*/ 953 PetscErrorCode DMSwarmGetSize(DM dm,PetscInt *n) 954 { 955 DM_Swarm *swarm = (DM_Swarm*)dm->data; 956 PetscErrorCode ierr; 957 PetscInt nlocal,ng; 958 959 PetscFunctionBegin; 960 ierr = DMSwarmDataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);CHKERRQ(ierr); 961 ierr = MPI_Allreduce(&nlocal,&ng,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRMPI(ierr); 962 if (n) { *n = ng; } 963 PetscFunctionReturn(0); 964 } 965 966 /*@C 967 DMSwarmRegisterPetscDatatypeField - Register a field to a DMSwarm with a native PETSc data type 968 969 Collective on dm 970 971 Input parameters: 972 + dm - a DMSwarm 973 . fieldname - the textual name to identify this field 974 . blocksize - the number of each data type 975 - type - a valid PETSc data type (PETSC_CHAR, PETSC_SHORT, PETSC_INT, PETSC_FLOAT, PETSC_REAL, PETSC_LONG) 976 977 Level: beginner 978 979 Notes: 980 The textual name for each registered field must be unique. 981 982 .seealso: DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField() 983 @*/ 984 PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm,const char fieldname[],PetscInt blocksize,PetscDataType type) 985 { 986 PetscErrorCode ierr; 987 DM_Swarm *swarm = (DM_Swarm*)dm->data; 988 size_t size; 989 990 PetscFunctionBegin; 991 if (!swarm->field_registration_initialized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmInitializeFieldRegister() first"); 992 if (swarm->field_registration_finalized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first"); 993 994 if (type == PETSC_OBJECT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 995 if (type == PETSC_FUNCTION) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 996 if (type == PETSC_STRING) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 997 if (type == PETSC_STRUCT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 998 if (type == PETSC_DATATYPE_UNKNOWN) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}"); 999 1000 ierr = PetscDataTypeGetSize(type, &size);CHKERRQ(ierr); 1001 /* Load a specific data type into data bucket, specifying textual name and its size in bytes */ 1002 ierr = DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterPetscDatatypeField",fieldname,blocksize*size,NULL);CHKERRQ(ierr); 1003 { 1004 DMSwarmDataField gfield; 1005 1006 ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr); 1007 ierr = DMSwarmDataFieldSetBlockSize(gfield,blocksize);CHKERRQ(ierr); 1008 } 1009 swarm->db->field[swarm->db->nfields-1]->petsc_type = type; 1010 PetscFunctionReturn(0); 1011 } 1012 1013 /*@C 1014 DMSwarmRegisterUserStructField - Register a user defined struct to a DMSwarm 1015 1016 Collective on dm 1017 1018 Input parameters: 1019 + dm - a DMSwarm 1020 . fieldname - the textual name to identify this field 1021 - size - the size in bytes of the user struct of each data type 1022 1023 Level: beginner 1024 1025 Notes: 1026 The textual name for each registered field must be unique. 1027 1028 .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmRegisterUserDatatypeField() 1029 @*/ 1030 PetscErrorCode DMSwarmRegisterUserStructField(DM dm,const char fieldname[],size_t size) 1031 { 1032 PetscErrorCode ierr; 1033 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1034 1035 PetscFunctionBegin; 1036 ierr = DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterUserStructField",fieldname,size,NULL);CHKERRQ(ierr); 1037 swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_STRUCT ; 1038 PetscFunctionReturn(0); 1039 } 1040 1041 /*@C 1042 DMSwarmRegisterUserDatatypeField - Register a user defined data type to a DMSwarm 1043 1044 Collective on dm 1045 1046 Input parameters: 1047 + dm - a DMSwarm 1048 . fieldname - the textual name to identify this field 1049 . size - the size in bytes of the user data type 1050 - blocksize - the number of each data type 1051 1052 Level: beginner 1053 1054 Notes: 1055 The textual name for each registered field must be unique. 1056 1057 .seealso: DMSwarmRegisterPetscDatatypeField(), DMSwarmRegisterUserStructField(), DMSwarmRegisterUserDatatypeField() 1058 @*/ 1059 PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm,const char fieldname[],size_t size,PetscInt blocksize) 1060 { 1061 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1062 PetscErrorCode ierr; 1063 1064 PetscFunctionBegin; 1065 ierr = DMSwarmDataBucketRegisterField(swarm->db,"DMSwarmRegisterUserDatatypeField",fieldname,blocksize*size,NULL);CHKERRQ(ierr); 1066 { 1067 DMSwarmDataField gfield; 1068 1069 ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr); 1070 ierr = DMSwarmDataFieldSetBlockSize(gfield,blocksize);CHKERRQ(ierr); 1071 } 1072 swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_DATATYPE_UNKNOWN; 1073 PetscFunctionReturn(0); 1074 } 1075 1076 /*@C 1077 DMSwarmGetField - Get access to the underlying array storing all entries associated with a registered field 1078 1079 Not collective 1080 1081 Input parameters: 1082 + dm - a DMSwarm 1083 - fieldname - the textual name to identify this field 1084 1085 Output parameters: 1086 + blocksize - the number of each data type 1087 . type - the data type 1088 - data - pointer to raw array 1089 1090 Level: beginner 1091 1092 Notes: 1093 The array must be returned using a matching call to DMSwarmRestoreField(). 1094 1095 .seealso: DMSwarmRestoreField() 1096 @*/ 1097 PetscErrorCode DMSwarmGetField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data) 1098 { 1099 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1100 DMSwarmDataField gfield; 1101 PetscErrorCode ierr; 1102 1103 PetscFunctionBegin; 1104 if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); } 1105 ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr); 1106 ierr = DMSwarmDataFieldGetAccess(gfield);CHKERRQ(ierr); 1107 ierr = DMSwarmDataFieldGetEntries(gfield,data);CHKERRQ(ierr); 1108 if (blocksize) {*blocksize = gfield->bs; } 1109 if (type) { *type = gfield->petsc_type; } 1110 PetscFunctionReturn(0); 1111 } 1112 1113 /*@C 1114 DMSwarmRestoreField - Restore access to the underlying array storing all entries associated with a registered field 1115 1116 Not collective 1117 1118 Input parameters: 1119 + dm - a DMSwarm 1120 - fieldname - the textual name to identify this field 1121 1122 Output parameters: 1123 + blocksize - the number of each data type 1124 . type - the data type 1125 - data - pointer to raw array 1126 1127 Level: beginner 1128 1129 Notes: 1130 The user must call DMSwarmGetField() prior to calling DMSwarmRestoreField(). 1131 1132 .seealso: DMSwarmGetField() 1133 @*/ 1134 PetscErrorCode DMSwarmRestoreField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data) 1135 { 1136 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1137 DMSwarmDataField gfield; 1138 PetscErrorCode ierr; 1139 1140 PetscFunctionBegin; 1141 ierr = DMSwarmDataBucketGetDMSwarmDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr); 1142 ierr = DMSwarmDataFieldRestoreAccess(gfield);CHKERRQ(ierr); 1143 if (data) *data = NULL; 1144 PetscFunctionReturn(0); 1145 } 1146 1147 /*@ 1148 DMSwarmAddPoint - Add space for one new point in the DMSwarm 1149 1150 Not collective 1151 1152 Input parameter: 1153 . dm - a DMSwarm 1154 1155 Level: beginner 1156 1157 Notes: 1158 The new point will have all fields initialized to zero. 1159 1160 .seealso: DMSwarmAddNPoints() 1161 @*/ 1162 PetscErrorCode DMSwarmAddPoint(DM dm) 1163 { 1164 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1165 PetscErrorCode ierr; 1166 1167 PetscFunctionBegin; 1168 if (!swarm->issetup) {ierr = DMSetUp(dm);CHKERRQ(ierr);} 1169 ierr = PetscLogEventBegin(DMSWARM_AddPoints,0,0,0,0);CHKERRQ(ierr); 1170 ierr = DMSwarmDataBucketAddPoint(swarm->db);CHKERRQ(ierr); 1171 ierr = PetscLogEventEnd(DMSWARM_AddPoints,0,0,0,0);CHKERRQ(ierr); 1172 PetscFunctionReturn(0); 1173 } 1174 1175 /*@ 1176 DMSwarmAddNPoints - Add space for a number of new points in the DMSwarm 1177 1178 Not collective 1179 1180 Input parameters: 1181 + dm - a DMSwarm 1182 - npoints - the number of new points to add 1183 1184 Level: beginner 1185 1186 Notes: 1187 The new point will have all fields initialized to zero. 1188 1189 .seealso: DMSwarmAddPoint() 1190 @*/ 1191 PetscErrorCode DMSwarmAddNPoints(DM dm,PetscInt npoints) 1192 { 1193 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1194 PetscErrorCode ierr; 1195 PetscInt nlocal; 1196 1197 PetscFunctionBegin; 1198 ierr = PetscLogEventBegin(DMSWARM_AddPoints,0,0,0,0);CHKERRQ(ierr); 1199 ierr = DMSwarmDataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);CHKERRQ(ierr); 1200 nlocal = nlocal + npoints; 1201 ierr = DMSwarmDataBucketSetSizes(swarm->db,nlocal,DMSWARM_DATA_BUCKET_BUFFER_DEFAULT);CHKERRQ(ierr); 1202 ierr = PetscLogEventEnd(DMSWARM_AddPoints,0,0,0,0);CHKERRQ(ierr); 1203 PetscFunctionReturn(0); 1204 } 1205 1206 /*@ 1207 DMSwarmRemovePoint - Remove the last point from the DMSwarm 1208 1209 Not collective 1210 1211 Input parameter: 1212 . dm - a DMSwarm 1213 1214 Level: beginner 1215 1216 .seealso: DMSwarmRemovePointAtIndex() 1217 @*/ 1218 PetscErrorCode DMSwarmRemovePoint(DM dm) 1219 { 1220 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1221 PetscErrorCode ierr; 1222 1223 PetscFunctionBegin; 1224 ierr = PetscLogEventBegin(DMSWARM_RemovePoints,0,0,0,0);CHKERRQ(ierr); 1225 ierr = DMSwarmDataBucketRemovePoint(swarm->db);CHKERRQ(ierr); 1226 ierr = PetscLogEventEnd(DMSWARM_RemovePoints,0,0,0,0);CHKERRQ(ierr); 1227 PetscFunctionReturn(0); 1228 } 1229 1230 /*@ 1231 DMSwarmRemovePointAtIndex - Removes a specific point from the DMSwarm 1232 1233 Not collective 1234 1235 Input parameters: 1236 + dm - a DMSwarm 1237 - idx - index of point to remove 1238 1239 Level: beginner 1240 1241 .seealso: DMSwarmRemovePoint() 1242 @*/ 1243 PetscErrorCode DMSwarmRemovePointAtIndex(DM dm,PetscInt idx) 1244 { 1245 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1246 PetscErrorCode ierr; 1247 1248 PetscFunctionBegin; 1249 ierr = PetscLogEventBegin(DMSWARM_RemovePoints,0,0,0,0);CHKERRQ(ierr); 1250 ierr = DMSwarmDataBucketRemovePointAtIndex(swarm->db,idx);CHKERRQ(ierr); 1251 ierr = PetscLogEventEnd(DMSWARM_RemovePoints,0,0,0,0);CHKERRQ(ierr); 1252 PetscFunctionReturn(0); 1253 } 1254 1255 /*@ 1256 DMSwarmCopyPoint - Copy point pj to point pi in the DMSwarm 1257 1258 Not collective 1259 1260 Input parameters: 1261 + dm - a DMSwarm 1262 . pi - the index of the point to copy 1263 - pj - the point index where the copy should be located 1264 1265 Level: beginner 1266 1267 .seealso: DMSwarmRemovePoint() 1268 @*/ 1269 PetscErrorCode DMSwarmCopyPoint(DM dm,PetscInt pi,PetscInt pj) 1270 { 1271 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1272 PetscErrorCode ierr; 1273 1274 PetscFunctionBegin; 1275 if (!swarm->issetup) {ierr = DMSetUp(dm);CHKERRQ(ierr);} 1276 ierr = DMSwarmDataBucketCopyPoint(swarm->db,pi,swarm->db,pj);CHKERRQ(ierr); 1277 PetscFunctionReturn(0); 1278 } 1279 1280 PetscErrorCode DMSwarmMigrate_Basic(DM dm,PetscBool remove_sent_points) 1281 { 1282 PetscErrorCode ierr; 1283 1284 PetscFunctionBegin; 1285 ierr = DMSwarmMigrate_Push_Basic(dm,remove_sent_points);CHKERRQ(ierr); 1286 PetscFunctionReturn(0); 1287 } 1288 1289 /*@ 1290 DMSwarmMigrate - Relocates points defined in the DMSwarm to other MPI-ranks 1291 1292 Collective on dm 1293 1294 Input parameters: 1295 + dm - the DMSwarm 1296 - remove_sent_points - flag indicating if sent points should be removed from the current MPI-rank 1297 1298 Notes: 1299 The DM will be modified to accommodate received points. 1300 If remove_sent_points = PETSC_TRUE, any points that were sent will be removed from the DM. 1301 Different styles of migration are supported. See DMSwarmSetMigrateType(). 1302 1303 Level: advanced 1304 1305 .seealso: DMSwarmSetMigrateType() 1306 @*/ 1307 PetscErrorCode DMSwarmMigrate(DM dm,PetscBool remove_sent_points) 1308 { 1309 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1310 PetscErrorCode ierr; 1311 1312 PetscFunctionBegin; 1313 ierr = PetscLogEventBegin(DMSWARM_Migrate,0,0,0,0);CHKERRQ(ierr); 1314 switch (swarm->migrate_type) { 1315 case DMSWARM_MIGRATE_BASIC: 1316 ierr = DMSwarmMigrate_Basic(dm,remove_sent_points);CHKERRQ(ierr); 1317 break; 1318 case DMSWARM_MIGRATE_DMCELLNSCATTER: 1319 ierr = DMSwarmMigrate_CellDMScatter(dm,remove_sent_points);CHKERRQ(ierr); 1320 break; 1321 case DMSWARM_MIGRATE_DMCELLEXACT: 1322 SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE_DMCELLEXACT not implemented"); 1323 case DMSWARM_MIGRATE_USER: 1324 SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE_USER not implemented"); 1325 default: 1326 SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE type unknown"); 1327 } 1328 ierr = PetscLogEventEnd(DMSWARM_Migrate,0,0,0,0);CHKERRQ(ierr); 1329 ierr = DMClearGlobalVectors(dm);CHKERRQ(ierr); 1330 PetscFunctionReturn(0); 1331 } 1332 1333 PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize); 1334 1335 /* 1336 DMSwarmCollectViewCreate 1337 1338 * Applies a collection method and gathers point neighbour points into dm 1339 1340 Notes: 1341 Users should call DMSwarmCollectViewDestroy() after 1342 they have finished computations associated with the collected points 1343 */ 1344 1345 /*@ 1346 DMSwarmCollectViewCreate - Applies a collection method and gathers points 1347 in neighbour MPI-ranks into the DMSwarm 1348 1349 Collective on dm 1350 1351 Input parameter: 1352 . dm - the DMSwarm 1353 1354 Notes: 1355 Users should call DMSwarmCollectViewDestroy() after 1356 they have finished computations associated with the collected points 1357 Different collect methods are supported. See DMSwarmSetCollectType(). 1358 1359 Level: advanced 1360 1361 .seealso: DMSwarmCollectViewDestroy(), DMSwarmSetCollectType() 1362 @*/ 1363 PetscErrorCode DMSwarmCollectViewCreate(DM dm) 1364 { 1365 PetscErrorCode ierr; 1366 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1367 PetscInt ng; 1368 1369 PetscFunctionBegin; 1370 if (swarm->collect_view_active) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView currently active"); 1371 ierr = DMSwarmGetLocalSize(dm,&ng);CHKERRQ(ierr); 1372 switch (swarm->collect_type) { 1373 1374 case DMSWARM_COLLECT_BASIC: 1375 ierr = DMSwarmMigrate_GlobalToLocal_Basic(dm,&ng);CHKERRQ(ierr); 1376 break; 1377 case DMSWARM_COLLECT_DMDABOUNDINGBOX: 1378 SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT_DMDABOUNDINGBOX not implemented"); 1379 case DMSWARM_COLLECT_GENERAL: 1380 SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT_GENERAL not implemented"); 1381 default: 1382 SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT type unknown"); 1383 } 1384 swarm->collect_view_active = PETSC_TRUE; 1385 swarm->collect_view_reset_nlocal = ng; 1386 PetscFunctionReturn(0); 1387 } 1388 1389 /*@ 1390 DMSwarmCollectViewDestroy - Resets the DMSwarm to the size prior to calling DMSwarmCollectViewCreate() 1391 1392 Collective on dm 1393 1394 Input parameters: 1395 . dm - the DMSwarm 1396 1397 Notes: 1398 Users should call DMSwarmCollectViewCreate() before this function is called. 1399 1400 Level: advanced 1401 1402 .seealso: DMSwarmCollectViewCreate(), DMSwarmSetCollectType() 1403 @*/ 1404 PetscErrorCode DMSwarmCollectViewDestroy(DM dm) 1405 { 1406 PetscErrorCode ierr; 1407 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1408 1409 PetscFunctionBegin; 1410 if (!swarm->collect_view_active) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView is currently not active"); 1411 ierr = DMSwarmSetLocalSizes(dm,swarm->collect_view_reset_nlocal,-1);CHKERRQ(ierr); 1412 swarm->collect_view_active = PETSC_FALSE; 1413 PetscFunctionReturn(0); 1414 } 1415 1416 PetscErrorCode DMSwarmSetUpPIC(DM dm) 1417 { 1418 PetscInt dim; 1419 PetscErrorCode ierr; 1420 1421 PetscFunctionBegin; 1422 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 1423 if (dim < 1) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Dimension must be 1,2,3 - found %D",dim); 1424 if (dim > 3) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Dimension must be 1,2,3 - found %D",dim); 1425 ierr = DMSwarmRegisterPetscDatatypeField(dm,DMSwarmPICField_coor,dim,PETSC_DOUBLE);CHKERRQ(ierr); 1426 ierr = DMSwarmRegisterPetscDatatypeField(dm,DMSwarmPICField_cellid,1,PETSC_INT);CHKERRQ(ierr); 1427 PetscFunctionReturn(0); 1428 } 1429 1430 /*@ 1431 DMSwarmSetPointCoordinatesRandom - Sets initial coordinates for particles in each cell 1432 1433 Collective on dm 1434 1435 Input parameters: 1436 + dm - the DMSwarm 1437 - Npc - The number of particles per cell in the cell DM 1438 1439 Notes: 1440 The user must use DMSwarmSetCellDM() to set the cell DM first. The particles are placed randomly inside each cell. If only 1441 one particle is in each cell, it is placed at the centroid. 1442 1443 Level: intermediate 1444 1445 .seealso: DMSwarmSetCellDM() 1446 @*/ 1447 PetscErrorCode DMSwarmSetPointCoordinatesRandom(DM dm, PetscInt Npc) 1448 { 1449 DM cdm; 1450 PetscRandom rnd; 1451 DMPolytopeType ct; 1452 PetscBool simplex; 1453 PetscReal *centroid, *coords, *xi0, *v0, *J, *invJ, detJ; 1454 PetscInt dim, d, cStart, cEnd, c, p; 1455 PetscErrorCode ierr; 1456 1457 PetscFunctionBeginUser; 1458 ierr = PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rnd);CHKERRQ(ierr); 1459 ierr = PetscRandomSetInterval(rnd, -1.0, 1.0);CHKERRQ(ierr); 1460 ierr = PetscRandomSetType(rnd, PETSCRAND48);CHKERRQ(ierr); 1461 1462 ierr = DMSwarmGetCellDM(dm, &cdm);CHKERRQ(ierr); 1463 ierr = DMGetDimension(cdm, &dim);CHKERRQ(ierr); 1464 ierr = DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1465 ierr = DMPlexGetCellType(cdm, cStart, &ct);CHKERRQ(ierr); 1466 simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 1467 1468 ierr = PetscMalloc5(dim, ¢roid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ);CHKERRQ(ierr); 1469 for (d = 0; d < dim; ++d) xi0[d] = -1.0; 1470 ierr = DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 1471 for (c = cStart; c < cEnd; ++c) { 1472 if (Npc == 1) { 1473 ierr = DMPlexComputeCellGeometryFVM(cdm, c, NULL, centroid, NULL);CHKERRQ(ierr); 1474 for (d = 0; d < dim; ++d) coords[c*dim+d] = centroid[d]; 1475 } else { 1476 ierr = DMPlexComputeCellGeometryFEM(cdm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); /* affine */ 1477 for (p = 0; p < Npc; ++p) { 1478 const PetscInt n = c*Npc + p; 1479 PetscReal sum = 0.0, refcoords[3]; 1480 1481 for (d = 0; d < dim; ++d) { 1482 ierr = PetscRandomGetValueReal(rnd, &refcoords[d]);CHKERRQ(ierr); 1483 sum += refcoords[d]; 1484 } 1485 if (simplex && sum > 0.0) for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim)*sum; 1486 CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n*dim]); 1487 } 1488 } 1489 } 1490 ierr = DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr); 1491 ierr = PetscFree5(centroid, xi0, v0, J, invJ);CHKERRQ(ierr); 1492 ierr = PetscRandomDestroy(&rnd);CHKERRQ(ierr); 1493 PetscFunctionReturn(0); 1494 } 1495 1496 /*@ 1497 DMSwarmSetType - Set particular flavor of DMSwarm 1498 1499 Collective on dm 1500 1501 Input parameters: 1502 + dm - the DMSwarm 1503 - stype - the DMSwarm type (e.g. DMSWARM_PIC) 1504 1505 Level: advanced 1506 1507 .seealso: DMSwarmSetMigrateType(), DMSwarmSetCollectType() 1508 @*/ 1509 PetscErrorCode DMSwarmSetType(DM dm,DMSwarmType stype) 1510 { 1511 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1512 PetscErrorCode ierr; 1513 1514 PetscFunctionBegin; 1515 swarm->swarm_type = stype; 1516 if (swarm->swarm_type == DMSWARM_PIC) { 1517 ierr = DMSwarmSetUpPIC(dm);CHKERRQ(ierr); 1518 } 1519 PetscFunctionReturn(0); 1520 } 1521 1522 PetscErrorCode DMSetup_Swarm(DM dm) 1523 { 1524 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1525 PetscErrorCode ierr; 1526 PetscMPIInt rank; 1527 PetscInt p,npoints,*rankval; 1528 1529 PetscFunctionBegin; 1530 if (swarm->issetup) PetscFunctionReturn(0); 1531 1532 swarm->issetup = PETSC_TRUE; 1533 1534 if (swarm->swarm_type == DMSWARM_PIC) { 1535 /* check dmcell exists */ 1536 if (!swarm->dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires you call DMSwarmSetCellDM"); 1537 1538 if (swarm->dmcell->ops->locatepointssubdomain) { 1539 /* check methods exists for exact ownership identificiation */ 1540 ierr = PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->ops->LocatePointsSubdomain\n");CHKERRQ(ierr); 1541 swarm->migrate_type = DMSWARM_MIGRATE_DMCELLEXACT; 1542 } else { 1543 /* check methods exist for point location AND rank neighbor identification */ 1544 if (swarm->dmcell->ops->locatepoints) { 1545 ierr = PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->LocatePoints\n");CHKERRQ(ierr); 1546 } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires the method CellDM->ops->locatepoints be defined"); 1547 1548 if (swarm->dmcell->ops->getneighbors) { 1549 ierr = PetscInfo(dm, "DMSWARM_PIC: Using method CellDM->GetNeigbors\n");CHKERRQ(ierr); 1550 } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires the method CellDM->ops->getneighbors be defined"); 1551 1552 swarm->migrate_type = DMSWARM_MIGRATE_DMCELLNSCATTER; 1553 } 1554 } 1555 1556 ierr = DMSwarmFinalizeFieldRegister(dm);CHKERRQ(ierr); 1557 1558 /* check some fields were registered */ 1559 if (swarm->db->nfields <= 2) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"At least one field user must be registered via DMSwarmRegisterXXX()"); 1560 1561 /* check local sizes were set */ 1562 if (swarm->db->L == -1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Local sizes must be set via DMSwarmSetLocalSizes()"); 1563 1564 /* initialize values in pid and rank placeholders */ 1565 /* TODO: [pid - use MPI_Scan] */ 1566 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRMPI(ierr); 1567 ierr = DMSwarmDataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr); 1568 ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 1569 for (p=0; p<npoints; p++) { 1570 rankval[p] = (PetscInt)rank; 1571 } 1572 ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr); 1573 PetscFunctionReturn(0); 1574 } 1575 1576 extern PetscErrorCode DMSwarmSortDestroy(DMSwarmSort *_ctx); 1577 1578 PetscErrorCode DMDestroy_Swarm(DM dm) 1579 { 1580 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1581 PetscErrorCode ierr; 1582 1583 PetscFunctionBegin; 1584 if (--swarm->refct > 0) PetscFunctionReturn(0); 1585 ierr = DMSwarmDataBucketDestroy(&swarm->db);CHKERRQ(ierr); 1586 if (swarm->sort_context) { 1587 ierr = DMSwarmSortDestroy(&swarm->sort_context);CHKERRQ(ierr); 1588 } 1589 ierr = PetscFree(swarm);CHKERRQ(ierr); 1590 PetscFunctionReturn(0); 1591 } 1592 1593 PetscErrorCode DMSwarmView_Draw(DM dm, PetscViewer viewer) 1594 { 1595 DM cdm; 1596 PetscDraw draw; 1597 PetscReal *coords, oldPause, radius = 0.01; 1598 PetscInt Np, p, bs; 1599 PetscErrorCode ierr; 1600 1601 PetscFunctionBegin; 1602 ierr = PetscOptionsGetReal(NULL, ((PetscObject) dm)->prefix, "-dm_view_swarm_radius", &radius, NULL);CHKERRQ(ierr); 1603 ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr); 1604 ierr = DMSwarmGetCellDM(dm, &cdm);CHKERRQ(ierr); 1605 ierr = PetscDrawGetPause(draw, &oldPause);CHKERRQ(ierr); 1606 ierr = PetscDrawSetPause(draw, 0.0);CHKERRQ(ierr); 1607 ierr = DMView(cdm, viewer);CHKERRQ(ierr); 1608 ierr = PetscDrawSetPause(draw, oldPause);CHKERRQ(ierr); 1609 1610 ierr = DMSwarmGetLocalSize(dm, &Np);CHKERRQ(ierr); 1611 ierr = DMSwarmGetField(dm, DMSwarmPICField_coor, &bs, NULL, (void **) &coords);CHKERRQ(ierr); 1612 for (p = 0; p < Np; ++p) { 1613 const PetscInt i = p*bs; 1614 1615 ierr = PetscDrawEllipse(draw, coords[i], coords[i+1], radius, radius, PETSC_DRAW_BLUE);CHKERRQ(ierr); 1616 } 1617 ierr = DMSwarmRestoreField(dm, DMSwarmPICField_coor, &bs, NULL, (void **) &coords);CHKERRQ(ierr); 1618 ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 1619 ierr = PetscDrawPause(draw);CHKERRQ(ierr); 1620 ierr = PetscDrawSave(draw);CHKERRQ(ierr); 1621 PetscFunctionReturn(0); 1622 } 1623 1624 PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer) 1625 { 1626 DM_Swarm *swarm = (DM_Swarm*)dm->data; 1627 PetscBool iascii,ibinary,ishdf5,isvtk,isdraw; 1628 PetscErrorCode ierr; 1629 1630 PetscFunctionBegin; 1631 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1632 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 1633 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 1634 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY,&ibinary);CHKERRQ(ierr); 1635 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk);CHKERRQ(ierr); 1636 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5);CHKERRQ(ierr); 1637 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw);CHKERRQ(ierr); 1638 if (iascii) { 1639 ierr = DMSwarmDataBucketView(PetscObjectComm((PetscObject)dm),swarm->db,NULL,DATABUCKET_VIEW_STDOUT);CHKERRQ(ierr); 1640 } else if (ibinary) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO Binary support"); 1641 #if defined(PETSC_HAVE_HDF5) 1642 else if (ishdf5) { 1643 ierr = DMSwarmView_HDF5(dm, viewer);CHKERRQ(ierr); 1644 } 1645 #else 1646 else if (ishdf5) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"HDF5 not supported. Please reconfigure using --download-hdf5"); 1647 #endif 1648 else if (isvtk) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support"); 1649 else if (isdraw) { 1650 ierr = DMSwarmView_Draw(dm, viewer);CHKERRQ(ierr); 1651 } 1652 PetscFunctionReturn(0); 1653 } 1654 1655 /*@C 1656 DMSwarmGetCellSwarm - Extracts a single cell from the DMSwarm object, returns it as a single cell DMSwarm. 1657 The cell DM is filtered for fields of that cell, and the filtered DM is used as the cell DM of the new swarm object. 1658 1659 Important: Changes to this cell of the swarm will be lost if they are made prior to restoring this cell. 1660 1661 Noncollective 1662 1663 Input parameters: 1664 + sw - the DMSwarm 1665 - cellID - the integer id of the cell to be extracted and filtered 1666 1667 Output parameters: 1668 . cellswarm - The new DMSwarm 1669 1670 Level: beginner 1671 1672 Note: This presently only supports DMSWARM_PIC type 1673 1674 .seealso: DMSwarmRestoreCellSwarm() 1675 @*/ 1676 PETSC_EXTERN PetscErrorCode DMSwarmGetCellSwarm(DM sw, PetscInt cellID, DM cellswarm) 1677 { 1678 DM_Swarm *original = (DM_Swarm*) sw->data; 1679 DMLabel label; 1680 DM dmc, subdmc; 1681 PetscInt *pids, particles, dim; 1682 PetscErrorCode ierr; 1683 1684 PetscFunctionBegin; 1685 /* Configure new swarm */ 1686 ierr = DMSetType(cellswarm, DMSWARM);CHKERRQ(ierr); 1687 ierr = DMGetDimension(sw, &dim);CHKERRQ(ierr); 1688 ierr = DMSetDimension(cellswarm, dim);CHKERRQ(ierr); 1689 ierr = DMSwarmSetType(cellswarm, DMSWARM_PIC);CHKERRQ(ierr); 1690 /* Destroy the unused, unconfigured data bucket to prevent stragglers in memory */ 1691 ierr = DMSwarmDataBucketDestroy(&((DM_Swarm*)cellswarm->data)->db);CHKERRQ(ierr); 1692 ierr = DMSwarmSortGetAccess(sw);CHKERRQ(ierr); 1693 ierr = DMSwarmSortGetNumberOfPointsPerCell(sw, cellID, &particles);CHKERRQ(ierr); 1694 ierr = DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids);CHKERRQ(ierr); 1695 ierr = DMSwarmDataBucketCreateFromSubset(original->db, particles, pids, &((DM_Swarm*)cellswarm->data)->db);CHKERRQ(ierr); 1696 ierr = DMSwarmSortRestoreAccess(sw);CHKERRQ(ierr); 1697 ierr = PetscFree(pids);CHKERRQ(ierr); 1698 ierr = DMSwarmGetCellDM(sw, &dmc);CHKERRQ(ierr); 1699 ierr = DMLabelCreate(PetscObjectComm((PetscObject)sw), "singlecell", &label);CHKERRQ(ierr); 1700 ierr = DMAddLabel(dmc, label);CHKERRQ(ierr); 1701 ierr = DMLabelSetValue(label, cellID, 1);CHKERRQ(ierr); 1702 ierr = DMPlexFilter(dmc, label, 1, &subdmc);CHKERRQ(ierr); 1703 ierr = DMSwarmSetCellDM(cellswarm, subdmc);CHKERRQ(ierr); 1704 ierr = DMLabelDestroy(&label);CHKERRQ(ierr); 1705 PetscFunctionReturn(0); 1706 } 1707 1708 /*@C 1709 DMSwarmRestoreCellSwarm - Restores a DMSwarm object obtained with DMSwarmGetCellSwarm. All fields are copied back into the parent swarm. 1710 1711 Noncollective 1712 1713 Input parameters: 1714 + sw - the parent DMSwarm 1715 . cellID - the integer id of the cell to be copied back into the parent swarm 1716 - cellswarm - the cell swarm object 1717 1718 Level: beginner 1719 1720 Note: This only supports DMSWARM_PIC types of DMSwarms 1721 1722 .seealso: DMSwarmGetCellSwarm() 1723 @*/ 1724 PETSC_EXTERN PetscErrorCode DMSwarmRestoreCellSwarm(DM sw, PetscInt cellID, DM cellswarm) 1725 { 1726 DM dmc; 1727 PetscInt *pids, particles, p; 1728 PetscErrorCode ierr; 1729 1730 PetscFunctionBegin; 1731 ierr = DMSwarmSortGetAccess(sw);CHKERRQ(ierr); 1732 ierr = DMSwarmSortGetPointsPerCell(sw, cellID, &particles, &pids);CHKERRQ(ierr); 1733 ierr = DMSwarmSortRestoreAccess(sw);CHKERRQ(ierr); 1734 /* Pointwise copy of each particle based on pid. The parent swarm may not be altered during this process. */ 1735 for (p=0; p<particles; ++p) { 1736 ierr = DMSwarmDataBucketCopyPoint(((DM_Swarm*)cellswarm->data)->db,pids[p],((DM_Swarm*)sw->data)->db,pids[p]);CHKERRQ(ierr); 1737 } 1738 /* Free memory, destroy cell dm */ 1739 ierr = DMSwarmGetCellDM(cellswarm, &dmc);CHKERRQ(ierr); 1740 ierr = DMDestroy(&dmc);CHKERRQ(ierr); 1741 ierr = PetscFree(pids);CHKERRQ(ierr); 1742 PetscFunctionReturn(0); 1743 } 1744 1745 PETSC_INTERN PetscErrorCode DMClone_Swarm(DM, DM *); 1746 1747 static PetscErrorCode DMInitialize_Swarm(DM sw) 1748 { 1749 PetscFunctionBegin; 1750 sw->dim = 0; 1751 sw->ops->view = DMView_Swarm; 1752 sw->ops->load = NULL; 1753 sw->ops->setfromoptions = NULL; 1754 sw->ops->clone = DMClone_Swarm; 1755 sw->ops->setup = DMSetup_Swarm; 1756 sw->ops->createlocalsection = NULL; 1757 sw->ops->createdefaultconstraints = NULL; 1758 sw->ops->createglobalvector = DMCreateGlobalVector_Swarm; 1759 sw->ops->createlocalvector = DMCreateLocalVector_Swarm; 1760 sw->ops->getlocaltoglobalmapping = NULL; 1761 sw->ops->createfieldis = NULL; 1762 sw->ops->createcoordinatedm = NULL; 1763 sw->ops->getcoloring = NULL; 1764 sw->ops->creatematrix = DMCreateMatrix_Swarm; 1765 sw->ops->createinterpolation = NULL; 1766 sw->ops->createinjection = NULL; 1767 sw->ops->createmassmatrix = DMCreateMassMatrix_Swarm; 1768 sw->ops->refine = NULL; 1769 sw->ops->coarsen = NULL; 1770 sw->ops->refinehierarchy = NULL; 1771 sw->ops->coarsenhierarchy = NULL; 1772 sw->ops->globaltolocalbegin = NULL; 1773 sw->ops->globaltolocalend = NULL; 1774 sw->ops->localtoglobalbegin = NULL; 1775 sw->ops->localtoglobalend = NULL; 1776 sw->ops->destroy = DMDestroy_Swarm; 1777 sw->ops->createsubdm = NULL; 1778 sw->ops->getdimpoints = NULL; 1779 sw->ops->locatepoints = NULL; 1780 PetscFunctionReturn(0); 1781 } 1782 1783 PETSC_INTERN PetscErrorCode DMClone_Swarm(DM dm, DM *newdm) 1784 { 1785 DM_Swarm *swarm = (DM_Swarm *) dm->data; 1786 PetscErrorCode ierr; 1787 1788 PetscFunctionBegin; 1789 swarm->refct++; 1790 (*newdm)->data = swarm; 1791 ierr = PetscObjectChangeTypeName((PetscObject) *newdm, DMSWARM);CHKERRQ(ierr); 1792 ierr = DMInitialize_Swarm(*newdm);CHKERRQ(ierr); 1793 (*newdm)->dim = dm->dim; 1794 PetscFunctionReturn(0); 1795 } 1796 1797 /*MC 1798 1799 DMSWARM = "swarm" - A DM object used to represent arrays of data (fields) of arbitrary data type. 1800 This implementation was designed for particle methods in which the underlying 1801 data required to be represented is both (i) dynamic in length, (ii) and of arbitrary data type. 1802 1803 User data can be represented by DMSwarm through a registering "fields". 1804 To register a field, the user must provide; 1805 (a) a unique name; 1806 (b) the data type (or size in bytes); 1807 (c) the block size of the data. 1808 1809 For example, suppose the application requires a unique id, energy, momentum and density to be stored 1810 on a set of particles. Then the following code could be used 1811 1812 $ DMSwarmInitializeFieldRegister(dm) 1813 $ DMSwarmRegisterPetscDatatypeField(dm,"uid",1,PETSC_LONG); 1814 $ DMSwarmRegisterPetscDatatypeField(dm,"energy",1,PETSC_REAL); 1815 $ DMSwarmRegisterPetscDatatypeField(dm,"momentum",3,PETSC_REAL); 1816 $ DMSwarmRegisterPetscDatatypeField(dm,"density",1,PETSC_FLOAT); 1817 $ DMSwarmFinalizeFieldRegister(dm) 1818 1819 The fields represented by DMSwarm are dynamic and can be re-sized at any time. 1820 The only restriction imposed by DMSwarm is that all fields contain the same number of points. 1821 1822 To support particle methods, "migration" techniques are provided. These methods migrate data 1823 between MPI-ranks. 1824 1825 DMSwarm supports the methods DMCreateGlobalVector() and DMCreateLocalVector(). 1826 As a DMSwarm may internally define and store values of different data types, 1827 before calling DMCreateGlobalVector() or DMCreateLocalVector(), the user must inform DMSwarm which 1828 fields should be used to define a Vec object via 1829 DMSwarmVectorDefineField() 1830 The specified field can be changed at any time - thereby permitting vectors 1831 compatible with different fields to be created. 1832 1833 A dual representation of fields in the DMSwarm and a Vec object is permitted via 1834 DMSwarmCreateGlobalVectorFromField() 1835 Here the data defining the field in the DMSwarm is shared with a Vec. 1836 This is inherently unsafe if you alter the size of the field at any time between 1837 calls to DMSwarmCreateGlobalVectorFromField() and DMSwarmDestroyGlobalVectorFromField(). 1838 If the local size of the DMSwarm does not match the local size of the global vector 1839 when DMSwarmDestroyGlobalVectorFromField() is called, an error is thrown. 1840 1841 Additional high-level support is provided for Particle-In-Cell methods. 1842 Please refer to the man page for DMSwarmSetType(). 1843 1844 Level: beginner 1845 1846 .seealso: DMType, DMCreate(), DMSetType() 1847 M*/ 1848 PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm) 1849 { 1850 DM_Swarm *swarm; 1851 PetscErrorCode ierr; 1852 1853 PetscFunctionBegin; 1854 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1855 ierr = PetscNewLog(dm,&swarm);CHKERRQ(ierr); 1856 dm->data = swarm; 1857 ierr = DMSwarmDataBucketCreate(&swarm->db);CHKERRQ(ierr); 1858 ierr = DMSwarmInitializeFieldRegister(dm);CHKERRQ(ierr); 1859 swarm->refct = 1; 1860 swarm->vec_field_set = PETSC_FALSE; 1861 swarm->issetup = PETSC_FALSE; 1862 swarm->swarm_type = DMSWARM_BASIC; 1863 swarm->migrate_type = DMSWARM_MIGRATE_BASIC; 1864 swarm->collect_type = DMSWARM_COLLECT_BASIC; 1865 swarm->migrate_error_on_missing_point = PETSC_FALSE; 1866 swarm->dmcell = NULL; 1867 swarm->collect_view_active = PETSC_FALSE; 1868 swarm->collect_view_reset_nlocal = -1; 1869 ierr = DMInitialize_Swarm(dm);CHKERRQ(ierr); 1870 PetscFunctionReturn(0); 1871 } 1872