xref: /petsc/include/petscdmceed.h (revision 66af8762ec03dbef0e079729eb2a1734a35ed7ff)
1 #pragma once
2 
3 #include <petscdm.h>
4 
5 #if defined(PETSC_HAVE_LIBCEED)
6   #include <ceed.h>
7 
8   #if defined(PETSC_CLANG_STATIC_ANALYZER)
9 void PetscCallCEED(CeedErrorType);
10   #else
11     #define PetscCallCEED(...) \
12       do { \
13         CeedErrorType ierr_ceed_ = __VA_ARGS__; \
14         PetscCheck(ierr_ceed_ == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "libCEED error: %s", CeedErrorTypes[ierr_ceed_]); \
15       } while (0)
16   #endif /* PETSC_CLANG_STATIC_ANALYZER */
17   #define CHKERRQ_CEED(...) PetscCallCEED(__VA_ARGS__)
18 
19 PETSC_EXTERN PetscErrorCode DMGetCeed(DM, Ceed *);
20 
21 PETSC_EXTERN PetscErrorCode VecGetCeedVector(Vec, Ceed, CeedVector *);
22 PETSC_EXTERN PetscErrorCode VecGetCeedVectorRead(Vec, Ceed, CeedVector *);
23 PETSC_EXTERN PetscErrorCode VecRestoreCeedVector(Vec, CeedVector *);
24 PETSC_EXTERN PetscErrorCode VecRestoreCeedVectorRead(Vec, CeedVector *);
25 PETSC_INTERN PetscErrorCode DMCeedCreate_Internal(DM, IS, PetscBool, CeedQFunctionUser, const char *, DMCeed *);
26 PETSC_EXTERN PetscErrorCode DMCeedCreate(DM, PetscBool, CeedQFunctionUser, const char *);
27 
28 struct _PETSc_DMCEED {
29   CeedBasis           basis;      // Basis for element function space
30   CeedElemRestriction er;         // Map from PETSc local vector to element vectors
31   CeedQFunctionUser   func;       // Plex Function for this operator
32   char               *funcSource; // Plex Function source as text
33   CeedQFunction       qf;         // QFunction expressing the operator action
34   CeedOperator        op;         // Operator action for this object
35   DMCeed              geom;       // Operator computing geometric data at quadrature points
36   CeedElemRestriction erq;        // Map from PETSc local vector to quadrature points
37   CeedVector          qd;         // Geometric data at quadrature points used in calculating the qfunction
38 };
39 
40 #else
41 
42 struct _PETSc_DMCEED {
43   PetscInt dummy;
44 };
45 
46 #endif
47 
48 PETSC_EXTERN PetscErrorCode DMCeedComputeGeometry(DM, DMCeed);
49 PETSC_EXTERN PetscErrorCode DMCeedDestroy(DMCeed *);
50