16524c165SJacob Faibussowitsch #ifndef PETSCDMSWARM_H 226bd1501SBarry Smith #define PETSCDMSWARM_H 33d5c7219SDave May 43d5c7219SDave May #include <petscdm.h> 535b38c60SMatthew G. Knepley #include <petscdt.h> 63d5c7219SDave May 7ac09b921SBarry Smith /* SUBMANSEC = DMSwarm */ 8ac09b921SBarry Smith 9e7af74c8SDave May /*E 1087497f52SBarry Smith DMSwarmType - Defines the type of `DMSWARM` 11e7af74c8SDave May 12e7af74c8SDave May DMSWARM_BASIC defines N entries of varied data-types which the user may register. 13e7af74c8SDave May 14ac09b921SBarry Smith DMSWARM_PIC is suitable for particle-in-cell methods. Configured as DMSWARM_PIC, the swarm will be aware of, another DM which serves as the background mesh. 15ac09b921SBarry Smith 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 16ac09b921SBarry Smith the owning rank. The background mesh will (by default) define the spatial decomposition of the points defined in the swarm. DMSWARM_PIC provides support 17ac09b921SBarry Smith for particle-in-cell operations such as defining initial point coordinates, communicating particles between sub-domains, projecting particle data fields on to the mesh. 18e7af74c8SDave May 19e7af74c8SDave May Level: beginner 20e7af74c8SDave May 21db781477SPatrick Sanan .seealso: `DMSwarmSetType()` 22e7af74c8SDave May E*/ 23480eef7bSDave May typedef enum { 24480eef7bSDave May DMSWARM_BASIC = 0, 25480eef7bSDave May DMSWARM_PIC 26480eef7bSDave May } DMSwarmType; 27480eef7bSDave May 28480eef7bSDave May typedef enum { 29480eef7bSDave May DMSWARM_MIGRATE_BASIC = 0, 30853ec3c6SDave May DMSWARM_MIGRATE_DMCELLNSCATTER, 31853ec3c6SDave May DMSWARM_MIGRATE_DMCELLEXACT, 32853ec3c6SDave May DMSWARM_MIGRATE_USER 33480eef7bSDave May } DMSwarmMigrateType; 34480eef7bSDave May 35480eef7bSDave May typedef enum { 36480eef7bSDave May DMSWARM_COLLECT_BASIC = 0, 37480eef7bSDave May DMSWARM_COLLECT_DMDABOUNDINGBOX, 38853ec3c6SDave May DMSWARM_COLLECT_GENERAL, 39853ec3c6SDave May DMSWARM_COLLECT_USER 40480eef7bSDave May } DMSwarmCollectType; 41480eef7bSDave May 42e7af74c8SDave May /*E 430e2ec84fSDave May DMSwarmPICLayoutType - Defines the method used to define particle coordinates within each cell. The layouts are constructured using the reference cell geometry 440e2ec84fSDave May 450e2ec84fSDave May DMSWARMPIC_LAYOUT_REGULAR defines points on a regular ijk mesh. 460e2ec84fSDave May When using DMSWARMPIC_LAYOUT_REGULAR, the fill_param defines the number of points in each spatial direction. 470e2ec84fSDave May 480e2ec84fSDave May DMSWARMPIC_LAYOUT_GAUSS defines points using an npoint Gauss-Legendre tensor product quadrature rule. 490e2ec84fSDave May When using DMSWARMPIC_LAYOUT_GAUSS, the fill_param defines the number of quadrature points in each spatial direction. 500e2ec84fSDave May 510e2ec84fSDave May DMSWARMPIC_LAYOUT_SUBDIVISION defines points on the centroid of a sub-divided reference cell. 520e2ec84fSDave May When using DMSWARMPIC_LAYOUT_SUBDIVISION, the fill_param defines the number times the reference cell is sub-divided. 530e2ec84fSDave May 540e2ec84fSDave May Level: beginner 550e2ec84fSDave May 56*a1cb98faSBarry Smith .seealso: `DMSWARM`, `DM`, `DMSwarmInsertPointsUsingCellDM()` 57e7af74c8SDave May E*/ 58e2d107dbSDave May typedef enum { 59e2d107dbSDave May DMSWARMPIC_LAYOUT_REGULAR = 0, 60e2d107dbSDave May DMSWARMPIC_LAYOUT_GAUSS, 61e0e61dcbSKarl Rupp DMSWARMPIC_LAYOUT_SUBDIVISION 62e2d107dbSDave May } DMSwarmPICLayoutType; 63e2d107dbSDave May 64853ec3c6SDave May PETSC_EXTERN const char *DMSwarmTypeNames[]; 65853ec3c6SDave May PETSC_EXTERN const char *DMSwarmMigrateTypeNames[]; 66853ec3c6SDave May PETSC_EXTERN const char *DMSwarmCollectTypeNames[]; 67853ec3c6SDave May 68853ec3c6SDave May PETSC_EXTERN const char DMSwarmField_pid[]; 69853ec3c6SDave May PETSC_EXTERN const char DMSwarmField_rank[]; 70853ec3c6SDave May PETSC_EXTERN const char DMSwarmPICField_coor[]; 71e2d107dbSDave May PETSC_EXTERN const char DMSwarmPICField_cellid[]; 72853ec3c6SDave May 73934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM, const char[], Vec *); 74934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM, const char[], Vec *); 75fb1bcc12SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmCreateLocalVectorFromField(DM, const char[], Vec *); 76fb1bcc12SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM, const char[], Vec *); 77853ec3c6SDave May 78934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmInitializeFieldRegister(DM); 79934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmFinalizeFieldRegister(DM); 80934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmSetLocalSizes(DM, PetscInt, PetscInt); 81934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM, const char[], PetscInt, PetscDataType); 82934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserStructField(DM, const char[], size_t); 83934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserDatatypeField(DM, const char[], size_t, PetscInt); 84934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmGetField(DM, const char[], PetscInt *, PetscDataType *, void **); 85934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmRestoreField(DM, const char[], PetscInt *, PetscDataType *, void **); 86934315b8SMatthew G. Knepley 87934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmVectorDefineField(DM, const char[]); 88934315b8SMatthew G. Knepley 89934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmAddPoint(DM); 90934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmAddNPoints(DM, PetscInt); 91934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmRemovePoint(DM); 92934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmRemovePointAtIndex(DM, PetscInt); 93ba4fc9c6SDave May PETSC_EXTERN PetscErrorCode DMSwarmCopyPoint(DM dm, PetscInt, PetscInt); 94934315b8SMatthew G. Knepley 95934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmGetLocalSize(DM, PetscInt *); 96934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmGetSize(DM, PetscInt *); 97934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmMigrate(DM, PetscBool); 98934315b8SMatthew G. Knepley 99934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmCollectViewCreate(DM); 100934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmCollectViewDestroy(DM); 101934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmSetCellDM(DM, DM); 102934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmGetCellDM(DM, DM *); 103934315b8SMatthew G. Knepley 104934315b8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmSetType(DM, DMSwarmType); 105934315b8SMatthew G. Knepley 1060e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmSetPointsUniformCoordinates(DM, PetscReal *, PetscReal *, PetscInt *, InsertMode); 10794f7d2dcSDave May PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinates(DM, PetscInt, PetscReal *, PetscBool, InsertMode); 1080e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmInsertPointsUsingCellDM(DM, DMSwarmPICLayoutType, PetscInt); 10992e40656SDave May PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinatesCellwise(DM, PetscInt, PetscReal *); 11031403186SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinatesRandom(DM, PetscInt); 1110e2ec84fSDave May PETSC_EXTERN PetscErrorCode DMSwarmViewFieldsXDMF(DM, const char *, PetscInt, const char **); 11294f7d2dcSDave May PETSC_EXTERN PetscErrorCode DMSwarmViewXDMF(DM, const char *); 1130e2ec84fSDave May 11494f7d2dcSDave May PETSC_EXTERN PetscErrorCode DMSwarmSortGetAccess(DM); 11594f7d2dcSDave May PETSC_EXTERN PetscErrorCode DMSwarmSortRestoreAccess(DM); 11694f7d2dcSDave May PETSC_EXTERN PetscErrorCode DMSwarmSortGetPointsPerCell(DM, PetscInt, PetscInt *, PetscInt **); 11794f7d2dcSDave May PETSC_EXTERN PetscErrorCode DMSwarmSortGetNumberOfPointsPerCell(DM, PetscInt, PetscInt *); 11894f7d2dcSDave May PETSC_EXTERN PetscErrorCode DMSwarmSortGetIsValid(DM, PetscBool *); 11994f7d2dcSDave May PETSC_EXTERN PetscErrorCode DMSwarmSortGetSizes(DM, PetscInt *, PetscInt *); 12094f7d2dcSDave May 12194f7d2dcSDave May PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM, PetscInt, const char **, Vec **, PetscBool); 1224cc7f7b2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmCreateMassMatrixSquare(DM, DM, Mat *); 123dc5f5ce9SDave May 12435b38c60SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmGetCellSwarm(DM, PetscInt, DM); 12535b38c60SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmRestoreCellSwarm(DM, PetscInt, DM); 12635b38c60SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmGetNumSpecies(DM, PetscInt *); 12735b38c60SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmSetNumSpecies(DM, PetscInt); 1286c5a79ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmGetCoordinateFunction(DM, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *)); 1296c5a79ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmSetCoordinateFunction(DM, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *)); 1306c5a79ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmGetVelocityFunction(DM, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *)); 1316c5a79ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmSetVelocityFunction(DM, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *)); 13235b38c60SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmComputeLocalSize(DM, PetscInt, PetscProbFunc); 13335b38c60SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmComputeLocalSizeFromOptions(DM); 13435b38c60SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmInitializeCoordinates(DM); 13535b38c60SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmInitializeVelocities(DM, PetscProbFunc, const PetscReal[]); 13635b38c60SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmInitializeVelocitiesFromOptions(DM, const PetscReal[]); 137d0c080abSJoseph Pusztay 1383d5c7219SDave May #endif 139