Actual source code: petscdmplexegads.h

  1: #pragma once

  3: #include <petscdmplex.h>

  5: #if !defined(PETSC_HAVE_EGADS)
  6:   #error "PETSc not configured for EGADS; reconfigrue --with-egads or --download-egads"
  7: #endif

  9: /* Declarations below provide an interface to the EGADS/EGADSlite libraries, that can be automatically installed
 10:    by PETSc. These functions enable the creation, interrogation, and manipulation of CAD geometries attached to
 11:    a DMPlex. */
 12: #include <egads.h>
 13: #include <egads_lite.h>
 14: typedef ego PetscGeom;

 16: #define PetscCallEGADS(func, args) \
 17:   do { \
 18:     int _status; \
 19:     PetscStackPushExternal(#func); \
 20:     _status = func args; \
 21:     PetscStackPop; \
 22:     PetscCheck(_status >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in EGADS call %s() Status %d", #func, (int)_status); \
 23:   } while (0)

 25: PETSC_EXTERN PetscErrorCode DMPlexGeomDataAndGrads(DM, PetscBool);
 26: PETSC_EXTERN PetscErrorCode DMPlexModifyGeomModel(DM, MPI_Comm, PetscScalar[], PetscScalar[], PetscBool, PetscBool, const char[]);
 27: PETSC_EXTERN PetscErrorCode DMPlexInflateToGeomModelUseXYZ(DM);
 28: PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelTUV(DM);
 29: PETSC_EXTERN PetscErrorCode DMPlexInflateToGeomModelUseTUV(DM);
 30: PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelBodies(DM, PetscGeom **, PetscInt *);
 31: PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelBodyShells(DM, PetscGeom, PetscGeom **, PetscInt *);
 32: PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelBodyFaces(DM, PetscGeom, PetscGeom **, PetscInt *);
 33: PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelBodyLoops(DM, PetscGeom, PetscGeom **, PetscInt *);
 34: PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelBodyEdges(DM, PetscGeom, PetscGeom **, PetscInt *);
 35: PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelBodyNodes(DM, PetscGeom, PetscGeom **, PetscInt *);
 36: PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelShellFaces(DM, PetscGeom, PetscGeom, PetscGeom **, PetscInt *);
 37: PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelFaceLoops(DM, PetscGeom, PetscGeom, PetscGeom **, PetscInt *);
 38: PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelFaceEdges(DM, PetscGeom, PetscGeom, PetscGeom **, PetscInt *);
 39: PETSC_EXTERN PetscErrorCode DMPlexGetGeomModelEdgeNodes(DM, PetscGeom, PetscGeom, PetscGeom **, PetscInt *);
 40: PETSC_EXTERN PetscErrorCode DMPlexGetGeomBodyMassProperties(DM, PetscGeom, PetscScalar *, PetscScalar *, PetscScalar **, PetscInt *, PetscScalar **, PetscInt *);
 41: PETSC_EXTERN PetscErrorCode DMPlexRestoreGeomBodyMassProperties(DM, PetscGeom, PetscScalar *, PetscScalar *, PetscScalar **, PetscInt *, PetscScalar **, PetscInt *);
 42: PETSC_EXTERN PetscErrorCode DMPlexGetGeomID(DM, PetscGeom, PetscGeom, PetscInt *);
 43: PETSC_EXTERN PetscErrorCode DMPlexGetGeomObject(DM, PetscGeom, PetscInt, PetscInt, PetscGeom *);
 44: PETSC_EXTERN PetscErrorCode DMPlexGetGeomFaceNumOfControlPoints(DM, PetscGeom, PetscInt *);
 45: PETSC_EXTERN PetscErrorCode DMPlexFreeGeomObject(DM, PetscGeom *);
 46: PETSC_EXTERN PetscErrorCode DMPlexGetGeomCntrlPntAndWeightData(DM, PetscHMapI *, PetscInt *, PetscScalar **, PetscInt *, Mat *, PetscHMapI *, PetscInt *, PetscScalar **);
 47: PETSC_EXTERN PetscErrorCode DMPlexRestoreGeomCntrlPntAndWeightData(DM, PetscHMapI *, PetscInt *, PetscScalar **, PetscInt *, Mat *, PetscHMapI *, PetscInt *, PetscScalar **);
 48: PETSC_EXTERN PetscErrorCode DMPlexGetGeomGradData(DM, PetscHMapI *, Mat *, PetscInt *, PetscScalar **, PetscScalar **, PetscInt *, PetscScalar **, PetscScalar **);
 49: PETSC_EXTERN PetscErrorCode DMPlexRestoreGeomGradData(DM, PetscHMapI *, Mat *, PetscInt *, PetscScalar **, PetscScalar **, PetscInt *, PetscScalar **, PetscScalar **);
 50: PETSC_EXTERN PetscErrorCode DMPlexGetGeomCntrlPntMaps(DM, PetscInt *, PetscInt **, PetscInt **, PetscInt **, PetscInt **, PetscInt **, PetscInt **);