Actual source code: petscdmceed.h
1: #pragma once
3: #include <petscdm.h>
5: #if defined(PETSC_HAVE_LIBCEED)
6: #include <ceed.h>
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__)
19: PETSC_EXTERN PetscErrorCode DMGetCeed(DM, Ceed *);
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: PETSC_EXTERN PetscErrorCode DMCeedCreateFVM(DM, PetscBool, CeedQFunctionUser, const char *, CeedQFunctionContext);
29: struct _PETSc_DMCEED {
30: CeedBasis basis; // Basis for element function space
31: CeedElemRestriction er; // Map from PETSc local vector to FE element vectors
32: CeedElemRestriction erL; // Map from PETSc local vector to FV left cell vectors
33: CeedElemRestriction erR; // Map from PETSc local vector to FV right cell vectors
34: CeedQFunctionUser func; // Plex Function for this operator
35: char *funcSource; // Plex Function source as text
36: CeedQFunction qf; // QFunction expressing the operator action
37: CeedOperator op; // Operator action for this object
38: DMCeed geom; // Operator computing geometric data at quadrature points
39: CeedElemRestriction erq; // Map from PETSc local vector to quadrature points
40: CeedVector qd; // Geometric data at quadrature points used in calculating the qfunction
41: DMCeed info; // Mesh information at quadrature points
42: CeedElemRestriction eri; // Map from PETSc local vector to quadrature points
43: CeedVector qi; // Mesh information at quadrature points
44: };
46: #else
48: struct _PETSc_DMCEED {
49: PetscInt dummy;
50: };
52: #endif
54: PETSC_EXTERN PetscErrorCode DMCeedComputeGeometry(DM, DMCeed);
55: PETSC_EXTERN PetscErrorCode DMCeedComputeInfo(DM, DMCeed);
56: PETSC_EXTERN PetscErrorCode DMCeedDestroy(DMCeed *);