xref: /petsc/src/binding/petsc4py/src/petsc4py/PETSc/petscdmswarm.pxi (revision 9764e2963186706b1a86270e84a7b36d242f0ed1)
1# --------------------------------------------------------------------
2
3cdef extern from * nogil:
4
5    ctypedef enum PetscDMSwarmType "DMSwarmType":
6        DMSWARM_BASIC
7        DMSWARM_PIC
8
9    ctypedef enum PetscDMSwarmMigrateType "DMSwarmMigrateType":
10        DMSWARM_MIGRATE_BASIC
11        DMSWARM_MIGRATE_DMCELLNSCATTER
12        DMSWARM_MIGRATE_DMCELLEXACT
13        DMSWARM_MIGRATE_USER
14
15    ctypedef enum PetscDMSwarmCollectType "DMSwarmCollectType":
16        DMSWARM_COLLECT_BASIC
17        DMSWARM_COLLECT_DMDABOUNDINGBOX
18        DMSWARM_COLLECT_GENERAL
19        DMSWARM_COLLECT_USER
20
21    ctypedef enum PetscDMSwarmPICLayoutType 'DMSwarmPICLayoutType':
22        DMSWARMPIC_LAYOUT_REGULAR
23        DMSWARMPIC_LAYOUT_GAUSS
24        DMSWARMPIC_LAYOUT_SUBDIVISION
25
26    PetscErrorCode DMSwarmCreateGlobalVectorFromField(PetscDM, const char[], PetscVec*)
27    PetscErrorCode DMSwarmDestroyGlobalVectorFromField(PetscDM, const char[], PetscVec*)
28    PetscErrorCode DMSwarmCreateLocalVectorFromField(PetscDM, const char[], PetscVec*)
29    PetscErrorCode DMSwarmDestroyLocalVectorFromField(PetscDM, const char[], PetscVec*)
30    PetscErrorCode DMSwarmCreateGlobalVectorFromFields(PetscDM, PetscInt, const char *[], PetscVec *)
31    PetscErrorCode DMSwarmDestroyGlobalVectorFromFields(PetscDM, PetscInt, const char *[], PetscVec *)
32    PetscErrorCode DMSwarmCreateLocalVectorFromFields(PetscDM, PetscInt, const char *[], PetscVec *)
33    PetscErrorCode DMSwarmDestroyLocalVectorFromFields(PetscDM, PetscInt, const char *[], PetscVec *)
34
35    PetscErrorCode DMSwarmInitializeFieldRegister(PetscDM)
36    PetscErrorCode DMSwarmFinalizeFieldRegister(PetscDM)
37    PetscErrorCode DMSwarmSetLocalSizes(PetscDM, PetscInt, PetscInt)
38    PetscErrorCode DMSwarmRegisterPetscDatatypeField(PetscDM, const char[], PetscInt, PetscDataType)
39#    PetscErrorCode DMSwarmRegisterUserStructField(PetscDM, const char[], size_t)
40#    PetscErrorCode DMSwarmRegisterUserDatatypeField(PetscDM, const char[], size_t, PetscInt)
41    PetscErrorCode DMSwarmGetField(PetscDM, const char[], PetscInt*, PetscDataType*, void**)
42    PetscErrorCode DMSwarmRestoreField(PetscDM, const char[], PetscInt*, PetscDataType*, void**)
43    PetscErrorCode DMSwarmGetFieldInfo(PetscDM, const char[], PetscInt *, PetscDataType *)
44
45    PetscErrorCode DMSwarmVectorDefineField(PetscDM, const char[])
46
47    PetscErrorCode DMSwarmAddPoint(PetscDM)
48    PetscErrorCode DMSwarmAddNPoints(PetscDM, PetscInt)
49    PetscErrorCode DMSwarmRemovePoint(PetscDM)
50    PetscErrorCode DMSwarmRemovePointAtIndex(PetscDM, PetscInt)
51    PetscErrorCode DMSwarmCopyPoint(PetscDM, PetscInt, PetscInt)
52
53    PetscErrorCode DMSwarmGetLocalSize(PetscDM, PetscInt*)
54    PetscErrorCode DMSwarmGetSize(PetscDM, PetscInt*)
55    PetscErrorCode DMSwarmMigrate(PetscDM, PetscBool)
56
57    PetscErrorCode DMSwarmCollectViewCreate(PetscDM)
58    PetscErrorCode DMSwarmCollectViewDestroy(PetscDM)
59    PetscErrorCode DMSwarmSetCellDM(PetscDM, PetscDM)
60    PetscErrorCode DMSwarmGetCellDM(PetscDM, PetscDM*)
61    PetscErrorCode DMSwarmGetCellDMByName(PetscDM, const char[], PetscDMSwarmCellDM *)
62    PetscErrorCode DMSwarmGetCellDMNames(PetscDM, PetscInt *, const char **[])
63    PetscErrorCode DMSwarmSetCellDMActive(PetscDM, const char[])
64    PetscErrorCode DMSwarmGetCellDMActive(PetscDM, PetscDMSwarmCellDM *)
65    PetscErrorCode DMSwarmAddCellDM(PetscDM, PetscDMSwarmCellDM)
66
67    PetscErrorCode DMSwarmSetType(PetscDM, PetscDMSwarmType)
68
69    PetscErrorCode DMSwarmSetPointsUniformCoordinates(PetscDM, PetscReal[], PetscReal[], PetscInt[], PetscInsertMode)
70    PetscErrorCode DMSwarmSetPointCoordinates(PetscDM, PetscInt, PetscReal*, PetscBool, PetscInsertMode)
71    PetscErrorCode DMSwarmInsertPointsUsingCellDM(PetscDM, PetscDMSwarmPICLayoutType, PetscInt)
72    PetscErrorCode DMSwarmSetPointCoordinatesCellwise(PetscDM, PetscInt, PetscReal*)
73    PetscErrorCode DMSwarmViewFieldsXDMF(PetscDM, const char*, PetscInt, const char**)
74    PetscErrorCode DMSwarmViewXDMF(PetscDM, const char*)
75
76    PetscErrorCode DMSwarmSortGetAccess(PetscDM)
77    PetscErrorCode DMSwarmSortRestoreAccess(PetscDM)
78    PetscErrorCode DMSwarmSortGetPointsPerCell(PetscDM, PetscInt, PetscInt*, PetscInt**)
79    PetscErrorCode DMSwarmSortGetNumberOfPointsPerCell(PetscDM, PetscInt, PetscInt*)
80    PetscErrorCode DMSwarmSortGetIsValid(PetscDM, PetscBool*)
81    PetscErrorCode DMSwarmSortGetSizes(PetscDM, PetscInt*, PetscInt*)
82
83    PetscErrorCode DMSwarmProjectFields(PetscDM, PetscDM, PetscInt, const char**, PetscVec*, PetscScatterMode)
84
85    PetscErrorCode DMSwarmComputeMoments(PetscDM, const char[], const char[], PetscReal[])
86
87    PetscErrorCode DMSwarmCellDMCreate(PetscDM, PetscInt, const char *[], PetscInt, const char *[], PetscDMSwarmCellDM *)
88    PetscErrorCode DMSwarmCellDMDestroy(PetscDMSwarmCellDM *)
89    PetscErrorCode DMSwarmCellDMView(PetscDMSwarmCellDM, PetscViewer)
90    PetscErrorCode DMSwarmCellDMGetDM(PetscDMSwarmCellDM, PetscDM *)
91    PetscErrorCode DMSwarmCellDMGetFields(PetscDMSwarmCellDM, PetscInt *, const char **[])
92    PetscErrorCode DMSwarmCellDMGetCoordinateFields(PetscDMSwarmCellDM, PetscInt *, const char **[])
93    PetscErrorCode DMSwarmCellDMGetCellID(PetscDMSwarmCellDM, const char *[])
94    PetscErrorCode DMSwarmCellDMGetBlockSize(PetscDMSwarmCellDM, PetscDM, PetscInt *)
95