1e1589f56SBarry Smith /* 2e1589f56SBarry Smith Objects to manage the interactions between the mesh data structures and the algebraic objects 3e1589f56SBarry Smith */ 43c48a1e8SJed Brown #if !defined(__PETSCDM_H) 53c48a1e8SJed Brown #define __PETSCDM_H 62c8e378dSBarry Smith #include <petscmat.h> 71e25c274SJed Brown #include <petscdmtypes.h> 8e1589f56SBarry Smith 9607a6623SBarry Smith PETSC_EXTERN PetscErrorCode DMInitializePackage(void); 10e1589f56SBarry Smith 11014dd563SJed Brown PETSC_EXTERN PetscClassId DM_CLASSID; 12e1589f56SBarry Smith 1376bdecfbSBarry Smith /*J 14*8f6c3df8SBarry Smith DMType - String with the name of a PETSc DM 15e1589f56SBarry Smith 16e1589f56SBarry Smith Level: beginner 17e1589f56SBarry Smith 18e1589f56SBarry Smith .seealso: DMSetType(), DM 1976bdecfbSBarry Smith J*/ 2019fd82e9SBarry Smith typedef const char* DMType; 21e1589f56SBarry Smith #define DMDA "da" 22e1589f56SBarry Smith #define DMADDA "adda" 23e1589f56SBarry Smith #define DMCOMPOSITE "composite" 24e1589f56SBarry Smith #define DMSLICED "sliced" 25fe1899a2SJed Brown #define DMSHELL "shell" 26b30b9b2eSMatthew G Knepley #define DMMESH "mesh" 27ab7f58a0SBarry Smith #define DMPLEX "plex" 28b30b9b2eSMatthew G Knepley #define DMCARTESIAN "cartesian" 298ac4e037SJed Brown #define DMREDUNDANT "redundant" 303a19ef87SMatthew G Knepley #define DMPATCH "patch" 311d72bce8STim Tautges #define DMMOAB "moab" 32e1589f56SBarry Smith 33140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList DMList; 34014dd563SJed Brown PETSC_EXTERN PetscBool DMRegisterAllCalled; 35014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreate(MPI_Comm,DM*); 3619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetType(DM, DMType); 3719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMGetType(DM, DMType *); 38bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode DMRegister(const char[],PetscErrorCode (*)(DM)); 39607a6623SBarry Smith PETSC_EXTERN PetscErrorCode DMRegisterAll(void); 40014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRegisterDestroy(void); 41e1589f56SBarry Smith 42014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMView(DM,PetscViewer); 43014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLoad(DM,PetscViewer); 44014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDestroy(DM*); 45014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM,Vec*); 46014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM,Vec*); 47014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM,Vec *); 48014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM,Vec *); 49014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM,Vec *); 50014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM,Vec *); 51014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM); 52014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM,const char*,Vec*); 53014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM,const char*,Vec*); 542348bcf4SPeter Brune PETSC_EXTERN PetscErrorCode DMGetNamedLocalVector(DM,const char*,Vec*); 552348bcf4SPeter Brune PETSC_EXTERN PetscErrorCode DMRestoreNamedLocalVector(DM,const char*,Vec*); 56014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM,ISLocalToGlobalMapping*); 57014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMappingBlock(DM,ISLocalToGlobalMapping*); 58014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM,PetscInt*,char***,IS**); 59014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM,PetscInt*); 6019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateColoring(DM,ISColoringType,MatType,ISColoring*); 6119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM,MatType,Mat*); 62014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM,PetscBool); 63014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM,DM,Mat*,Vec*); 64014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateInjection(DM,DM,VecScatter*); 65aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM,PetscInt,PetscDataType,void*); 66aa1993deSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM,PetscInt,PetscDataType,void*); 67014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefine(DM,MPI_Comm,DM*); 68014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsen(DM,MPI_Comm,DM*); 69014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM,PetscInt,DM[]); 70014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM,PetscInt,DM[]); 71014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,Vec,Mat,DM,void*),void*); 72014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*); 73014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestrict(DM,Mat,Vec,Mat,DM); 74014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMInterpolate(DM,Mat,DM); 75014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM); 76146574abSBarry Smith PETSC_EXTERN PetscErrorCode DMViewFromOptions(DM,const char[], const char[]); 77ca266f36SBarry Smith 78014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetUp(DM); 79014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM,DM,Mat,Vec*); 80014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateAggregates(DM,DM,Mat*); 81baf369e7SPeter Brune PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*); 82014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM,Vec,InsertMode,Vec); 83014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM,Vec,InsertMode,Vec); 84014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM,Vec,InsertMode,Vec); 85014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM,Vec,InsertMode,Vec); 8619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMConvert(DM,DMType,DM*); 87e1589f56SBarry Smith 886636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM,DM*); 896636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM,Vec*); 906636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM,Vec); 916636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM,Vec*); 926636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM,Vec); 93e87bb0d3SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMLocatePoints(DM,Vec,IS*); 946636e97aSMatthew G Knepley 955dbd56e3SPeter Brune /* block hook interface */ 96be081cd6SPeter Brune PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*); 97be081cd6SPeter Brune PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM,VecScatter,VecScatter,DM); 985dbd56e3SPeter Brune 99014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM,const char []); 10019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetVecType(DM,VecType); 10119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetMatType(DM,MatType); 102014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM,void*); 103014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM,PetscErrorCode (*)(void**)); 104014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM,void*); 105014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM,PetscErrorCode (*)(DM,Vec,Vec)); 106014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM,PetscBool *); 107b0ae01b7SPeter Brune PETSC_EXTERN PetscErrorCode DMHasColoring(DM,PetscBool *); 108014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM,Vec,Vec); 10993d92d96SBarry Smith 11081d26defSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, PetscInt[], IS *, DM *); 11116621825SDmitry Karpeev PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM,PetscInt*,char***,IS**,DM**); 1128d4ac253SDmitry Karpeev PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM,PetscInt*,char***,IS**,IS**,DM**); 113e30e807fSPeter Brune PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**); 114e7c4fc90SDmitry Karpeev 115014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM,PetscInt*); 116014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM,PetscInt*); 117014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMFinalizePackage(void); 118e1589f56SBarry Smith 1195f1ad066SMatthew G Knepley PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM*); 1205f1ad066SMatthew G Knepley PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM); 121c688c046SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM*); 122c688c046SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM); 1235f1ad066SMatthew G Knepley 124e1589f56SBarry Smith typedef struct NLF_DAAD* NLF; 125e1589f56SBarry Smith 126bc2bf880SBarry Smith #define DM_FILE_CLASSID 1211221 1277da65231SMatthew G Knepley 1287da65231SMatthew G Knepley /* FEM support */ 129014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char [], PetscInt, const PetscScalar []); 130014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char [], PetscInt, PetscInt, const PetscScalar []); 1317da65231SMatthew G Knepley 132014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetDefaultSection(DM, PetscSection *); 133014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetDefaultSection(DM, PetscSection); 134014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetDefaultGlobalSection(DM, PetscSection *); 135eaf8d80aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetDefaultGlobalSection(DM, PetscSection); 136014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetDefaultSF(DM, PetscSF *); 137014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetDefaultSF(DM, PetscSF); 138014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateDefaultSF(DM, PetscSection, PetscSection); 139b21d0597SMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *); 140057b4bcdSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF); 14188ed4aceSMatthew G Knepley 142af122d2aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *); 143af122d2aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt); 144af122d2aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, PetscObject *); 145af122d2aSMatthew G Knepley 1467da65231SMatthew G Knepley typedef struct { 1477da65231SMatthew G Knepley PetscInt numQuadPoints; /* The number of quadrature points on an element */ 1487da65231SMatthew G Knepley const PetscReal *quadPoints; /* The quadrature point coordinates */ 1497da65231SMatthew G Knepley const PetscReal *quadWeights; /* The quadrature weights */ 1507da65231SMatthew G Knepley PetscInt numBasisFuncs; /* The number of finite element basis functions on an element */ 1517da65231SMatthew G Knepley PetscInt numComponents; /* The number of components for each basis function */ 1527da65231SMatthew G Knepley const PetscReal *basis; /* The basis functions tabulated at the quadrature points */ 1537da65231SMatthew G Knepley const PetscReal *basisDer; /* The basis function derivatives tabulated at the quadrature points */ 1547da65231SMatthew G Knepley } PetscQuadrature; 15583e5549fSMatthew G Knepley 15692fc3d87SMatthew G Knepley typedef struct { 15792fc3d87SMatthew G Knepley PetscQuadrature *quad; 15892fc3d87SMatthew G Knepley void (**f0Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */ 15992fc3d87SMatthew G Knepley void (**f1Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */ 16092fc3d87SMatthew G Knepley void (**g0Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_0 functions for each field pair */ 16192fc3d87SMatthew G Knepley void (**g1Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_1 functions for each field pair */ 16292fc3d87SMatthew G Knepley void (**g2Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_2 functions for each field pair */ 16392fc3d87SMatthew G Knepley void (**g3Funcs)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_3 functions for each field pair */ 164bef01ed6SMatthew G. Knepley void (**bcFuncs)(const PetscReal[], PetscScalar *); /* The boundary condition function for each field component */ 165652b88e8SMatthew G. Knepley PetscQuadrature *quadBd; 166652b88e8SMatthew G. Knepley void (**f0BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */ 167652b88e8SMatthew G. Knepley void (**f1BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */ 16892fc3d87SMatthew G Knepley } PetscFEM; 16992fc3d87SMatthew G Knepley 17083e5549fSMatthew G Knepley typedef enum {PETSC_UNIT_LENGTH, PETSC_UNIT_MASS, PETSC_UNIT_TIME, PETSC_UNIT_CURRENT, PETSC_UNIT_TEMPERATURE, PETSC_UNIT_AMOUNT, PETSC_UNIT_LUMINOSITY, NUM_PETSC_UNITS} PetscUnit; 171e87bb0d3SMatthew G Knepley 172e87bb0d3SMatthew G Knepley struct _DMInterpolationInfo { 173e87bb0d3SMatthew G Knepley MPI_Comm comm; 174e87bb0d3SMatthew G Knepley PetscInt dim; /*1 The spatial dimension of points */ 175e87bb0d3SMatthew G Knepley PetscInt nInput; /* The number of input points */ 176e87bb0d3SMatthew G Knepley PetscReal *points; /* The input point coordinates */ 177e87bb0d3SMatthew G Knepley PetscInt *cells; /* The cell containing each point */ 178e87bb0d3SMatthew G Knepley PetscInt n; /* The number of local points */ 179e87bb0d3SMatthew G Knepley Vec coords; /* The point coordinates */ 180e87bb0d3SMatthew G Knepley PetscInt dof; /* The number of components to interpolate */ 181e87bb0d3SMatthew G Knepley }; 182e87bb0d3SMatthew G Knepley typedef struct _DMInterpolationInfo *DMInterpolationInfo; 183e87bb0d3SMatthew G Knepley 184e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *); 185e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt); 186e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *); 187e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt); 188e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *); 189e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]); 190e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool); 191e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *); 192e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *); 193e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *); 194e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec); 195e87bb0d3SMatthew G Knepley PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *); 196e1589f56SBarry Smith #endif 197