xref: /petsc/include/petscdmplexegads.h (revision 5552b385b77b214b234683fbe1b434751e6350f0)
1*5552b385SBrandon #pragma once
2*5552b385SBrandon 
3*5552b385SBrandon #include <petscdmplex.h>
4*5552b385SBrandon 
5*5552b385SBrandon #if !defined(PETSC_HAVE_EGADS)
6*5552b385SBrandon   #error "PETSc not configured for EGADS; reconfigrue --with-egads or --download-egads"
7*5552b385SBrandon #endif
8*5552b385SBrandon 
9*5552b385SBrandon /* Declarations below provide an interface to the EGADS/EGADSlite libraries, that can be automatically installed
10*5552b385SBrandon    by PETSc. These functions enable the creation, interrogation, and manipulation of CAD geometries attached to
11*5552b385SBrandon    a DMPlex. */
12*5552b385SBrandon #include <egads.h>
13*5552b385SBrandon #include <egads_lite.h>
14*5552b385SBrandon typedef ego PetscGeom;
15*5552b385SBrandon 
16*5552b385SBrandon #define PetscCallEGADS(func, args) \
17*5552b385SBrandon   do { \
18*5552b385SBrandon     int _status; \
19*5552b385SBrandon     PetscStackPushExternal(#func); \
20*5552b385SBrandon     _status = func args; \
21*5552b385SBrandon     PetscStackPop; \
22*5552b385SBrandon     PetscCheck(_status >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in EGADS call %s() Status %d", #func, (int)_status); \
23*5552b385SBrandon   } while (0)
24*5552b385SBrandon 
25*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexGeomDataAndGrads(DM, PetscBool);
26*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexModifyGeomModel(DM, MPI_Comm, PetscScalar[], PetscScalar[], PetscBool, PetscBool, const char[]);
27*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexInflateToGeomModelUseXYZ(DM);
28*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelTUV(DM);
29*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexInflateToGeomModelUseTUV(DM);
30*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelBodies(DM, PetscGeom **, PetscInt *);
31*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelBodyShells(DM, PetscGeom, PetscGeom **, PetscInt *);
32*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelBodyFaces(DM, PetscGeom, PetscGeom **, PetscInt *);
33*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelBodyLoops(DM, PetscGeom, PetscGeom **, PetscInt *);
34*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelBodyEdges(DM, PetscGeom, PetscGeom **, PetscInt *);
35*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelBodyNodes(DM, PetscGeom, PetscGeom **, PetscInt *);
36*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelShellFaces(DM, PetscGeom, PetscGeom, PetscGeom **, PetscInt *);
37*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelFaceLoops(DM, PetscGeom, PetscGeom, PetscGeom **, PetscInt *);
38*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelFaceEdges(DM, PetscGeom, PetscGeom, PetscGeom **, PetscInt *);
39*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelEdgeNodes(DM, PetscGeom, PetscGeom, PetscGeom **, PetscInt *);
40*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexGetGeomBodyMassProperties(DM, PetscGeom, PetscScalar *, PetscScalar *, PetscScalar **, PetscInt *, PetscScalar **, PetscInt *);
41*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexRestoreGeomBodyMassProperties(DM, PetscGeom, PetscScalar *, PetscScalar *, PetscScalar **, PetscInt *, PetscScalar **, PetscInt *);
42*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexGetGeomID(DM, PetscGeom, PetscGeom, PetscInt *);
43*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexGetGeomObject(DM, PetscGeom, PetscInt, PetscInt, PetscGeom *);
44*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexGetGeomFaceNumOfControlPoints(DM, PetscGeom, PetscInt *);
45*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexFreeGeomObject(DM, PetscGeom *);
46*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexGetGeomCntrlPntAndWeightData(DM, PetscHMapI *, PetscInt *, PetscScalar **, PetscInt *, Mat *, PetscHMapI *, PetscInt *, PetscScalar **);
47*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexRestoreGeomCntrlPntAndWeightData(DM, PetscHMapI *, PetscInt *, PetscScalar **, PetscInt *, Mat *, PetscHMapI *, PetscInt *, PetscScalar **);
48*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexGetGeomGradData(DM, PetscHMapI *, Mat *, PetscInt *, PetscScalar **, PetscScalar **, PetscInt *, PetscScalar **, PetscScalar **);
49*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexRestoreGeomGradData(DM, PetscHMapI *, Mat *, PetscInt *, PetscScalar **, PetscScalar **, PetscInt *, PetscScalar **, PetscScalar **);
50*5552b385SBrandon PETSC_EXTERN PetscErrorCode DMPlexGetGeomCntrlPntMaps(DM, PetscInt *, PetscInt **, PetscInt **, PetscInt **, PetscInt **, PetscInt **, PetscInt **);
51