xref: /petsc/include/petscds.h (revision 2764a2aaaec701d8975de44069275628a72629e7)
1*2764a2aaSMatthew G. Knepley /*
2*2764a2aaSMatthew G. Knepley       Objects which encapsulate discretizations+continuum residuals
3*2764a2aaSMatthew G. Knepley */
4*2764a2aaSMatthew G. Knepley #if !defined(__PETSCDS_H)
5*2764a2aaSMatthew G. Knepley #define __PETSCDS_H
6*2764a2aaSMatthew G. Knepley #include <petscfe.h>
7*2764a2aaSMatthew G. Knepley #include <petscfv.h>
8*2764a2aaSMatthew G. Knepley #include <petscdstypes.h>
9*2764a2aaSMatthew G. Knepley 
10*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSInitializePackage(void);
11*2764a2aaSMatthew G. Knepley 
12*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCDS_CLASSID;
13*2764a2aaSMatthew G. Knepley 
14*2764a2aaSMatthew G. Knepley /*J
15*2764a2aaSMatthew G. Knepley   PetscDSType - String with the name of a PETSc discrete system
16*2764a2aaSMatthew G. Knepley 
17*2764a2aaSMatthew G. Knepley   Level: beginner
18*2764a2aaSMatthew G. Knepley 
19*2764a2aaSMatthew G. Knepley .seealso: PetscDSSetType(), PetscDS
20*2764a2aaSMatthew G. Knepley J*/
21*2764a2aaSMatthew G. Knepley typedef const char *PetscDSType;
22*2764a2aaSMatthew G. Knepley #define PETSCDSBASIC "basic"
23*2764a2aaSMatthew G. Knepley 
24*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscDSList;
25*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscBool         PetscDSRegisterAllCalled;
26*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSCreate(MPI_Comm, PetscDS *);
27*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSDestroy(PetscDS *);
28*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetType(PetscDS, PetscDSType);
29*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetType(PetscDS, PetscDSType *);
30*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetUp(PetscDS);
31*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetFromOptions(PetscDS);
32*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSViewFromOptions(PetscDS,const char[],const char[]);
33*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSView(PetscDS,PetscViewer);
34*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSRegister(const char [], PetscErrorCode (*)(PetscDS));
35*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSRegisterAll(void);
36*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSRegisterDestroy(void);
37*2764a2aaSMatthew G. Knepley 
38*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetSpatialDimension(PetscDS, PetscInt *);
39*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetNumFields(PetscDS, PetscInt *);
40*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetTotalDimension(PetscDS, PetscInt *);
41*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetTotalBdDimension(PetscDS, PetscInt *);
42*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetTotalComponents(PetscDS, PetscInt *);
43*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetFieldOffset(PetscDS, PetscInt, PetscInt *);
44*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetBdFieldOffset(PetscDS, PetscInt, PetscInt *);
45*2764a2aaSMatthew G. Knepley 
46*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetDiscretization(PetscDS, PetscInt, PetscObject *);
47*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetDiscretization(PetscDS, PetscInt, PetscObject);
48*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSAddDiscretization(PetscDS, PetscObject);
49*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetBdDiscretization(PetscDS, PetscInt, PetscObject *);
50*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetBdDiscretization(PetscDS, PetscInt, PetscObject);
51*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSAddBdDiscretization(PetscDS, PetscObject);
52*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetObjective(PetscDS, PetscInt, void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]));
53*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetObjective(PetscDS, PetscInt, void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]));
54*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetResidual(PetscDS, PetscInt,
55*2764a2aaSMatthew G. Knepley                                                     void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
56*2764a2aaSMatthew G. Knepley                                                     void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]));
57*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetResidual(PetscDS, PetscInt,
58*2764a2aaSMatthew G. Knepley                                                     void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
59*2764a2aaSMatthew G. Knepley                                                     void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]));
60*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetJacobian(PetscDS, PetscInt, PetscInt,
61*2764a2aaSMatthew G. Knepley                                                     void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
62*2764a2aaSMatthew G. Knepley                                                     void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
63*2764a2aaSMatthew G. Knepley                                                     void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
64*2764a2aaSMatthew G. Knepley                                                     void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]));
65*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetJacobian(PetscDS, PetscInt, PetscInt,
66*2764a2aaSMatthew G. Knepley                                                     void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
67*2764a2aaSMatthew G. Knepley                                                     void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
68*2764a2aaSMatthew G. Knepley                                                     void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
69*2764a2aaSMatthew G. Knepley                                                     void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]));
70*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetBdResidual(PetscDS, PetscInt,
71*2764a2aaSMatthew G. Knepley                                                       void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
72*2764a2aaSMatthew G. Knepley                                                       void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]));
73*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetBdResidual(PetscDS, PetscInt,
74*2764a2aaSMatthew G. Knepley                                                       void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
75*2764a2aaSMatthew G. Knepley                                                       void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]));
76*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetBdJacobian(PetscDS, PetscInt, PetscInt,
77*2764a2aaSMatthew G. Knepley                                                       void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
78*2764a2aaSMatthew G. Knepley                                                       void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
79*2764a2aaSMatthew G. Knepley                                                       void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
80*2764a2aaSMatthew G. Knepley                                                       void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]));
81*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetBdJacobian(PetscDS, PetscInt, PetscInt,
82*2764a2aaSMatthew G. Knepley                                                       void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
83*2764a2aaSMatthew G. Knepley                                                       void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
84*2764a2aaSMatthew G. Knepley                                                       void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
85*2764a2aaSMatthew G. Knepley                                                       void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]));
86*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetTabulation(PetscDS, PetscReal ***, PetscReal ***);
87*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetBdTabulation(PetscDS, PetscReal ***, PetscReal ***);
88*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetEvaluationArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **);
89*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetWeakFormArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **);
90*2764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetRefCoordArrays(PetscDS, PetscReal **, PetscScalar **);
91*2764a2aaSMatthew G. Knepley 
92*2764a2aaSMatthew G. Knepley #endif
93