xref: /petsc/include/petscdstypes.h (revision 2192575e8c1428e47f7eac7b348fe8e57717183e)
1a4963045SJacob Faibussowitsch #pragma once
22764a2aaSMatthew G. Knepley 
36528b96dSMatthew G. Knepley #include <petscdmlabel.h>
46528b96dSMatthew G. Knepley 
5ce78bad3SBarry Smith /* MANSEC = DM */
6ac09b921SBarry Smith /* SUBMANSEC = DT */
7ac09b921SBarry Smith 
82764a2aaSMatthew G. Knepley /*S
987497f52SBarry Smith   PetscDS - PETSc object that manages a discrete system, which is a set of discretizations + continuum equations from a `PetscWeakForm`
102764a2aaSMatthew G. Knepley 
112764a2aaSMatthew G. Knepley   Level: intermediate
122764a2aaSMatthew G. Knepley 
13db781477SPatrick Sanan .seealso: `PetscDSCreate()`, `PetscDSSetType()`, `PetscDSType`, `PetscWeakForm`, `PetscFECreate()`, `PetscFVCreate()`
142764a2aaSMatthew G. Knepley S*/
152764a2aaSMatthew G. Knepley typedef struct _p_PetscDS *PetscDS;
162764a2aaSMatthew G. Knepley 
176528b96dSMatthew G. Knepley /*S
186528b96dSMatthew G. Knepley   PetscWeakForm - PETSc object that manages a sets of pointwise functions defining a system of equations
196528b96dSMatthew G. Knepley 
206528b96dSMatthew G. Knepley   Level: intermediate
216528b96dSMatthew G. Knepley 
22db781477SPatrick Sanan .seealso: `PetscWeakFormCreate()`, `PetscDS`, `PetscFECreate()`, `PetscFVCreate()`
236528b96dSMatthew G. Knepley S*/
246528b96dSMatthew G. Knepley typedef struct _p_PetscWeakForm *PetscWeakForm;
256528b96dSMatthew G. Knepley 
2606ad1575SMatthew G. Knepley /*S
2706ad1575SMatthew G. Knepley   PetscFormKey - This key indicates how to use a set of pointwise functions defining part of a system of equations
2806ad1575SMatthew G. Knepley 
29a5b23f4aSJose E. Roman   The subdomain on which to integrate is specified by (label, value), the test function field by (field), and the
3006ad1575SMatthew G. Knepley   piece of the equation by (part). For example, LHS = 0 and RHS = 1 in IMEX methods. More pieces can be present for
3106ad1575SMatthew G. Knepley   operator splitting methods.
3206ad1575SMatthew G. Knepley 
3306ad1575SMatthew G. Knepley   Level: intermediate
3406ad1575SMatthew G. Knepley 
355d83a8b1SBarry Smith   Note:
365d83a8b1SBarry Smith   This is a struct, not a `PetscObject`
375d83a8b1SBarry Smith 
38db781477SPatrick Sanan .seealso: `DMPlexSNESComputeResidualFEM()`, `DMPlexSNESComputeJacobianFEM()`, `DMPlexSNESComputeBoundaryFEM()`
3906ad1575SMatthew G. Knepley S*/
40ce78bad3SBarry Smith typedef struct {
4106ad1575SMatthew G. Knepley   DMLabel  label; /* The (label, value) select a subdomain */
426528b96dSMatthew G. Knepley   PetscInt value;
4306ad1575SMatthew G. Knepley   PetscInt field; /* Selects the field for the test function */
4406ad1575SMatthew G. Knepley   PetscInt part;  /* Selects the equation part. For example, LHS = 0 and RHS = 1 in IMEX methods. More pieces can be present for operator splitting methods. */
4506ad1575SMatthew G. Knepley } PetscFormKey;
4606ad1575SMatthew G. Knepley 
4706ad1575SMatthew G. Knepley /*E
4806ad1575SMatthew G. Knepley   PetscWeakFormKind - The kind of weak form. The specific forms are given in the documentation for the integraton functions.
4906ad1575SMatthew G. Knepley 
5016a05f60SBarry Smith   Values:
5116a05f60SBarry Smith + OBJECTIVE                  - Objective form
5216a05f60SBarry Smith . F0, F1                     - Residual forms
5316a05f60SBarry Smith . G0, G1, G2, G3             - Jacobian forms
5416a05f60SBarry Smith . GP0, GP1, GP2, GP3         - Jacobian preconditioner matrix forms
5516a05f60SBarry Smith . GT0, GT1, GT2, GT3         - Dynamic Jacobian matrix forms
5616a05f60SBarry Smith . BDF0, BDF1                 - Boundary Residual forms
5716a05f60SBarry Smith . BDG0, BDG1, BDG2, BDG3     - Jacobian forms
5816a05f60SBarry Smith . BDGP0, BDGP1, BDGP2, BDGP3 - Jacobian preconditioner matrix forms
59d2b2dc1eSMatthew G. Knepley . R                          - Riemann solver
60d2b2dc1eSMatthew G. Knepley - CEED                       - libCEED QFunction
6106ad1575SMatthew G. Knepley 
6206ad1575SMatthew G. Knepley   Level: beginner
6306ad1575SMatthew G. Knepley 
6416a05f60SBarry Smith .seealso: `PetscWeakForm`, `PetscFEIntegrateResidual()`, `PetscFEIntegrateJacobian()`, `PetscFEIntegrateBdResidual()`, `PetscFEIntegrateBdJacobian()`,
6516a05f60SBarry Smith           `PetscFVIntegrateRHSFunction()`, `PetscWeakFormSetIndexResidual()`, `PetscWeakFormClearIndex()`
6606ad1575SMatthew G. Knepley E*/
679371c9d4SSatish Balay typedef enum {
689371c9d4SSatish Balay   PETSC_WF_OBJECTIVE,
699371c9d4SSatish Balay   PETSC_WF_F0,
709371c9d4SSatish Balay   PETSC_WF_F1,
719371c9d4SSatish Balay   PETSC_WF_G0,
729371c9d4SSatish Balay   PETSC_WF_G1,
739371c9d4SSatish Balay   PETSC_WF_G2,
749371c9d4SSatish Balay   PETSC_WF_G3,
759371c9d4SSatish Balay   PETSC_WF_GP0,
769371c9d4SSatish Balay   PETSC_WF_GP1,
779371c9d4SSatish Balay   PETSC_WF_GP2,
789371c9d4SSatish Balay   PETSC_WF_GP3,
799371c9d4SSatish Balay   PETSC_WF_GT0,
809371c9d4SSatish Balay   PETSC_WF_GT1,
819371c9d4SSatish Balay   PETSC_WF_GT2,
829371c9d4SSatish Balay   PETSC_WF_GT3,
839371c9d4SSatish Balay   PETSC_WF_BDF0,
849371c9d4SSatish Balay   PETSC_WF_BDF1,
859371c9d4SSatish Balay   PETSC_WF_BDG0,
869371c9d4SSatish Balay   PETSC_WF_BDG1,
879371c9d4SSatish Balay   PETSC_WF_BDG2,
889371c9d4SSatish Balay   PETSC_WF_BDG3,
899371c9d4SSatish Balay   PETSC_WF_BDGP0,
909371c9d4SSatish Balay   PETSC_WF_BDGP1,
919371c9d4SSatish Balay   PETSC_WF_BDGP2,
929371c9d4SSatish Balay   PETSC_WF_BDGP3,
939371c9d4SSatish Balay   PETSC_WF_R,
94d2b2dc1eSMatthew G. Knepley   PETSC_WF_CEED,
959371c9d4SSatish Balay   PETSC_NUM_WF
969371c9d4SSatish Balay } PetscWeakFormKind;
9706ad1575SMatthew G. Knepley PETSC_EXTERN const char *const PetscWeakFormKinds[];
98509b31aaSMatthew G. Knepley 
99*2192575eSBarry Smith /*S
100*2192575eSBarry Smith   PetscPointFn - A prototype of a pointwise function that can be passed to, for example, `PetscDSSetObjective()`
101*2192575eSBarry Smith 
102*2192575eSBarry Smith   Calling Sequence:
103*2192575eSBarry Smith + dim          - the coordinate dimension
104*2192575eSBarry Smith . Nf           - the number of fields
105*2192575eSBarry Smith . NfAux        - the number of auxiliary fields
106*2192575eSBarry Smith . uOff         - the offset into `u`[] and `u_t`[] for each field
107*2192575eSBarry Smith . uOff_x       - the offset into `u_x`[] for each field
108*2192575eSBarry Smith . u            - each field evaluated at the current point
109*2192575eSBarry Smith . u_t          - the time derivative of each field evaluated at the current point
110*2192575eSBarry Smith . u_x          - the gradient of each field evaluated at the current point
111*2192575eSBarry Smith . aOff         - the offset into `a`[] and `a_t`[] for each auxiliary field
112*2192575eSBarry Smith . aOff_x       - the offset into `a_x`[] for each auxiliary field
113*2192575eSBarry Smith . a            - each auxiliary field evaluated at the current point
114*2192575eSBarry Smith . a_t          - the time derivative of each auxiliary field evaluated at the current point
115*2192575eSBarry Smith . a_x          - the gradient of auxiliary each field evaluated at the current point
116*2192575eSBarry Smith . t            - current time
117*2192575eSBarry Smith . x            - coordinates of the current point
118*2192575eSBarry Smith . numConstants - number of constant parameters
119*2192575eSBarry Smith . constants    - constant parameters
120*2192575eSBarry Smith - obj          - output values at the current point
121*2192575eSBarry Smith 
122*2192575eSBarry Smith   Level: beginner
123*2192575eSBarry Smith 
124*2192575eSBarry Smith .seealso: `PetscPointFn`, `PetscDSSetObjective()`, `PetscDSGetObjective()`, PetscDSGetResidual()`, `PetscDSSetResidual()`,
125*2192575eSBarry Smith           `PetscDSGetRHSResidual()`, `PetscDSGetRHSResidual()`, `PetscDSSetUpdate()`, `PetscDSGetUpdate()`, `DMPlexSetCoordinateMap()`
126*2192575eSBarry Smith S*/
127*2192575eSBarry Smith typedef void PetscPointFn(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar result[]);
128*2192575eSBarry Smith 
129*2192575eSBarry Smith /*S
130*2192575eSBarry Smith   PetscPointJacFn - A prototype of a pointwise function that can be passed to, for example, `PetscDSSetJacobian()` for computing Jacobians
131*2192575eSBarry Smith 
132*2192575eSBarry Smith   Calling Sequence:
133*2192575eSBarry Smith + dim          - the coordinate dimension
134*2192575eSBarry Smith . Nf           - the number of fields
135*2192575eSBarry Smith . NfAux        - the number of auxiliary fields
136*2192575eSBarry Smith . uOff         - the offset into `u`[] and `u_t`[] for each field
137*2192575eSBarry Smith . uOff_x       - the offset into `u_x`[] for each field
138*2192575eSBarry Smith . u            - each field evaluated at the current point
139*2192575eSBarry Smith . u_t          - the time derivative of each field evaluated at the current point
140*2192575eSBarry Smith . u_x          - the gradient of each field evaluated at the current point
141*2192575eSBarry Smith . aOff         - the offset into `a`[] and `a_t`[] for each auxiliary field
142*2192575eSBarry Smith . aOff_x       - the offset into a_`x`[] for each auxiliary field
143*2192575eSBarry Smith . a            - each auxiliary field evaluated at the current point
144*2192575eSBarry Smith . a_t          - the time derivative of each auxiliary field evaluated at the current point
145*2192575eSBarry Smith . a_x          - the gradient of auxiliary each field evaluated at the current point
146*2192575eSBarry Smith . t            - current time
147*2192575eSBarry Smith . u_tShift     - the multiplier `a` for $dF/dU_t$
148*2192575eSBarry Smith . x            - coordinates of the current point
149*2192575eSBarry Smith . numConstants - number of constant parameters
150*2192575eSBarry Smith . constants    - constant parameters
151*2192575eSBarry Smith - g            - output values at the current point
152*2192575eSBarry Smith 
153*2192575eSBarry Smith   Level: beginner
154*2192575eSBarry Smith 
155*2192575eSBarry Smith .seealso: `PetscPointFn`, `PetscDSSetJacobian()`, `PetscDSGetJacobian()`, PetscDSSetJacobianPreconditioner()`, `PetscDSGetJacobianPreconditioner()`,
156*2192575eSBarry Smith           `PetscDSSetDynamicJacobian()`, `PetscDSGetDynamicJacobian()`
157*2192575eSBarry Smith S*/
158*2192575eSBarry Smith typedef void PetscPointJacFn(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g[]);
159*2192575eSBarry Smith 
160*2192575eSBarry Smith /*S
161*2192575eSBarry Smith   PetscBdPointFn - A prototype of a pointwise boundary function that can be passed to, for example, `PetscDSSetBdResidual()`
162*2192575eSBarry Smith 
163*2192575eSBarry Smith   Calling Sequence:
164*2192575eSBarry Smith + dim          - the coordinate dimension
165*2192575eSBarry Smith . Nf           - the number of fields
166*2192575eSBarry Smith . NfAux        - the number of auxiliary fields
167*2192575eSBarry Smith . uOff         - the offset into `u`[] and `u_t`[] for each field
168*2192575eSBarry Smith . uOff_x       - the offset into `u_x`[] for each field
169*2192575eSBarry Smith . u            - each field evaluated at the current point
170*2192575eSBarry Smith . u_t          - the time derivative of each field evaluated at the current point
171*2192575eSBarry Smith . u_x          - the gradient of each field evaluated at the current point
172*2192575eSBarry Smith . aOff         - the offset into `a`[] and `a_t`[] for each auxiliary field
173*2192575eSBarry Smith . aOff_x       - the offset into `a_x`[] for each auxiliary field
174*2192575eSBarry Smith . a            - each auxiliary field evaluated at the current point
175*2192575eSBarry Smith . a_t          - the time derivative of each auxiliary field evaluated at the current point
176*2192575eSBarry Smith . a_x          - the gradient of auxiliary each field evaluated at the current point
177*2192575eSBarry Smith . t            - current time
178*2192575eSBarry Smith . x            - coordinates of the current point
179*2192575eSBarry Smith . n            - unit normal at the current point
180*2192575eSBarry Smith . numConstants - number of constant parameters
181*2192575eSBarry Smith . constants    - constant parameters
182*2192575eSBarry Smith - f            - output values at the current point
183*2192575eSBarry Smith 
184*2192575eSBarry Smith   Level: beginner
185*2192575eSBarry Smith 
186*2192575eSBarry Smith .seealso: `PetscPointFn`, `PetscDSSetBdResidual()`, `PetscDSGetBdResidual()`, `PetscDSSetObjective()`, `PetscDSGetObjective()`, PetscDSGetResidual()`,
187*2192575eSBarry Smith           `PetscDSGetRHSResidual()`, `PetscDSGetRHSResidual()`, `PetscDSSetUpdate()`, `PetscDSGetUpdate()`, `DMPlexSetCoordinateMap()`,
188*2192575eSBarry Smith           `PetscDSSetResidual()`, `PetscPointJacFn`
189*2192575eSBarry Smith S*/
190*2192575eSBarry Smith typedef void PetscBdPointFn(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]);
191*2192575eSBarry Smith 
192*2192575eSBarry Smith /*S
193*2192575eSBarry Smith   PetscBdPointJacFn - A prototype of a pointwise boundary function that can be passed to, for example, `PetscDSSetBdJacobian()`
194*2192575eSBarry Smith 
195*2192575eSBarry Smith   Calling Sequence:
196*2192575eSBarry Smith + dim          - the coordinate dimension
197*2192575eSBarry Smith . Nf           - the number of fields
198*2192575eSBarry Smith . NfAux        - the number of auxiliary fields
199*2192575eSBarry Smith . uOff         - the offset into `u`[] and `u_t`[] for each field
200*2192575eSBarry Smith . uOff_x       - the offset into `u_x`[] for each field
201*2192575eSBarry Smith . u            - each field evaluated at the current point
202*2192575eSBarry Smith . u_t          - the time derivative of each field evaluated at the current point
203*2192575eSBarry Smith . u_x          - the gradient of each field evaluated at the current point
204*2192575eSBarry Smith . aOff         - the offset into `a`[] and `a_t`[] for each auxiliary field
205*2192575eSBarry Smith . aOff_x       - the offset into `a_x`[] for each auxiliary field
206*2192575eSBarry Smith . a            - each auxiliary field evaluated at the current point
207*2192575eSBarry Smith . a_t          - the time derivative of each auxiliary field evaluated at the current point
208*2192575eSBarry Smith . a_x          - the gradient of auxiliary each field evaluated at the current point
209*2192575eSBarry Smith . t            - current time
210*2192575eSBarry Smith . u_tShift     - the multiplier `a` for $dF/dU_t$
211*2192575eSBarry Smith . x            - coordinates of the current point
212*2192575eSBarry Smith . n            - normal at the current point
213*2192575eSBarry Smith . numConstants - number of constant parameters
214*2192575eSBarry Smith . constants    - constant parameters
215*2192575eSBarry Smith - g            - output values at the current point
216*2192575eSBarry Smith 
217*2192575eSBarry Smith   Level: beginner
218*2192575eSBarry Smith 
219*2192575eSBarry Smith .seealso: `PetscPointFn`, `PetscDSSetBdJacobian()`, PetscDSGetBdJacobian()`, `PetscDSSetBdJacobianPreconditioner()`, `PetscDSGetBdJacobianPreconditioner()`,
220*2192575eSBarry Smith           `PetscDSSetBdResidual()`, `PetscDSGetBdResidual()`, `PetscDSSetObjective()`, `PetscDSGetObjective()`, PetscDSGetResidual()`,
221*2192575eSBarry Smith           `PetscDSGetRHSResidual()`, `PetscDSGetRHSResidual()`, `PetscDSSetUpdate()`, `PetscDSGetUpdate()`, `DMPlexSetCoordinateMap()`,
222*2192575eSBarry Smith           `PetscDSSetResidual()`, `PetscPointJacFn`
223*2192575eSBarry Smith S*/
224*2192575eSBarry Smith typedef void PetscBdPointJacFn(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]);
225*2192575eSBarry Smith 
226*2192575eSBarry Smith /*S
227*2192575eSBarry Smith   PetscPointExactSolutionFn - A prototype of a pointwise function that computes the exact solution to a PDE. Used with, for example,
228*2192575eSBarry Smith   `PetscDSSetExactSolution()`
229*2192575eSBarry Smith 
230*2192575eSBarry Smith   Calling Sequence:
231*2192575eSBarry Smith + dim - the coordinate dimension
232*2192575eSBarry Smith . t   - current time
233*2192575eSBarry Smith . x   - coordinates of the current point
234*2192575eSBarry Smith . Nc  - the number of field components
235*2192575eSBarry Smith . u   - the solution field evaluated at the current point
236*2192575eSBarry Smith - ctx - a user context, set with `PetscDSSetExactSolution()` or `PetscDSSetExactSolutionTimeDerivative()`
237*2192575eSBarry Smith 
238*2192575eSBarry Smith   Level: beginner
239*2192575eSBarry Smith 
240*2192575eSBarry Smith .seealso: `PetscPointFn`, `PetscDSSetExactSolution()`, `PetscDSGetExactSolution()`, `PetscDSSetExactSolutionTimeDerivative()`, `PetscDSGetExactSolutionTimeDerivative()`
241*2192575eSBarry Smith S*/
242*2192575eSBarry Smith typedef PetscErrorCode PetscPointExactSolutionFn(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx);
243*2192575eSBarry Smith 
244*2192575eSBarry Smith /*S
245*2192575eSBarry Smith   PetscRiemannFn - A prototype of a pointwise function that can be passed to, for example, `PetscDSSetRiemannSolver()`
246*2192575eSBarry Smith 
247*2192575eSBarry Smith   Calling Sequence:
248*2192575eSBarry Smith + dim          - the coordinate dimension
249*2192575eSBarry Smith . Nf           - The number of fields
250*2192575eSBarry Smith . x            - The coordinates at a point on the interface
251*2192575eSBarry Smith . n            - The normal vector to the interface
252*2192575eSBarry Smith . uL           - The state vector to the left of the interface
253*2192575eSBarry Smith . uR           - The state vector to the right of the interface
254*2192575eSBarry Smith . numConstants - number of constant parameters
255*2192575eSBarry Smith . constants    - constant parameters
256*2192575eSBarry Smith . flux         - output array of flux through the interface
257*2192575eSBarry Smith - ctx          - optional user context
258*2192575eSBarry Smith 
259*2192575eSBarry Smith   Level: beginner
260*2192575eSBarry Smith 
261*2192575eSBarry Smith .seealso: `PetscPointFn`, `PetscDSSetRiemannSolver()`, `PetscDSGetRiemannSolver()`
262*2192575eSBarry Smith S*/
263*2192575eSBarry Smith typedef void PetscRiemannFn(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscInt numConstants, const PetscScalar constants[], PetscScalar flux[], void *ctx);
264509b31aaSMatthew G. Knepley 
265509b31aaSMatthew G. Knepley /*S
266509b31aaSMatthew G. Knepley   PetscSimplePointFn - A prototype of a simple pointwise function that can be passed to, for example, `DMPlexTransformExtrudeSetNormalFunction()`
267509b31aaSMatthew G. Knepley 
268509b31aaSMatthew G. Knepley   Calling Sequence:
269509b31aaSMatthew G. Knepley + dim  - The coordinate dimension of the original mesh (usually a surface)
270509b31aaSMatthew G. Knepley . time - The current time, or 0.
271509b31aaSMatthew G. Knepley . x    - The location of the current normal, in the coordinate space of the original mesh
272509b31aaSMatthew G. Knepley . r    - The layer number of this point
273509b31aaSMatthew G. Knepley . u    - The user provides the computed normal on output
274*2192575eSBarry Smith - ctx  - An optional user context, this context may be obtained by the calling code with `DMGetApplicationContext()`
275509b31aaSMatthew G. Knepley 
276509b31aaSMatthew G. Knepley   Level: beginner
277509b31aaSMatthew G. Knepley 
278*2192575eSBarry Smith   Developer Note:
279*2192575eSBarry Smith   The handling of `ctx` in the use of such functions may not be ideal since the context is not provided when the function pointer is provided with, for example, `DMSwarmSetCoordinateFunction()`
280509b31aaSMatthew G. Knepley 
281*2192575eSBarry Smith .seealso: `PetscPointFn`, `DMPlexTransformExtrudeSetNormalFunction()`, `DMSwarmSetCoordinateFunction()`
282509b31aaSMatthew G. Knepley S*/
283*2192575eSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PetscSimplePointFn(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt r, PetscScalar u[], void *ctx);
284509b31aaSMatthew G. Knepley 
285*2192575eSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscSimplePointFn *PetscSimplePointFunc PETSC_DEPRECATED_TYPEDEF(3, 24, 0, "PetscSimplePointFn*", );
286*2192575eSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscPointFn       *PetscPointFunc PETSC_DEPRECATED_TYPEDEF(3, 24, 0, "PetscPointFn*", );
287*2192575eSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscPointJacFn    *PetscPointJac PETSC_DEPRECATED_TYPEDEF(3, 24, 0, "PetscPointJacFn*", );
288*2192575eSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscBdPointFn     *PetscBdPointFunc PETSC_DEPRECATED_TYPEDEF(3, 24, 0, "PetscBdPointFn*", );
289*2192575eSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscBdPointJacFn  *PetscBdPointJac PETSC_DEPRECATED_TYPEDEF(3, 24, 0, "PetscBdPointJacFn*", );
290*2192575eSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscRiemannFn     *PetscRiemannFunc PETSC_DEPRECATED_TYPEDEF(3, 24, 0, "PetscRiemannFn*", );
291*2192575eSBarry Smith 
292*2192575eSBarry Smith /*S
293*2192575eSBarry Smith   PetscPointBoundFn - A prototype of a pointwise function that can be passed to, for example, `PetscDSSetLowerBound()`
294*2192575eSBarry Smith 
295*2192575eSBarry Smith   Calling Sequence:
296*2192575eSBarry Smith + dim - the coordinate dimension
297*2192575eSBarry Smith . t   - current time
298*2192575eSBarry Smith . x   - coordinates of the current point
299*2192575eSBarry Smith . Nc  - the number of field components
300*2192575eSBarry Smith . u   - the lower bound evaluated at the current point
301*2192575eSBarry Smith - ctx - a user context, passed in with, for example, `PetscDSSetLowerBound()`
302*2192575eSBarry Smith 
303*2192575eSBarry Smith   Level: beginner
304*2192575eSBarry Smith 
305*2192575eSBarry Smith .seealso: `PetscPointFn`, `PetscDSSetLowerBound()`, `PetscDSSetUpperBound()`
306*2192575eSBarry Smith S*/
307*2192575eSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PetscPointBoundFn(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], void *ctx);
308