xref: /petsc/include/petscds.h (revision fbdf3bed89ca4550907b2a3c2c22a34a69e8be89)
1 /*
2       Objects which encapsulate discretizations+continuum residuals
3 */
4 #if !defined(__PETSCDS_H)
5 #define __PETSCDS_H
6 #include <petscfe.h>
7 #include <petscfv.h>
8 #include <petscdstypes.h>
9 
10 PETSC_EXTERN PetscErrorCode PetscDSInitializePackage(void);
11 
12 PETSC_EXTERN PetscClassId PETSCDS_CLASSID;
13 
14 /*J
15   PetscDSType - String with the name of a PETSc discrete system
16 
17   Level: beginner
18 
19 .seealso: PetscDSSetType(), PetscDS
20 J*/
21 typedef const char *PetscDSType;
22 #define PETSCDSBASIC "basic"
23 
24 typedef void (*PetscPointFunc)(PetscInt, PetscInt,
25                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
26                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
27                                PetscReal, const PetscReal[], PetscScalar[]);
28 typedef void (*PetscBdPointFunc)(PetscInt, PetscInt,
29                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
30                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
31                                  PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]);
32 typedef void (*PetscRiemannFunc)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *);
33 
34 
35 PETSC_EXTERN PetscFunctionList PetscDSList;
36 PETSC_EXTERN PetscErrorCode PetscDSCreate(MPI_Comm, PetscDS *);
37 PETSC_EXTERN PetscErrorCode PetscDSDestroy(PetscDS *);
38 PETSC_EXTERN PetscErrorCode PetscDSSetType(PetscDS, PetscDSType);
39 PETSC_EXTERN PetscErrorCode PetscDSGetType(PetscDS, PetscDSType *);
40 PETSC_EXTERN PetscErrorCode PetscDSSetUp(PetscDS);
41 PETSC_EXTERN PetscErrorCode PetscDSSetFromOptions(PetscDS);
42 PETSC_EXTERN PetscErrorCode PetscDSViewFromOptions(PetscDS,const char[],const char[]);
43 PETSC_EXTERN PetscErrorCode PetscDSView(PetscDS,PetscViewer);
44 PETSC_EXTERN PetscErrorCode PetscDSRegister(const char [], PetscErrorCode (*)(PetscDS));
45 PETSC_EXTERN PetscErrorCode PetscDSRegisterDestroy(void);
46 
47 PETSC_EXTERN PetscErrorCode PetscDSGetSpatialDimension(PetscDS, PetscInt *);
48 PETSC_EXTERN PetscErrorCode PetscDSGetNumFields(PetscDS, PetscInt *);
49 PETSC_EXTERN PetscErrorCode PetscDSGetTotalDimension(PetscDS, PetscInt *);
50 PETSC_EXTERN PetscErrorCode PetscDSGetTotalBdDimension(PetscDS, PetscInt *);
51 PETSC_EXTERN PetscErrorCode PetscDSGetTotalComponents(PetscDS, PetscInt *);
52 PETSC_EXTERN PetscErrorCode PetscDSGetFieldOffset(PetscDS, PetscInt, PetscInt *);
53 PETSC_EXTERN PetscErrorCode PetscDSGetBdFieldOffset(PetscDS, PetscInt, PetscInt *);
54 PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffset(PetscDS, PetscInt, PetscInt *);
55 PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffsets(PetscDS, PetscInt *[]);
56 PETSC_EXTERN PetscErrorCode PetscDSGetComponentBdOffsets(PetscDS, PetscInt *[]);
57 PETSC_EXTERN PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS, PetscInt *[]);
58 PETSC_EXTERN PetscErrorCode PetscDSGetComponentBdDerivativeOffsets(PetscDS, PetscInt *[]);
59 
60 PETSC_EXTERN PetscErrorCode PetscDSGetDiscretization(PetscDS, PetscInt, PetscObject *);
61 PETSC_EXTERN PetscErrorCode PetscDSSetDiscretization(PetscDS, PetscInt, PetscObject);
62 PETSC_EXTERN PetscErrorCode PetscDSAddDiscretization(PetscDS, PetscObject);
63 PETSC_EXTERN PetscErrorCode PetscDSGetBdDiscretization(PetscDS, PetscInt, PetscObject *);
64 PETSC_EXTERN PetscErrorCode PetscDSSetBdDiscretization(PetscDS, PetscInt, PetscObject);
65 PETSC_EXTERN PetscErrorCode PetscDSAddBdDiscretization(PetscDS, PetscObject);
66 PETSC_EXTERN PetscErrorCode PetscDSGetImplicit(PetscDS, PetscInt, PetscBool*);
67 PETSC_EXTERN PetscErrorCode PetscDSSetImplicit(PetscDS, PetscInt, PetscBool);
68 PETSC_EXTERN PetscErrorCode PetscDSGetAdjacency(PetscDS, PetscInt, PetscBool*, PetscBool*);
69 PETSC_EXTERN PetscErrorCode PetscDSSetAdjacency(PetscDS, PetscInt, PetscBool,  PetscBool);
70 PETSC_EXTERN PetscErrorCode PetscDSGetObjective(PetscDS, PetscInt,
71                                                 void (**)(PetscInt, PetscInt,
72                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
73                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
74                                                           const PetscReal, const PetscReal[], PetscScalar[]));
75 PETSC_EXTERN PetscErrorCode PetscDSSetObjective(PetscDS, PetscInt,
76                                                 void (*)(PetscInt, PetscInt,
77                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
78                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
79                                                          const PetscReal, const PetscReal[], PetscScalar[]));
80 PETSC_EXTERN PetscErrorCode PetscDSGetResidual(PetscDS, PetscInt,
81                                                void (**)(PetscInt, PetscInt,
82                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
83                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
84                                                          const PetscReal, const PetscReal[], PetscScalar[]),
85                                                void (**)(PetscInt, PetscInt,
86                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
87                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
88                                                          const PetscReal, const PetscReal[], PetscScalar[]));
89 PETSC_EXTERN PetscErrorCode PetscDSSetResidual(PetscDS, PetscInt,
90                                                void (*)(PetscInt, PetscInt,
91                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
92                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
93                                                         const PetscReal, const PetscReal[], PetscScalar[]),
94                                                void (*)(PetscInt, PetscInt,
95                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
96                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
97                                                         const PetscReal, const PetscReal[], PetscScalar[]));
98 PETSC_EXTERN PetscErrorCode PetscDSGetJacobian(PetscDS, PetscInt, PetscInt,
99                                                void (**)(PetscInt, PetscInt,
100                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
101                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
102                                                          const PetscReal, const PetscReal[], PetscScalar[]),
103                                                void (**)(PetscInt, PetscInt,
104                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
105                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
106                                                          const PetscReal, const PetscReal[], PetscScalar[]),
107                                                void (**)(PetscInt, PetscInt,
108                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
109                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
110                                                          const PetscReal, const PetscReal[], PetscScalar[]),
111                                                void (**)(PetscInt, PetscInt,
112                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
113                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
114                                                          const PetscReal, const PetscReal[], PetscScalar[]));
115 PETSC_EXTERN PetscErrorCode PetscDSSetJacobian(PetscDS, PetscInt, PetscInt,
116                                                void (*)(PetscInt, PetscInt,
117                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
118                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
119                                                         const PetscReal, const PetscReal[], PetscScalar[]),
120                                                void (*)(PetscInt, PetscInt,
121                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
122                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
123                                                         const PetscReal, const PetscReal[], PetscScalar[]),
124                                                void (*)(PetscInt, PetscInt,
125                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
126                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
127                                                         const PetscReal, const PetscReal[], PetscScalar[]),
128                                                void (*)(PetscInt, PetscInt,
129                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
130                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
131                                                         const PetscReal, const PetscReal[], PetscScalar[]));
132 PETSC_EXTERN PetscErrorCode PetscDSGetRiemannSolver(PetscDS, PetscInt,
133                                                     void (**)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *));
134 PETSC_EXTERN PetscErrorCode PetscDSSetRiemannSolver(PetscDS, PetscInt,
135                                                     void (*)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *));
136 PETSC_EXTERN PetscErrorCode PetscDSGetContext(PetscDS, PetscInt, void **);
137 PETSC_EXTERN PetscErrorCode PetscDSSetContext(PetscDS, PetscInt, void *);
138 PETSC_EXTERN PetscErrorCode PetscDSGetBdResidual(PetscDS, PetscInt,
139                                                  void (**)(PetscInt, PetscInt,
140                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
141                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
142                                                            const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
143                                                  void (**)(PetscInt, PetscInt,
144                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
145                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
146                                                            const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
147 PETSC_EXTERN PetscErrorCode PetscDSSetBdResidual(PetscDS, PetscInt,
148                                                  void (*)(PetscInt, PetscInt,
149                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
150                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
151                                                           const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
152                                                  void (*)(PetscInt, PetscInt,
153                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
154                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
155                                                           const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
156 PETSC_EXTERN PetscErrorCode PetscDSGetBdJacobian(PetscDS, PetscInt, PetscInt,
157                                                  void (**)(PetscInt, PetscInt,
158                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
159                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
160                                                            const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
161                                                  void (**)(PetscInt, PetscInt,
162                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
163                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
164                                                            const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
165                                                  void (**)(PetscInt, PetscInt,
166                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
167                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
168                                                            const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
169                                                  void (**)(PetscInt, PetscInt,
170                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
171                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
172                                                            const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
173 PETSC_EXTERN PetscErrorCode PetscDSSetBdJacobian(PetscDS, PetscInt, PetscInt,
174                                                  void (*)(PetscInt, PetscInt,
175                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
176                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
177                                                           const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
178                                                  void (*)(PetscInt, PetscInt,
179                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
180                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
181                                                           const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
182                                                  void (*)(PetscInt, PetscInt,
183                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
184                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
185                                                           const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
186                                                  void (*)(PetscInt, PetscInt,
187                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
188                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
189                                                           const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
190 PETSC_EXTERN PetscErrorCode PetscDSGetTabulation(PetscDS, PetscReal ***, PetscReal ***);
191 PETSC_EXTERN PetscErrorCode PetscDSGetBdTabulation(PetscDS, PetscReal ***, PetscReal ***);
192 PETSC_EXTERN PetscErrorCode PetscDSGetEvaluationArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **);
193 PETSC_EXTERN PetscErrorCode PetscDSGetWeakFormArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **);
194 PETSC_EXTERN PetscErrorCode PetscDSGetRefCoordArrays(PetscDS, PetscReal **, PetscScalar **);
195 
196 #endif
197