1 #pragma once 2 3 #include <petscdm.h> 4 #include <petscdt.h> 5 6 typedef struct _p_DMSwarmDataField *DMSwarmDataField; 7 typedef struct _p_DMSwarmDataBucket *DMSwarmDataBucket; 8 typedef struct _p_DMSwarmSort *DMSwarmSort; 9 10 /* SUBMANSEC = DMSwarm */ 11 12 /*E 13 DMSwarmType - Defines the type of `DMSWARM` 14 15 Values: 16 + `DMSWARM_BASIC` - defines N entries of varied data-types which the user may register. 17 - `DMSWARM_PIC` - suitable for particle-in-cell methods. Configured as `DMSWARM_PIC`, the swarm will be aware of, another `DM` which serves as the background mesh. 18 19 Fields specific to particle-in-cell methods are registered by default. These include spatial coordinates, a unique identifier, a cell index and an index for 20 the owning rank. The background mesh will (by default) define the spatial decomposition of the points defined in the swarm. `DMSWARM_PIC` provides support 21 for particle-in-cell operations such as defining initial point coordinates, communicating particles between sub-domains, projecting particle data fields on to the mesh. 22 23 Level: beginner 24 25 .seealso: [](ch_dmbase), `DMSWARM`, `DMSwarmSetType()` 26 E*/ 27 typedef enum { 28 DMSWARM_BASIC = 0, 29 DMSWARM_PIC 30 } DMSwarmType; 31 32 typedef enum { 33 DMSWARM_MIGRATE_BASIC = 0, 34 DMSWARM_MIGRATE_DMCELLNSCATTER, 35 DMSWARM_MIGRATE_DMCELLEXACT, 36 DMSWARM_MIGRATE_USER 37 } DMSwarmMigrateType; 38 39 typedef enum { 40 DMSWARM_COLLECT_BASIC = 0, 41 DMSWARM_COLLECT_DMDABOUNDINGBOX, 42 DMSWARM_COLLECT_GENERAL, 43 DMSWARM_COLLECT_USER 44 } DMSwarmCollectType; 45 46 /*E 47 DMSwarmPICLayoutType - Defines the method used to define particle coordinates within each cell. The layouts are constructured using the reference cell geometry 48 49 Values: 50 + `DMSWARMPIC_LAYOUT_REGULAR` - defines points on a regular ijk mesh. In this case 51 the `fill_param` argument of `DMSwarmInsertPointsUsingCellDM()` defines the number of points in each spatial direction. 52 . `DMSWARMPIC_LAYOUT_GAUSS` - defines points using an npoint Gauss-Legendre tensor product quadrature rule. In this case 53 the `fill_param` argument of `DMSwarmInsertPointsUsingCellDM()` defines the number of quadrature points in each spatial direction. 54 - `DMSWARMPIC_LAYOUT_SUBDIVISION` - defines points on the centroid of a sub-divided reference cell. In this case 55 the `fill_param` argument of `DMSwarmInsertPointsUsingCellDM()` defines the number times the reference cell is sub-divided. 56 57 Level: beginner 58 59 .seealso: [](ch_dmbase), `DMSWARM`, `DM`, `DMSwarmInsertPointsUsingCellDM()` 60 E*/ 61 typedef enum { 62 DMSWARMPIC_LAYOUT_REGULAR = 0, 63 DMSWARMPIC_LAYOUT_GAUSS, 64 DMSWARMPIC_LAYOUT_SUBDIVISION 65 } DMSwarmPICLayoutType; 66 67 PETSC_EXTERN const char *DMSwarmTypeNames[]; 68 PETSC_EXTERN const char *DMSwarmMigrateTypeNames[]; 69 PETSC_EXTERN const char *DMSwarmCollectTypeNames[]; 70 71 PETSC_EXTERN const char DMSwarmField_pid[]; 72 PETSC_EXTERN const char DMSwarmField_rank[]; 73 PETSC_EXTERN const char DMSwarmPICField_coor[]; 74 PETSC_EXTERN const char DMSwarmPICField_cellid[]; 75 76 PETSC_EXTERN PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM, const char[], Vec *); 77 PETSC_EXTERN PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM, const char[], Vec *); 78 PETSC_EXTERN PetscErrorCode DMSwarmCreateLocalVectorFromField(DM, const char[], Vec *); 79 PETSC_EXTERN PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM, const char[], Vec *); 80 81 PETSC_EXTERN PetscErrorCode DMSwarmInitializeFieldRegister(DM); 82 PETSC_EXTERN PetscErrorCode DMSwarmFinalizeFieldRegister(DM); 83 PETSC_EXTERN PetscErrorCode DMSwarmSetLocalSizes(DM, PetscInt, PetscInt); 84 PETSC_EXTERN PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM, const char[], PetscInt, PetscDataType); 85 PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserStructField(DM, const char[], size_t); 86 PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserDatatypeField(DM, const char[], size_t, PetscInt); 87 PETSC_EXTERN PetscErrorCode DMSwarmGetField(DM, const char[], PetscInt *, PetscDataType *, void **); 88 PETSC_EXTERN PetscErrorCode DMSwarmRestoreField(DM, const char[], PetscInt *, PetscDataType *, void **); 89 90 PETSC_EXTERN PetscErrorCode DMSwarmVectorDefineField(DM, const char[]); 91 92 PETSC_EXTERN PetscErrorCode DMSwarmAddPoint(DM); 93 PETSC_EXTERN PetscErrorCode DMSwarmAddNPoints(DM, PetscInt); 94 PETSC_EXTERN PetscErrorCode DMSwarmRemovePoint(DM); 95 PETSC_EXTERN PetscErrorCode DMSwarmRemovePointAtIndex(DM, PetscInt); 96 PETSC_EXTERN PetscErrorCode DMSwarmCopyPoint(DM dm, PetscInt, PetscInt); 97 98 PETSC_EXTERN PetscErrorCode DMSwarmGetLocalSize(DM, PetscInt *); 99 PETSC_EXTERN PetscErrorCode DMSwarmGetSize(DM, PetscInt *); 100 PETSC_EXTERN PetscErrorCode DMSwarmGetMigrateType(DM, DMSwarmMigrateType *); 101 PETSC_EXTERN PetscErrorCode DMSwarmSetMigrateType(DM, DMSwarmMigrateType); 102 PETSC_EXTERN PetscErrorCode DMSwarmMigrate(DM, PetscBool); 103 104 PETSC_EXTERN PetscErrorCode DMSwarmCollectViewCreate(DM); 105 PETSC_EXTERN PetscErrorCode DMSwarmCollectViewDestroy(DM); 106 PETSC_EXTERN PetscErrorCode DMSwarmSetCellDM(DM, DM); 107 PETSC_EXTERN PetscErrorCode DMSwarmGetCellDM(DM, DM *); 108 109 PETSC_EXTERN PetscErrorCode DMSwarmSetType(DM, DMSwarmType); 110 111 PETSC_EXTERN PetscErrorCode DMSwarmSetPointsUniformCoordinates(DM, PetscReal *, PetscReal *, PetscInt *, InsertMode); 112 PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinates(DM, PetscInt, PetscReal *, PetscBool, InsertMode); 113 PETSC_EXTERN PetscErrorCode DMSwarmInsertPointsUsingCellDM(DM, DMSwarmPICLayoutType, PetscInt); 114 PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinatesCellwise(DM, PetscInt, PetscReal *); 115 PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinatesRandom(DM, PetscInt); 116 PETSC_EXTERN PetscErrorCode DMSwarmViewFieldsXDMF(DM, const char *, PetscInt, const char **); 117 PETSC_EXTERN PetscErrorCode DMSwarmViewXDMF(DM, const char *); 118 119 PETSC_EXTERN PetscErrorCode DMSwarmSortGetAccess(DM); 120 PETSC_EXTERN PetscErrorCode DMSwarmSortRestoreAccess(DM); 121 PETSC_EXTERN PetscErrorCode DMSwarmSortGetPointsPerCell(DM, PetscInt, PetscInt *, PetscInt **); 122 PETSC_EXTERN PetscErrorCode DMSwarmSortGetNumberOfPointsPerCell(DM, PetscInt, PetscInt *); 123 PETSC_EXTERN PetscErrorCode DMSwarmSortGetIsValid(DM, PetscBool *); 124 PETSC_EXTERN PetscErrorCode DMSwarmSortGetSizes(DM, PetscInt *, PetscInt *); 125 126 PETSC_EXTERN PetscErrorCode DMSwarmCreateMassMatrixSquare(DM, DM, Mat *); 127 128 PETSC_EXTERN PetscErrorCode DMSwarmGetCellSwarm(DM, PetscInt, DM); 129 PETSC_EXTERN PetscErrorCode DMSwarmRestoreCellSwarm(DM, PetscInt, DM); 130 PETSC_EXTERN PetscErrorCode DMSwarmGetNumSpecies(DM, PetscInt *); 131 PETSC_EXTERN PetscErrorCode DMSwarmSetNumSpecies(DM, PetscInt); 132 PETSC_EXTERN PetscErrorCode DMSwarmGetCoordinateFunction(DM, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *)); 133 PETSC_EXTERN PetscErrorCode DMSwarmSetCoordinateFunction(DM, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *)); 134 PETSC_EXTERN PetscErrorCode DMSwarmGetVelocityFunction(DM, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *)); 135 PETSC_EXTERN PetscErrorCode DMSwarmSetVelocityFunction(DM, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *)); 136 PETSC_EXTERN PetscErrorCode DMSwarmComputeLocalSize(DM, PetscInt, PetscProbFunc); 137 PETSC_EXTERN PetscErrorCode DMSwarmComputeLocalSizeFromOptions(DM); 138 PETSC_EXTERN PetscErrorCode DMSwarmInitializeCoordinates(DM); 139 PETSC_EXTERN PetscErrorCode DMSwarmInitializeVelocities(DM, PetscProbFunc, const PetscReal[]); 140 PETSC_EXTERN PetscErrorCode DMSwarmInitializeVelocitiesFromOptions(DM, const PetscReal[]); 141 142 // Interface to internal storage 143 PETSC_EXTERN PetscErrorCode DMSwarmDataFieldGetEntries(const DMSwarmDataField, void **); 144 PETSC_EXTERN PetscErrorCode DMSwarmDataFieldRestoreEntries(const DMSwarmDataField, void **); 145 PETSC_EXTERN PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFieldByName(DMSwarmDataBucket, const char[], DMSwarmDataField *); 146 PETSC_EXTERN PetscErrorCode DMSwarmDataBucketGetDMSwarmDataFieldIdByName(DMSwarmDataBucket, const char[], PetscInt *); 147 PETSC_EXTERN PetscErrorCode DMSwarmDataBucketQueryDMSwarmDataFieldByName(DMSwarmDataBucket, const char[], PetscBool *); 148