xref: /petsc/include/petscdmswarm.h (revision a4f09dd6775dcf4381066c075a763b56863c40b8)
1 #if !defined(__PETSCDMSWARM_H)
2 #define __PETSCDMSWARM_H
3 
4 #include <petscdm.h>
5 
6 typedef enum {
7   DMSWARM_BASIC=0,
8   DMSWARM_PIC
9 } DMSwarmType;
10 
11 typedef enum {
12   DMSWARM_MIGRATE_BASIC=0,
13   DMSWARM_MIGRATE_DMCELLNSCATTER,
14   DMSWARM_MIGRATE_DMCELLEXACT,
15   DMSWARM_MIGRATE_USER
16 } DMSwarmMigrateType;
17 
18 typedef enum {
19   DMSWARM_COLLECT_BASIC=0,
20   DMSWARM_COLLECT_DMDABOUNDINGBOX,
21   DMSWARM_COLLECT_GENERAL,
22   DMSWARM_COLLECT_USER
23 } DMSwarmCollectType;
24 
25 /*J
26    DMSwarmPICLayoutType - Defines the method used to define particle coordinates within each cell. The layouts are constructured using the reference cell geometry
27 
28    DMSWARMPIC_LAYOUT_REGULAR defines points on a regular ijk mesh.
29    When using DMSWARMPIC_LAYOUT_REGULAR, the fill_param defines the number of points in each spatial direction.
30 
31    DMSWARMPIC_LAYOUT_GAUSS defines points using an npoint Gauss-Legendre tensor product quadrature rule.
32    When using DMSWARMPIC_LAYOUT_GAUSS, the fill_param defines the number of quadrature points in each spatial direction.
33 
34    DMSWARMPIC_LAYOUT_SUBDIVISION defines points on the centroid of a sub-divided reference cell.
35    When using DMSWARMPIC_LAYOUT_SUBDIVISION, the fill_param defines the number times the reference cell is sub-divided.
36 
37    Level: beginner
38 
39 .seealso DMSwarmInsertPointsUsingCellDM()
40 J*/
41 typedef enum {
42   DMSWARMPIC_LAYOUT_REGULAR=0,
43   DMSWARMPIC_LAYOUT_GAUSS,
44   DMSWARMPIC_LAYOUT_SUBDIVISION,
45 } DMSwarmPICLayoutType;
46 
47 PETSC_EXTERN const char* DMSwarmTypeNames[];
48 PETSC_EXTERN const char* DMSwarmMigrateTypeNames[];
49 PETSC_EXTERN const char* DMSwarmCollectTypeNames[];
50 
51 PETSC_EXTERN const char DMSwarmField_pid[];
52 PETSC_EXTERN const char DMSwarmField_rank[];
53 PETSC_EXTERN const char DMSwarmPICField_coor[];
54 PETSC_EXTERN const char DMSwarmPICField_cellid[];
55 
56 PETSC_EXTERN PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM,const char[],Vec*);
57 PETSC_EXTERN PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM,const char[],Vec*);
58 PETSC_EXTERN PetscErrorCode DMSwarmCreateLocalVectorFromField(DM,const char[],Vec*);
59 PETSC_EXTERN PetscErrorCode DMSwarmDestroyLocalVectorFromField(DM,const char[],Vec*);
60 
61 PETSC_EXTERN PetscErrorCode DMSwarmInitializeFieldRegister(DM);
62 PETSC_EXTERN PetscErrorCode DMSwarmFinalizeFieldRegister(DM);
63 PETSC_EXTERN PetscErrorCode DMSwarmSetLocalSizes(DM,PetscInt,PetscInt);
64 PETSC_EXTERN PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM,const char[],PetscInt,PetscDataType);
65 PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserStructField(DM,const char[],size_t);
66 PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserDatatypeField(DM,const char[],size_t,PetscInt);
67 PETSC_EXTERN PetscErrorCode DMSwarmGetField(DM,const char[],PetscInt*,PetscDataType*,void**);
68 PETSC_EXTERN PetscErrorCode DMSwarmRestoreField(DM,const char[],PetscInt*,PetscDataType*,void**);
69 
70 PETSC_EXTERN PetscErrorCode DMSwarmVectorDefineField(DM,const char[]);
71 
72 PETSC_EXTERN PetscErrorCode DMSwarmAddPoint(DM);
73 PETSC_EXTERN PetscErrorCode DMSwarmAddNPoints(DM,PetscInt);
74 PETSC_EXTERN PetscErrorCode DMSwarmRemovePoint(DM);
75 PETSC_EXTERN PetscErrorCode DMSwarmRemovePointAtIndex(DM,PetscInt);
76 
77 PETSC_EXTERN PetscErrorCode DMSwarmGetLocalSize(DM,PetscInt*);
78 PETSC_EXTERN PetscErrorCode DMSwarmGetSize(DM,PetscInt*);
79 PETSC_EXTERN PetscErrorCode DMSwarmMigrate(DM,PetscBool);
80 
81 PETSC_EXTERN PetscErrorCode DMSwarmCollectViewCreate(DM);
82 PETSC_EXTERN PetscErrorCode DMSwarmCollectViewDestroy(DM);
83 PETSC_EXTERN PetscErrorCode DMSwarmSetCellDM(DM,DM);
84 PETSC_EXTERN PetscErrorCode DMSwarmGetCellDM(DM,DM*);
85 
86 PETSC_EXTERN PetscErrorCode DMSwarmSetType(DM,DMSwarmType);
87 
88 PETSC_EXTERN PetscErrorCode DMSwarmSetPointsUniformCoordinates(DM,PetscReal*,PetscReal*,PetscInt*,InsertMode);
89 PETSC_EXTERN PetscErrorCode DMSwarmSetPointCoordinates(DM dm,PetscInt npoints,PetscReal coor[],PetscBool redundant,InsertMode mode);
90 PETSC_EXTERN PetscErrorCode DMSwarmInsertPointsUsingCellDM(DM,DMSwarmPICLayoutType,PetscInt);
91 PETSC_EXTERN PetscErrorCode DMSwarmViewFieldsXDMF(DM,const char*,PetscInt,const char**);
92 PETSC_EXTERN PetscErrorCode DMSwarmViewXDMF(DM dm,const char filename[]);
93 
94 PETSC_EXTERN PetscErrorCode DMSwarmSortGetAccess(DM dm);
95 PETSC_EXTERN PetscErrorCode DMSwarmSortRestoreAccess(DM dm);
96 PETSC_EXTERN PetscErrorCode DMSwarmSortGetPointsPerCell(DM dm,PetscInt e,PetscInt *npoints,PetscInt **pidlist);
97 PETSC_EXTERN PetscErrorCode DMSwarmSortGetNumberOfPointsPerCell(DM dm,PetscInt cell_idx,PetscInt *npoints);
98 PETSC_EXTERN PetscErrorCode DMSwarmSortGetIsValid(DM dm,PetscBool *isvalid);
99 PETSC_EXTERN PetscErrorCode DMSwarmSortGetSizes(DM dm,PetscInt *ncells,PetscInt *npoints);
100 
101 #endif
102 
103