xref: /petsc/include/petscdstypes.h (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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
547addb90fSBarry Smith . GP0, GP1, GP2, GP3         - Jacobian forms used to construct the preconditioner
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
587addb90fSBarry Smith . BDGP0, BDGP1, BDGP2, BDGP3 - Jacobian forms used to construct the preconditioner
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 
992192575eSBarry Smith /*S
1002192575eSBarry Smith   PetscPointFn - A prototype of a pointwise function that can be passed to, for example, `PetscDSSetObjective()`
1012192575eSBarry Smith 
1022192575eSBarry Smith   Calling Sequence:
1032192575eSBarry Smith + dim          - the coordinate dimension
1042192575eSBarry Smith . Nf           - the number of fields
1052192575eSBarry Smith . NfAux        - the number of auxiliary fields
1062192575eSBarry Smith . uOff         - the offset into `u`[] and `u_t`[] for each field
1072192575eSBarry Smith . uOff_x       - the offset into `u_x`[] for each field
1082192575eSBarry Smith . u            - each field evaluated at the current point
1092192575eSBarry Smith . u_t          - the time derivative of each field evaluated at the current point
1102192575eSBarry Smith . u_x          - the gradient of each field evaluated at the current point
1112192575eSBarry Smith . aOff         - the offset into `a`[] and `a_t`[] for each auxiliary field
1122192575eSBarry Smith . aOff_x       - the offset into `a_x`[] for each auxiliary field
1132192575eSBarry Smith . a            - each auxiliary field evaluated at the current point
1142192575eSBarry Smith . a_t          - the time derivative of each auxiliary field evaluated at the current point
1152192575eSBarry Smith . a_x          - the gradient of auxiliary each field evaluated at the current point
1162192575eSBarry Smith . t            - current time
1172192575eSBarry Smith . x            - coordinates of the current point
1182192575eSBarry Smith . numConstants - number of constant parameters
1192192575eSBarry Smith . constants    - constant parameters
1202192575eSBarry Smith - obj          - output values at the current point
1212192575eSBarry Smith 
1222192575eSBarry Smith   Level: beginner
1232192575eSBarry Smith 
1242192575eSBarry Smith .seealso: `PetscPointFn`, `PetscDSSetObjective()`, `PetscDSGetObjective()`, PetscDSGetResidual()`, `PetscDSSetResidual()`,
1252192575eSBarry Smith           `PetscDSGetRHSResidual()`, `PetscDSGetRHSResidual()`, `PetscDSSetUpdate()`, `PetscDSGetUpdate()`, `DMPlexSetCoordinateMap()`
1262192575eSBarry Smith S*/
1272192575eSBarry 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[]);
1282192575eSBarry Smith 
1292192575eSBarry Smith /*S
1302192575eSBarry Smith   PetscPointJacFn - A prototype of a pointwise function that can be passed to, for example, `PetscDSSetJacobian()` for computing Jacobians
1312192575eSBarry Smith 
1322192575eSBarry Smith   Calling Sequence:
1332192575eSBarry Smith + dim          - the coordinate dimension
1342192575eSBarry Smith . Nf           - the number of fields
1352192575eSBarry Smith . NfAux        - the number of auxiliary fields
1362192575eSBarry Smith . uOff         - the offset into `u`[] and `u_t`[] for each field
1372192575eSBarry Smith . uOff_x       - the offset into `u_x`[] for each field
1382192575eSBarry Smith . u            - each field evaluated at the current point
1392192575eSBarry Smith . u_t          - the time derivative of each field evaluated at the current point
1402192575eSBarry Smith . u_x          - the gradient of each field evaluated at the current point
1412192575eSBarry Smith . aOff         - the offset into `a`[] and `a_t`[] for each auxiliary field
1422192575eSBarry Smith . aOff_x       - the offset into a_`x`[] for each auxiliary field
1432192575eSBarry Smith . a            - each auxiliary field evaluated at the current point
1442192575eSBarry Smith . a_t          - the time derivative of each auxiliary field evaluated at the current point
1452192575eSBarry Smith . a_x          - the gradient of auxiliary each field evaluated at the current point
1462192575eSBarry Smith . t            - current time
1472192575eSBarry Smith . u_tShift     - the multiplier `a` for $dF/dU_t$
1482192575eSBarry Smith . x            - coordinates of the current point
1492192575eSBarry Smith . numConstants - number of constant parameters
1502192575eSBarry Smith . constants    - constant parameters
1512192575eSBarry Smith - g            - output values at the current point
1522192575eSBarry Smith 
1532192575eSBarry Smith   Level: beginner
1542192575eSBarry Smith 
1552192575eSBarry Smith .seealso: `PetscPointFn`, `PetscDSSetJacobian()`, `PetscDSGetJacobian()`, PetscDSSetJacobianPreconditioner()`, `PetscDSGetJacobianPreconditioner()`,
1562192575eSBarry Smith           `PetscDSSetDynamicJacobian()`, `PetscDSGetDynamicJacobian()`
1572192575eSBarry Smith S*/
1582192575eSBarry 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[]);
1592192575eSBarry Smith 
1602192575eSBarry Smith /*S
1612192575eSBarry Smith   PetscBdPointFn - A prototype of a pointwise boundary function that can be passed to, for example, `PetscDSSetBdResidual()`
1622192575eSBarry Smith 
1632192575eSBarry Smith   Calling Sequence:
1642192575eSBarry Smith + dim          - the coordinate dimension
1652192575eSBarry Smith . Nf           - the number of fields
1662192575eSBarry Smith . NfAux        - the number of auxiliary fields
1672192575eSBarry Smith . uOff         - the offset into `u`[] and `u_t`[] for each field
1682192575eSBarry Smith . uOff_x       - the offset into `u_x`[] for each field
1692192575eSBarry Smith . u            - each field evaluated at the current point
1702192575eSBarry Smith . u_t          - the time derivative of each field evaluated at the current point
1712192575eSBarry Smith . u_x          - the gradient of each field evaluated at the current point
1722192575eSBarry Smith . aOff         - the offset into `a`[] and `a_t`[] for each auxiliary field
1732192575eSBarry Smith . aOff_x       - the offset into `a_x`[] for each auxiliary field
1742192575eSBarry Smith . a            - each auxiliary field evaluated at the current point
1752192575eSBarry Smith . a_t          - the time derivative of each auxiliary field evaluated at the current point
1762192575eSBarry Smith . a_x          - the gradient of auxiliary each field evaluated at the current point
1772192575eSBarry Smith . t            - current time
1782192575eSBarry Smith . x            - coordinates of the current point
1792192575eSBarry Smith . n            - unit normal at the current point
1802192575eSBarry Smith . numConstants - number of constant parameters
1812192575eSBarry Smith . constants    - constant parameters
1822192575eSBarry Smith - f            - output values at the current point
1832192575eSBarry Smith 
1842192575eSBarry Smith   Level: beginner
1852192575eSBarry Smith 
1862192575eSBarry Smith .seealso: `PetscPointFn`, `PetscDSSetBdResidual()`, `PetscDSGetBdResidual()`, `PetscDSSetObjective()`, `PetscDSGetObjective()`, PetscDSGetResidual()`,
1872192575eSBarry Smith           `PetscDSGetRHSResidual()`, `PetscDSGetRHSResidual()`, `PetscDSSetUpdate()`, `PetscDSGetUpdate()`, `DMPlexSetCoordinateMap()`,
1882192575eSBarry Smith           `PetscDSSetResidual()`, `PetscPointJacFn`
1892192575eSBarry Smith S*/
1902192575eSBarry 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[]);
1912192575eSBarry Smith 
1922192575eSBarry Smith /*S
1932192575eSBarry Smith   PetscBdPointJacFn - A prototype of a pointwise boundary function that can be passed to, for example, `PetscDSSetBdJacobian()`
1942192575eSBarry Smith 
1952192575eSBarry Smith   Calling Sequence:
1962192575eSBarry Smith + dim          - the coordinate dimension
1972192575eSBarry Smith . Nf           - the number of fields
1982192575eSBarry Smith . NfAux        - the number of auxiliary fields
1992192575eSBarry Smith . uOff         - the offset into `u`[] and `u_t`[] for each field
2002192575eSBarry Smith . uOff_x       - the offset into `u_x`[] for each field
2012192575eSBarry Smith . u            - each field evaluated at the current point
2022192575eSBarry Smith . u_t          - the time derivative of each field evaluated at the current point
2032192575eSBarry Smith . u_x          - the gradient of each field evaluated at the current point
2042192575eSBarry Smith . aOff         - the offset into `a`[] and `a_t`[] for each auxiliary field
2052192575eSBarry Smith . aOff_x       - the offset into `a_x`[] for each auxiliary field
2062192575eSBarry Smith . a            - each auxiliary field evaluated at the current point
2072192575eSBarry Smith . a_t          - the time derivative of each auxiliary field evaluated at the current point
2082192575eSBarry Smith . a_x          - the gradient of auxiliary each field evaluated at the current point
2092192575eSBarry Smith . t            - current time
2102192575eSBarry Smith . u_tShift     - the multiplier `a` for $dF/dU_t$
2112192575eSBarry Smith . x            - coordinates of the current point
2122192575eSBarry Smith . n            - normal at the current point
2132192575eSBarry Smith . numConstants - number of constant parameters
2142192575eSBarry Smith . constants    - constant parameters
2152192575eSBarry Smith - g            - output values at the current point
2162192575eSBarry Smith 
2172192575eSBarry Smith   Level: beginner
2182192575eSBarry Smith 
2192192575eSBarry Smith .seealso: `PetscPointFn`, `PetscDSSetBdJacobian()`, PetscDSGetBdJacobian()`, `PetscDSSetBdJacobianPreconditioner()`, `PetscDSGetBdJacobianPreconditioner()`,
2202192575eSBarry Smith           `PetscDSSetBdResidual()`, `PetscDSGetBdResidual()`, `PetscDSSetObjective()`, `PetscDSGetObjective()`, PetscDSGetResidual()`,
2212192575eSBarry Smith           `PetscDSGetRHSResidual()`, `PetscDSGetRHSResidual()`, `PetscDSSetUpdate()`, `PetscDSGetUpdate()`, `DMPlexSetCoordinateMap()`,
2222192575eSBarry Smith           `PetscDSSetResidual()`, `PetscPointJacFn`
2232192575eSBarry Smith S*/
2242192575eSBarry 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[]);
2252192575eSBarry Smith 
2262192575eSBarry Smith /*S
2272192575eSBarry Smith   PetscPointExactSolutionFn - A prototype of a pointwise function that computes the exact solution to a PDE. Used with, for example,
2282192575eSBarry Smith   `PetscDSSetExactSolution()`
2292192575eSBarry Smith 
2302192575eSBarry Smith   Calling Sequence:
2312192575eSBarry Smith + dim - the coordinate dimension
2322192575eSBarry Smith . t   - current time
2332192575eSBarry Smith . x   - coordinates of the current point
2342192575eSBarry Smith . Nc  - the number of field components
2352192575eSBarry Smith . u   - the solution field evaluated at the current point
2362192575eSBarry Smith - ctx - a user context, set with `PetscDSSetExactSolution()` or `PetscDSSetExactSolutionTimeDerivative()`
2372192575eSBarry Smith 
2382192575eSBarry Smith   Level: beginner
2392192575eSBarry Smith 
2402192575eSBarry Smith .seealso: `PetscPointFn`, `PetscDSSetExactSolution()`, `PetscDSGetExactSolution()`, `PetscDSSetExactSolutionTimeDerivative()`, `PetscDSGetExactSolutionTimeDerivative()`
2412192575eSBarry Smith S*/
242*2a8381b2SBarry Smith typedef PetscErrorCode PetscPointExactSolutionFn(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], PetscCtx ctx);
2432192575eSBarry Smith 
2442192575eSBarry Smith /*S
2452192575eSBarry Smith   PetscRiemannFn - A prototype of a pointwise function that can be passed to, for example, `PetscDSSetRiemannSolver()`
2462192575eSBarry Smith 
2472192575eSBarry Smith   Calling Sequence:
2482192575eSBarry Smith + dim          - the coordinate dimension
2492192575eSBarry Smith . Nf           - The number of fields
2502192575eSBarry Smith . x            - The coordinates at a point on the interface
2512192575eSBarry Smith . n            - The normal vector to the interface
2522192575eSBarry Smith . uL           - The state vector to the left of the interface
2532192575eSBarry Smith . uR           - The state vector to the right of the interface
2542192575eSBarry Smith . numConstants - number of constant parameters
2552192575eSBarry Smith . constants    - constant parameters
2562192575eSBarry Smith . flux         - output array of flux through the interface
2572192575eSBarry Smith - ctx          - optional user context
2582192575eSBarry Smith 
2592192575eSBarry Smith   Level: beginner
2602192575eSBarry Smith 
2612192575eSBarry Smith .seealso: `PetscPointFn`, `PetscDSSetRiemannSolver()`, `PetscDSGetRiemannSolver()`
2622192575eSBarry Smith S*/
263*2a8381b2SBarry 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[], PetscCtx 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
2742192575eSBarry 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 
2782192575eSBarry Smith   Developer Note:
2792192575eSBarry 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 
2812192575eSBarry Smith .seealso: `PetscPointFn`, `DMPlexTransformExtrudeSetNormalFunction()`, `DMSwarmSetCoordinateFunction()`
282509b31aaSMatthew G. Knepley S*/
283*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PetscSimplePointFn(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt r, PetscScalar u[], PetscCtx ctx);
284509b31aaSMatthew G. Knepley 
2852192575eSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscSimplePointFn *PetscSimplePointFunc PETSC_DEPRECATED_TYPEDEF(3, 24, 0, "PetscSimplePointFn*", );
2862192575eSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscPointFn       *PetscPointFunc PETSC_DEPRECATED_TYPEDEF(3, 24, 0, "PetscPointFn*", );
2872192575eSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscPointJacFn    *PetscPointJac PETSC_DEPRECATED_TYPEDEF(3, 24, 0, "PetscPointJacFn*", );
2882192575eSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscBdPointFn     *PetscBdPointFunc PETSC_DEPRECATED_TYPEDEF(3, 24, 0, "PetscBdPointFn*", );
2892192575eSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscBdPointJacFn  *PetscBdPointJac PETSC_DEPRECATED_TYPEDEF(3, 24, 0, "PetscBdPointJacFn*", );
2902192575eSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscRiemannFn     *PetscRiemannFunc PETSC_DEPRECATED_TYPEDEF(3, 24, 0, "PetscRiemannFn*", );
2912192575eSBarry Smith 
2922192575eSBarry Smith /*S
2932192575eSBarry Smith   PetscPointBoundFn - A prototype of a pointwise function that can be passed to, for example, `PetscDSSetLowerBound()`
2942192575eSBarry Smith 
2952192575eSBarry Smith   Calling Sequence:
2962192575eSBarry Smith + dim - the coordinate dimension
2972192575eSBarry Smith . t   - current time
2982192575eSBarry Smith . x   - coordinates of the current point
2992192575eSBarry Smith . Nc  - the number of field components
3002192575eSBarry Smith . u   - the lower bound evaluated at the current point
3012192575eSBarry Smith - ctx - a user context, passed in with, for example, `PetscDSSetLowerBound()`
3022192575eSBarry Smith 
3032192575eSBarry Smith   Level: beginner
3042192575eSBarry Smith 
3052192575eSBarry Smith .seealso: `PetscPointFn`, `PetscDSSetLowerBound()`, `PetscDSSetUpperBound()`
3062192575eSBarry Smith S*/
307*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PetscPointBoundFn(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nc, PetscScalar u[], PetscCtx ctx);
308