xref: /petsc/include/petscds.h (revision 249df28457db0a897b34333bdf9e5494258ac347)
12764a2aaSMatthew G. Knepley /*
22764a2aaSMatthew G. Knepley       Objects which encapsulate discretizations+continuum residuals
32764a2aaSMatthew G. Knepley */
42764a2aaSMatthew G. Knepley #if !defined(__PETSCDS_H)
52764a2aaSMatthew G. Knepley #define __PETSCDS_H
62764a2aaSMatthew G. Knepley #include <petscfe.h>
72764a2aaSMatthew G. Knepley #include <petscfv.h>
82764a2aaSMatthew G. Knepley #include <petscdstypes.h>
92764a2aaSMatthew G. Knepley 
102764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSInitializePackage(void);
112764a2aaSMatthew G. Knepley 
122764a2aaSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCDS_CLASSID;
132764a2aaSMatthew G. Knepley 
142764a2aaSMatthew G. Knepley /*J
152764a2aaSMatthew G. Knepley   PetscDSType - String with the name of a PETSc discrete system
162764a2aaSMatthew G. Knepley 
172764a2aaSMatthew G. Knepley   Level: beginner
182764a2aaSMatthew G. Knepley 
192764a2aaSMatthew G. Knepley .seealso: PetscDSSetType(), PetscDS
202764a2aaSMatthew G. Knepley J*/
212764a2aaSMatthew G. Knepley typedef const char *PetscDSType;
222764a2aaSMatthew G. Knepley #define PETSCDSBASIC "basic"
232764a2aaSMatthew G. Knepley 
242764a2aaSMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscDSList;
252764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSCreate(MPI_Comm, PetscDS *);
262764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSDestroy(PetscDS *);
272764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetType(PetscDS, PetscDSType);
282764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetType(PetscDS, PetscDSType *);
292764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetUp(PetscDS);
302764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetFromOptions(PetscDS);
312764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSViewFromOptions(PetscDS,const char[],const char[]);
322764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSView(PetscDS,PetscViewer);
332764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSRegister(const char [], PetscErrorCode (*)(PetscDS));
342764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSRegisterDestroy(void);
352764a2aaSMatthew G. Knepley 
362764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetSpatialDimension(PetscDS, PetscInt *);
372764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetNumFields(PetscDS, PetscInt *);
382764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetTotalDimension(PetscDS, PetscInt *);
392764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetTotalBdDimension(PetscDS, PetscInt *);
402764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetTotalComponents(PetscDS, PetscInt *);
412764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetFieldOffset(PetscDS, PetscInt, PetscInt *);
422764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetBdFieldOffset(PetscDS, PetscInt, PetscInt *);
436ce16762SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffset(PetscDS, PetscInt, PetscInt *);
442764a2aaSMatthew G. Knepley 
452764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetDiscretization(PetscDS, PetscInt, PetscObject *);
462764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetDiscretization(PetscDS, PetscInt, PetscObject);
472764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSAddDiscretization(PetscDS, PetscObject);
482764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetBdDiscretization(PetscDS, PetscInt, PetscObject *);
492764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetBdDiscretization(PetscDS, PetscInt, PetscObject);
502764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSAddBdDiscretization(PetscDS, PetscObject);
51*249df284SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetImplicit(PetscDS, PetscInt, PetscBool*);
52*249df284SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetImplicit(PetscDS, PetscInt, PetscBool);
532764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetObjective(PetscDS, PetscInt, void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]));
542764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetObjective(PetscDS, PetscInt, void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]));
552764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetResidual(PetscDS, PetscInt,
562764a2aaSMatthew G. Knepley                                                     void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
572764a2aaSMatthew G. Knepley                                                     void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]));
582764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetResidual(PetscDS, PetscInt,
592764a2aaSMatthew G. Knepley                                                     void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
602764a2aaSMatthew G. Knepley                                                     void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]));
612764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetJacobian(PetscDS, PetscInt, PetscInt,
622764a2aaSMatthew G. Knepley                                                     void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
632764a2aaSMatthew G. Knepley                                                     void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
642764a2aaSMatthew G. Knepley                                                     void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
652764a2aaSMatthew G. Knepley                                                     void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]));
662764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetJacobian(PetscDS, PetscInt, PetscInt,
672764a2aaSMatthew G. Knepley                                                     void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
682764a2aaSMatthew G. Knepley                                                     void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
692764a2aaSMatthew G. Knepley                                                     void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
702764a2aaSMatthew G. Knepley                                                     void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]));
710c2f2876SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetRiemannSolver(PetscDS, PetscInt,
720c2f2876SMatthew G. Knepley                                                     void (**)(const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *));
730c2f2876SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetRiemannSolver(PetscDS, PetscInt,
740c2f2876SMatthew G. Knepley                                                     void (*)(const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *));
750c2f2876SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetContext(PetscDS, PetscInt, void **);
760c2f2876SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetContext(PetscDS, PetscInt, void *);
772764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetBdResidual(PetscDS, PetscInt,
782764a2aaSMatthew G. Knepley                                                       void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
792764a2aaSMatthew G. Knepley                                                       void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]));
802764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetBdResidual(PetscDS, PetscInt,
812764a2aaSMatthew G. Knepley                                                       void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
822764a2aaSMatthew G. Knepley                                                       void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]));
832764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetBdJacobian(PetscDS, PetscInt, PetscInt,
842764a2aaSMatthew G. Knepley                                                       void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
852764a2aaSMatthew G. Knepley                                                       void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
862764a2aaSMatthew G. Knepley                                                       void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
872764a2aaSMatthew G. Knepley                                                       void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]));
882764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSSetBdJacobian(PetscDS, PetscInt, PetscInt,
892764a2aaSMatthew G. Knepley                                                       void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
902764a2aaSMatthew G. Knepley                                                       void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
912764a2aaSMatthew G. Knepley                                                       void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
922764a2aaSMatthew G. Knepley                                                       void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]));
932764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetTabulation(PetscDS, PetscReal ***, PetscReal ***);
942764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetBdTabulation(PetscDS, PetscReal ***, PetscReal ***);
952764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetEvaluationArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **);
962764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetWeakFormArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **);
972764a2aaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDSGetRefCoordArrays(PetscDS, PetscReal **, PetscScalar **);
982764a2aaSMatthew G. Knepley 
992764a2aaSMatthew G. Knepley #endif
100